Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 219

Historique | Voir | Annoter | Télécharger (39,34 ko)

1
#!/usr/bin/env python3
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
import time,string
30
from numpy.random import randint as nprnd
31
import sys
32
import getopt
33
import matplotlib.pyplot as plt
34
from socket import gethostname
35

    
36
def DictionariesAPI():
37
    PhysicsList={'Einstein':0,'Newton':1}
38
    return(PhysicsList)
39

    
40
BlobOpenCL= """
41

42
#define PI (float)3.14159265359
43
#define nbr 256
44

45
#define EINSTEIN 0
46
#define NEWTON 1
47

48
#define TRACKPOINTS 2048
49

50
float atanp(float x,float y)
51
{
52
  float angle;
53

54
  angle=atan2(y,x);
55

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

61
  return angle;
62
}
63

64
float f(float v)
65
{
66
  return v;
67
}
68

69
#if PHYSICS == NEWTON
70
float g(float u,float m,float b)
71
{
72
  return (-u);
73
}
74
#else
75
float g(float u,float m,float b)
76
{
77
  return (3.e0f*m/b*pow(u,2)-u);
78
}
79
#endif
80

81
void calcul(float *us,float *vs,float up,float vp,
82
            float h,float m,float b)
83
{
84
  float c0,c1,c2,c3,d0,d1,d2,d3;
85

86
  c0=h*f(vp);
87
  c1=h*f(vp+c0/2.);
88
  c2=h*f(vp+c1/2.);
89
  c3=h*f(vp+c2);
90
  d0=h*g(up,m,b);
91
  d1=h*g(up+d0/2.,m,b);
92
  d2=h*g(up+d1/2.,m,b);
93
  d3=h*g(up+d2,m,b);
94

95
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
96
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
97
}
98

99
void rungekutta(float *ps,float *us,float *vs,
100
                float pp,float up,float vp,
101
                float h,float m,float b)
102
{
103
  calcul(us,vs,up,vp,h,m,b);
104
  *ps=pp+h;
105
}
106

107
float decalage_spectral(float r,float b,float phi,
108
                         float tho,float m)
109
{
110
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
111
}
112

113
float spectre(float rf,int q,float b,float db,
114
             float h,float r,float m,float bss)
115
{
116
  float flx;
117

118
  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
119
  return(flx);
120
}
121

122
float spectre_cn(float rf32,float b32,float db32,
123
                 float h32,float r32,float m32,float bss32)
124
{
125
  
126
#define MYFLOAT float
127

128
  MYFLOAT rf=(MYFLOAT)(rf32);
129
  MYFLOAT b=(MYFLOAT)(b32);
130
  MYFLOAT db=(MYFLOAT)(db32);
131
  MYFLOAT h=(MYFLOAT)(h32);
132
  MYFLOAT r=(MYFLOAT)(r32);
133
  MYFLOAT m=(MYFLOAT)(m32);
134
  MYFLOAT bss=(MYFLOAT)(bss32);
135

136
  MYFLOAT flx;
137
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
138
  int fi,posfreq;
139

140
#define planck 6.62e-34
141
#define k 1.38e-23
142
#define c2 9.e16
143
#define temp 3.e7
144
#define m_point 1.
145

146
#define lplanck (log(6.62)-34.*log(10.))
147
#define lk (log(1.38)-23.*log(10.))
148
#define lc2 (log(9.)+16.*log(10.))
149

150
  MYFLOAT v=1.-3./r;
151

152
  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 ));
153

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

156
  flux_int=0.;
157
  flx=0.;
158

159
  for (fi=0;fi<nbr;fi++)
160
    {
161
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
162
      nu_rec=nu_em*rf;
163
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
164
      if ((posfreq>0)&&(posfreq<nbr))
165
        {
166
          // Initial version
167
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
168
          // Version with log used
169
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
170
          // flux_int*=pow(rf,3)*b*db*h;
171
          //flux_int*=exp(3.*log(rf))*b*db*h;
172
          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;
173

174
          flx+=flux_int;
175
        }
176
    }
177

178
  return((float)(flx));
179
}
180

181
void impact(float phi,float r,float b,float tho,float m,
182
            float *zp,float *fp,
183
            int q,float db,
184
            float h,int raie)
185
{
186
  float flx,rf,bss;
187

188
  rf=decalage_spectral(r,b,phi,tho,m);
189

190
  if (raie==0)
191
    {
192
      bss=1.e19;
193
      flx=spectre_cn(rf,b,db,h,r,m,bss);
194
    }
195
  else
196
    { 
197
      bss=2.;
198
      flx=spectre(rf,q,b,db,h,r,m,bss);
199
    }
200
  
201
  *zp=1./rf;
202
  *fp=flx;
203

204
}
205

206
__kernel void EachPixel(__global float *zImage,__global float *fImage,
207
                        float Mass,float InternalRadius,
208
                        float ExternalRadius,float Angle,
209
                        int Line)
210
{
211
   uint xi=(uint)get_global_id(0);
212
   uint yi=(uint)get_global_id(1);
213
   uint sizex=(uint)get_global_size(0);
214
   uint sizey=(uint)get_global_size(1);
215

216
   // Perform trajectory for each pixel, exit on hit
217

218
  private float m,rs,ri,re,tho;
219
  private int q,raie;
220

221
  m=Mass;
222
  rs=2.*m;
223
  ri=InternalRadius;
224
  re=ExternalRadius;
225
  tho=Angle;
226
  q=-2;
227
  raie=Line;
228

229
  private float d,bmx,db,b,h;
230
  private float rp0,rpp,rps;
231
  private float phi,thi,phd,php,nr,r;
232
  private int nh;
233
  private float zp,fp;
234

235
  // Autosize for image
236
  bmx=1.25*re;
237
  b=0.;
238

239
  h=4.e0f*PI/(float)TRACKPOINTS;
240

241
  // set origin as center of image
242
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
243
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
244
  // angle extracted from cylindric symmetry
245
  phi=atanp(x,y);
246
  phd=atanp(cos(phi)*sin(tho),cos(tho));
247

248
  float up,vp,pp,us,vs,ps;
249

250
  // impact parameter
251
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
252
  // step of impact parameter;
253
  db=bmx/(float)(sizex);
254

255
  up=0.;
256
  vp=1.;
257
      
258
  pp=0.;
259
  nh=0;
260

261
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
262
    
263
  rps=fabs(b/us);
264
  rp0=rps;
265

266
  int ExitOnImpact=0;
267

268
  do
269
  {
270
     nh++;
271
     pp=ps;
272
     up=us;
273
     vp=vs;
274
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
275
     rpp=rps;
276
     rps=fabs(b/us);
277
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;       
278

279
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
280

281
  if (ExitOnImpact==1) {
282
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
283
  }
284
  else
285
  {
286
     zp=0.;
287
     fp=0.;
288
  }
289

290
  barrier(CLK_GLOBAL_MEM_FENCE);
291

292
  zImage[yi+sizex*xi]=(float)zp;
293
  fImage[yi+sizex*xi]=(float)fp;  
294
}
295

296
__kernel void Pixel(__global float *zImage,__global float *fImage,
297
                    __global float *Trajectories,__global int *IdLast,
298
                    uint ImpactParameter,uint TrackPoints,
299
                    float Mass,float InternalRadius,
300
                    float ExternalRadius,float Angle,
301
                    int Line)
302
{
303
   uint xi=(uint)get_global_id(0);
304
   uint yi=(uint)get_global_id(1);
305
   uint sizex=(uint)get_global_size(0);
306
   uint sizey=(uint)get_global_size(1);
307

308
   // Perform trajectory for each pixel
309

310
  float m,rs,ri,re,tho;
311
  int q,raie;
312

313
  m=Mass;
314
  rs=2.*m;
315
  ri=InternalRadius;
316
  re=ExternalRadius;
317
  tho=Angle;
318
  q=-2;
319
  raie=Line;
320

321
  float d,bmx,db,b,h;
322
  float phi,thi,phd,php,nr,r;
323
  int nh;
324
  float zp=0,fp=0;
325

326
  // Autosize for image, 25% greater than external radius
327
  bmx=1.25*re;
328

329
  // Angular step of integration
330
  h=4.e0f*PI/(float)TrackPoints;
331

332
  // Step of Impact Parameter
333
  db=bmx/(2.e0*(float)ImpactParameter);
334

335
  // set origin as center of image
336
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
337
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
338

339
  // angle extracted from cylindric symmetry
340
  phi=atanp(x,y);
341
  phd=atanp(cos(phi)*sin(tho),cos(tho));
342

343
  // Real Impact Parameter
344
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
345

346
  // Integer Impact Parameter
347
  uint bi=(uint)sqrt(x*x+y*y);
348

349
  int HalfLap=0,ExitOnImpact=0,ni;
350

351
  if (bi<ImpactParameter)
352
  {
353
    do
354
    {
355
      php=phd+(float)HalfLap*PI;
356
      nr=php/h;
357
      ni=(int)nr;
358

359
      if (ni<IdLast[bi])
360
      {
361
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
362
      }
363
      else
364
      {
365
        r=Trajectories[bi*TrackPoints+ni];
366
      }
367
           
368
      if ((r<=re)&&(r>=ri))
369
      {
370
        ExitOnImpact=1;
371
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
372
      }
373
              
374
      HalfLap++;
375
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
376

377
  }
378

379
  barrier(CLK_GLOBAL_MEM_FENCE);
380

381
  zImage[yi+sizex*xi]=zp;
382
  fImage[yi+sizex*xi]=fp;
383
}
384

385
__kernel void Circle(__global float *Trajectories,__global int *IdLast,
386
                     __global float *zImage,__global float *fImage,
387
                     int TrackPoints,
388
                     float Mass,float InternalRadius,
389
                     float ExternalRadius,float Angle,
390
                     int Line)
391
{
392
   // Integer Impact Parameter ID 
393
   int bi=get_global_id(0);
394
   // Integer points on circle
395
   int i=get_global_id(1);
396
   // Integer Impact Parameter Size (half of image)
397
   int bmaxi=get_global_size(0);
398
   // Integer Points on circle
399
   int imx=get_global_size(1);
400

401
   // Perform trajectory for each pixel
402

403
  float m,rs,ri,re,tho;
404
  int q,raie;
405

406
  m=Mass;
407
  rs=2.*m;
408
  ri=InternalRadius;
409
  re=ExternalRadius;
410
  tho=Angle;
411
  raie=Line;
412

413
  float bmx,db,b,h;
414
  float phi,thi,phd;
415
  int nh;
416
  float zp=0,fp=0;
417

418
  // Autosize for image
419
  bmx=1.25*re;
420

421
  // Angular step of integration
422
  h=4.e0f*PI/(float)TrackPoints;
423

424
  // impact parameter
425
  b=(float)bi/(float)bmaxi*bmx;
426
  db=bmx/(2.e0*(float)bmaxi);
427

428
  phi=2.*PI/(float)imx*(float)i;
429
  phd=atanp(cos(phi)*sin(tho),cos(tho));
430
  int yi=(int)((float)bi*sin(phi))+bmaxi;
431
  int xi=(int)((float)bi*cos(phi))+bmaxi;
432

433
  int HalfLap=0,ExitOnImpact=0,ni;
434
  float php,nr,r;
435

436
  do
437
  {
438
     php=phd+(float)HalfLap*PI;
439
     nr=php/h;
440
     ni=(int)nr;
441

442
     if (ni<IdLast[bi])
443
     {
444
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
445
     }
446
     else
447
     {
448
        r=Trajectories[bi*TrackPoints+ni];
449
     }
450
           
451
     if ((r<=re)&&(r>=ri))
452
     {
453
        ExitOnImpact=1;
454
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
455
     }
456
              
457
     HalfLap++;
458
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
459
  
460
  zImage[yi+2*bmaxi*xi]=zp;
461
  fImage[yi+2*bmaxi*xi]=fp;  
462

463
  barrier(CLK_GLOBAL_MEM_FENCE);
464

465
}
466

467
__kernel void Trajectory(__global float *Trajectories,
468
                         __global int *IdLast,int TrackPoints,
469
                         float Mass,float InternalRadius,
470
                         float ExternalRadius,float Angle,
471
                         int Line)
472
{
473
  // Integer Impact Parameter ID 
474
  int bi=get_global_id(0);
475
  // Integer Impact Parameter Size (half of image)
476
  int bmaxi=get_global_size(0);
477

478
  // Perform trajectory for each pixel
479

480
  float m,rs,ri,re,tho;
481
  int raie,q;
482

483
  m=Mass;
484
  rs=2.*m;
485
  ri=InternalRadius;
486
  re=ExternalRadius;
487
  tho=Angle;
488
  q=-2;
489
  raie=Line;
490

491
  float d,bmx,db,b,h;
492
  float phi,thi,phd,php,nr,r;
493
  int nh;
494
  float zp,fp;
495

496
  // Autosize for image
497
  bmx=1.25*re;
498

499
  // Angular step of integration
500
  h=4.e0f*PI/(float)TrackPoints;
501

502
  // impact parameter
503
  b=(float)bi/(float)bmaxi*bmx;
504

505
  float up,vp,pp,us,vs,ps;
506

507
  up=0.;
508
  vp=1.;
509
      
510
  pp=0.;
511
  nh=0;
512

513
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
514
  
515
  // b versus us
516
  float bvus=fabs(b/us);
517
  float bvus0=bvus;
518
  Trajectories[bi*TrackPoints+nh]=bvus;
519

520
  do
521
  {
522
     nh++;
523
     pp=ps;
524
     up=us;
525
     vp=vs;
526
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
527
     bvus=fabs(b/us);
528
     Trajectories[bi*TrackPoints+nh]=bvus;
529

530
  } while ((bvus>=rs)&&(bvus<=bvus0));
531

532
  IdLast[bi]=nh;
533

534
  barrier(CLK_GLOBAL_MEM_FENCE);
535
 
536
}
537
"""
538

    
539
def KernelCodeCuda():
540
    BlobCUDA= """
541

542
#define PI (float)3.14159265359
543
#define nbr 256
544

545
#define EINSTEIN 0
546
#define NEWTON 1
547

548
#define TRACKPOINTS 2048
549

550
__device__ float nothing(float x)
551
{
552
  return(x);
553
}
554

555
__device__ float atanp(float x,float y)
556
{
557
  float angle;
558

559
  angle=atan2(y,x);
560

561
  if (angle<0.e0f)
562
    {
563
      angle+=(float)2.e0f*PI;
564
    }
565

566
  return(angle);
567
}
568

569
__device__ float f(float v)
570
{
571
  return(v);
572
}
573

574
#if PHYSICS == NEWTON
575
__device__ float g(float u,float m,float b)
576
{
577
  return (-u);
578
}
579
#else
580
__device__ float g(float u,float m,float b)
581
{
582
  return (3.e0f*m/b*pow(u,2)-u);
583
}
584
#endif
585

586
__device__ void calcul(float *us,float *vs,float up,float vp,
587
                       float h,float m,float b)
588
{
589
  float c0,c1,c2,c3,d0,d1,d2,d3;
590

591
  c0=h*f(vp);
592
  c1=h*f(vp+c0/2.);
593
  c2=h*f(vp+c1/2.);
594
  c3=h*f(vp+c2);
595
  d0=h*g(up,m,b);
596
  d1=h*g(up+d0/2.,m,b);
597
  d2=h*g(up+d1/2.,m,b);
598
  d3=h*g(up+d2,m,b);
599

600
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
601
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
602
}
603

604
__device__ void rungekutta(float *ps,float *us,float *vs,
605
                           float pp,float up,float vp,
606
                           float h,float m,float b)
607
{
608
  calcul(us,vs,up,vp,h,m,b);
609
  *ps=pp+h;
610
}
611

612
__device__ float decalage_spectral(float r,float b,float phi,
613
                                   float tho,float m)
614
{
615
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
616
}
617

618
__device__ float spectre(float rf,int q,float b,float db,
619
                         float h,float r,float m,float bss)
620
{
621
  float flx;
622

623
  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
624
  return(flx);
625
}
626

627
__device__ float spectre_cn(float rf32,float b32,float db32,
628
                            float h32,float r32,float m32,float bss32)
629
{
630
  
631
#define MYFLOAT float
632

633
  MYFLOAT rf=(MYFLOAT)(rf32);
634
  MYFLOAT b=(MYFLOAT)(b32);
635
  MYFLOAT db=(MYFLOAT)(db32);
636
  MYFLOAT h=(MYFLOAT)(h32);
637
  MYFLOAT r=(MYFLOAT)(r32);
638
  MYFLOAT m=(MYFLOAT)(m32);
639
  MYFLOAT bss=(MYFLOAT)(bss32);
640

641
  MYFLOAT flx;
642
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
643
  int fi,posfreq;
644

645
#define planck 6.62e-34
646
#define k 1.38e-23
647
#define c2 9.e16
648
#define temp 3.e7
649
#define m_point 1.
650

651
#define lplanck (log(6.62)-34.*log(10.))
652
#define lk (log(1.38)-23.*log(10.))
653
#define lc2 (log(9.)+16.*log(10.))
654

655
  MYFLOAT v=1.-3./r;
656

657
  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 ));
658

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

661
  flux_int=0.;
662
  flx=0.;
663

664
  for (fi=0;fi<nbr;fi++)
665
    {
666
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
667
      nu_rec=nu_em*rf;
668
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
669
      if ((posfreq>0)&&(posfreq<nbr))
670
        {
671
          // Initial version
672
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
673
          // Version with log used
674
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
675
          // flux_int*=pow(rf,3)*b*db*h;
676
          //flux_int*=exp(3.*log(rf))*b*db*h;
677
          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;
678

679
          flx+=flux_int;
680
        }
681
    }
682

683
  return((float)(flx));
684
}
685

686
__device__ void impact(float phi,float r,float b,float tho,float m,
687
                       float *zp,float *fp,
688
                       int q,float db,
689
                       float h,int raie)
690
{
691
  float flx,rf,bss;
692

693
  rf=decalage_spectral(r,b,phi,tho,m);
694

695
  if (raie==0)
696
    {
697
      bss=1.e19;
698
      flx=spectre_cn(rf,b,db,h,r,m,bss);
699
    }
700
  else
701
    { 
702
      bss=2.;
703
      flx=spectre(rf,q,b,db,h,r,m,bss);
704
    }
705
  
706
  *zp=1./rf;
707
  *fp=flx;
708

709
}
710

711
__global__ void EachPixel(float *zImage,float *fImage,
712
                          float Mass,float InternalRadius,
713
                          float ExternalRadius,float Angle,
714
                          int Line)
715
{
716
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
717
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
718
   uint sizex=(uint)gridDim.x*blockDim.x;
719
   uint sizey=(uint)gridDim.y*blockDim.y;
720

721
   // Perform trajectory for each pixel, exit on hit
722

723
  float m,rs,ri,re,tho;
724
  int q,raie;
725

726
  m=Mass;
727
  rs=2.*m;
728
  ri=InternalRadius;
729
  re=ExternalRadius;
730
  tho=Angle;
731
  q=-2;
732
  raie=Line;
733

734
  float d,bmx,db,b,h;
735
  float rp0,rpp,rps;
736
  float phi,thi,phd,php,nr,r;
737
  int nh;
738
  float zp,fp;
739

740
  // Autosize for image
741
  bmx=1.25*re;
742
  b=0.;
743

744
  h=4.e0f*PI/(float)TRACKPOINTS;
745

746
  // set origin as center of image
747
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
748
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
749
  // angle extracted from cylindric symmetry
750
  phi=atanp(x,y);
751
  phd=atanp(cos(phi)*sin(tho),cos(tho));
752

753
  float up,vp,pp,us,vs,ps;
754

755
  // impact parameter
756
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
757
  // step of impact parameter;
758
//  db=bmx/(float)(sizex/2);
759
  db=bmx/(float)(sizex);
760

761
  up=0.;
762
  vp=1.;
763
      
764
  pp=0.;
765
  nh=0;
766

767
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
768
    
769
  rps=fabs(b/us);
770
  rp0=rps;
771

772
  int ExitOnImpact=0;
773

774
  do
775
  {
776
     nh++;
777
     pp=ps;
778
     up=us;
779
     vp=vs;
780
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
781
     rpp=rps;          
782
     rps=fabs(b/us);
783
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;          
784

785
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
786

787
  if (ExitOnImpact==1) {
788
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
789
  }
790
  else
791
  {
792
     zp=0.;
793
     fp=0.;
794
  }
795

796
  __syncthreads();
797

798
  zImage[yi+sizex*xi]=(float)zp;
799
  fImage[yi+sizex*xi]=(float)fp;
800
}
801

802
__global__ void Pixel(float *zImage,float *fImage,
803
                      float *Trajectories,int *IdLast,
804
                      uint ImpactParameter,uint TrackPoints,
805
                      float Mass,float InternalRadius,
806
                      float ExternalRadius,float Angle,
807
                      int Line)
808
{
809
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
810
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
811
   uint sizex=(uint)gridDim.x*blockDim.x;
812
   uint sizey=(uint)gridDim.y*blockDim.y;
813
//  uint xi=(uint)blockIdx.x;
814
//  uint yi=(uint)blockIdx.y;
815
//  uint sizex=(uint)gridDim.x*blockDim.x;
816
//  uint sizey=(uint)gridDim.y*blockDim.y;
817

818
  // Perform trajectory for each pixel
819

820
  float m,rs,ri,re,tho;
821
  int q,raie;
822

823
  m=Mass;
824
  rs=2.*m;
825
  ri=InternalRadius;
826
  re=ExternalRadius;
827
  tho=Angle;
828
  q=-2;
829
  raie=Line;
830

831
  float d,bmx,db,b,h;
832
  float phi,thi,phd,php,nr,r;
833
  int nh;
834
  float zp=0,fp=0;
835
  // Autosize for image, 25% greater than external radius
836
  bmx=1.25*re;
837

838
  // Angular step of integration
839
  h=4.e0f*PI/(float)TrackPoints;
840

841
  // Step of Impact Parameter
842
  db=bmx/(2.e0*(float)ImpactParameter);
843

844
  // set origin as center of image
845
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
846
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
847
  // angle extracted from cylindric symmetry
848
  phi=atanp(x,y);
849
  phd=atanp(cos(phi)*sin(tho),cos(tho));
850

851
  // Real Impact Parameter
852
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
853

854
  // Integer Impact Parameter
855
  uint bi=(uint)sqrt(x*x+y*y);
856

857
  int HalfLap=0,ExitOnImpact=0,ni;
858

859
  if (bi<ImpactParameter)
860
  {
861
    do
862
    {
863
      php=phd+(float)HalfLap*PI;
864
      nr=php/h;
865
      ni=(int)nr;
866

867
      if (ni<IdLast[bi])
868
      {
869
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
870
      }
871
      else
872
      {
873
        r=Trajectories[bi*TrackPoints+ni];
874
      }
875
           
876
      if ((r<=re)&&(r>=ri))
877
      {
878
        ExitOnImpact=1;
879
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
880
      }
881
              
882
      HalfLap++;
883
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
884

885
  }
886

887
  zImage[yi+sizex*xi]=zp;
888
  fImage[yi+sizex*xi]=fp;
889
}
890

891
__global__ void Circle(float *Trajectories,int *IdLast,
892
                       float *zImage,float *fImage,
893
                       int TrackPoints,
894
                       float Mass,float InternalRadius,
895
                       float ExternalRadius,float Angle,
896
                       int Line)
897
{
898
   // Integer Impact Parameter ID 
899
   int bi=blockIdx.x*blockDim.x+threadIdx.x;
900
   // Integer points on circle
901
   int i=blockIdx.y*blockDim.y+threadIdx.y;
902
   // Integer Impact Parameter Size (half of image)
903
   int bmaxi=gridDim.x*blockDim.x;
904
   // Integer Points on circle
905
   int imx=gridDim.y*blockDim.y;
906

907

908
   // Integer Impact Parameter ID 
909
//   int bi=blockIdx.x;
910
   // Integer points on circle
911
//   int i=blockIdx.y;
912
   // Integer Impact Parameter Size (half of image)
913
//   int bmaxi=gridDim.x*blockDim.x;
914
   // Integer Points on circle
915
//   int imx=gridDim.y*blockDim.y;
916

917
   // Perform trajectory for each pixel
918

919
  float m,rs,ri,re,tho;
920
  int q,raie;
921

922
  m=Mass;
923
  rs=2.*m;
924
  ri=InternalRadius;
925
  re=ExternalRadius;
926
  tho=Angle;
927
  raie=Line;
928

929
  float bmx,db,b,h;
930
  float phi,thi,phd;
931
  int nh;
932
  float zp=0,fp=0;
933

934
  // Autosize for image
935
  bmx=1.25*re;
936

937
  // Angular step of integration
938
  h=4.e0f*PI/(float)TrackPoints;
939

940
  // impact parameter
941
  b=(float)bi/(float)bmaxi*bmx;
942
  db=bmx/(2.e0*(float)bmaxi);
943

944
  phi=2.*PI/(float)imx*(float)i;
945
  phd=atanp(cos(phi)*sin(tho),cos(tho));
946
  int yi=(int)((float)bi*sin(phi))+bmaxi;
947
  int xi=(int)((float)bi*cos(phi))+bmaxi;
948

949
  int HalfLap=0,ExitOnImpact=0,ni;
950
  float php,nr,r;
951

952
  do
953
  {
954
     php=phd+(float)HalfLap*PI;
955
     nr=php/h;
956
     ni=(int)nr;
957

958
     if (ni<IdLast[bi])
959
     {
960
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
961
     }
962
     else
963
     {
964
        r=Trajectories[bi*TrackPoints+ni];
965
     }
966
           
967
     if ((r<=re)&&(r>=ri))
968
     {
969
        ExitOnImpact=1;
970
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
971
     }
972
              
973
     HalfLap++;
974
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
975
  
976
  zImage[yi+2*bmaxi*xi]=zp;
977
  fImage[yi+2*bmaxi*xi]=fp;  
978

979
}
980

981
__global__ void Trajectory(float *Trajectories,
982
                           int *IdLast,int TrackPoints,
983
                           float Mass,float InternalRadius,
984
                           float ExternalRadius,float Angle,
985
                           int Line)
986
{
987
  // Integer Impact Parameter ID 
988
  //  int bi=blockIdx.x;
989
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
990
  // Integer Impact Parameter Size (half of image)
991
  int bmaxi=gridDim.x*blockDim.x;
992

993
  // Perform trajectory for each pixel
994

995
  float m,rs,ri,re,tho;
996
  int raie,q;
997

998
  m=Mass;
999
  rs=2.*m;
1000
  ri=InternalRadius;
1001
  re=ExternalRadius;
1002
  tho=Angle;
1003
  q=-2;
1004
  raie=Line;
1005

1006
  float d,bmx,db,b,h;
1007
  float phi,thi,phd,php,nr,r;
1008
  int nh;
1009
  float zp,fp;
1010

1011
  // Autosize for image
1012
  bmx=1.25*re;
1013

1014
  // Angular step of integration
1015
  h=4.e0f*PI/(float)TrackPoints;
1016

1017
  // impact parameter
1018
  b=(float)bi/(float)bmaxi*bmx;
1019

1020
  float up,vp,pp,us,vs,ps;
1021

1022
  up=0.;
1023
  vp=1.;
1024
      
1025
  pp=0.;
1026
  nh=0;
1027

1028
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1029
  
1030
  // b versus us
1031
  float bvus=fabs(b/us);
1032
  float bvus0=bvus;
1033
  Trajectories[bi*TrackPoints+nh]=bvus;
1034

1035
  do
1036
  {
1037
     nh++;
1038
     pp=ps;
1039
     up=us;
1040
     vp=vs;
1041
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1042
     bvus=fabs(b/us);
1043
     Trajectories[bi*TrackPoints+nh]=bvus;
1044

1045
  } while ((bvus>=rs)&&(bvus<=bvus0));
1046

1047
  IdLast[bi]=nh;
1048

1049
}
1050
"""
1051
    return(BlobCUDA)
