Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 210

Historique | Voir | Annoter | Télécharger (21,97 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
#
18
# Thanks to :
19
# 
20
# - Herve Aussel for his part of code of black body spectrum
21
# - Didier Pelat for his help to perform this work
22
# - Jean-Pierre Luminet for his article published in 1979
23
# - Numerical Recipies for Runge Kutta recipies
24
# - Luc Blanchet for his disponibility about my questions in General Relativity
25
# - Pierre Lena for his passion about science and vulgarisation
26

    
27
import pyopencl as cl
28
import numpy
29
from PIL import Image
30
import time,string
31
from numpy.random import randint as nprnd
32
import sys
33
import getopt
34
import matplotlib.pyplot as plt
35

    
36

    
37

    
38

    
39

    
40

    
41

    
42

    
43
KERNEL_CODE=string.Template("""
44

45
#define PI (float)3.14159265359
46
#define nbr 256
47

48
float atanp(float x,float y)
49
{
50
  float angle;
51

52
  angle=atan2(y,x);
53

54
  if (angle<0.e0f)
55
    {
56
      angle+=(float)2.e0f*PI;
57
    }
58

59
  return angle;
60
}
61

62
float f(float v)
63
{
64
  return v;
65
}
66

67
float g(float u,float m,float b)
68
{
69
  return (3.e0f*m/b*pow(u,2)-u);
70
}
71

72

73
void calcul(float *us,float *vs,float up,float vp,
74
            float h,float m,float b)
75
{
76
  float c0,c1,c2,c3,d0,d1,d2,d3;
77

78
  c0=h*f(vp);
79
  c1=h*f(vp+c0/2.);
80
  c2=h*f(vp+c1/2.);
81
  c3=h*f(vp+c2);
82
  d0=h*g(up,m,b);
83
  d1=h*g(up+d0/2.,m,b);
84
  d2=h*g(up+d1/2.,m,b);
85
  d3=h*g(up+d2,m,b);
86

87
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
88
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
89
}
90

91
void rungekutta(float *ps,float *us,float *vs,
92
                float pp,float up,float vp,
93
                float h,float m,float b)
94
{
95
  calcul(us,vs,up,vp,h,m,b);
96
  *ps=pp+h;
97
}
98

99
float decalage_spectral(float r,float b,float phi,
100
                         float tho,float m)
101
{
102
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
103
}
104

105
float spectre(float rf,int q,float b,float db,
106
             float h,float r,float m,float bss)
107
{
108
  float flx;
109

110
  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
111
  return(flx);
112
}
113

114
float spectre_cn(float rf32,float b32,float db32,
115
                 float h32,float r32,float m32,float bss32)
116
{
117
  
118
#define MYFLOAT double
119

120
  MYFLOAT rf=(MYFLOAT)(rf32);
121
  MYFLOAT b=(MYFLOAT)(b32);
122
  MYFLOAT db=(MYFLOAT)(db32);
123
  MYFLOAT h=(MYFLOAT)(h32);
124
  MYFLOAT r=(MYFLOAT)(r32);
125
  MYFLOAT m=(MYFLOAT)(m32);
126
  MYFLOAT bss=(MYFLOAT)(bss32);
127

128
  MYFLOAT flx;
129
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
130
  int fi,posfreq;
131

132
#define planck 6.62e-34
133
#define k 1.38e-23
134
#define c2 9.e16
135
#define temp 3.e7
136
#define m_point 1.
137

138
#define lplanck (log(6.62)-34.*log(10.))
139
#define lk (log(1.38)-23.*log(10.))
140
#define lc2 (log(9.)+16.*log(10.))
141

142
  MYFLOAT v=1.-3./r;
143

144
  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 ));
145

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

148
  flux_int=0.;
149
  flx=0.;
150

151
  for (fi=0;fi<nbr;fi++)
152
    {
153
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
154
      nu_rec=nu_em*rf;
155
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
156
      if ((posfreq>0)&&(posfreq<nbr))
157
        {
158
          // Initial version
159
          flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
160
          // Version with log used
161
          flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
162
          flux_int*=pow(rf,3)*b*db*h;
163
          flx+=flux_int;
164
        }
165
    }
166

167
  return((float)(flx));
168
}
169

170
void impact(float phi,float r,float b,float tho,float m,
171
            float *zp,float *fp,
172
            int q,float db,
173
            float h,int raie)
174
{
175
  float flx,rf,bss;
176

177
  rf=decalage_spectral(r,b,phi,tho,m);
178

179
  if (raie==0)
180
    {
181
      bss=1.e19;
182
      flx=spectre_cn(rf,b,db,h,r,m,bss);
183
    }
184
  else
185
    { 
186
      bss=2.;
187
      flx=spectre(rf,q,b,db,h,r,m,bss);
188
    }
189
  
190
  *zp=1./rf;
191
  *fp=flx;
192

193
}
194

195
__kernel void EachPixel(__global float *zImage,__global float *fImage,
196
                        float Mass,float InternalRadius,
197
                        float ExternalRadius,float Angle,
198
                        int Line)
199
{
200
   uint xi=(uint)get_global_id(0);
201
   uint yi=(uint)get_global_id(1);
202
   uint sizex=(uint)get_global_size(0);
203
   uint sizey=(uint)get_global_size(1);
204

205
   // Perform trajectory for each pixel, exit on hit
206

207
  float m,rs,ri,re,tho;
208
  int q,raie;
209

210
  m=Mass;
211
  rs=2.*m;
212
  ri=InternalRadius;
213
  re=ExternalRadius;
214
  tho=Angle;
215
  q=-2;
216
  raie=Line;
217

218
  float d,bmx,db,b,h;
219
  float rp[2048];
220
  float phi,thi,phd,php,nr,r;
221
  int nh;
222
  float zp,fp;
223

224
  // Autosize for image
225
  bmx=1.25*re;
226
  b=0.;
227

228
  h=4.e0f*PI/(float)2048;
229

230
  // set origin as center of image
231
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
232
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
233
  // angle extracted from cylindric symmetry
234
  phi=atanp(x,y);
235
  phd=atanp(cos(phi)*sin(tho),cos(tho));
236

237
  float up,vp,pp,us,vs,ps;
238

239
  // impact parameter
240
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
241
  // step of impact parameter;
242
//  db=bmx/(float)(sizex/2);
243
  db=bmx/(float)(sizex);
244

245
  up=0.;
246
  vp=1.;
247
      
248
  pp=0.;
249
  nh=0;
250

251
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
252
    
253
  rp[nh]=fabs(b/us);
254

255
  int ExitOnImpact=0;
256

257
  do
258
  {
259
     nh++;
260
     pp=ps;
261
     up=us;
262
     vp=vs;
263
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);          
264
     rp[nh]=fabs(b/us);
265
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rp[nh]>ri)&&(rp[nh]<re)?1:0;          
266

267
  } while ((rp[nh]>=rs)&&(rp[nh]<=rp[0])&&(ExitOnImpact==0));
268

269
  if (ExitOnImpact==1) {
270
     impact(phi,rp[nh-1],b,tho,m,&zp,&fp,q,db,h,raie);
271
  }
272
  else
273
  {
274
     zp=0.;
275
     fp=0.;
276
  }
277

278
  barrier(CLK_GLOBAL_MEM_FENCE);
279

280
  zImage[yi+sizex*xi]=(float)zp;
281
  fImage[yi+sizex*xi]=(float)fp;  
282
}
283

284
__kernel void Pixel(__global float *zImage,__global float *fImage,
285
                    __global float *Trajectories,__global int *IdLast,
286
                    uint ImpactParameter,uint TrackPoints,
287
                    float Mass,float InternalRadius,
288
                    float ExternalRadius,float Angle,
289
                    int Line)
290
{
291
   uint xi=(uint)get_global_id(0);
292
   uint yi=(uint)get_global_id(1);
293
   uint sizex=(uint)get_global_size(0);
294
   uint sizey=(uint)get_global_size(1);
295

296
   // Perform trajectory for each pixel
297

298
  float m,rs,ri,re,tho;
299
  int q,raie;
300

301
  m=Mass;
302
  rs=2.*m;
303
  ri=InternalRadius;
304
  re=ExternalRadius;
305
  tho=Angle;
306
  q=-2;
307
  raie=Line;
308

309
  float d,bmx,db,b,h;
310
  float phi,thi,phd,php,nr,r;
311
  int nh;
312
  float zp=0,fp=0;
313

314
  // Autosize for image, 25% greater than external radius
315
  bmx=1.25*re;
316

317
  // Angular step of integration
318
  h=4.e0f*PI/(float)TrackPoints;
319

320
  // Step of Impact Parameter
321
  db=bmx/(2.e0*(float)ImpactParameter);
322

323
  // set origin as center of image
324
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
325
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
326

327
  // angle extracted from cylindric symmetry
328
  phi=atanp(x,y);
329
  phd=atanp(cos(phi)*sin(tho),cos(tho));
330

331
  // Real Impact Parameter
332
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
333

334
  // Integer Impact Parameter
335
  uint bi=(uint)sqrt(x*x+y*y);
336

337
  int HalfLap=0,ExitOnImpact=0,ni;
338

339
  if (bi<ImpactParameter)
340
  {
341
    do
342
    {
343
      php=phd+(float)HalfLap*PI;
344
      nr=php/h;
345
      ni=(int)nr;
346

347
      if (ni<IdLast[bi])
348
      {
349
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
350
      }
351
      else
352
      {
353
        r=Trajectories[bi*TrackPoints+ni];
354
      }
355
           
356
      if ((r<=re)&&(r>=ri))
357
      {
358
        ExitOnImpact=1;
359
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
360
      }
361
              
362
      HalfLap++;
363
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
364

365
  }
366

367
  barrier(CLK_GLOBAL_MEM_FENCE);
368

369
  zImage[yi+sizex*xi]=zp;
370
  fImage[yi+sizex*xi]=fp;
371
}
372

373
__kernel void Circle(__global float *Trajectories,__global int *IdLast,
374
                     __global float *zImage,__global float *fImage,
375
                     int TrackPoints,
376
                     float Mass,float InternalRadius,
377
                     float ExternalRadius,float Angle,
378
                     int Line)
379
{
380
   // Integer Impact Parameter ID 
381
   int bi=get_global_id(0);
382
   // Integer points on circle
383
   int i=get_global_id(1);
384
   // Integer Impact Parameter Size (half of image)
385
   int bmaxi=get_global_size(0);
386
   // Integer Points on circle
387
   int imx=get_global_size(1);
388

389
   // Perform trajectory for each pixel
390

391
  float m,rs,ri,re,tho;
392
  int q,raie;
393

394
  m=Mass;
395
  rs=2.*m;
396
  ri=InternalRadius;
397
  re=ExternalRadius;
398
  tho=Angle;
399
  raie=Line;
400

401
  float bmx,db,b,h;
402
  float phi,thi,phd;
403
  int nh;
404
  float zp=0,fp=0;
405

406
  // Autosize for image
407
  bmx=1.25*re;
408

409
  // Angular step of integration
410
  h=4.e0f*PI/(float)TrackPoints;
411

412
  // impact parameter
413
  b=(float)bi/(float)bmaxi*bmx;
414
  db=bmx/(2.e0*(float)bmaxi);
415

416
  phi=2.*PI/(float)imx*(float)i;
417
  phd=atanp(cos(phi)*sin(tho),cos(tho));
418
  int yi=(int)((float)bi*sin(phi))+bmaxi;
419
  int xi=(int)((float)bi*cos(phi))+bmaxi;
420

421
  int HalfLap=0,ExitOnImpact=0,ni;
422
  float php,nr,r;
423

424
  do
425
  {
426
     php=phd+(float)HalfLap*PI;
427
     nr=php/h;
428
     ni=(int)nr;
429

430
     if (ni<IdLast[bi])
431
     {
432
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
433
     }
434
     else
435
     {
436
        r=Trajectories[bi*TrackPoints+ni];
437
     }
438
           
439
     if ((r<=re)&&(r>=ri))
440
     {
441
        ExitOnImpact=1;
442
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
443
     }
444
              
445
     HalfLap++;
446
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
447
  
448
  zImage[yi+2*bmaxi*xi]=zp;
449
  fImage[yi+2*bmaxi*xi]=fp;  
450

451
  barrier(CLK_GLOBAL_MEM_FENCE);
452

453
}
454

455
__kernel void Trajectory(__global float *Trajectories,
456
                         __global int *IdLast,int TrackPoints,
457
                         float Mass,float InternalRadius,
458
                         float ExternalRadius,float Angle,
459
                         int Line)
460
{
461
  // Integer Impact Parameter ID 
462
  int bi=get_global_id(0);
463
  // Integer Impact Parameter Size (half of image)
464
  int bmaxi=get_global_size(0);
465

466
  // Perform trajectory for each pixel
467

468
  float m,rs,ri,re,tho;
469
  int raie,q;
470

471
  m=Mass;
472
  rs=2.*m;
473
  ri=InternalRadius;
474
  re=ExternalRadius;
475
  tho=Angle;
476
  q=-2;
477
  raie=Line;
478

479
  float d,bmx,db,b,h;
480
  float phi,thi,phd,php,nr,r;
481
  int nh;
482
  float zp,fp;
483

484
  // Autosize for image
485
  bmx=1.25*re;
486

487
  // Angular step of integration
488
  h=4.e0f*PI/(float)TrackPoints;
489

490
  // impact parameter
491
  b=(float)bi/(float)bmaxi*bmx;
492

493
  float up,vp,pp,us,vs,ps;
494

495
  up=0.;
496
  vp=1.;
497
      
498
  pp=0.;
499
  nh=0;
500

501
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
502
  
503
  // b versus us
504
  float bvus=fabs(b/us);
505
  float bvus0=bvus;
506
  Trajectories[bi*TrackPoints+nh]=bvus;
507

508
  do
509
  {
510
     nh++;
511
     pp=ps;
512
     up=us;
513
     vp=vs;
514
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
515
     bvus=fabs(b/us);
516
     Trajectories[bi*TrackPoints+nh]=bvus;
517

518
  } while ((bvus>=rs)&&(bvus<=bvus0));
519

520
  IdLast[bi]=nh;
521

522
  barrier(CLK_GLOBAL_MEM_FENCE);
523
 
524
}
525
""")
526

    
527
def ImageOutput(sigma,prefix):
528
    Max=sigma.max()
529
    Min=sigma.min()
530
    
531
    # Normalize value as 8bits Integer
532
    SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
533
    image = Image.fromarray(SigmaInt)
534
    image.save("%s.jpg" % prefix)
535
    
536
def BlackHoleCL(zImage,fImage,InputCL):
537

    
538
    kernel_params = {}
539

    
540
    print(InputCL)
541

    
542
    Device=InputCL['Device']
543
    GpuStyle=InputCL['GpuStyle']
544
    VariableType=InputCL['VariableType']
545
    Size=InputCL['Size']
546
    Mass=InputCL['Mass']
547
    InternalRadius=InputCL['InternalRadius']
548
    ExternalRadius=InputCL['ExternalRadius']
549
    Angle=InputCL['Angle']
550
    Method=InputCL['Method']
551
    TrackPoints=InputCL['TrackPoints']
552

    
553
    if InputCL['BlackBody']:
554
        # Spectrum is Black Body one
555
        Line=0
556
    else:
557
        # Spectrum is Monochromatic Line one
558
        Line=1
559

    
560
    Trajectories=numpy.zeros((InputCL['Size']/2,InputCL['TrackPoints']),
561
                              dtype=numpy.float32)
562
    IdLast=numpy.zeros(InputCL['Size']/2,dtype=numpy.int32)
563

    
564
    # Je detecte un peripherique GPU dans la liste des peripheriques
565
    Id=0
566
    HasXPU=False
567
    for platform in cl.get_platforms():
568
        for device in platform.get_devices():
569
            if Id==Device:
570
                XPU=device
571
                print "CPU/GPU selected: ",device.name.lstrip()
572
                HasXPU=True
573
            Id+=1
574

    
575
    if HasXPU==False:
576
        print "No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1)
577
        sys.exit()                
578
        
579
    ctx = cl.Context([XPU])
580
    queue = cl.CommandQueue(ctx,
581
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
582

    
583
    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
584
    
585
    # Je recupere les flag possibles pour les buffers
586
    mf = cl.mem_flags
587

    
588
    print(zImage)
589
    
590
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
591
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
592
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
593
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
594
    
595
    start_time=time.time() 
596

    
597
    print(Trajectories.shape)
598
    if Method=='EachPixel':
599
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
600
                                       None,zImageCL,fImageCL,
601
                                       numpy.float32(Mass),
602
                                       numpy.float32(InternalRadius),
603
                                       numpy.float32(ExternalRadius),
604
                                       numpy.float32(Angle),
605
                                       numpy.int32(Line))
606
        CLLaunch.wait()
607
    elif Method=='TrajectoCircle':
608
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
609
                                        None,TrajectoriesCL,IdLastCL,
610
                                        numpy.uint32(Trajectories.shape[1]),
611
                                        numpy.float32(Mass),
612
                                        numpy.float32(InternalRadius),
613
                                        numpy.float32(ExternalRadius),
614
                                        numpy.float32(Angle),
615
                                        numpy.int32(Line))
