Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 230

Historique | Voir | Annoter | Télécharger (50 ko)

1
#!/usr/bin/env python
2
#
3
# TrouNoir model using PyOpenCL 
4
#
5
# CC BY-NC-SA 2019 : <emmanuel.quemener@ens-lyon.fr> 
6
#
7
# Thanks to Andreas Klockner for PyOpenCL:
8
# http://mathema.tician.de/software/pyopencl
9
# 
10
# Original code programmed in Fortran 77 in mars 1994
11
# for Practical Work of Numerical Simulation
12
# DEA (old Master2) in astrophysics and spatial techniques in Meudon
13
# by Herve Aussel & Emmanuel Quemener
14
#
15
# Conversion in C done by Emmanuel Quemener in august 1997
16
# GPUfication in OpenCL under Python in july 2019
17
# GPUfication in CUDA under Python in august 2019
18
#
19
# Thanks to :
20
# 
21
# - Herve Aussel for his part of code of black body spectrum
22
# - Didier Pelat for his help to perform this work
23
# - Jean-Pierre Luminet for his article published in 1979
24
# - Numerical Recipies for Runge Kutta recipies
25
# - Luc Blanchet for his disponibility about my questions in General Relativity
26
# - Pierre Lena for his passion about science and vulgarisation
27

    
28
# If crash on OpenCL Intel implementation, add following options and force
29
#export PYOPENCL_COMPILER_OUTPUT=1
30
#export CL_CONFIG_USE_VECTORIZER=True
31
#export CL_CONFIG_CPU_VECTORIZER_MODE=16
32

    
33
import pyopencl as cl
34
import numpy
35
import time,string
36
from numpy.random import randint as nprnd
37
import sys
38
import getopt
39
import matplotlib.pyplot as plt
40
from socket import gethostname
41

    
42
def DictionariesAPI():
43
    PhysicsList={'Einstein':0,'Newton':1}
44
    return(PhysicsList)