1052
    
1053
# def ImageOutput(sigma,prefix):
1054
#     from PIL import Image
1055
#     Max=sigma.max()
1056
#     Min=sigma.min()
1057
    
1058
#     # Normalize value as 8bits Integer
1059
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1060
#     image = Image.fromarray(SigmaInt)
1061
#     image.save("%s.jpg" % prefix)
1062

    
1063
def ImageOutput(sigma,prefix,Colors):
1064
    import matplotlib.pyplot as plt
1065
    if Colors == 'Red2Yellow':
1066
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1067
    else:
1068
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1069

    
1070
def BlackHoleCL(zImage,fImage,InputCL):
1071

    
1072
    kernel_params = {}
1073

    
1074
    print(InputCL)
1075

    
1076
    Device=InputCL['Device']
1077
    GpuStyle=InputCL['GpuStyle']
1078
    VariableType=InputCL['VariableType']
1079
    Size=InputCL['Size']
1080
    Mass=InputCL['Mass']
1081
    InternalRadius=InputCL['InternalRadius']
1082
    ExternalRadius=InputCL['ExternalRadius']
1083
    Angle=InputCL['Angle']
1084
    Method=InputCL['Method']
1085
    TrackPoints=InputCL['TrackPoints']
1086
    Physics=InputCL['Physics']