616
        
617
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
618
                                           zImage.shape[0]*4),None,
619
                                    TrajectoriesCL,IdLastCL,
620
                                    zImageCL,fImageCL,
621
                                    numpy.uint32(Trajectories.shape[1]),
622
                                    numpy.float32(Mass),
623
                                    numpy.float32(InternalRadius),
624
                                    numpy.float32(ExternalRadius),
625
                                    numpy.float32(Angle),
626
                                    numpy.int32(Line))
627
        CLLaunch.wait()
628
    else:
629
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
630
                                        None,TrajectoriesCL,IdLastCL,
631
                                        numpy.uint32(Trajectories.shape[1]),
632
                                        numpy.float32(Mass),
633
                                        numpy.float32(InternalRadius),
634
                                        numpy.float32(ExternalRadius),
635
                                        numpy.float32(Angle),
636
                                        numpy.int32(Line))
637
        
638
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],
639
                                                  zImage.shape[1]),None,
640
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
641
                                   numpy.uint32(Trajectories.shape[0]),
642
                                   numpy.uint32(Trajectories.shape[1]),
643
                                        numpy.float32(Mass),
644
                                        numpy.float32(InternalRadius),
645
                                        numpy.float32(ExternalRadius),
646
                                        numpy.float32(Angle),
647
                                        numpy.int32(Line))