45

    
46
#
47
# Blank space below to simplify debugging on OpenCL code
48
#
49

    
50

    
51

    
52

    
53

    
54

    
55

    
56

    
57

    
58

    
59

    
60

    
61

    
62

    
63

    
64

    
65

    
66

    
67

    
68

    
69

    
70

    
71

    
72

    
73

    
74

    
75

    
76

    
77

    
78

    
79

    
80

    
81

    
82

    
83

    
84

    
85

    
86

    
87

    
88

    
89

    
90

    
91

    
92

    
93

    
94

    
95

    
96

    
97

    
98

    
99

    
100

    
101

    
102

    
103
BlobOpenCL= """
104

105
#define PI (float)3.14159265359
106
#define nbr 256
107

108
#define EINSTEIN 0
109
#define NEWTON 1
110

111
#ifdef SETTRACKPOINTS
112
#define TRACKPOINTS SETTRACKPOINTS
113
#else
114
#define TRACKPOINTS 2048
115
#endif
116

117
float atanp(float x,float y)
118
{
119
  float angle;
120

121
  angle=atan2(y,x);
122

123
  if (angle<0.e0f)
124
    {
125
      angle+=(float)2.e0f*PI;
126
    }
127

128
  return angle;
129
}
130

131
float f(float v)
132
{
133
  return v;
134
}
135

136
#if PHYSICS == NEWTON
137
float g(float u,float m,float b)
138
{
139
  return (-u);
140
}
141
#else
142
float g(float u,float m,float b)
143
{
144
  return (3.e0f*m/b*pow(u,2)-u);
145
}
146
#endif
147

148
void calcul(float *us,float *vs,float up,float vp,
149
            float h,float m,float b)
150
{
151
  float c0,c1,c2,c3,d0,d1,d2,d3;
152

153
  c0=h*f(vp);
154
  c1=h*f(vp+c0/2.);
155
  c2=h*f(vp+c1/2.);
156
  c3=h*f(vp+c2);
157
  d0=h*g(up,m,b);
158
  d1=h*g(up+d0/2.,m,b);
159
  d2=h*g(up+d1/2.,m,b);
160
  d3=h*g(up+d2,m,b);
161

162
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
163
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
164
}
165

166
void rungekutta(float *ps,float *us,float *vs,
167
                float pp,float up,float vp,
168
                float h,float m,float b)
169
{
170
  calcul(us,vs,up,vp,h,m,b);
171
  *ps=pp+h;
172
}
173

174
float decalage_spectral(float r,float b,float phi,
175
                         float tho,float m)
176
{
177
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
178
}
179

180
float spectre(float rf,int q,float b,float db,
181
             float h,float r,float m,float bss)
182
{
183
  float flx;
184

185
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
186
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
187
  return(flx);
188
}
189

190
float spectre_cn(float rf32,float b32,float db32,
191
                 float h32,float r32,float m32,float bss32)
192
{
193
  
194
#define MYFLOAT float
195

196
  MYFLOAT rf=(MYFLOAT)(rf32);
197
  MYFLOAT b=(MYFLOAT)(b32);
198
  MYFLOAT db=(MYFLOAT)(db32);
199
  MYFLOAT h=(MYFLOAT)(h32);
200
  MYFLOAT r=(MYFLOAT)(r32);
201
  MYFLOAT m=(MYFLOAT)(m32);
202
  MYFLOAT bss=(MYFLOAT)(bss32);
203

204
  MYFLOAT flx;
205
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
206
  int fi,posfreq;
207

208
#define planck 6.62e-34
209
#define k 1.38e-23
210
#define c2 9.e16
211
#define temp 3.e7
212
#define m_point 1.
213

214
#define lplanck (log(6.62)-34.*log(10.))
215
#define lk (log(1.38)-23.*log(10.))
216
#define lc2 (log(9.)+16.*log(10.))
217

218
  MYFLOAT v=1.-3./r;
219

220
  qu=1./sqrt((1.-3./r)*r)*(sqrt(r)-sqrt(6.)+sqrt(3.)/2.*log((sqrt(r)+sqrt(3.))/(sqrt(r)-sqrt(3.))* 0.17157287525380988 ));
221

222
  temp_em=temp*sqrt(m)*exp(0.25*log(m_point)-0.75*log(r)-0.125*log(v)+0.25*log(fabs(qu)));
223

224
  flux_int=0.;
225
  flx=0.;
226

227
  for (fi=0;fi<nbr;fi++)
228
    {
229
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
230
      nu_rec=nu_em*rf;
231
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
232
      if ((posfreq>0)&&(posfreq<nbr))
233
        {
234
          // Initial version
235
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
236
          // Version with log used
237
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
238
          // flux_int*=pow(rf,3)*b*db*h;
239
          //flux_int*=exp(3.*log(rf))*b*db*h;
240
          flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.)*exp(3.*log(rf))*b*db*h;
241

242
          flx+=flux_int;
243
        }
244
    }
245

246
  return((float)(flx));
247
}
248

249
void impact(float phi,float r,float b,float tho,float m,
250
            float *zp,float *fp,
251
            int q,float db,
252
            float h,int raie)
253
{
254
  float flx,rf,bss;
255

256
  rf=decalage_spectral(r,b,phi,tho,m);
257

258
  if (raie==0)
259
    {
260
      bss=1.e19;
261
      flx=spectre_cn(rf,b,db,h,r,m,bss);
262
    }
263
  else
264
    { 
265
      bss=2.;
266
      flx=spectre(rf,q,b,db,h,r,m,bss);
267
    }
268
  
269
  *zp=1./rf;
270
  *fp=flx;
271

272
}
273

274
__kernel void EachPixel(__global float *zImage,__global float *fImage,
275
                        float Mass,float InternalRadius,
276
                        float ExternalRadius,float Angle,
277
                        int Line)
278
{
279
   uint xi=(uint)get_global_id(0);
280
   uint yi=(uint)get_global_id(1);
281
   uint sizex=(uint)get_global_size(0);
282
   uint sizey=(uint)get_global_size(1);
283

284
   // Perform trajectory for each pixel, exit on hit
285

286
  float m,rs,ri,re,tho;
287
  int q,raie;
288

289
  m=Mass;
290
  rs=2.e0f*m;
291
  ri=InternalRadius;
292
  re=ExternalRadius;
293
  tho=Angle;
294
  q=-2;
295
  raie=Line;
296

297
  float bmx,db,b,h;
298
  float rp0,rpp,rps;
299
  float phi,phd;
300
  uint nh=0;
301
  float zp=0.e0f,fp=0.e0f;
302

303
  // Autosize for image
304
  bmx=1.25e0f*re;
305

306
  h=4.e0f*PI/(float)TRACKPOINTS;
307

308
  // set origin as center of image
309
  float x=(float)xi-(float)(sizex/2)+(float)5.e-1f;
310
  float y=(float)yi-(float)(sizey/2)+(float)5.e-1f;
311
  // angle extracted from cylindric symmetry
312
  phi=atanp(x,y);
313
  phd=atanp(cos(phi)*sin(tho),cos(tho));
314

315

316
  float up,vp,pp,us,vs,ps;
317

318
  // impact parameter
319
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
320
  // step of impact parameter;
321
  db=bmx/(float)(sizex);
322

323
  up=0.e0f;
324
  vp=1.e0f;
325
  pp=0.e0f;
326

327
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
328

329
  rps=fabs(b/us);
330
  rp0=rps;
331

332
  int ExitOnImpact=0;
333

334
  do
335
  {
336
     nh++;
337
     pp=ps;
338
     up=us;
339
     vp=vs;
340
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
341
     rps=fabs(b/us);
342
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>=ri)&&(rps<=re)?1:0;       
343

344
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0)&&(nh<TRACKPOINTS));
345

346

347
  if (ExitOnImpact==1) {
348
     impact(phi,rps,b,tho,m,&zp,&fp,q,db,h,raie);
349
  }
350
  else
351
  {
352
     zp=0.e0f;
353
     fp=0.e0f;
354
  }
355

356
  barrier(CLK_GLOBAL_MEM_FENCE);
357

358
  zImage[yi+sizex*xi]=(float)zp;
359
  fImage[yi+sizex*xi]=(float)fp;  
360
}
361

362
__kernel void Pixel(__global float *zImage,__global float *fImage,
363
                    __global float *Trajectories,__global int *IdLast,
364
                    uint ImpactParameter,
365
                    float Mass,float InternalRadius,
366
                    float ExternalRadius,float Angle,
367
                    int Line)
368
{
369
   uint xi=(uint)get_global_id(0);
370
   uint yi=(uint)get_global_id(1);
371
   uint sizex=(uint)get_global_size(0);
372
   uint sizey=(uint)get_global_size(1);
373

374
   // Perform trajectory for each pixel
375

376
  float m,ri,re,tho;
377
  int q,raie;
378

379
  m=Mass;
380
  ri=InternalRadius;
381
  re=ExternalRadius;
382
  tho=Angle;
383
  q=-2;
384
  raie=Line;
385

386
  float bmx,db,b,h;
387
  float phi,phd,php,nr,r;
388
  float zp=0.,fp=0.;
389

390
  // Autosize for image, 25% greater than external radius
391
  bmx=1.25*re;
392

393
  // Angular step of integration
394
  h=4.e0f*PI/(float)TRACKPOINTS;
395

396
  // Step of Impact Parameter
397
  db=bmx/(2.e0*(float)ImpactParameter);
398

399
  // set origin as center of image
400
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
401
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
402

403
  // angle extracted from cylindric symmetry
404
  phi=atanp(x,y);
405
  phd=atanp(cos(phi)*sin(tho),cos(tho));
406

407
  // Real Impact Parameter
408
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
409

410
  // Integer Impact Parameter
411
  uint bi=(uint)sqrt(x*x+y*y);
412

413
  int HalfLap=0,ExitOnImpact=0,ni;
414

415
  if (bi<ImpactParameter)
416
  {
417
    do
418
    {
419
      php=phd+(float)HalfLap*PI;
420
      nr=php/h;
421
      ni=(int)nr;
422

423
      if (ni<IdLast[bi])
424
      {
425
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
426
      }
427
      else
428
      {
429
        r=Trajectories[bi*TRACKPOINTS+ni];
430
      }
431
           
432
      if ((r<=re)&&(r>=ri))
433
      {
434
        ExitOnImpact=1;
435
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
436
      }
437
              
438
      HalfLap++;
439
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
440

441
  }
442

443
  barrier(CLK_GLOBAL_MEM_FENCE);
444

445
  zImage[yi+sizex*xi]=zp;
446
  fImage[yi+sizex*xi]=fp;
447
}
448

449
__kernel void Circle(__global float *Trajectories,__global int *IdLast,
450
                     __global float *zImage,__global float *fImage,
451
                     float Mass,float InternalRadius,
452
                     float ExternalRadius,float Angle,
453
                     int Line)
454
{
455
   // Integer Impact Parameter ID 
456
   int bi=get_global_id(0);
457
   // Integer points on circle
458
   int i=get_global_id(1);
459
   // Integer Impact Parameter Size (half of image)
460
   int bmaxi=get_global_size(0);
461
   // Integer Points on circle
462
   int imx=get_global_size(1);
463

464
   // Perform trajectory for each pixel
465

466
  float m,ri,re,tho;
467
  int q,raie;
468

469
  m=Mass;
470
  ri=InternalRadius;
471
  re=ExternalRadius;
472
  tho=Angle;
473
  raie=Line;
474

475
  float bmx,db,b,h;
476
  float phi,phd;
477
  float zp=0,fp=0;
478

479
  // Autosize for image
480
  bmx=1.25*re;
481

482
  // Angular step of integration
483
  h=4.e0f*PI/(float)TRACKPOINTS;
484

485
  // impact parameter
486
  b=(float)bi/(float)bmaxi*bmx;
487
  db=bmx/(2.e0f*(float)bmaxi);
488

489
  phi=2.*PI/(float)imx*(float)i;
490
  phd=atanp(cos(phi)*sin(tho),cos(tho));
491
  int yi=(int)((float)bi*sin(phi))+bmaxi;
492
  int xi=(int)((float)bi*cos(phi))+bmaxi;
493

494
  int HalfLap=0,ExitOnImpact=0,ni;
495
  float php,nr,r;
496

497
  do
498
  {
499
     php=phd+(float)HalfLap*PI;
500
     nr=php/h;
501
     ni=(int)nr;
502

503
     if (ni<IdLast[bi])
504
     {
505
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
506
     }
507
     else
508
     {
509
        r=Trajectories[bi*TRACKPOINTS+ni];
510
     }
511
           
512
     if ((r<=re)&&(r>=ri))
513
     {
514
        ExitOnImpact=1;
515
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
516
     }
517
              
518
     HalfLap++;
519
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
520
  
521
  zImage[yi+2*bmaxi*xi]=zp;
522
  fImage[yi+2*bmaxi*xi]=fp;  
523

524
  barrier(CLK_GLOBAL_MEM_FENCE);
525

526
}
527

528
__kernel void Trajectory(__global float *Trajectories,__global int *IdLast,
529
                         float Mass,float InternalRadius,
530
                         float ExternalRadius,float Angle,
531
                         int Line)
532
{
533
  // Integer Impact Parameter ID 
534
  int bi=get_global_id(0);
535
  // Integer Impact Parameter Size (half of image)
536
  int bmaxi=get_global_size(0);
537

538
  // Perform trajectory for each pixel
539

540
  float m,rs,re;
541

542
  m=Mass;
543
  rs=2.*m;
544
  re=ExternalRadius;
545

546
  float bmx,b,h;
547
  int nh;
548

549
  // Autosize for image
550
  bmx=1.25*re;
551

552
  // Angular step of integration
553
  h=4.e0f*PI/(float)TRACKPOINTS;
554

555
  // impact parameter
556
  b=(float)bi/(float)bmaxi*bmx;
557

558
  float up,vp,pp,us,vs,ps;
559

560
  up=0.;
561
  vp=1.;
562
      
563
  pp=0.;
564
  nh=0;
565

566
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
567
  
568
  // b versus us
569
  float bvus=fabs(b/us);
570
  float bvus0=bvus;
571
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
572

573
  do
574
  {
575
     nh++;
576
     pp=ps;
577
     up=us;
578
     vp=vs;
579
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
580
     bvus=fabs(b/us);
581
     Trajectories[bi*TRACKPOINTS+nh]=bvus;
582

583
  } while ((bvus>=rs)&&(bvus<=bvus0));
584

585
  IdLast[bi]=nh;
586

587
  barrier(CLK_GLOBAL_MEM_FENCE);
588
 
589
}
590

591
__kernel void EachCircle(__global float *zImage,__global float *fImage,
592
                         float Mass,float InternalRadius,
593
                         float ExternalRadius,float Angle,
594
                         int Line)
595
{
596
  // Integer Impact Parameter ID 
597
  int bi=get_global_id(0);
598
  // Integer Impact Parameter Size (half of image)
599
  int bmaxi=get_global_size(0);
600

601
  float Trajectory[TRACKPOINTS];
602

603
  float m,rs,ri,re,tho;
604
  int raie,q;
605

606
  m=Mass;
607
  rs=2.*m;
608
  ri=InternalRadius;
609
  re=ExternalRadius;
610
  tho=Angle;
611
  q=-2;
612
  raie=Line;
613

614
  float bmx,db,b,h;
615
  int nh;
616

617
  // Autosize for image
618
  bmx=1.25e0f*re;
619

620
  // Angular step of integration
621
  h=4.e0f*PI/(float)TRACKPOINTS;
622

623
  // impact parameter
624
  b=(float)bi/(float)bmaxi*bmx;
625
  db=bmx/(2.e0f*(float)bmaxi);
626

627
  float up,vp,pp,us,vs,ps;
628

629
  up=0.;
630
  vp=1.;
631
      
632
  pp=0.;
633
  nh=0;
634

635
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
636
  
637
  // b versus us
638
  float bvus=fabs(b/us);
639
  float bvus0=bvus;
640
  Trajectory[nh]=bvus;
641

642
  do
643
  {
644
     nh++;
645
     pp=ps;
646
     up=us;
647
     vp=vs;
648
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
649
     bvus=fabs(b/us);
650
     Trajectory[nh]=bvus;
651

652
  } while ((bvus>=rs)&&(bvus<=bvus0));
653

654
  for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
655
     Trajectory[i]=0.e0f;
656
  }
657

658
  int imx=(int)(16*bi);
659

660
  for (int i=0;i<imx;i++)
661
  {
662
     float zp=0.,fp=0.;
663
     float phi=2.*PI/(float)imx*(float)i;
664
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
665
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
666
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
667

668
     int HalfLap=0,ExitOnImpact=0,ni;
669
     float php,nr,r;
670

671
     do
672
     {
673
        php=phd+(float)HalfLap*PI;
674
        nr=php/h;
675
        ni=(int)nr;
676

677
        if (ni<nh)
678
        {
679
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
680
        }
681
        else
682
        {
683
           r=Trajectory[ni];
684
        }
685
           
686
        if ((r<=re)&&(r>=ri))
687
        {
688
           ExitOnImpact=1;
689
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
690
        }
691
              
692
        HalfLap++;
693

694
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
695
 
696
     zImage[yi+2*bmaxi*xi]=zp;
697
     fImage[yi+2*bmaxi*xi]=fp;
698

699
  }
700

701
  barrier(CLK_GLOBAL_MEM_FENCE);
702
 
703
}
704

705
__kernel void Original(__global float *zImage,__global float *fImage,
706
                       uint Size,float Mass,float InternalRadius,
707
                       float ExternalRadius,float Angle,
708
                       int Line)
709
{
710
   // Integer Impact Parameter Size (half of image)
711
   uint bmaxi=(uint)Size;
712

713
   float Trajectory[TRACKPOINTS];
714

715
   // Perform trajectory for each pixel
716

717
   float m,rs,ri,re,tho;
718
   int raie,q;
719

720
   m=Mass;
721
   rs=2.*m;
722
   ri=InternalRadius;
723
   re=ExternalRadius;
724
   tho=Angle;
725
   q=-2;
726
   raie=Line;
727

728
   float bmx,db,b,h;
729
   int nh;
730

731
   // Autosize for image
732
   bmx=1.25e0f*re;
733

734
   // Angular step of integration
735
   h=4.e0f*PI/(float)TRACKPOINTS;
736

737
   // Integer Impact Parameter ID 
738
   for (int bi=0;bi<bmaxi;bi++)
739
   {
740
      // impact parameter
741
      b=(float)bi/(float)bmaxi*bmx;
742
      db=bmx/(2.e0f*(float)bmaxi);
743

744
      float up,vp,pp,us,vs,ps;
745

746
      up=0.;
747
      vp=1.;
748
      
749
      pp=0.;
750
      nh=0;
751

752
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
753
  
754
      // b versus us
755
      float bvus=fabs(b/us);
756
      float bvus0=bvus;
757
      Trajectory[nh]=bvus;
758

759
      do
760
      {
761
         nh++;
762
         pp=ps;
763
         up=us;
764
         vp=vs;
765
         rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
766
         bvus=fabs(b/us);
767
         Trajectory[nh]=bvus;
768

769
      } while ((bvus>=rs)&&(bvus<=bvus0));
770

771
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
772
         Trajectory[i]=0.e0f;
773
      }
774

775
      int imx=(int)(16*bi);
776

777
      for (int i=0;i<imx;i++)
778
      {
779
         float zp=0,fp=0;
780
         float phi=2.e0f*PI/(float)imx*(float)i;
781
         float phd=atanp(cos(phi)*sin(tho),cos(tho));
782
         uint yi=(uint)((float)bi*sin(phi)+bmaxi);
783
         uint xi=(uint)((float)bi*cos(phi)+bmaxi);
784

785
         int HalfLap=0,ExitOnImpact=0,ni;
786
         float php,nr,r;
787

788
         do
789
         {
790
            php=phd+(float)HalfLap*PI;
791
            nr=php/h;
792
            ni=(int)nr;
793

794
            if (ni<nh)
795
            {
796
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
797
            }
798
            else
799
            {
800
               r=Trajectory[ni];
801
            }
802
           
803
            if ((r<=re)&&(r>=ri))
804
            {
805
               ExitOnImpact=1;
806
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
807
            }
808
               
809
            HalfLap++;
810
 
811
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
812
  
813
         zImage[yi+2*bmaxi*xi]=zp;
814
         fImage[yi+2*bmaxi*xi]=fp;
815
  
816
      }
817

818
   } 
819
   
820
   barrier(CLK_GLOBAL_MEM_FENCE);
821
 
822
}
823
"""
824

    
825
def KernelCodeCuda():
826
    BlobCUDA= """
827

828
#define PI (float)3.14159265359
829
#define nbr 256
830

831
#define EINSTEIN 0
832
#define NEWTON 1
833

834
#ifdef SETTRACKPOINTS
835
#define TRACKPOINTS SETTRACKPOINTS
836
#else
837
#define TRACKPOINTS
838
#endif
839

840
__device__ float nothing(float x)
841
{
842
  return(x);
843
}
844

845
__device__ float atanp(float x,float y)
846
{
847
  float angle;
848

849
  angle=atan2(y,x);
850

851
  if (angle<0.e0f)
852
    {
853
      angle+=(float)2.e0f*PI;
854
    }
855

856
  return(angle);
857
}
858

859
__device__ float f(float v)
860
{
861
  return(v);
862
}
863

864
#if PHYSICS == NEWTON
865
__device__ float g(float u,float m,float b)
866
{
867
  return (-u);
868
}
869
#else
870
__device__ float g(float u,float m,float b)
871
{
872
  return (3.e0f*m/b*pow(u,2)-u);
873
}
874
#endif
875

876
__device__ void calcul(float *us,float *vs,float up,float vp,
877
                       float h,float m,float b)
878
{
879
  float c0,c1,c2,c3,d0,d1,d2,d3;
880

881
  c0=h*f(vp);
882
  c1=h*f(vp+c0/2.);
883
  c2=h*f(vp+c1/2.);
884
  c3=h*f(vp+c2);
885
  d0=h*g(up,m,b);
886
  d1=h*g(up+d0/2.,m,b);
887
  d2=h*g(up+d1/2.,m,b);
888
  d3=h*g(up+d2,m,b);
889

890
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
891
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
892
}
893

894
__device__ void rungekutta(float *ps,float *us,float *vs,
895
                           float pp,float up,float vp,
896
                           float h,float m,float b)
897
{
898
  calcul(us,vs,up,vp,h,m,b);
899
  *ps=pp+h;
900
}
901

902
__device__ float decalage_spectral(float r,float b,float phi,
903
                                   float tho,float m)
904
{
905
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
906
}
907

908
__device__ float spectre(float rf,int q,float b,float db,
909
                         float h,float r,float m,float bss)
910
{
911
  float flx;
912

913
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
914
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
915
  return(flx);
916
}
917

918
__device__ float spectre_cn(float rf32,float b32,float db32,
919
                            float h32,float r32,float m32,float bss32)
920
{
921
  
922
#define MYFLOAT float
923

924
  MYFLOAT rf=(MYFLOAT)(rf32);
925
  MYFLOAT b=(MYFLOAT)(b32);
926
  MYFLOAT db=(MYFLOAT)(db32);
927
  MYFLOAT h=(MYFLOAT)(h32);
928
  MYFLOAT r=(MYFLOAT)(r32);
929
  MYFLOAT m=(MYFLOAT)(m32);
930
  MYFLOAT bss=(MYFLOAT)(bss32);
931

932
  MYFLOAT flx;
933
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
934
  int fi,posfreq;
935

936
#define planck 6.62e-34
937
#define k 1.38e-23
938
#define c2 9.e16
939
#define temp 3.e7
940
#define m_point 1.
941

942
#define lplanck (log(6.62)-34.*log(10.))
943
#define lk (log(1.38)-23.*log(10.))
944
#define lc2 (log(9.)+16.*log(10.))
945

946
  MYFLOAT v=1.-3./r;
947

948
  qu=1./sqrt((1.-3./r)*r)*(sqrt(r)-sqrt(6.)+sqrt(3.)/2.*log((sqrt(r)+sqrt(3.))/(sqrt(r)-sqrt(3.))* 0.17157287525380988 ));
949

950
  temp_em=temp*sqrt(m)*exp(0.25*log(m_point)-0.75*log(r)-0.125*log(v)+0.25*log(fabs(qu)));
951

952
  flux_int=0.;
953
  flx=0.;
954

955
  for (fi=0;fi<nbr;fi++)
956
    {
957
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
958
      nu_rec=nu_em*rf;
959
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
960
      if ((posfreq>0)&&(posfreq<nbr))
961
        {
962
          // Initial version
963
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
964
          // Version with log used
965
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
966
          // flux_int*=pow(rf,3)*b*db*h;
967
          //flux_int*=exp(3.*log(rf))*b*db*h;
968
          flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.)*exp(3.*log(rf))*b*db*h;
969

970
          flx+=flux_int;
971
        }
972
    }
973

974
  return((float)(flx));
975
}
976

977
__device__ void impact(float phi,float r,float b,float tho,float m,
978
                       float *zp,float *fp,
979
                       int q,float db,
980
                       float h,int raie)
981
{
982
  float flx,rf,bss;
983

984
  rf=decalage_spectral(r,b,phi,tho,m);
985

986
  if (raie==0)
987
    {
988
      bss=1.e19;
989
      flx=spectre_cn(rf,b,db,h,r,m,bss);
990
    }
991
  else
992
    { 
993
      bss=2.;
994
      flx=spectre(rf,q,b,db,h,r,m,bss);
995
    }
996
  
997
  *zp=1./rf;
998
  *fp=flx;
999

1000
}
1001

1002
__global__ void EachPixel(float *zImage,float *fImage,
1003
                          float Mass,float InternalRadius,
1004
                          float ExternalRadius,float Angle,
1005
                          int Line)
1006
{
1007
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
1008
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
1009
   uint sizex=(uint)gridDim.x*blockDim.x;
1010
   uint sizey=(uint)gridDim.y*blockDim.y;
1011

1012
   // Perform trajectory for each pixel, exit on hit
1013

1014
  float m,rs,ri,re,tho;
1015
  int q,raie;
1016

1017
  m=Mass;
1018
  rs=2.*m;
1019
  ri=InternalRadius;
1020
  re=ExternalRadius;
1021
  tho=Angle;
1022
  q=-2;
1023
  raie=Line;
1024

1025
  float d,bmx,db,b,h;
1026
  float rp0,rpp,rps;
1027
  float phi,thi,phd,php,nr,r;
1028
  int nh;
1029
  float zp,fp;
1030

1031
  // Autosize for image
1032
  bmx=1.25*re;
1033
  b=0.;
1034

1035
  h=4.e0f*PI/(float)TRACKPOINTS;
1036

1037
  // set origin as center of image
1038
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
1039
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
1040
  // angle extracted from cylindric symmetry
1041
  phi=atanp(x,y);
1042
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1043

1044
  float up,vp,pp,us,vs,ps;
1045

1046
  // impact parameter
1047
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
1048
  // step of impact parameter;
1049
//  db=bmx/(float)(sizex/2);
1050
  db=bmx/(float)(sizex);
1051

1052
  up=0.;
1053
  vp=1.;
1054
      
1055
  pp=0.;
1056
  nh=0;
1057

1058
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1059
    
1060
  rps=fabs(b/us);
1061
  rp0=rps;
1062

1063
  int ExitOnImpact=0;
1064

1065
  do
1066
  {
1067
     nh++;
1068
     pp=ps;
1069
     up=us;
1070
     vp=vs;
1071
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1072
     rpp=rps;          
1073
     rps=fabs(b/us);
1074
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;          
1075

1076
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
1077

1078
  if (ExitOnImpact==1) {
1079
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
1080
  }
1081
  else
1082
  {
1083
     zp=0.;
1084
     fp=0.;
1085
  }
1086

1087
  __syncthreads();
1088

1089
  zImage[yi+sizex*xi]=(float)zp;
1090
  fImage[yi+sizex*xi]=(float)fp;
1091
}
1092

1093
__global__ void Pixel(float *zImage,float *fImage,
1094
                      float *Trajectories,int *IdLast,
1095
                      uint ImpactParameter,
1096
                      float Mass,float InternalRadius,
1097
                      float ExternalRadius,float Angle,
1098
                      int Line)
1099
{
1100
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
1101
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
1102
   uint sizex=(uint)gridDim.x*blockDim.x;
1103
   uint sizey=(uint)gridDim.y*blockDim.y;
1104

1105
  // Perform trajectory for each pixel
1106

1107
  float m,rs,ri,re,tho;
1108
  int q,raie;
1109

1110
  m=Mass;
1111
  rs=2.*m;
1112
  ri=InternalRadius;
1113
  re=ExternalRadius;
1114
  tho=Angle;
1115
  q=-2;
1116
  raie=Line;
1117

1118
  float d,bmx,db,b,h;
1119
  float phi,thi,phd,php,nr,r;
1120
  int nh;
1121
  float zp=0,fp=0;
1122
  // Autosize for image, 25% greater than external radius
1123
  bmx=1.25*re;
1124

1125
  // Angular step of integration
1126
  h=4.e0f*PI/(float)TRACKPOINTS;
1127

1128
  // Step of Impact Parameter
1129
  db=bmx/(2.e0*(float)ImpactParameter);
1130

1131
  // set origin as center of image
1132
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
1133
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
1134
  // angle extracted from cylindric symmetry
1135
  phi=atanp(x,y);
1136
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1137

1138
  // Real Impact Parameter
1139
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
1140

1141
  // Integer Impact Parameter
1142
  uint bi=(uint)sqrt(x*x+y*y);
1143

1144
  int HalfLap=0,ExitOnImpact=0,ni;
1145

1146
  if (bi<ImpactParameter)
1147
  {
1148
    do
1149
    {
1150
      php=phd+(float)HalfLap*PI;
1151
      nr=php/h;
1152
      ni=(int)nr;
1153

1154
      if (ni<IdLast[bi])
1155
      {
1156
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1157
      }
1158
      else
1159
      {
1160
        r=Trajectories[bi*TRACKPOINTS+ni];
1161
      }
1162
           
1163
      if ((r<=re)&&(r>=ri))
1164
      {
1165
        ExitOnImpact=1;
1166
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1167
      }
1168
              
1169
      HalfLap++;
1170
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1171

1172
  }
1173

1174
  zImage[yi+sizex*xi]=zp;
1175
  fImage[yi+sizex*xi]=fp;
1176
}
1177

1178
__global__ void Circle(float *Trajectories,int *IdLast,
1179
                       float *zImage,float *fImage,
1180
                       float Mass,float InternalRadius,
1181
                       float ExternalRadius,float Angle,
1182
                       int Line)
1183
{
1184
   // Integer Impact Parameter ID 
1185
   int bi=blockIdx.x*blockDim.x+threadIdx.x;
1186
   // Integer points on circle
1187
   int i=blockIdx.y*blockDim.y+threadIdx.y;
1188
   // Integer Impact Parameter Size (half of image)
1189
   int bmaxi=gridDim.x*blockDim.x;
1190
   // Integer Points on circle
1191
   int imx=gridDim.y*blockDim.y;
1192

1193
   // Perform trajectory for each pixel
1194

1195
  float m,rs,ri,re,tho;
1196
  int q,raie;
1197

1198
  m=Mass;
1199
  rs=2.*m;
1200
  ri=InternalRadius;
1201
  re=ExternalRadius;
1202
  tho=Angle;
1203
  raie=Line;
1204

1205
  float bmx,db,b,h;
1206
  float phi,thi,phd;
1207
  int nh;
1208
  float zp=0,fp=0;
1209

1210
  // Autosize for image
1211
  bmx=1.25*re;
1212

1213
  // Angular step of integration
1214
  h=4.e0f*PI/(float)TRACKPOINTS;
1215

1216
  // impact parameter
1217
  b=(float)bi/(float)bmaxi*bmx;
1218
  db=bmx/(2.e0*(float)bmaxi);
1219

1220
  phi=2.*PI/(float)imx*(float)i;
1221
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1222
  int yi=(int)((float)bi*sin(phi))+bmaxi;
1223
  int xi=(int)((float)bi*cos(phi))+bmaxi;
1224

1225
  int HalfLap=0,ExitOnImpact=0,ni;
1226
  float php,nr,r;
1227

1228
  do
1229
  {
1230
     php=phd+(float)HalfLap*PI;
1231
     nr=php/h;
1232
     ni=(int)nr;
1233

1234
     if (ni<IdLast[bi])
1235
     {
1236
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1237
     }
1238
     else
1239
     {
1240
        r=Trajectories[bi*TRACKPOINTS+ni];
1241
     }
1242
           
1243
     if ((r<=re)&&(r>=ri))
1244
     {
1245
        ExitOnImpact=1;
1246
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1247
     }
1248
              
1249
     HalfLap++;
1250
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1251
  
1252
  zImage[yi+2*bmaxi*xi]=zp;
1253
  fImage[yi+2*bmaxi*xi]=fp;  
1254

1255
}
1256

1257
__global__ void Trajectory(float *Trajectories,int *IdLast,
1258
                           float Mass,float InternalRadius,
1259
                           float ExternalRadius,float Angle,
1260
                           int Line)
1261
{
1262
  // Integer Impact Parameter ID 
1263
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1264
  // Integer Impact Parameter Size (half of image)
1265
  int bmaxi=gridDim.x*blockDim.x;
1266

1267
  // Perform trajectory for each pixel
1268

1269
  float m,rs,ri,re,tho;
1270
  int raie,q;
1271

1272
  m=Mass;
1273
  rs=2.*m;
1274
  ri=InternalRadius;
1275
  re=ExternalRadius;
1276
  tho=Angle;
1277
  q=-2;
1278
  raie=Line;
1279

1280
  float d,bmx,db,b,h;
1281
  float phi,thi,phd,php,nr,r;
1282
  int nh;
1283
  float zp,fp;
1284

1285
  // Autosize for image
1286
  bmx=1.25*re;
1287

1288
  // Angular step of integration
1289
  h=4.e0f*PI/(float)TRACKPOINTS;
1290

1291
  // impact parameter
1292
  b=(float)bi/(float)bmaxi*bmx;
1293

1294
  float up,vp,pp,us,vs,ps;
1295

1296
  up=0.;
1297
  vp=1.;
1298
      
1299
  pp=0.;
1300
  nh=0;
1301

1302
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1303
  
1304
  // b versus us
1305
  float bvus=fabs(b/us);
1306
  float bvus0=bvus;
1307
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
1308

1309
  do
1310
  {
1311
     nh++;
1312
     pp=ps;
1313
     up=us;
1314
     vp=vs;
1315
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1316
     bvus=fabs(b/us);
1317
     Trajectories[bi*TRACKPOINTS+nh]=bvus;
1318

1319
  } while ((bvus>=rs)&&(bvus<=bvus0));
1320

1321
  IdLast[bi]=nh;
1322

1323
}
1324

1325
__global__ void EachCircle(float *zImage,float *fImage,
1326
                           float Mass,float InternalRadius,
1327
                           float ExternalRadius,float Angle,
1328
                           int Line)
1329
{
1330
  // Integer Impact Parameter ID 
1331
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1332

1333
  // Integer Impact Parameter Size (half of image)
1334
  int bmaxi=gridDim.x*blockDim.x;
1335

1336
  float Trajectory[2048];
1337

1338
  // Perform trajectory for each pixel
1339

1340
  float m,rs,ri,re,tho;
1341
  int raie,q;
1342

1343
  m=Mass;
1344
  rs=2.*m;
1345
  ri=InternalRadius;
1346
  re=ExternalRadius;
1347
  tho=Angle;
1348
  q=-2;
1349
  raie=Line;
1350

1351
  float bmx,db,b,h;
1352
  int nh;
1353

1354
  // Autosize for image
1355
  bmx=1.25e0f*re;
1356

1357
  // Angular step of integration
1358
  h=4.e0f*PI/(float)TRACKPOINTS;
1359

1360
  // impact parameter
1361
  b=(float)bi/(float)bmaxi*bmx;
1362
  db=bmx/(2.e0f*(float)bmaxi);
1363

1364
  float up,vp,pp,us,vs,ps;
1365

1366
  up=0.;
1367
  vp=1.;
1368
      
1369
  pp=0.;
1370
  nh=0;
1371

1372
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1373
  
1374
  // b versus us
1375
  float bvus=fabs(b/us);
1376
  float bvus0=bvus;
1377
  Trajectory[nh]=bvus;
1378

1379
  do
1380
  {
1381
     nh++;
1382
     pp=ps;
1383
     up=us;
1384
     vp=vs;
1385
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1386
     bvus=fabs(b/us);
1387
     Trajectory[nh]=bvus;
1388

1389
  } while ((bvus>=rs)&&(bvus<=bvus0));
1390

1391
  int imx=(int)(16*bi);
1392

1393
  for (int i=0;i<imx;i++)
1394
  {
1395
     float zp=0,fp=0;
1396
     float phi=2.*PI/(float)imx*(float)i;
1397
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
1398
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1399
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1400

1401
     int HalfLap=0,ExitOnImpact=0,ni;
1402
     float php,nr,r;
1403

1404
     do
1405
     {
1406
        php=phd+(float)HalfLap*PI;
1407
        nr=php/h;
1408
        ni=(int)nr;
1409

1410
        if (ni<nh)
1411
        {
1412
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1413
        }
1414
        else
1415
        {
1416
           r=Trajectory[ni];
1417
        }
1418
           
1419
        if ((r<=re)&&(r>=ri))
1420
        {
1421
           ExitOnImpact=1;
1422
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1423
        }
1424
              
1425
        HalfLap++;
1426

1427
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
1428

1429
   __syncthreads();
1430

1431
   zImage[yi+2*bmaxi*xi]=zp;
1432
   fImage[yi+2*bmaxi*xi]=fp;
1433

1434
  }
1435
 
1436
}
1437

1438
__global__ void Original(float *zImage,float *fImage,
1439
                         uint Size,float Mass,float InternalRadius,
1440
                         float ExternalRadius,float Angle,
1441
                         int Line)
1442
{
1443
   // Integer Impact Parameter Size (half of image)
1444
   uint bmaxi=(uint)Size;
1445

1446
   float Trajectory[TRACKPOINTS];
1447

1448
   // Perform trajectory for each pixel
1449

1450
   float m,rs,ri,re,tho;
1451
   int raie,q;
1452

1453
   m=Mass;
1454
   rs=2.*m;
1455
   ri=InternalRadius;
1456
   re=ExternalRadius;
1457
   tho=Angle;
1458
   q=-2;
1459
   raie=Line;
1460

1461
   float bmx,db,b,h;
1462
   int nh;
1463

1464
   // Autosize for image
1465
   bmx=1.25e0f*re;
1466

1467
   // Angular step of integration
1468
   h=4.e0f*PI/(float)TRACKPOINTS;
1469

1470
   // Integer Impact Parameter ID 
1471
   for (int bi=0;bi<bmaxi;bi++)
1472
   {
1473
      // impact parameter
1474
      b=(float)bi/(float)bmaxi*bmx;
1475
      db=bmx/(2.e0f*(float)bmaxi);
1476

1477
      float up,vp,pp,us,vs,ps;
1478

1479
      up=0.;
1480
      vp=1.;
1481
      
1482
      pp=0.;
1483
      nh=0;
1484

1485
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1486
  
1487
      // b versus us
1488
      float bvus=fabs(b/us);
1489
      float bvus0=bvus;
1490
      Trajectory[nh]=bvus;
1491

1492
      do
1493
      {
1494
         nh++;
1495
         pp=ps;
1496
         up=us;
1497
         vp=vs;
1498
         rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1499
         bvus=fabs(b/us);
1500
         Trajectory[nh]=bvus;
1501

1502
      } while ((bvus>=rs)&&(bvus<=bvus0));
1503

1504
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
1505
         Trajectory[i]=0.e0f;
1506
      }
1507

1508
      int imx=(int)(16*bi);
1509

1510
      for (int i=0;i<imx;i++)
1511
      {
1512
         float zp=0,fp=0;
1513
         float phi=2.e0f*PI/(float)imx*(float)i;
1514
         float phd=atanp(cos(phi)*sin(tho),cos(tho));
1515
         uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1516
         uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1517

1518
         int HalfLap=0,ExitOnImpact=0,ni;
1519
         float php,nr,r;
1520

1521
         do
1522
         {
1523
            php=phd+(float)HalfLap*PI;
1524
            nr=php/h;
1525
            ni=(int)nr;
1526

1527
            if (ni<nh)
1528
            {
1529
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1530
            }
1531
            else
1532
            {
1533
               r=Trajectory[ni];
1534
            }
1535
           
1536
            if ((r<=re)&&(r>=ri))
1537
            {
1538
               ExitOnImpact=1;
1539
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1540
            }
1541
               
1542
            HalfLap++;
1543
 
1544
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
1545
  
1546
         zImage[yi+2*bmaxi*xi]=zp;
1547
         fImage[yi+2*bmaxi*xi]=fp;
1548
  
1549
      }
1550

1551
   } 
1552
   
1553
}
1554
"""
1555
    return(BlobCUDA)