1087

    
1088
    PhysicsList=DictionariesAPI()
1089
    
1090
    if InputCL['BlackBody']:
1091
        # Spectrum is Black Body one
1092
        Line=0
1093
    else:
1094
        # Spectrum is Monochromatic Line one
1095
        Line=1
1096

    
1097
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1098
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1099

    
1100
    # Je detecte un peripherique GPU dans la liste des peripheriques
1101
    Id=0
1102
    HasXPU=False
1103
    for platform in cl.get_platforms():
1104
        for device in platform.get_devices():
1105
            if Id==Device:
1106
                XPU=device
1107
                print("CPU/GPU selected: ",device.name.lstrip())
1108
                HasXPU=True
1109
            Id+=1
1110

    
1111
    if HasXPU==False:
1112
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1113
        sys.exit()                
1114
        
1115
    ctx = cl.Context([XPU])
1116
    queue = cl.CommandQueue(ctx,
1117
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1118

    
1119
    
1120
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1121

    
1122
    BuildOptions="-cl-mad-enable -DPHYSICS=%i " % (PhysicsList[Physics])
1123

    
1124
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1125
    
1126
    # Je recupere les flag possibles pour les buffers
1127
    mf = cl.mem_flags
1128

    
1129
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
1130
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1131
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1132
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1133
    
1134
    start_time=time.time() 
1135

    
1136
    if Method=='EachPixel':
1137
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1138
                                       None,zImageCL,fImageCL,
1139
                                       numpy.float32(Mass),
1140
                                       numpy.float32(InternalRadius),
1141
                                       numpy.float32(ExternalRadius),
1142
                                       numpy.float32(Angle),
1143
                                       numpy.int32(Line))