648
        CLLaunch.wait()
649
    
650
    elapsed = time.time()-start_time
651
    print("Elapsed %f: " % elapsed)
652
    
653
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
654
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
655
    cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
656
    cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
657
    
658
    print(zImage.max(),numpy.where(zImage[:,:]==zImage.max()))
659
    print(fImage.max(),numpy.where(fImage[:,:]==fImage.max()))
660
    zImageCL.release()
661
    fImageCL.release()
662

    
663
    TrajectoriesCL.release()
664
    IdLastCL.release()
665

    
666
    return(elapsed)
667

    
668
if __name__=='__main__':
669
        
670
    GpuStyle = 'OpenCL'
671
    Mass = 1.
672
    # Internal Radius 3 times de Schwarzschild Radius
673
    InternalRadius=6.*Mass
674
    #
675
    ExternalRadius=12.
676
    #
677
    # Angle with normal to disc 10 degrees
678
    Angle = numpy.pi/180.*(90.-10.)
679
    # Radiation of disc : BlackBody or Monochromatic
680
    BlackBody = False
681
    # Size of image
682
    Size=256
683
    # Variable Type
684
    VariableType='FP32'
685
    # ?
686
    q=-2
687
    # Method of resolution
688
    Method='TrajectoPixel'
689
          
690
    HowToUse='%s -h [Help] -b [BlackBodyEmission] -s <SizeInPixels> -m <Mass> -i <DiscInternalRadius> -x <DiscExternalRadius> -a <AngleAboveDisc> -d <DeviceId> -g <CUDA/OpenCL> -t <EachPixel/TrajectoCircle/TrajectoPixel> -v <FP32/FP64>'