1556
    
1557
# def ImageOutput(sigma,prefix):
1558
#     from PIL import Image
1559
#     Max=sigma.max()
1560
#     Min=sigma.min()
1561
    
1562
#     # Normalize value as 8bits Integer
1563
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1564
#     image = Image.fromarray(SigmaInt)
1565
#     image.save("%s.jpg" % prefix)
1566

    
1567
def ImageOutput(sigma,prefix,Colors):
1568
    import matplotlib.pyplot as plt
1569
    start_time=time.time() 
1570
    if Colors == 'Red2Yellow':
1571
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1572
    else:
1573
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1574
    save_time = time.time()-start_time
1575
    print("Save Time : %f" % save_time)
1576

    
1577
def BlackHoleCL(zImage,fImage,InputCL):
1578

    
1579
    kernel_params = {}
1580

    
1581
    print(InputCL)
1582

    
1583
    Device=InputCL['Device']
1584
    GpuStyle=InputCL['GpuStyle']
1585
    VariableType=InputCL['VariableType']
1586
    Size=InputCL['Size']
1587
    Mass=InputCL['Mass']
1588
    InternalRadius=InputCL['InternalRadius']
1589
    ExternalRadius=InputCL['ExternalRadius']
1590
    Angle=InputCL['Angle']
1591
    Method=InputCL['Method']