1144
        CLLaunch.wait()
1145
    elif Method=='TrajectoCircle':
1146
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1147
                                        None,TrajectoriesCL,IdLastCL,
1148
                                        numpy.uint32(Trajectories.shape[1]),
1149
                                        numpy.float32(Mass),
1150
                                        numpy.float32(InternalRadius),
1151
                                        numpy.float32(ExternalRadius),
1152
                                        numpy.float32(Angle),
1153
                                        numpy.int32(Line))
1154
        
1155
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1156
                                           zImage.shape[0]*4),None,
1157
                                    TrajectoriesCL,IdLastCL,
1158
                                    zImageCL,fImageCL,
1159
                                    numpy.uint32(Trajectories.shape[1]),
1160
                                    numpy.float32(Mass),
1161
                                    numpy.float32(InternalRadius),
1162
                                    numpy.float32(ExternalRadius),
1163
                                    numpy.float32(Angle),
1164
                                    numpy.int32(Line))
1165
        CLLaunch.wait()
1166
    else:
1167
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1168
                                        None,TrajectoriesCL,IdLastCL,
1169
                                        numpy.uint32(Trajectories.shape[1]),
1170
                                        numpy.float32(Mass),
1171
                                        numpy.float32(InternalRadius),
1172
                                        numpy.float32(ExternalRadius),
