Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 222

Historique | Voir | Annoter | Télécharger (39,17 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
import pyopencl as cl
29
import numpy
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
from socket import gethostname
36

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

    
41
BlobOpenCL= """
42

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

46
#define EINSTEIN 0
47
#define NEWTON 1
48

49
#define TRACKPOINTS 2048
50

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

55
  angle=atan2(y,x);
56

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

62
  return angle;
63
}
64

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

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

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

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

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

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

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

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

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

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

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

138
  MYFLOAT flx;
139
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
140
  int fi,posfreq;
141

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

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

152
  MYFLOAT v=1.-3./r;
153

154
  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 ));
155

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

158
  flux_int=0.;
159
  flx=0.;
160

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

176
          flx+=flux_int;
177
        }
178
    }
179

180
  return((float)(flx));
181
}
182

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

190
  rf=decalage_spectral(r,b,phi,tho,m);
191

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

206
}
207

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

218
   // Perform trajectory for each pixel, exit on hit
219

220
  private float m,rs,ri,re,tho;
221
  private int q,raie;
222

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

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

237
  // Autosize for image
238
  bmx=1.25*re;
239
  b=0.;
240

241
  h=4.e0f*PI/(float)TRACKPOINTS;
242

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

250
  float up,vp,pp,us,vs,ps;
251

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

257
  up=0.;
258
  vp=1.;
259
      
260
  pp=0.;
261
  nh=0;
262

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

268
  int ExitOnImpact=0;
269

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

281
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
282

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

292
  barrier(CLK_GLOBAL_MEM_FENCE);
293

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

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

310
   // Perform trajectory for each pixel
311

312
  float m,rs,ri,re,tho;
313
  int q,raie;
314

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

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

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

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

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

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

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

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

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

351
  int HalfLap=0,ExitOnImpact=0,ni;
352

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

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

379
  }
380

381
  barrier(CLK_GLOBAL_MEM_FENCE);
382

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

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

403
   // Perform trajectory for each pixel
404

405
  float m,rs,ri,re,tho;
406
  int q,raie;
407

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

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

420
  // Autosize for image
421
  bmx=1.25*re;
422

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

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

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

435
  int HalfLap=0,ExitOnImpact=0,ni;
436
  float php,nr,r;
437

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

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

465
  barrier(CLK_GLOBAL_MEM_FENCE);
466

467
}
468

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

480
  // Perform trajectory for each pixel
481

482
  float m,rs,ri,re,tho;
483
  int raie,q;
484

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

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

498
  // Autosize for image
499
  bmx=1.25*re;
500

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

504
  // impact parameter
505
  b=(float)bi/(float)bmaxi*bmx;
506

507
  float up,vp,pp,us,vs,ps;
508

509
  up=0.;
510
  vp=1.;
511
      
512
  pp=0.;
513
  nh=0;
514

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

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

532
  } while ((bvus>=rs)&&(bvus<=bvus0));
533

534
  IdLast[bi]=nh;
535

536
  barrier(CLK_GLOBAL_MEM_FENCE);
537
 
538
}
539
"""
540

    
541
def KernelCodeCuda():
542
    BlobCUDA= """
543

544
#define PI (float)3.14159265359
545
#define nbr 256
546

547
#define EINSTEIN 0
548
#define NEWTON 1
549

550
#define TRACKPOINTS 2048
551

552
__device__ float nothing(float x)
553
{
554
  return(x);
555
}
556

557
__device__ float atanp(float x,float y)
558
{
559
  float angle;
560

561
  angle=atan2(y,x);
562

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

568
  return(angle);
569
}
570

571
__device__ float f(float v)
572
{
573
  return(v);
574
}
575

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

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

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

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

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

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

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

625
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
626
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
627
  return(flx);
628
}
629

630
__device__ float spectre_cn(float rf32,float b32,float db32,
631
                            float h32,float r32,float m32,float bss32)
632
{
633
  
634
#define MYFLOAT float
635

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

644
  MYFLOAT flx;
645
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
646
  int fi,posfreq;
647

648
#define planck 6.62e-34
649
#define k 1.38e-23
650
#define c2 9.e16
651
#define temp 3.e7
652
#define m_point 1.
653

654
#define lplanck (log(6.62)-34.*log(10.))
655
#define lk (log(1.38)-23.*log(10.))
656
#define lc2 (log(9.)+16.*log(10.))
657

658
  MYFLOAT v=1.-3./r;
659

660
  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 ));
661

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

664
  flux_int=0.;
665
  flx=0.;
666

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

682
          flx+=flux_int;
683
        }
684
    }
685

686
  return((float)(flx));
687
}
688

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

696
  rf=decalage_spectral(r,b,phi,tho,m);
697

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

712
}
713

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

724
   // Perform trajectory for each pixel, exit on hit
725

726
  float m,rs,ri,re,tho;
727
  int q,raie;
728

729
  m=Mass;
730
  rs=2.*m;
731
  ri=InternalRadius;
732
  re=ExternalRadius;
733
  tho=Angle;
734
  q=-2;
735
  raie=Line;
736

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

743
  // Autosize for image
744
  bmx=1.25*re;
745
  b=0.;
746

747
  h=4.e0f*PI/(float)TRACKPOINTS;
748

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

756
  float up,vp,pp,us,vs,ps;
757

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

764
  up=0.;
765
  vp=1.;
766
      
767
  pp=0.;
768
  nh=0;
769

770
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
771
    
772
  rps=fabs(b/us);
773
  rp0=rps;
774

775
  int ExitOnImpact=0;
776

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

788
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
789

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

799
  __syncthreads();
800

801
  zImage[yi+sizex*xi]=(float)zp;
802
  fImage[yi+sizex*xi]=(float)fp;
803
}
804

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

817
  // Perform trajectory for each pixel
818

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

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

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

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

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

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

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

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

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

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

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

884
  }
885

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

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

906
   // Perform trajectory for each pixel
907

908
  float m,rs,ri,re,tho;
909
  int q,raie;
910

911
  m=Mass;
912
  rs=2.*m;
913
  ri=InternalRadius;
914
  re=ExternalRadius;
915
  tho=Angle;
916
  raie=Line;
917

918
  float bmx,db,b,h;
919
  float phi,thi,phd;
920
  int nh;
921
  float zp=0,fp=0;
922

923
  // Autosize for image
924
  bmx=1.25*re;
925

926
  // Angular step of integration
927
  h=4.e0f*PI/(float)TrackPoints;
928

929
  // impact parameter
930
  b=(float)bi/(float)bmaxi*bmx;
931
  db=bmx/(2.e0*(float)bmaxi);
932

933
  phi=2.*PI/(float)imx*(float)i;
934
  phd=atanp(cos(phi)*sin(tho),cos(tho));
935
  int yi=(int)((float)bi*sin(phi))+bmaxi;
936
  int xi=(int)((float)bi*cos(phi))+bmaxi;
937

938
  int HalfLap=0,ExitOnImpact=0,ni;
939
  float php,nr,r;
940

941
  do
942
  {
943
     php=phd+(float)HalfLap*PI;
944
     nr=php/h;
945
     ni=(int)nr;
946

947
     if (ni<IdLast[bi])
948
     {
949
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
950
     }
951
     else
952
     {
953
        r=Trajectories[bi*TrackPoints+ni];
954
     }
955
           
956
     if ((r<=re)&&(r>=ri))
957
     {
958
        ExitOnImpact=1;
959
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
960
     }
961
              
962
     HalfLap++;
963
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
964
  
965
  zImage[yi+2*bmaxi*xi]=zp;
966
  fImage[yi+2*bmaxi*xi]=fp;  
967

968
}
969

970
__global__ void Trajectory(float *Trajectories,
971
                           int *IdLast,int TrackPoints,
972
                           float Mass,float InternalRadius,
973
                           float ExternalRadius,float Angle,
974
                           int Line)
975
{
976
  // Integer Impact Parameter ID 
977
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
978
  // Integer Impact Parameter Size (half of image)
979
  int bmaxi=gridDim.x*blockDim.x;
980

981
  // Perform trajectory for each pixel
982

983
  float m,rs,ri,re,tho;
984
  int raie,q;
985

986
  m=Mass;
987
  rs=2.*m;
988
  ri=InternalRadius;
989
  re=ExternalRadius;
990
  tho=Angle;
991
  q=-2;
992
  raie=Line;
993

994
  float d,bmx,db,b,h;
995
  float phi,thi,phd,php,nr,r;
996
  int nh;
997
  float zp,fp;
998

999
  // Autosize for image
1000
  bmx=1.25*re;
1001

1002
  // Angular step of integration
1003
  h=4.e0f*PI/(float)TrackPoints;
1004

1005
  // impact parameter
1006
  b=(float)bi/(float)bmaxi*bmx;
1007

1008
  float up,vp,pp,us,vs,ps;
1009

1010
  up=0.;
1011
  vp=1.;
1012
      
1013
  pp=0.;
1014
  nh=0;
1015

1016
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1017
  
1018
  // b versus us
1019
  float bvus=fabs(b/us);
1020
  float bvus0=bvus;
1021
  Trajectories[bi*TrackPoints+nh]=bvus;
1022

1023
  do
1024
  {
1025
     nh++;
1026
     pp=ps;
1027
     up=us;
1028
     vp=vs;
1029
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1030
     bvus=fabs(b/us);
1031
     Trajectories[bi*TrackPoints+nh]=bvus;
1032

1033
  } while ((bvus>=rs)&&(bvus<=bvus0));
1034

1035
  IdLast[bi]=nh;
1036

1037
}
1038
"""
1039
    return(BlobCUDA)
1040
    
1041
# def ImageOutput(sigma,prefix):
1042
#     from PIL import Image
1043
#     Max=sigma.max()
1044
#     Min=sigma.min()
1045
    
1046
#     # Normalize value as 8bits Integer
1047
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1048
#     image = Image.fromarray(SigmaInt)
1049
#     image.save("%s.jpg" % prefix)
1050

    
1051
def ImageOutput(sigma,prefix,Colors):
1052
    import matplotlib.pyplot as plt
1053
    start_time=time.time() 
1054
    if Colors == 'Red2Yellow':
1055
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1056
    else:
1057
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1058
    save_time = time.time()-start_time
1059
    print("Save Time : %f" % save_time)
1060

    
1061
def BlackHoleCL(zImage,fImage,InputCL):
1062

    
1063
    kernel_params = {}
1064

    
1065
    print(InputCL)
1066

    
1067
    Device=InputCL['Device']
1068
    GpuStyle=InputCL['GpuStyle']
1069
    VariableType=InputCL['VariableType']
1070
    Size=InputCL['Size']
1071
    Mass=InputCL['Mass']
1072
    InternalRadius=InputCL['InternalRadius']
1073
    ExternalRadius=InputCL['ExternalRadius']
1074
    Angle=InputCL['Angle']
1075
    Method=InputCL['Method']
1076
    TrackPoints=InputCL['TrackPoints']
1077
    Physics=InputCL['Physics']
1078

    
1079
    PhysicsList=DictionariesAPI()
1080
    
1081
    if InputCL['BlackBody']:
1082
        # Spectrum is Black Body one
1083
        Line=0
1084
    else:
1085
        # Spectrum is Monochromatic Line one
1086
        Line=1
1087

    
1088
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1089
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1090

    
1091
    # Je detecte un peripherique GPU dans la liste des peripheriques
1092
    Id=0
1093
    HasXPU=False
1094
    for platform in cl.get_platforms():
1095
        for device in platform.get_devices():
1096
            if Id==Device:
1097
                XPU=device
1098
                print("CPU/GPU selected: ",device.name.lstrip())
1099
                HasXPU=True
1100
            Id+=1
1101

    
1102
    if HasXPU==False:
1103
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1104
        sys.exit()                
1105
        
1106
    ctx = cl.Context([XPU])
1107
    queue = cl.CommandQueue(ctx,
1108
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1109

    
1110
    
1111
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1112

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

    
1115
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1116
    
1117
    # Je recupere les flag possibles pour les buffers
1118
    mf = cl.mem_flags
1119

    
1120
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
1121
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1122
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1123
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1124
    
1125
    start_time=time.time() 
1126

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

    
1189
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1190
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1191
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1192
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1193
    zImageCL.release()
1194
    fImageCL.release()
1195

    
1196
    TrajectoriesCL.release()
1197
    IdLastCL.release()
1198

    
1199
    return(elapsed)
1200

    
1201
def BlackHoleCUDA(zImage,fImage,InputCL):
1202

    
1203
    kernel_params = {}
1204

    
1205
    print(InputCL)
1206

    
1207
    Device=InputCL['Device']
1208
    GpuStyle=InputCL['GpuStyle']
1209
    VariableType=InputCL['VariableType']
1210
    Size=InputCL['Size']
1211
    Mass=InputCL['Mass']
1212
    InternalRadius=InputCL['InternalRadius']
1213
    ExternalRadius=InputCL['ExternalRadius']
1214
    Angle=InputCL['Angle']
1215
    Method=InputCL['Method']
1216
    TrackPoints=InputCL['TrackPoints']
1217
    Physics=InputCL['Physics']
1218

    
1219
    PhysicsList=DictionariesAPI()
1220
    
1221
    if InputCL['BlackBody']:
1222
        # Spectrum is Black Body one
1223
        Line=0
1224
    else:
1225
        # Spectrum is Monochromatic Line one
1226
        Line=1
1227

    
1228
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1229
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1230

    
1231
    try:
1232
        # For PyCUDA import
1233
        import pycuda.driver as cuda
1234
        from pycuda.compiler import SourceModule
1235
        
1236
        cuda.init()
1237
        for Id in range(cuda.Device.count()):
1238
            if Id==Device:
1239
                XPU=cuda.Device(Id)
1240
                print("GPU selected %s" % XPU.name())
1241
        print
1242

    
1243
    except ImportError:
1244
        print("Platform does not seem to support CUDA")
1245

    
1246
    Context=XPU.make_context()
1247

    
1248
    try:
1249
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i' % (PhysicsList[Physics])])
1250
        print("Compilation seems to be OK")
1251
    except:
1252
        print("Compilation seems to break")
1253

    
1254
    EachPixelCU=mod.get_function("EachPixel")
1255
    TrajectoryCU=mod.get_function("Trajectory")
1256
    PixelCU=mod.get_function("Pixel")
1257
    CircleCU=mod.get_function("Circle")
1258
    
1259
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1260
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1261
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1262
    cuda.memcpy_htod(zImageCU, zImage)
1263
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1264
    cuda.memcpy_htod(zImageCU, fImage)
1265
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1266
    cuda.memcpy_htod(IdLastCU, IdLast)
1267

    
1268
    start_time=time.time() 
1269

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

    
1321
    Context.synchronize()
1322

    
1323
    compute = time.time()-start_time
1324

    
1325
    cuda.memcpy_dtoh(zImage,zImageCU)
1326
    cuda.memcpy_dtoh(fImage,fImageCU)
1327
    elapsed = time.time()-start_time
1328
    print("\nCompute Time : %f" % compute)
1329
    print("Elapsed Time : %f\n" % elapsed)
1330

    
1331
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1332
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1333
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1334
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1335

    
1336
    
1337
    Context.pop()
1338

    
1339
    Context.detach()
1340

    
1341
    return(elapsed)
1342

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

    
1373
    try:
1374
        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="])
1375
    except getopt.GetoptError:
1376
        print(HowToUse % sys.argv[0])
1377
        sys.exit(2)
1378

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

    
1444
    print("Device Identification selected : %s" % Device)
1445
    print("GpuStyle used : %s" % GpuStyle)
1446
    print("VariableType : %s" % VariableType)
1447
    print("Size : %i" % Size)
1448
    print("Mass : %f" % Mass)
1449
    print("Internal Radius : %f" % InternalRadius)
1450
    print("External Radius : %f" % ExternalRadius)
1451
    print("Angle with normal of (in radians) : %f" % Angle)
1452
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
1453
    print("Method of resolution : %s" % Method)
1454
    print("Colors for output images : %s" % Colors)
1455
    print("Physics used for Trajectories : %s" % Physics)
1456

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

    
1473
    if GpuStyle=='OpenCL':
1474
        print("\nSelection of OpenCL device") 
1475
        try:
1476
            # For PyOpenCL import
1477
            import pyopencl as cl
1478
            Id=0
1479
            for platform in cl.get_platforms():
1480
                for device in platform.get_devices():
1481
                    #deviceType=cl.device_type.to_string(device.type)
1482
                    deviceType="xPU"
1483
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
1484

    
1485
                    if Id in Devices:
1486
                    # Set the Alu as detected Device Type
1487
                        Alu[Id]=deviceType
1488
                    Id=Id+1
1489
        except ImportError:
1490
            print("Platform does not seem to support OpenCL")
1491

    
1492
#    print(Devices,Alu)
1493

    
1494
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
1495
    TrackPoints=2048
1496
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1497
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1498

    
1499
    InputCL={}
1500
    InputCL['Device']=Device
1501
    InputCL['GpuStyle']=GpuStyle
1502
    InputCL['VariableType']=VariableType
1503
    InputCL['Size']=Size
1504
    InputCL['Mass']=Mass
1505
    InputCL['InternalRadius']=InternalRadius
1506
    InputCL['ExternalRadius']=ExternalRadius
1507
    InputCL['Angle']=Angle
1508
    InputCL['BlackBody']=BlackBody
1509
    InputCL['Method']=Method
1510
    InputCL['TrackPoints']=TrackPoints
1511
    InputCL['Physics']=Physics
1512

    
1513
    if GpuStyle=='OpenCL':
1514
        duration=BlackHoleCL(zImage,fImage,InputCL)
1515
    else:
1516
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
1517
        
1518
    Hostname=gethostname()
1519
    Date=time.strftime("%Y%m%d_%H%M%S")
1520
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1521
    
1522
    
1523
    if not NoImage:
1524
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
1525
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)