1592
    TrackPoints=InputCL['TrackPoints']
1593
    Physics=InputCL['Physics']
1594

    
1595
    PhysicsList=DictionariesAPI()
1596
    
1597
    if InputCL['BlackBody']:
1598
        # Spectrum is Black Body one
1599
        Line=0
1600
    else:
1601
        # Spectrum is Monochromatic Line one
1602
        Line=1
1603

    
1604
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1605
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1606

    
1607
    # Je detecte un peripherique GPU dans la liste des peripheriques
1608
    Id=0
1609
    HasXPU=False
1610
    for platform in cl.get_platforms():
1611
        for device in platform.get_devices():
1612
            if Id==Device:
1613
                XPU=device
1614
                print("CPU/GPU selected: ",device.name.lstrip())
1615
                HasXPU=True
1616
            Id+=1
1617

    
1618
    if HasXPU==False:
1619
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1620
        sys.exit()                
1621
        
1622
    ctx = cl.Context([XPU])
1623
    queue = cl.CommandQueue(ctx,
1624
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1625

    
1626
    
1627
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1628

    
1629
    BuildOptions="-DPHYSICS=%i -DSETTRACKPOINTS=%i " % (PhysicsList[Physics],InputCL['TrackPoints'])
1630

    
1631
    print('My Platform is ',platform.name)
1632
    
1633
    if 'Intel' in platform.name or 'Experimental' in platform.name or 'Clover' in platform.name or 'Portable' in platform.name :
1634
        print('No extra options for Intel and Clover!')
1635
    else:
1636
        BuildOptions = BuildOptions+" -cl-mad-enable"
1637

    
1638
    print(BuildOptions)
1639
    
1640
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1641
    
1642
    # Je recupere les flag possibles pour les buffers
1643
    mf = cl.mem_flags
1644

    
1645
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1646
        TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories) 