1173
                                        numpy.float32(Angle),
1174
                                        numpy.int32(Line))
1175
        
1176
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],
1177
                                          zImage.shape[1]),None,
1178
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1179
                                   numpy.uint32(Trajectories.shape[0]),
1180
                                   numpy.uint32(Trajectories.shape[1]),
1181
                                        numpy.float32(Mass),
1182
                                        numpy.float32(InternalRadius),
1183
                                        numpy.float32(ExternalRadius),
1184
                                        numpy.float32(Angle),
1185
                                        numpy.int32(Line))
1186
        CLLaunch.wait()
1187
    
1188
    compute = time.time()-start_time
1189
    
1190
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1191
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1192
    cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1193
    cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1194
    elapsed = time.time()-start_time
1195
    print("\nCompute Time : %f" % compute)
1196
    print("Elapsed Time : %f\n" % elapsed)
1197

    
1198
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1199
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1200
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1201
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1202
    zImageCL.release()
1203
    fImageCL.release()
1204

    
1205
    TrajectoriesCL.release()
1206
    IdLastCL.release()
1207

    
1208
    return(elapsed)
1209

    
1210
def BlackHoleCUDA(zImage,fImage,InputCL):
1211

    
1212
    kernel_params = {}
1213

    
1214
    print(InputCL)