691

    
692
    try:
693
        opts, args = getopt.getopt(sys.argv[1:],"hbs:m:i:x:a:d:g:v:t:",["blackbody","camera","size=","mass=","internal=","external=","angle=","device=","gpustyle=","variabletype=","method="])
694
    except getopt.GetoptError:
695
        print(HowToUse % sys.argv[0])
696
        sys.exit(2)
697

    
698
    # List of Devices
699
    Devices=[]
700
    Alu={}
701
        
702
    for opt, arg in opts:
703
        if opt == '-h':
704
            print(HowToUse % sys.argv[0])
705
            
706
            print("\nInformations about devices detected under OpenCL API:")
707
            # For PyOpenCL import
708
            try:
709
                import pyopencl as cl
710
                Id=0
711
                for platform in cl.get_platforms():
712
                    for device in platform.get_devices():
713
                        #deviceType=cl.device_type.to_string(device.type)
714
                        deviceType="xPU"
715
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
716
                        Id=Id+1
717
                        
718
            except:
719
                print("Your platform does not seem to support OpenCL")
720
                
721
            print("\nInformations about devices detected under CUDA API:")
722
                # For PyCUDA import
723
            try:
724
                import pycuda.driver as cuda
725
                cuda.init()
726
                for Id in range(cuda.Device.count()):