1647
        IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1648

    
1649
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1650
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1651
    
1652
    start_time=time.time() 
1653

    
1654
    if Method=='EachPixel':
1655
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1656
                                       None,zImageCL,fImageCL,
1657
                                       numpy.float32(Mass),
1658
                                       numpy.float32(InternalRadius),
1659
                                       numpy.float32(ExternalRadius),
1660
                                       numpy.float32(Angle),
1661
                                       numpy.int32(Line))
1662
        CLLaunch.wait()
1663
    elif Method=='Original':
1664
        CLLaunch=BlackHoleCL.Original(queue,(1,),
1665
                                      None,zImageCL,fImageCL,
1666
                                      numpy.uint32(zImage.shape[0]/2),
1667
                                      numpy.float32(Mass),
1668
                                      numpy.float32(InternalRadius),
1669
                                      numpy.float32(ExternalRadius),
1670
                                      numpy.float32(Angle),
1671
                                      numpy.int32(Line))
1672
        CLLaunch.wait()
1673
    elif Method=='EachCircle':
1674
        CLLaunch=BlackHoleCL.EachCircle(queue,(zImage.shape[0]/2,),
1675
                                        None,zImageCL,fImageCL,
1676
                                        numpy.float32(Mass),
1677
                                        numpy.float32(InternalRadius),
1678
                                        numpy.float32(ExternalRadius),
1679
                                        numpy.float32(Angle),
1680
                                        numpy.int32(Line))