1215

    
1216
    Device=InputCL['Device']
1217
    GpuStyle=InputCL['GpuStyle']
1218
    VariableType=InputCL['VariableType']
1219
    Size=InputCL['Size']
1220
    Mass=InputCL['Mass']
1221
    InternalRadius=InputCL['InternalRadius']
1222
    ExternalRadius=InputCL['ExternalRadius']
1223
    Angle=InputCL['Angle']
1224
    Method=InputCL['Method']
1225
    TrackPoints=InputCL['TrackPoints']
1226
    Physics=InputCL['Physics']
1227

    
1228
    PhysicsList=DictionariesAPI()
1229
    
1230
    if InputCL['BlackBody']:
1231
        # Spectrum is Black Body one
1232
        Line=0
1233
    else:
1234
        # Spectrum is Monochromatic Line one
1235
        Line=1
1236

    
1237
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1238
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1239

    
1240
    try:
1241
        # For PyCUDA import
1242
        import pycuda.driver as cuda
1243
        from pycuda.compiler import SourceModule
1244
        
1245
        cuda.init()
1246
        for Id in range(cuda.Device.count()):
1247
            if Id==Device:
1248
                XPU=cuda.Device(Id)
1249
                print("GPU selected %s" % XPU.name())
1250
        print
1251

    
1252
    except ImportError:
1253
        print("Platform does not seem to support CUDA")
1254

    
1255
    Context=XPU.make_context()
1256

    
1257
    try:
1258
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i' % (PhysicsList[Physics])])
1259
        print("Compilation seems to be OK")
1260
    except:
1261
        print("Compilation seems to break")
1262

    
1263
    EachPixelCU=mod.get_function("EachPixel")
1264
    TrajectoryCU=mod.get_function("Trajectory")
1265
    PixelCU=mod.get_function("Pixel")
1266
    CircleCU=mod.get_function("Circle")
1267
    
1268
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1269
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1270
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1271
    cuda.memcpy_htod(zImageCU, zImage)
1272
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1273
    cuda.memcpy_htod(zImageCU, fImage)
1274
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1275
    cuda.memcpy_htod(IdLastCU, IdLast)
1276

    
1277
    start_time=time.time() 
1278

    
1279
    if Method=='EachPixel':
1280
        EachPixelCU(zImageCU,fImageCU,
1281
                    numpy.float32(Mass),
1282
                    numpy.float32(InternalRadius),
1283
                    numpy.float32(ExternalRadius),
1284
                    numpy.float32(Angle),
1285
                    numpy.int32(Line),
1286
                    grid=(zImage.shape[0]/32,zImage.shape[1]/32),
1287
                    block=(32,32,1))
1288
    elif Method=='TrajectoCircle':
1289
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1290
                     numpy.uint32(Trajectories.shape[1]),
1291
                     numpy.float32(Mass),
1292
                     numpy.float32(InternalRadius),
1293
                     numpy.float32(ExternalRadius),
1294
                     numpy.float32(Angle),
1295
                     numpy.int32(Line),
1296
                     grid=(Trajectories.shape[0]/32,1),
1297
                     block=(32,1,1))
1298
        
1299
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1300
                 numpy.uint32(Trajectories.shape[1]),
1301
                 numpy.float32(Mass),
1302
                 numpy.float32(InternalRadius),
1303
                 numpy.float32(ExternalRadius),
1304
                 numpy.float32(Angle),
1305
                 numpy.int32(Line),
1306
                 grid=(Trajectories.shape[0]/32,zImage.shape[0]*4/32),
1307
                 block=(32,32,1))
1308
    else:
1309
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1310
                     numpy.uint32(Trajectories.shape[1]),
1311
                     numpy.float32(Mass),
1312
                     numpy.float32(InternalRadius),
1313
                     numpy.float32(ExternalRadius),
1314
                     numpy.float32(Angle),
1315
                     numpy.int32(Line),
1316
                     grid=(Trajectories.shape[0]/32,1),
1317
                     block=(32,1,1))
1318
        
1319
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1320
                numpy.uint32(Trajectories.shape[0]),
1321
                numpy.uint32(Trajectories.shape[1]),
1322
                numpy.float32(Mass),
1323
                numpy.float32(InternalRadius),
1324
                numpy.float32(ExternalRadius),
1325
                numpy.float32(Angle),
1326
                numpy.int32(Line),
1327
                grid=(zImage.shape[0]/32,zImage.shape[1]/32,1),
1328
                block=(32,32,1))
1329
    
1330

    
1331
    compute = time.time()-start_time
1332

    
1333
    cuda.memcpy_dtoh(zImage,zImageCU)
1334
    cuda.memcpy_dtoh(fImage,fImageCU)
1335
    elapsed = time.time()-start_time
1336
    print("\nCompute Time : %f" % compute)
1337
    print("Elapsed Time : %f\n" % elapsed)
1338

    
1339
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1340
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1341
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1342
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1343

    
1344
    Context.pop()
1345
    
1346
    Context.detach()
1347

    
1348
    return(elapsed)
1349

    
1350
if __name__=='__main__':
1351
        
1352
    GpuStyle = 'OpenCL'
1353
    Mass = 1.
1354
    # Internal Radius 3 times de Schwarzschild Radius
1355
    InternalRadius=6.*Mass
1356
    #
1357
    ExternalRadius=12.
1358
    #
1359
    # Angle with normal to disc 10 degrees
1360
    Angle = numpy.pi/180.*(90.-10.)
1361
    # Radiation of disc : BlackBody or Monochromatic
1362
    BlackBody = False
1363
    # Size of image
1364
    Size=256
1365
    # Variable Type
1366
    VariableType='FP32'
1367
    # ?
1368
    q=-2
1369
    # Method of resolution
1370
    Method='TrajectoPixel'
1371
    # Colors for output image
1372
    Colors='Greyscale'
1373
    # Physics
1374
    Physics='Einstein'
1375
    # No output as image