727
                    device=cuda.Device(Id)
728
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
729
                print
730
            except:
731
                print("Your platform does not seem to support CUDA")
732
        
733
            sys.exit()
734
        
735
        elif opt in ("-d", "--device"):
736
#            Devices.append(int(arg))
737
            Device=int(arg)
738
        elif opt in ("-g", "--gpustyle"):
739
            GpuStyle = arg
740
        elif opt in ("-v", "--variabletype"):
741
            VariableType = arg
742
        elif opt in ("-s", "--size"):
743
            Size = int(arg)
744
        elif opt in ("-m", "--mass"):
745
            Mass = float(arg)
746
        elif opt in ("-i", "--internal"):
747
            InternalRadius = float(arg)
748
        elif opt in ("-e", "--external"):
749
            ExternalRadius = float(arg)
750
        elif opt in ("-a", "--angle"):
751
            Angle = numpy.pi/180.*(90.-float(arg))
752
        elif opt in ("-b", "--blackbody"):
753
            BlackBody = True
754
        elif opt in ("-t", "--method"):
755
            Method = arg
756

    
757
    print("Device Identification selected : %s" % Device)
758
    print("GpuStyle used : %s" % GpuStyle)
759
    print("VariableType : %s" % VariableType)
760
    print("Size : %i" % Size)