1681
        CLLaunch.wait()
1682
    elif Method=='TrajectoCircle':
1683
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1684
                                        None,TrajectoriesCL,IdLastCL,
1685
                                        numpy.float32(Mass),
1686
                                        numpy.float32(InternalRadius),
1687
                                        numpy.float32(ExternalRadius),
1688
                                        numpy.float32(Angle),
1689
                                        numpy.int32(Line))
1690
        
1691
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1692
                                           zImage.shape[0]*4),None,
1693
                                    TrajectoriesCL,IdLastCL,
1694
                                    zImageCL,fImageCL,
1695
                                    numpy.float32(Mass),
1696
                                    numpy.float32(InternalRadius),
1697
                                    numpy.float32(ExternalRadius),
1698
                                    numpy.float32(Angle),
1699
                                    numpy.int32(Line))
1700
        CLLaunch.wait()
1701
    else:
1702
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1703
                                        None,TrajectoriesCL,IdLastCL,
1704
                                        numpy.float32(Mass),
1705
                                        numpy.float32(InternalRadius),
1706
                                        numpy.float32(ExternalRadius),
1707
                                        numpy.float32(Angle),
1708
                                        numpy.int32(Line))
1709
        
1710
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],zImage.shape[1]),None,
1711
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1712
                                   numpy.uint32(Trajectories.shape[0]),