1376
    NoImage = False
1377
    
1378
    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> -v <FP32/FP64>'
1379

    
1380
    try:
1381
        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="])
1382
    except getopt.GetoptError:
1383
        print(HowToUse % sys.argv[0])
1384
        sys.exit(2)
1385

    
1386
    # List of Devices
1387
    Devices=[]
1388
    Alu={}
1389
        
1390
    for opt, arg in opts:
1391
        if opt == '-h':
1392
            print(HowToUse % sys.argv[0])
1393
            
1394
            print("\nInformations about devices detected under OpenCL API:")
1395
            # For PyOpenCL import
1396
            try:
1397
                import pyopencl as cl
1398
                Id=0
1399
                for platform in cl.get_platforms():
1400
                    for device in platform.get_devices():
1401
                        #deviceType=cl.device_type.to_string(device.type)
1402
                        deviceType="xPU"
1403
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
1404
                        Id=Id+1
1405
                        
1406
            except:
1407
                print("Your platform does not seem to support OpenCL")
1408
                
1409
            print("\nInformations about devices detected under CUDA API:")
1410
                # For PyCUDA import
1411
            try:
1412
                import pycuda.driver as cuda
1413
                cuda.init()
1414
                for Id in range(cuda.Device.count()):
1415
                    device=cuda.Device(Id)
1416
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
1417
                print
1418
            except:
1419
                print("Your platform does not seem to support CUDA")
1420
        
1421
            sys.exit()
1422
        
1423
        elif opt in ("-d", "--device"):
1424
#            Devices.append(int(arg))
1425
            Device=int(arg)
1426
        elif opt in ("-g", "--gpustyle"):
1427
            GpuStyle = arg
1428
        elif opt in ("-v", "--variabletype"):
1429
            VariableType = arg
1430
        elif opt in ("-s", "--size"):
1431
            Size = int(arg)
1432
        elif opt in ("-m", "--mass"):
1433
            Mass = float(arg)
1434
        elif opt in ("-i", "--internal"):
1435
            InternalRadius = float(arg)
1436
        elif opt in ("-e", "--external"):
1437
            ExternalRadius = float(arg)
1438
        elif opt in ("-a", "--angle"):
1439
            Angle = numpy.pi/180.*(90.-float(arg))
1440
        elif opt in ("-b", "--blackbody"):
1441
            BlackBody = True
1442
        elif opt in ("-n", "--noimage"):
1443
            NoImage = True
1444
        elif opt in ("-t", "--method"):
1445
            Method = arg
1446
        elif opt in ("-c", "--colors"):
1447
            Colors = arg
1448
        elif opt in ("-p", "--physics"):
1449
            Physics = arg
1450

    
1451
    print("Device Identification selected : %s" % Device)
1452
    print("GpuStyle used : %s" % GpuStyle)
1453
    print("VariableType : %s" % VariableType)
1454
    print("Size : %i" % Size)
1455
    print("Mass : %f" % Mass)
1456
    print("Internal Radius : %f" % InternalRadius)
1457
    print("External Radius : %f" % ExternalRadius)
1458
    print("Angle with normal of (in radians) : %f" % Angle)
1459
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
1460
    print("Method of resolution : %s" % Method)
1461
    print("Colors for output images : %s" % Colors)
1462
    print("Physics used for Trajectories : %s" % Physics)
1463

    
1464
    if GpuStyle=='CUDA':
1465
        print("\nSelection of CUDA device") 
1466
        try:
1467
            # For PyCUDA import
1468
            import pycuda.driver as cuda
1469
            
1470
            cuda.init()
1471
            for Id in range(cuda.Device.count()):
1472
                device=cuda.Device(Id)
1473
                print("Device #%i of type GPU : %s" % (Id,device.name()))
1474
                if Id in Devices:
1475
                    Alu[Id]='GPU'
1476
            
1477
        except ImportError:
1478
            print("Platform does not seem to support CUDA")
1479

    
1480
    if GpuStyle=='OpenCL':
1481
        print("\nSelection of OpenCL device") 
1482
        try:
1483
            # For PyOpenCL import
1484
            import pyopencl as cl
1485
            Id=0
1486
            for platform in cl.get_platforms():
1487
                for device in platform.get_devices():
1488
                    #deviceType=cl.device_type.to_string(device.type)
1489
                    deviceType="xPU"
1490
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
1491

    
1492
                    if Id in Devices:
1493
                    # Set the Alu as detected Device Type
1494
                        Alu[Id]=deviceType
1495
                    Id=Id+1
1496
        except ImportError:
1497
            print("Platform does not seem to support OpenCL")
1498

    
1499
#    print(Devices,Alu)
1500

    
1501
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
1502
    TrackPoints=2048
1503
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1504
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1505

    
1506
    InputCL={}
1507
    InputCL['Device']=Device
1508
    InputCL['GpuStyle']=GpuStyle
1509
    InputCL['VariableType']=VariableType
1510
    InputCL['Size']=Size
1511
    InputCL['Mass']=Mass
1512
    InputCL['InternalRadius']=InternalRadius
1513
    InputCL['ExternalRadius']=ExternalRadius
1514
    InputCL['Angle']=Angle
1515
    InputCL['BlackBody']=BlackBody
1516
    InputCL['Method']=Method
1517
    InputCL['TrackPoints']=TrackPoints
1518
    InputCL['Physics']=Physics
1519

    
1520
    if GpuStyle=='OpenCL':
1521
        duration=BlackHoleCL(zImage,fImage,InputCL)
1522
    else:
1523
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
1524
        
1525
    Hostname=gethostname()
1526
    Date=time.strftime("%Y%m%d_%H%M%S")
1527
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1528
    
1529
    
1530
    if not NoImage:
1531
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
1532
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)