761
    print("Mass : %f" % Mass)
762
    print("Internal Radius : %f" % InternalRadius)
763
    print("External Radius : %f" % ExternalRadius)
764
    print("Angle with normal of (in radians) : %f" % Angle)
765
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
766
    print("Method of resolution : %s" % Method)
767

    
768
    if GpuStyle=='CUDA':
769
        print("\nSelection of CUDA device") 
770
        try:
771
            # For PyCUDA import
772
            import pycuda.driver as cuda
773
            
774
            cuda.init()
775
            for Id in range(cuda.Device.count()):
776
                device=cuda.Device(Id)
777
                print("Device #%i of type GPU : %s" % (Id,device.name()))
778
                if Id in Devices:
779
                    Alu[Id]='GPU'
780
            
781
        except ImportError:
782
            print("Platform does not seem to support CUDA")
783

    
784
    if GpuStyle=='OpenCL':
785
        print("\nSelection of OpenCL device") 
786
        try:
787
            # For PyOpenCL import
788
            import pyopencl as cl
789
            Id=0
790
            for platform in cl.get_platforms():
791
                for device in platform.get_devices():
792
                    #deviceType=cl.device_type.to_string(device.type)
793
                    deviceType="xPU"
794
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
795

    
796
                    if Id in Devices:
797
                    # Set the Alu as detected Device Type
798
                        Alu[Id]=deviceType
799
                    Id=Id+1
800
        except ImportError:
801
            print("Platform does not seem to support OpenCL")
802

    
803
#    print(Devices,Alu)
804

    
805
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
806
    TrackPoints=2048
807
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
808
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
809

    
810
    InputCL={}
811
    InputCL['Device']=Device
812
    InputCL['GpuStyle']=GpuStyle
813
    InputCL['VariableType']=VariableType
814
    InputCL['Size']=Size
815
    InputCL['Mass']=Mass
816
    InputCL['InternalRadius']=InternalRadius
817
    InputCL['ExternalRadius']=ExternalRadius
818
    InputCL['Angle']=Angle
819
    InputCL['BlackBody']=BlackBody
820
    InputCL['Method']=Method
821
    InputCL['TrackPoints']=TrackPoints
822
    
823
    duration=BlackHoleCL(zImage,fImage,InputCL)
824

    
825
    ImageOutput(zImage,"TrouNoirZ_%s" % Method)
826
    ImageOutput(fImage,"TrouNoirF_%s" % Method)