1713
                                   numpy.float32(Mass),
1714
                                   numpy.float32(InternalRadius),
1715
                                   numpy.float32(ExternalRadius),
1716
                                   numpy.float32(Angle),
1717
                                   numpy.int32(Line))
1718
        CLLaunch.wait()
1719
    
1720
    compute = time.time()-start_time
1721
    
1722
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1723
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1724
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1725
        cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1726
        cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1727
    elapsed = time.time()-start_time
1728
    print("\nCompute Time : %f" % compute)
1729
    print("Elapsed Time : %f\n" % elapsed)
1730

    
1731
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1732
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1733
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1734
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1735
    zImageCL.release()
1736
    fImageCL.release()
1737

    
1738
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1739
        TrajectoriesCL.release()
1740
        IdLastCL.release()
1741

    
1742
    return(elapsed)
1743

    
1744
def BlackHoleCUDA(zImage,fImage,InputCL):
1745

    
1746
    kernel_params = {}
1747

    
1748
    print(InputCL)
1749

    
1750
    Device=InputCL['Device']
1751
    GpuStyle=InputCL['GpuStyle']
1752
    VariableType=InputCL['VariableType']
1753
    Size=InputCL['Size']
1754
    Mass=InputCL['Mass']
1755
    InternalRadius=InputCL['InternalRadius']
1756
    ExternalRadius=InputCL['ExternalRadius']
1757
    Angle=InputCL['Angle']
1758
    Method=InputCL['Method']
1759
    TrackPoints=InputCL['TrackPoints']
1760
    Physics=InputCL['Physics']
1761

    
1762
    PhysicsList=DictionariesAPI()
1763
    
1764
    if InputCL['BlackBody']:
1765
        # Spectrum is Black Body one
1766
        Line=0
1767
    else:
1768
        # Spectrum is Monochromatic Line one
1769
        Line=1
1770

    
1771
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1772
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1773

    
1774
    try:
1775
        # For PyCUDA import
1776
        import pycuda.driver as cuda
1777
        from pycuda.compiler import SourceModule
1778
        
1779
        cuda.init()
1780
        for Id in range(cuda.Device.count()):
1781
            if Id==Device:
1782
                XPU=cuda.Device(Id)
1783
                print("GPU selected %s" % XPU.name())
1784
        print
1785

    
1786
    except ImportError:
1787
        print("Platform does not seem to support CUDA")
1788

    
1789
    Context=XPU.make_context()
1790

    
1791
    try:
1792
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1793
        print("Compilation seems to be OK")
1794
    except:
1795
        print("Compilation seems to break")
1796

    
1797
    EachPixelCU=mod.get_function("EachPixel")
1798
    OriginalCU=mod.get_function("Original")
1799
    EachCircleCU=mod.get_function("EachCircle")
1800
    TrajectoryCU=mod.get_function("Trajectory")
1801
    PixelCU=mod.get_function("Pixel")
1802
    CircleCU=mod.get_function("Circle")
1803
    
1804
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1805
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1806
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1807
    cuda.memcpy_htod(zImageCU, zImage)
1808
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1809
    cuda.memcpy_htod(zImageCU, fImage)
1810
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1811
    cuda.memcpy_htod(IdLastCU, IdLast)
1812

    
1813
    start_time=time.time() 
1814

    
1815
    if Method=='EachPixel':
1816
        EachPixelCU(zImageCU,fImageCU,
1817
                    numpy.float32(Mass),
1818
                    numpy.float32(InternalRadius),
1819
                    numpy.float32(ExternalRadius),
1820
                    numpy.float32(Angle),
1821
                    numpy.int32(Line),
1822
                    grid=(zImage.shape[0]/32,zImage.shape[1]/32),
1823
                    block=(32,32,1))
1824
    elif Method=='EachCircle':
1825
        EachCircleCU(zImageCU,fImageCU,
1826
                   numpy.float32(Mass),
1827
                   numpy.float32(InternalRadius),
1828
                   numpy.float32(ExternalRadius),
1829
                   numpy.float32(Angle),
1830
                   numpy.int32(Line),
1831
                   grid=(zImage.shape[0]/32/2,1),
1832
                   block=(32,1,1))
1833
    elif Method=='Original':
1834
        OriginalCU(zImageCU,fImageCU,
1835
                   numpy.uint32(zImage.shape[0]/2),
1836
                   numpy.float32(Mass),
1837
                   numpy.float32(InternalRadius),
1838
                   numpy.float32(ExternalRadius),
1839
                   numpy.float32(Angle),
1840
                   numpy.int32(Line),
1841
                   grid=(1,1),
1842
                   block=(1,1,1))
1843
    elif Method=='TrajectoCircle':
1844
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1845
                     numpy.float32(Mass),
1846
                     numpy.float32(InternalRadius),
1847
                     numpy.float32(ExternalRadius),
1848
                     numpy.float32(Angle),
1849
                     numpy.int32(Line),
1850
                     grid=(Trajectories.shape[0]/32,1),
1851
                     block=(32,1,1))
1852
        
1853
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1854
                 numpy.float32(Mass),
1855
                 numpy.float32(InternalRadius),
1856
                 numpy.float32(ExternalRadius),
1857
                 numpy.float32(Angle),
1858
                 numpy.int32(Line),
1859
                 grid=(Trajectories.shape[0]/32,zImage.shape[0]*4/32),
1860
                 block=(32,32,1))
1861
    else:
1862
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1863
                     numpy.float32(Mass),
1864
                     numpy.float32(InternalRadius),
1865
                     numpy.float32(ExternalRadius),
1866
                     numpy.float32(Angle),
1867
                     numpy.int32(Line),
1868
                     grid=(Trajectories.shape[0]/32,1),
1869
                     block=(32,1,1))
1870
        
1871
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1872
                numpy.uint32(Trajectories.shape[0]),
1873
                numpy.float32(Mass),
1874
                numpy.float32(InternalRadius),
1875
                numpy.float32(ExternalRadius),
1876
                numpy.float32(Angle),
1877
                numpy.int32(Line),
1878
                grid=(zImage.shape[0]/32,zImage.shape[1]/32,1),
1879
                block=(32,32,1))
1880

    
1881
    Context.synchronize()
1882

    
1883
    compute = time.time()-start_time
1884

    
1885
    cuda.memcpy_dtoh(zImage,zImageCU)
1886
    cuda.memcpy_dtoh(fImage,fImageCU)
1887
    elapsed = time.time()-start_time
1888
    print("\nCompute Time : %f" % compute)
1889
    print("Elapsed Time : %f\n" % elapsed)
1890

    
1891
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1892
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1893
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1894
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1895

    
1896
    
1897
    Context.pop()
1898

    
1899
    Context.detach()
1900

    
1901
    return(elapsed)
1902

    
1903
if __name__=='__main__':
1904
        
1905
    GpuStyle = 'OpenCL'
1906
    Mass = 1.
1907
    # Internal Radius 3 times de Schwarzschild Radius
1908
    InternalRadius=6.*Mass
1909
    #
1910
    ExternalRadius=12.
1911
    #
1912
    # Angle with normal to disc 10 degrees
1913
    Angle = numpy.pi/180.*(90.-10.)
1914
    # Radiation of disc : BlackBody or Monochromatic
1915
    BlackBody = False
1916
    # Size of image
1917
    Size=256
1918
    # Variable Type
1919
    VariableType='FP32'
1920
    # ?
1921
    q=-2
1922
    # Method of resolution
1923
    Method='TrajectoPixel'
1924
    # Colors for output image
1925
    Colors='Greyscale'
1926
    # Physics
1927
    Physics='Einstein'
1928
    # No output as image
1929
    NoImage = False
1930
    
1931
    HowToUse='%s -h [Help] -b [BlackBodyEmission] -n [NoImage] -p <Einstein/Newton> -s <SizeInPixels> -m <Mass> -i <DiscInternalRadius> -x <DiscExternalRadius> -a <AngleAboveDisc> -d <DeviceId> -c <Greyscale/Red2Yellow> -g <CUDA/OpenCL> -t <EachPixel/TrajectoCircle/TrajectoPixel/EachCircle/Original> -v <FP32/FP64>'
1932

    
1933
    try:
1934
        opts, args = getopt.getopt(sys.argv[1:],"hbns:m:i:x:a:d:g:v:t:c:p:",["blackbody","noimage","camera","size=","mass=","internal=","external=","angle=","device=","gpustyle=","variabletype=","method=","colors=","physics="])
1935
    except getopt.GetoptError:
1936
        print(HowToUse % sys.argv[0])
1937
        sys.exit(2)
1938

    
1939
    # List of Devices
1940
    Devices=[]
1941
    Alu={}
1942
        
1943
    for opt, arg in opts:
1944
        if opt == '-h':
1945
            print(HowToUse % sys.argv[0])
1946
            
1947
            print("\nInformations about devices detected under OpenCL API:")
1948
            # For PyOpenCL import
1949
            try:
1950
                import pyopencl as cl
1951
                Id=0
1952
                for platform in cl.get_platforms():
1953
                    for device in platform.get_devices():
1954
                        #deviceType=cl.device_type.to_string(device.type)
1955
                        deviceType="xPU"
1956
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
1957
                        Id=Id+1
1958
                        
1959
            except:
1960
                print("Your platform does not seem to support OpenCL")
1961
                
1962
            print("\nInformations about devices detected under CUDA API:")
1963
                # For PyCUDA import
1964
            try:
1965
                import pycuda.driver as cuda
1966
                cuda.init()
1967
                for Id in range(cuda.Device.count()):
1968
                    device=cuda.Device(Id)
1969
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
1970
                print
1971
            except:
1972
                print("Your platform does not seem to support CUDA")
1973
        
1974
            sys.exit()
1975
        
1976
        elif opt in ("-d", "--device"):
1977
#            Devices.append(int(arg))
1978
            Device=int(arg)
1979
        elif opt in ("-g", "--gpustyle"):
1980
            GpuStyle = arg
1981
        elif opt in ("-v", "--variabletype"):
1982
            VariableType = arg
1983
        elif opt in ("-s", "--size"):
1984
            Size = int(arg)
1985
        elif opt in ("-m", "--mass"):
1986
            Mass = float(arg)
1987
        elif opt in ("-i", "--internal"):
1988
            InternalRadius = float(arg)
1989
        elif opt in ("-e", "--external"):
1990
            ExternalRadius = float(arg)
1991
        elif opt in ("-a", "--angle"):
1992
            Angle = numpy.pi/180.*(90.-float(arg))
1993
        elif opt in ("-b", "--blackbody"):
1994
            BlackBody = True
1995
        elif opt in ("-n", "--noimage"):
1996
            NoImage = True
1997
        elif opt in ("-t", "--method"):
1998
            Method = arg
1999
        elif opt in ("-c", "--colors"):
2000
            Colors = arg
2001
        elif opt in ("-p", "--physics"):
2002
            Physics = arg
2003

    
2004
    print("Device Identification selected : %s" % Device)
2005
    print("GpuStyle used : %s" % GpuStyle)
2006
    print("VariableType : %s" % VariableType)
2007
    print("Size : %i" % Size)
2008
    print("Mass : %f" % Mass)
2009
    print("Internal Radius : %f" % InternalRadius)
2010
    print("External Radius : %f" % ExternalRadius)
2011
    print("Angle with normal of (in radians) : %f" % Angle)
2012
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
2013
    print("Method of resolution : %s" % Method)
2014
    print("Colors for output images : %s" % Colors)
2015
    print("Physics used for Trajectories : %s" % Physics)
2016

    
2017
    if GpuStyle=='CUDA':
2018
        print("\nSelection of CUDA device") 
2019
        try:
2020
            # For PyCUDA import
2021
            import pycuda.driver as cuda
2022
            
2023
            cuda.init()
2024
            for Id in range(cuda.Device.count()):
2025
                device=cuda.Device(Id)
2026
                print("Device #%i of type GPU : %s" % (Id,device.name()))
2027
                if Id in Devices:
2028
                    Alu[Id]='GPU'
2029
            
2030
        except ImportError:
2031
            print("Platform does not seem to support CUDA")
2032

    
2033
    if GpuStyle=='OpenCL':
2034
        print("\nSelection of OpenCL device") 
2035
        try:
2036
            # For PyOpenCL import
2037
            import pyopencl as cl
2038
            Id=0
2039
            for platform in cl.get_platforms():
2040
                for device in platform.get_devices():
2041
                    #deviceType=cl.device_type.to_string(device.type)
2042
                    deviceType="xPU"
2043
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
2044

    
2045
                    if Id in Devices:
2046
                    # Set the Alu as detected Device Type
2047
                        Alu[Id]=deviceType
2048
                    Id=Id+1
2049
        except ImportError:
2050
            print("Platform does not seem to support OpenCL")
2051

    
2052
#    print(Devices,Alu)
2053

    
2054
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
2055
    TrackPoints=2048
2056
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2057
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2058

    
2059
    InputCL={}
2060
    InputCL['Device']=Device
2061
    InputCL['GpuStyle']=GpuStyle
2062
    InputCL['VariableType']=VariableType
2063
    InputCL['Size']=Size
2064
    InputCL['Mass']=Mass
2065
    InputCL['InternalRadius']=InternalRadius
2066
    InputCL['ExternalRadius']=ExternalRadius
2067
    InputCL['Angle']=Angle
2068
    InputCL['BlackBody']=BlackBody
2069
    InputCL['Method']=Method
2070
    InputCL['TrackPoints']=TrackPoints
2071
    InputCL['Physics']=Physics
2072

    
2073
    if GpuStyle=='OpenCL':
2074
        duration=BlackHoleCL(zImage,fImage,InputCL)
2075
    else:
2076
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
2077
        
2078
    Hostname=gethostname()
2079
    Date=time.strftime("%Y%m%d_%H%M%S")
2080
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2081
    
2082
    
2083
    if not NoImage:
2084
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
2085
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)