Statistiques
| Révision :

root / TrouNoir / TrouNoir.py

Historique | Voir | Annoter | Télécharger (53,99 ko)

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

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

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

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

    
48
#
49
# Blank space below to simplify debugging on OpenCL code
50
#
51

    
52

    
53

    
54

    
55

    
56

    
57

    
58

    
59

    
60

    
61

    
62

    
63

    
64

    
65

    
66

    
67

    
68

    
69

    
70

    
71

    
72

    
73

    
74

    
75

    
76

    
77

    
78

    
79

    
80

    
81

    
82

    
83

    
84

    
85

    
86

    
87

    
88

    
89

    
90

    
91

    
92

    
93

    
94

    
95

    
96

    
97

    
98

    
99

    
100

    
101

    
102

    
103

    
104

    
105
BlobOpenCL= """
106

107
#define PI (float)3.14159265359e0f
108
#define nbr 256
109

110
#define EINSTEIN 0
111
#define NEWTON 1
112

113
#ifdef SETTRACKPOINTS
114
#define TRACKPOINTS SETTRACKPOINTS
115
#else
116
#define TRACKPOINTS 2048
117
#endif
118

119
float atanp(float x,float y)
120
{
121
  float angle;
122

123
  angle=atan2(y,x);
124

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

130
  return angle;
131
}
132

133
float f(float v)
134
{
135
  return v;
136
}
137

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

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

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

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

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

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

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

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

192
float spectre_cn(float rf32,float b32,float db32,
193
                 float h32,float r32,float m32,float bss32)
194
{
195

196
#define MYFLOAT float
197

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

206
  MYFLOAT flx;
207
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
208
  int fi,posfreq;
209

210
#define planck 6.62e-34f
211
#define k 1.38e-23f
212
#define c2 9.e16f
213
#define temp 3.e7f
214
#define m_point 1.e0f
215

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

220
  MYFLOAT v=1.e0f-3.e0f/r;
221

222
  qu=1.e0f/sqrt((1.e0f-3.e0f/r)*r)*(sqrt(r)-sqrt(6.e0f)+sqrt(3.e0f)/2.e0f*log((sqrt(r)+sqrt(3.e0f))/(sqrt(r)-sqrt(3.e0f))* 0.17157287525380988e0f ));
223

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

226
  flux_int=0.e0f;
227
  flx=0.e0f;
228

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

244
          flx+=flux_int;
245
        }
246
    }
247

248
  return((float)(flx));
249
}
250

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

258
  rf=decalage_spectral(r,b,phi,tho,m);
259

260
  if (raie==0)
261
    {
262
      bss=1.e19f;
263
      flx=spectre_cn(rf,b,db,h,r,m,bss);
264
    }
265
  else
266
    {
267
      bss=2.e0f;
268
      flx=spectre(rf,q,b,db,h,r,m,bss);
269
    }
270

271
  *zp=1.e0f/rf;
272
  *fp=flx;
273

274
}
275

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

286
   // Perform trajectory for each pixel, exit on hit
287

288
  float m,rs,ri,re,tho;
289
  int q,raie;
290

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

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

305
  // Autosize for image
306
  bmx=1.25e0f*re;
307

308
  h=4.e0f*PI/(float)TRACKPOINTS;
309

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

317

318
  float up,vp,pp,us,vs,ps;
319

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

325
  up=0.e0f;
326
  vp=1.e0f;
327
  pp=0.e0f;
328

329
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
330

331
  rps=fabs(b/us);
332
  rp0=rps;
333

334
  int ExitOnImpact=0;
335

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

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

348

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

358
  barrier(CLK_GLOBAL_MEM_FENCE);
359

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

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

376
   // Perform trajectory for each pixel
377

378
  float m,ri,re,tho;
379
  int q,raie;
380

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

388
  float bmx,db,b,h;
389
  float phi,phd,php,nr,r;
390
  float zp=0.e0f,fp=0.e0f;
391

392
  // Autosize for image, 25% greater than external radius
393
  bmx=1.25e0f*re;
394

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

398
  // Step of Impact Parameter
399
  db=bmx/(2.e0f*(float)ImpactParameter);
400

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

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

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

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

415
  int HalfLap=0,ExitOnImpact=0,ni;
416

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

425
      if (ni<IdLast[bi])
426
      {
427
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.e0f)+Trajectories[bi*TRACKPOINTS+ni];
428
      }
429
      else
430
      {
431
        r=Trajectories[bi*TRACKPOINTS+ni];
432
      }
433

434
      if ((r<=re)&&(r>=ri))
435
      {
436
        ExitOnImpact=1;
437
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
438
      }
439

440
      HalfLap++;
441
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
442

443
  }
444

445
  barrier(CLK_GLOBAL_MEM_FENCE);
446

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

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

466
   // Perform trajectory for each pixel
467

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

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

477
  float bmx,db,b,h;
478
  float phi,phd;
479
  float zp=0.e0f,fp=0.e0f;
480

481
  // Autosize for image
482
  bmx=1.25e0f*re;
483

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

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

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

496
  int HalfLap=0,ExitOnImpact=0,ni;
497
  float php,nr,r;
498

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

505
     if (ni<IdLast[bi])
506
     {
507
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.e0f)+Trajectories[bi*TRACKPOINTS+ni];
508
     }
509
     else
510
     {
511
        r=Trajectories[bi*TRACKPOINTS+ni];
512
     }
513

514
     if ((r<=re)&&(r>=ri))
515
     {
516
        ExitOnImpact=1;
517
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
518
     }
519

520
     HalfLap++;
521
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
522

523
  zImage[yi+2*bmaxi*xi]=zp;
524
  fImage[yi+2*bmaxi*xi]=fp;
525

526
  barrier(CLK_GLOBAL_MEM_FENCE);
527

528
}
529

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

540
  // Perform trajectory for each pixel
541

542
  float m,rs,re;
543

544
  m=Mass;
545
  rs=2.e0f*m;
546
  re=ExternalRadius;
547

548
  float bmx,b,h;
549
  int nh;
550

551
  // Autosize for image
552
  bmx=1.25e0f*re;
553

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

557
  // impact parameter
558
  b=(float)bi/(float)bmaxi*bmx;
559

560
  float up,vp,pp,us,vs,ps;
561

562
  up=0.e0f;
563
  vp=1.e0f;
564

565
  pp=0.e0f;
566
  nh=0;
567

568
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
569

570
  // b versus us
571
  float bvus=fabs(b/us);
572
  float bvus0=bvus;
573
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
574

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

585
  } while ((bvus>=rs)&&(bvus<=bvus0));
586

587
  IdLast[bi]=nh;
588

589
  barrier(CLK_GLOBAL_MEM_FENCE);
590

591
}
592

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

603
  private float Trajectory[TRACKPOINTS];
604

605
  float m,rs,ri,re,tho;
606
  int raie,q;
607

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

616
  float bmx,db,b,h;
617
  uint nh;
618

619

620
  // Autosize for image
621
  bmx=1.25e0f*re;
622

623
  // Angular step of integration
624
  h=4.e0f*PI/(float)TRACKPOINTS;
625

626
  // impact parameter
627
  b=(float)bi/(float)bmaxi*bmx;
628
  db=bmx/(2.e0f*(float)bmaxi);
629

630
  float up,vp,pp,us,vs,ps;
631

632
  up=0.e0f;
633
  vp=1.e0f;
634

635
  pp=0.e0f;
636
  nh=0;
637

638
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
639

640
  // b versus us
641
  float bvus=fabs(b/us);
642
  float bvus0=bvus;
643
  Trajectory[nh]=bvus;
644

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

655
  } while ((bvus>=rs)&&(bvus<=bvus0));
656

657

658
  for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
659
     Trajectory[i]=0.e0f;
660
  }
661

662

663
  uint imx=(uint)(16*bi);
664

665
  for (uint i=0;i<imx;i++)
666
  {
667
     float zp=0.e0f,fp=0.e0f;
668
     float phi=2.e0f*PI/(float)imx*(float)i;
669
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
670
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
671
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
672

673
     uint HalfLap=0,ExitOnImpact=0,ni;
674
     float php,nr,r;
675

676
     do
677
     {
678
        php=phd+(float)HalfLap*PI;
679
        nr=php/h;
680
        ni=(int)nr;
681

682
        if (ni<nh)
683
        {
684
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.e0f)+Trajectory[ni];
685
        }
686
        else
687
        {
688
           r=Trajectory[ni];
689
        }
690

691
        if ((r<=re)&&(r>=ri))
692
        {
693
           ExitOnImpact=1;
694
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
695
        }
696

697
        HalfLap++;
698

699
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
700

701
     zImage[yi+2*bmaxi*xi]=zp;
702
     fImage[yi+2*bmaxi*xi]=fp;
703

704
  }
705

706
  barrier(CLK_GLOBAL_MEM_FENCE);
707

708
}
709

710
__kernel void Original(__global float *zImage,__global float *fImage,
711
                       uint Size,float Mass,float InternalRadius,
712
                       float ExternalRadius,float Angle,
713
                       int Line)
714
{
715
   // Integer Impact Parameter Size (half of image)
716
   uint bmaxi=(uint)Size;
717

718
   float Trajectory[TRACKPOINTS];
719

720
   // Perform trajectory for each pixel
721

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

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

733
   float bmx,db,b,h;
734
   uint nh;
735

736
   // Autosize for image
737
   bmx=1.25e0f*re;
738

739
   // Angular step of integration
740
   h=4.e0f*PI/(float)TRACKPOINTS;
741

742
   // Integer Impact Parameter ID
743
   for (int bi=0;bi<bmaxi;bi++)
744
   {
745
      // impact parameter
746
      b=(float)bi/(float)bmaxi*bmx;
747
      db=bmx/(2.e0f*(float)bmaxi);
748

749
      float up,vp,pp,us,vs,ps;
750

751
      up=0.e0f;
752
      vp=1.e0f;
753

754
      pp=0.e0f;
755
      nh=0;
756

757
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
758

759
      // b versus us
760
      float bvus=fabs(b/us);
761
      float bvus0=bvus;
762
      Trajectory[nh]=bvus;
763

764
      do
765
      {
766
         nh++;
767
         pp=ps;
768
         up=us;
769
         vp=vs;
770
         rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
771
         bvus=fabs(b/us);
772
         Trajectory[nh]=bvus;
773

774
      } while ((bvus>=rs)&&(bvus<=bvus0));
775

776
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
777
         Trajectory[i]=0.e0f;
778
      }
779

780
      int imx=(int)(16*bi);
781

782
      for (int i=0;i<imx;i++)
783
      {
784
         float zp=0.e0f,fp=0.e0f;
785
         float phi=2.e0f*PI/(float)imx*(float)i;
786
         float phd=atanp(cos(phi)*sin(tho),cos(tho));
787
         uint yi=(uint)((float)bi*sin(phi)+bmaxi);
788
         uint xi=(uint)((float)bi*cos(phi)+bmaxi);
789

790
         uint HalfLap=0,ExitOnImpact=0,ni;
791
         float php,nr,r;
792

793
         do
794
         {
795
            php=phd+(float)HalfLap*PI;
796
            nr=php/h;
797
            ni=(int)nr;
798

799
            if (ni<nh)
800
            {
801
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.e0f)+Trajectory[ni];
802
            }
803
            else
804
            {
805
               r=Trajectory[ni];
806
            }
807

808
            if ((r<=re)&&(r>=ri))
809
            {
810
               ExitOnImpact=1;
811
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
812
            }
813

814
            HalfLap++;
815

816
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
817

818
         zImage[yi+2*bmaxi*xi]=zp;
819
         fImage[yi+2*bmaxi*xi]=fp;
820

821
      }
822

823
   }
824

825
   barrier(CLK_GLOBAL_MEM_FENCE);
826

827
}
828
"""
829

    
830

    
831

    
832

    
833

    
834

    
835

    
836

    
837

    
838

    
839

    
840

    
841

    
842

    
843

    
844

    
845

    
846

    
847

    
848

    
849

    
850

    
851

    
852

    
853

    
854

    
855

    
856

    
857

    
858

    
859

    
860

    
861

    
862

    
863

    
864

    
865

    
866

    
867

    
868

    
869

    
870

    
871

    
872

    
873

    
874

    
875

    
876

    
877

    
878

    
879

    
880

    
881

    
882

    
883

    
884

    
885

    
886

    
887

    
888

    
889

    
890

    
891

    
892

    
893

    
894

    
895

    
896

    
897

    
898

    
899

    
900

    
901
def KernelCodeCuda():
902
    BlobCUDA= """
903

904
#define PI (float)3.14159265359
905
#define nbr 256
906

907
#define EINSTEIN 0
908
#define NEWTON 1
909

910
#ifdef SETTRACKPOINTS
911
#define TRACKPOINTS SETTRACKPOINTS
912
#else
913
#define TRACKPOINTS
914
#endif
915
__device__ float nothing(float x)
916
{
917
  return(x);
918
}
919

920
__device__ float atanp(float x,float y)
921
{
922
  float angle;
923

924
  angle=atan2(y,x);
925

926
  if (angle<0.e0f)
927
    {
928
      angle+=(float)2.e0f*PI;
929
    }
930

931
  return(angle);
932
}
933

934
__device__ float f(float v)
935
{
936
  return(v);
937
}
938

939
#if PHYSICS == NEWTON
940
__device__ float g(float u,float m,float b)
941
{
942
  return (-u);
943
}
944
#else
945
__device__ float g(float u,float m,float b)
946
{
947
  return (3.e0f*m/b*pow(u,2)-u);
948
}
949
#endif
950

951
__device__ void calcul(float *us,float *vs,float up,float vp,
952
                       float h,float m,float b)
953
{
954
  float c0,c1,c2,c3,d0,d1,d2,d3;
955

956
  c0=h*f(vp);
957
  c1=h*f(vp+c0/2.);
958
  c2=h*f(vp+c1/2.);
959
  c3=h*f(vp+c2);
960
  d0=h*g(up,m,b);
961
  d1=h*g(up+d0/2.,m,b);
962
  d2=h*g(up+d1/2.,m,b);
963
  d3=h*g(up+d2,m,b);
964

965
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
966
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
967
}
968

969
__device__ void rungekutta(float *ps,float *us,float *vs,
970
                           float pp,float up,float vp,
971
                           float h,float m,float b)
972
{
973
  calcul(us,vs,up,vp,h,m,b);
974
  *ps=pp+h;
975
}
976

977
__device__ float decalage_spectral(float r,float b,float phi,
978
                                   float tho,float m)
979
{
980
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
981
}
982

983
__device__ float spectre(float rf,int q,float b,float db,
984
                         float h,float r,float m,float bss)
985
{
986
  float flx;
987

988
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
989
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
990
  return(flx);
991
}
992

993
__device__ float spectre_cn(float rf32,float b32,float db32,
994
                            float h32,float r32,float m32,float bss32)
995
{
996

997
#define MYFLOAT float
998

999
  MYFLOAT rf=(MYFLOAT)(rf32);
1000
  MYFLOAT b=(MYFLOAT)(b32);
1001
  MYFLOAT db=(MYFLOAT)(db32);
1002
  MYFLOAT h=(MYFLOAT)(h32);
1003
  MYFLOAT r=(MYFLOAT)(r32);
1004
  MYFLOAT m=(MYFLOAT)(m32);
1005
  MYFLOAT bss=(MYFLOAT)(bss32);
1006

1007
  MYFLOAT flx;
1008
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
1009
  int fi,posfreq;
1010

1011
#define planck 6.62e-34
1012
#define k 1.38e-23
1013
#define c2 9.e16
1014
#define temp 3.e7
1015
#define m_point 1.
1016

1017
#define lplanck (log(6.62)-34.*log(10.))
1018
#define lk (log(1.38)-23.*log(10.))
1019
#define lc2 (log(9.)+16.*log(10.))
1020

1021
  MYFLOAT v=1.-3./r;
1022

1023
  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 ));
1024

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

1027
  flux_int=0.;
1028
  flx=0.;
1029

1030
  for (fi=0;fi<nbr;fi++)
1031
    {
1032
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
1033
      nu_rec=nu_em*rf;
1034
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
1035
      if ((posfreq>0)&&(posfreq<nbr))
1036
        {
1037
          // Initial version
1038
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
1039
          // Version with log used
1040
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
1041
          // flux_int*=pow(rf,3)*b*db*h;
1042
          //flux_int*=exp(3.*log(rf))*b*db*h;
1043
          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;
1044

1045
          flx+=flux_int;
1046
        }
1047
    }
1048

1049
  return((float)(flx));
1050
}
1051

1052
__device__ void impact(float phi,float r,float b,float tho,float m,
1053
                       float *zp,float *fp,
1054
                       int q,float db,
1055
                       float h,int raie)
1056
{
1057
  float flx,rf,bss;
1058

1059
  rf=decalage_spectral(r,b,phi,tho,m);
1060

1061
  if (raie==0)
1062
    {
1063
      bss=1.e19;
1064
      flx=spectre_cn(rf,b,db,h,r,m,bss);
1065
    }
1066
  else
1067
    {
1068
      bss=2.;
1069
      flx=spectre(rf,q,b,db,h,r,m,bss);
1070
    }
1071

1072
  *zp=1./rf;
1073
  *fp=flx;
1074

1075
}
1076

1077
__global__ void EachPixel(float *zImage,float *fImage,
1078
                          float Mass,float InternalRadius,
1079
                          float ExternalRadius,float Angle,
1080
                          int Line)
1081
{
1082
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
1083
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
1084
   uint sizex=(uint)gridDim.x*blockDim.x;
1085
   uint sizey=(uint)gridDim.y*blockDim.y;
1086

1087

1088
   // Perform trajectory for each pixel, exit on hit
1089

1090
  float m,rs,ri,re,tho;
1091
  int q,raie;
1092

1093
  m=Mass;
1094
  rs=2.*m;
1095
  ri=InternalRadius;
1096
  re=ExternalRadius;
1097
  tho=Angle;
1098
  q=-2;
1099
  raie=Line;
1100

1101
  float bmx,db,b,h;
1102
  float rp0,rpp,rps;
1103
  float phi,phd;
1104
  int nh;
1105
  float zp,fp;
1106

1107
  // Autosize for image
1108
  bmx=1.25*re;
1109
  b=0.;
1110

1111
  h=4.e0f*PI/(float)TRACKPOINTS;
1112

1113
  // set origin as center of image
1114
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
1115
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
1116
  // angle extracted from cylindric symmetry
1117
  phi=atanp(x,y);
1118
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1119

1120
  float up,vp,pp,us,vs,ps;
1121

1122
  // impact parameter
1123
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
1124
  // step of impact parameter;
1125
//  db=bmx/(float)(sizex/2);
1126
  db=bmx/(float)(sizex);
1127

1128
  up=0.;
1129
  vp=1.;
1130
  pp=0.;
1131
  nh=0;
1132

1133
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1134

1135
  rps=fabs(b/us);
1136
  rp0=rps;
1137

1138
  int ExitOnImpact=0;
1139

1140
  do
1141
  {
1142
     nh++;
1143
     pp=ps;
1144
     up=us;
1145
     vp=vs;
1146
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1147
     rpp=rps;
1148
     rps=fabs(b/us);
1149
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;          
1150

1151
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
1152

1153
  if (ExitOnImpact==1) {
1154
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
1155
  }
1156
  else
1157
  {
1158
     zp=0.e0f;
1159
     fp=0.e0f;
1160
  }
1161

1162
  __syncthreads();
1163

1164
  zImage[yi+sizex*xi]=(float)zp;
1165
  fImage[yi+sizex*xi]=(float)fp;
1166
}
1167

1168
__global__ void Pixel(float *zImage,float *fImage,
1169
                      float *Trajectories,int *IdLast,
1170
                      uint ImpactParameter,
1171
                      float Mass,float InternalRadius,
1172
                      float ExternalRadius,float Angle,
1173
                      int Line)
1174
{
1175
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
1176
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
1177
   uint sizex=(uint)gridDim.x*blockDim.x;
1178
   uint sizey=(uint)gridDim.y*blockDim.y;
1179

1180
  // Perform trajectory for each pixel
1181

1182
  float m,ri,re,tho;
1183
  int q,raie;
1184

1185
  m=Mass;
1186
  ri=InternalRadius;
1187
  re=ExternalRadius;
1188
  tho=Angle;
1189
  q=-2;
1190
  raie=Line;
1191

1192
  float bmx,db,b,h;
1193
  float phi,phd,php,nr,r;
1194
  float zp=0,fp=0;
1195
  // Autosize for image, 25% greater than external radius
1196
  bmx=1.25e0f*re;
1197

1198
  // Angular step of integration
1199
  h=4.e0f*PI/(float)TRACKPOINTS;
1200

1201
  // Step of Impact Parameter
1202
  db=bmx/(2.e0f*(float)ImpactParameter);
1203

1204
  // set origin as center of image
1205
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
1206
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
1207
  // angle extracted from cylindric symmetry
1208
  phi=atanp(x,y);
1209
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1210

1211
  // Real Impact Parameter
1212
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
1213

1214
  // Integer Impact Parameter
1215
  uint bi=(uint)sqrt(x*x+y*y);
1216

1217
  int HalfLap=0,ExitOnImpact=0,ni;
1218

1219
  if (bi<ImpactParameter)
1220
  {
1221
    do
1222
    {
1223
      php=phd+(float)HalfLap*PI;
1224
      nr=php/h;
1225
      ni=(int)nr;
1226

1227
      if (ni<IdLast[bi])
1228
      {
1229
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.e0f)+Trajectories[bi*TRACKPOINTS+ni];
1230
      }
1231
      else
1232
      {
1233
        r=Trajectories[bi*TRACKPOINTS+ni];
1234
      }
1235

1236
      if ((r<=re)&&(r>=ri))
1237
      {
1238
        ExitOnImpact=1;
1239
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1240
      }
1241

1242
      HalfLap++;
1243
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1244

1245
  }
1246

1247
  zImage[yi+sizex*xi]=zp;
1248
  fImage[yi+sizex*xi]=fp;
1249
}
1250

1251
__global__ void Circle(float *Trajectories,int *IdLast,
1252
                       float *zImage,float *fImage,
1253
                       float Mass,float InternalRadius,
1254
                       float ExternalRadius,float Angle,
1255
                       int Line)
1256
{
1257
   // Integer Impact Parameter ID
1258
   int bi=blockIdx.x*blockDim.x+threadIdx.x;
1259
   // Integer points on circle
1260
   int i=blockIdx.y*blockDim.y+threadIdx.y;
1261
   // Integer Impact Parameter Size (half of image)
1262
   int bmaxi=gridDim.x*blockDim.x;
1263
   // Integer Points on circle
1264
   int imx=gridDim.y*blockDim.y;
1265

1266
   // Perform trajectory for each pixel
1267

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

1271
  m=Mass;
1272
  ri=InternalRadius;
1273
  re=ExternalRadius;
1274
  tho=Angle;
1275
  raie=Line;
1276

1277
  float bmx,db,b,h;
1278
  float phi,phd;
1279
  float zp=0,fp=0;
1280

1281
  // Autosize for image
1282
  bmx=1.25e0f*re;
1283

1284
  // Angular step of integration
1285
  h=4.e0f*PI/(float)TRACKPOINTS;
1286

1287
  // impact parameter
1288
  b=(float)bi/(float)bmaxi*bmx;
1289
  db=bmx/(2.e0f*(float)bmaxi);
1290

1291
  phi=2.e0f*PI/(float)imx*(float)i;
1292
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1293
  int yi=(int)((float)bi*sin(phi))+bmaxi;
1294
  int xi=(int)((float)bi*cos(phi))+bmaxi;
1295

1296
  int HalfLap=0,ExitOnImpact=0,ni;
1297
  float php,nr,r;
1298

1299
  do
1300
  {
1301
     php=phd+(float)HalfLap*PI;
1302
     nr=php/h;
1303
     ni=(int)nr;
1304

1305
     if (ni<IdLast[bi])
1306
     {
1307
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.e0f)+Trajectories[bi*TRACKPOINTS+ni];
1308
     }
1309
     else
1310
     {
1311
        r=Trajectories[bi*TRACKPOINTS+ni];
1312
     }
1313

1314
     if ((r<=re)&&(r>=ri))
1315
     {
1316
        ExitOnImpact=1;
1317
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1318
     }
1319

1320
     HalfLap++;
1321
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1322

1323
  zImage[yi+2*bmaxi*xi]=zp;
1324
  fImage[yi+2*bmaxi*xi]=fp;
1325

1326
}
1327

1328
__global__ void Trajectory(float *Trajectories,int *IdLast,
1329
                           float Mass,float InternalRadius,
1330
                           float ExternalRadius,float Angle,
1331
                           int Line)
1332
{
1333
  // Integer Impact Parameter ID
1334
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1335
  // Integer Impact Parameter Size (half of image)
1336
  int bmaxi=gridDim.x*blockDim.x;
1337

1338
  // Perform trajectory for each pixel
1339

1340
  float m,rs,re;
1341

1342
  m=Mass;
1343
  rs=2.e0f*m;
1344
  re=ExternalRadius;
1345

1346
  float bmx,b,h;
1347
  int nh;
1348

1349
  // Autosize for image
1350
  bmx=1.25e0f*re;
1351

1352
  // Angular step of integration
1353
  h=4.e0f*PI/(float)TRACKPOINTS;
1354

1355
  // impact parameter
1356
  b=(float)bi/(float)bmaxi*bmx;
1357

1358
  float up,vp,pp,us,vs,ps;
1359

1360
  up=0.e0f;
1361
  vp=1.e0f;
1362
  pp=0.e0f;
1363
  nh=0;
1364

1365
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1366

1367
  // b versus us
1368
  float bvus=fabs(b/us);
1369
  float bvus0=bvus;
1370
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
1371

1372
  do
1373
  {
1374
     nh++;
1375
     pp=ps;
1376
     up=us;
1377
     vp=vs;
1378
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1379
     bvus=fabs(b/us);
1380
     Trajectories[bi*TRACKPOINTS+nh]=bvus;
1381

1382
  } while ((bvus>=rs)&&(bvus<=bvus0));
1383

1384
  IdLast[bi]=nh;
1385

1386
}
1387

1388
__global__ void EachCircle(float *zImage,float *fImage,
1389
                           float Mass,float InternalRadius,
1390
                           float ExternalRadius,float Angle,
1391
                           int Line)
1392
{
1393
  // Integer Impact Parameter ID
1394
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1395

1396
  // Integer Impact Parameter Size (half of image)
1397
  int bmaxi=gridDim.x*blockDim.x;
1398

1399
  float Trajectory[TRACKPOINTS];
1400

1401
  // Perform trajectory for each pixel
1402

1403
  float m,rs,ri,re,tho;
1404
  int raie,q;
1405

1406
  m=Mass;
1407
  rs=2.*m;
1408
  ri=InternalRadius;
1409
  re=ExternalRadius;
1410
  tho=Angle;
1411
  q=-2;
1412
  raie=Line;
1413

1414
  float bmx,db,b,h;
1415
  int nh;
1416

1417
  // Autosize for image
1418
  bmx=1.25e0f*re;
1419

1420
  // Angular step of integration
1421
  h=4.e0f*PI/(float)TRACKPOINTS;
1422

1423
  // impact parameter
1424
  b=(float)bi/(float)bmaxi*bmx;
1425
  db=bmx/(2.e0f*(float)bmaxi);
1426

1427
  float up,vp,pp,us,vs,ps;
1428

1429
  up=0.;
1430
  vp=1.;
1431
  pp=0.;
1432
  nh=0;
1433

1434
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1435

1436
  // b versus us
1437
  float bvus=fabs(b/us);
1438
  float bvus0=bvus;
1439
  Trajectory[nh]=bvus;
1440

1441
  do
1442
  {
1443
     nh++;
1444
     pp=ps;
1445
     up=us;
1446
     vp=vs;
1447
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1448
     bvus=fabs(b/us);
1449
     Trajectory[nh]=bvus;
1450

1451
  } while ((bvus>=rs)&&(bvus<=bvus0));
1452

1453
  for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
1454
     Trajectory[i]=0.e0f;
1455
  }
1456

1457
  int imx=(int)(16*bi);
1458

1459
  for (int i=0;i<imx;i++)
1460
  {
1461
     float zp=0,fp=0;
1462
     float phi=2.*PI/(float)imx*(float)i;
1463
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
1464
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1465
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1466

1467
     int HalfLap=0,ExitOnImpact=0,ni;
1468
     float php,nr,r;
1469

1470
     do
1471
     {
1472
        php=phd+(float)HalfLap*PI;
1473
        nr=php/h;
1474
        ni=(int)nr;
1475

1476
        if (ni<nh)
1477
        {
1478
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1479
        }
1480
        else
1481
        {
1482
           r=Trajectory[ni];
1483
        }
1484

1485
        if ((r<=re)&&(r>=ri))
1486
        {
1487
           ExitOnImpact=1;
1488
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1489
        }
1490

1491
        HalfLap++;
1492

1493
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
1494

1495
   __syncthreads();
1496

1497
   zImage[yi+2*bmaxi*xi]=zp;
1498
   fImage[yi+2*bmaxi*xi]=fp;
1499

1500
  }
1501

1502
}
1503

1504
__global__ void Original(float *zImage,float *fImage,
1505
                         uint Size,float Mass,float InternalRadius,
1506
                         float ExternalRadius,float Angle,
1507
                         int Line)
1508
{
1509
   // Integer Impact Parameter Size (half of image)
1510
   uint bmaxi=(uint)Size;
1511

1512
   float Trajectory[TRACKPOINTS];
1513

1514
   // Perform trajectory for each pixel
1515

1516
   float m,rs,ri,re,tho;
1517
   int raie,q;
1518

1519
   m=Mass;
1520
   rs=2.e0f*m;
1521
   ri=InternalRadius;
1522
   re=ExternalRadius;
1523
   tho=Angle;
1524
   q=-2;
1525
   raie=Line;
1526

1527
   float bmx,db,b,h;
1528
   int nh;
1529

1530
   // Autosize for image
1531
   bmx=1.25e0f*re;
1532

1533
   // Angular step of integration
1534
   h=4.e0f*PI/(float)TRACKPOINTS;
1535

1536
   // Integer Impact Parameter ID
1537
   for (int bi=0;bi<bmaxi;bi++)
1538
   {
1539
      // impact parameter
1540
      b=(float)bi/(float)bmaxi*bmx;
1541
      db=bmx/(2.e0f*(float)bmaxi);
1542

1543
      float up,vp,pp,us,vs,ps;
1544

1545
      up=0.;
1546
      vp=1.;
1547
      pp=0.;
1548
      nh=0;
1549

1550
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1551

1552
      // b versus us
1553
      float bvus=fabs(b/us);
1554
      float bvus0=bvus;
1555
      Trajectory[nh]=bvus;
1556

1557
      do
1558
      {
1559
         nh++;
1560
         pp=ps;
1561
         up=us;
1562
         vp=vs;
1563
         rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1564
         bvus=fabs(b/us);
1565
         Trajectory[nh]=bvus;
1566

1567
      } while ((bvus>=rs)&&(bvus<=bvus0));
1568

1569
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
1570
         Trajectory[i]=0.e0f;
1571
      }
1572

1573
      int imx=(int)(16*bi);
1574

1575
      for (int i=0;i<imx;i++)
1576
      {
1577
         float zp=0,fp=0;
1578
         float phi=2.e0f*PI/(float)imx*(float)i;
1579
         float phd=atanp(cos(phi)*sin(tho),cos(tho));
1580
         uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1581
         uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1582

1583
         int HalfLap=0,ExitOnImpact=0,ni;
1584
         float php,nr,r;
1585

1586
         do
1587
         {
1588
            php=phd+(float)HalfLap*PI;
1589
            nr=php/h;
1590
            ni=(int)nr;
1591

1592
            if (ni<nh)
1593
            {
1594
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1595
            }
1596
            else
1597
            {
1598
               r=Trajectory[ni];
1599
            }
1600

1601
            if ((r<=re)&&(r>=ri))
1602
            {
1603
               ExitOnImpact=1;
1604
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1605
            }
1606

1607
            HalfLap++;
1608

1609
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
1610

1611
         zImage[yi+2*bmaxi*xi]=zp;
1612
         fImage[yi+2*bmaxi*xi]=fp;
1613

1614
      }
1615

1616
   }
1617

1618
}
1619
"""
1620
    return(BlobCUDA)
1621

    
1622
# def ImageOutput(sigma,prefix):
1623
#     from PIL import Image
1624
#     Max=sigma.max()
1625
#     Min=sigma.min()
1626
#     # Normalize value as 8bits Integer
1627
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1628
#     image = Image.fromarray(SigmaInt)
1629
#     image.save("%s.jpg" % prefix)
1630

    
1631
def ImageOutput(sigma,prefix,Colors):
1632
    try:
1633
        import matplotlib.pyplot as plt
1634
        start_time=time.time()
1635
        if Colors == 'Red2Yellow':
1636
            plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1637
        else:
1638
            plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1639
        save_time = time.time()-start_time
1640
        print("Save image as %s.png file" % prefix)
1641
        print("Save Time : %f" % save_time)
1642
    except:
1643
        from PIL import Image
1644
        start_time=time.time()
1645
        Max = sigma.max()
1646
        Min = sigma.min()
1647
        # Normalize value as 8bits Integer
1648
        SigmaInt = (255 * (sigma - Min) / (Max - Min)).astype("uint8")
1649
        image = Image.fromarray(SigmaInt)
1650
        image.save("%s.jpg" % prefix)
1651
        save_time = time.time()-start_time
1652
        print("Save image as %s.png file" % prefix)
1653
        print("Save Time : %f" % save_time)
1654
            
1655
def BlackHoleCL(zImage,fImage,InputCL):
1656

    
1657
    kernel_params = {}
1658

    
1659
    Device=InputCL['Device']
1660
    GpuStyle=InputCL['GpuStyle']
1661
    VariableType=InputCL['VariableType']
1662
    Size=InputCL['Size']
1663
    Mass=InputCL['Mass']
1664
    InternalRadius=InputCL['InternalRadius']
1665
    ExternalRadius=InputCL['ExternalRadius']
1666
    Angle=InputCL['Angle']
1667
    Method=InputCL['Method']
1668
    TrackPoints=InputCL['TrackPoints']
1669
    Physics=InputCL['Physics']
1670
    NoImage=InputCL['NoImage']
1671
    TrackSave=InputCL['TrackSave']
1672

    
1673
    PhysicsList=DictionariesAPI()
1674
    
1675
    if InputCL['BlackBody']:
1676
        # Spectrum is Black Body one
1677
        Line=0
1678
    else:
1679
        # Spectrum is Monochromatic Line one
1680
        Line=1
1681

    
1682
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1683
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1684

    
1685
    # Je detecte un peripherique GPU dans la liste des peripheriques
1686
    Id=0
1687
    HasXPU=False
1688
    for platform in cl.get_platforms():
1689
        for device in platform.get_devices():
1690
            if Id==Device:
1691
                PF4XPU=platform.name
1692
                XPU=device
1693
                print("CPU/GPU selected: ",device.name.lstrip())
1694
                HasXPU=True
1695
            Id+=1
1696

    
1697
    if HasXPU==False:
1698
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1699
        sys.exit()
1700

    
1701
    ctx = cl.Context([XPU])
1702
    queue = cl.CommandQueue(ctx,
1703
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1704

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

    
1707
    print('My Platform is ',PF4XPU)
1708
    
1709
    if 'Intel' in PF4XPU or 'Experimental' in PF4XPU or 'Clover' in PF4XPU or 'Portable' in PF4XPU :
1710
        print('No extra options for Intel and Clover!')
1711
    else:
1712
        BuildOptions = BuildOptions+" -cl-mad-enable"
1713

    
1714
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1715
    
1716
    # Je recupere les flag possibles pour les buffers
1717
    mf = cl.mem_flags
1718

    
1719
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1720
    #     TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories) 
1721
    #     IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1722

    
1723
    # zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1724
    # fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1725

    
1726
        TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY, Trajectories.nbytes) 
1727
        IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY, IdLast.nbytes)
1728

    
1729
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY, zImage.nbytes)
1730
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY, fImage.nbytes)
1731
    
1732
    start_time=time.time()
1733

    
1734
    if Method=='EachPixel':
1735
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1736
                                       None,zImageCL,fImageCL,
1737
                                       numpy.float32(Mass),
1738
                                       numpy.float32(InternalRadius),
1739
                                       numpy.float32(ExternalRadius),
1740
                                       numpy.float32(Angle),
1741
                                       numpy.int32(Line))
1742
        CLLaunch.wait()
1743
    elif Method=='Original':
1744
        CLLaunch=BlackHoleCL.Original(queue,(1,),
1745
                                      None,zImageCL,fImageCL,
1746
                                      numpy.uint32(zImage.shape[0]/2),
1747
                                      numpy.float32(Mass),
1748
                                      numpy.float32(InternalRadius),
1749
                                      numpy.float32(ExternalRadius),
1750
                                      numpy.float32(Angle),
1751
                                      numpy.int32(Line))
1752
        CLLaunch.wait()
1753
    elif Method=='EachCircle':
1754
        CLLaunch=BlackHoleCL.EachCircle(queue,(int(zImage.shape[0]/2),),
1755
                                        None,zImageCL,fImageCL,
1756
                                        numpy.float32(Mass),
1757
                                        numpy.float32(InternalRadius),
1758
                                        numpy.float32(ExternalRadius),
1759
                                        numpy.float32(Angle),
1760
                                        numpy.int32(Line))
1761
        CLLaunch.wait()
1762
    elif Method=='TrajectoCircle':
1763
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1764
                                        None,TrajectoriesCL,IdLastCL,
1765
                                        numpy.float32(Mass),
1766
                                        numpy.float32(InternalRadius),
1767
                                        numpy.float32(ExternalRadius),
1768
                                        numpy.float32(Angle),
1769
                                        numpy.int32(Line))
1770
        
1771
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1772
                                           int(zImage.shape[0]*4)),None,
1773
                                    TrajectoriesCL,IdLastCL,
1774
                                    zImageCL,fImageCL,
1775
                                    numpy.float32(Mass),
1776
                                    numpy.float32(InternalRadius),
1777
                                    numpy.float32(ExternalRadius),
1778
                                    numpy.float32(Angle),
1779
                                    numpy.int32(Line))
1780
        CLLaunch.wait()
1781
    else:
1782
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1783
                                        None,TrajectoriesCL,IdLastCL,
1784
                                        numpy.float32(Mass),
1785
                                        numpy.float32(InternalRadius),
1786
                                        numpy.float32(ExternalRadius),
1787
                                        numpy.float32(Angle),
1788
                                        numpy.int32(Line))
1789
        
1790
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],zImage.shape[1]),None,
1791
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1792
                                   numpy.uint32(Trajectories.shape[0]),
1793
                                   numpy.float32(Mass),
1794
                                   numpy.float32(InternalRadius),
1795
                                   numpy.float32(ExternalRadius),
1796
                                   numpy.float32(Angle),
1797
                                   numpy.int32(Line))
1798
        CLLaunch.wait()
1799
    
1800
    compute = time.time()-start_time
1801
    
1802
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1803
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1804
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1805
        cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1806
        cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1807
    elapsed = time.time()-start_time
1808
    print("\nCompute Time : %f" % compute)
1809
    print("Elapsed Time : %f\n" % elapsed)
1810

    
1811
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1812
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1813
    print("Z max @(%f,%f) : %f" % ((1.*zMaxPosition[1][0]/zImage.shape[1]-0.5,1.*zMaxPosition[0][0]/zImage.shape[0]-0.5,zImage.max())))
1814
    print("Flux max @(%f,%f) : %f" % ((1.*fMaxPosition[1][0]/fImage.shape[1]-0.5,1.*fMaxPosition[0][0]/fImage.shape[0]-0.5,fImage.max())))
1815
    zImageCL.release()
1816
    fImageCL.release()
1817

    
1818
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1819
        if not NoImage:
1820
            AngleStep=4*numpy.pi/TrackPoints
1821
            Angles=numpy.arange(0.,4*numpy.pi,AngleStep)
1822
            Angles.shape=(1,TrackPoints)
1823
            Hostname=gethostname()
1824
            Date=time.strftime("%Y%m%d_%H%M%S")
1825
            ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1826
            
1827
            if TrackSave:
1828
                # numpy.savetxt("TrouNoirTrajectories_%s.csv" % ImageInfo,
1829
                #               numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1830
                #               delimiter=' ', fmt='%.2e')
1831
                numpy.savetxt("TrouNoirTrajectories.csv",
1832
                            numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),delimiter=' ', fmt='%.2e')
1833
        
1834
        TrajectoriesCL.release()
1835
        IdLastCL.release()
1836
        
1837
    return(elapsed)
1838

    
1839
def BlackHoleCUDA(zImage,fImage,InputCL):
1840

    
1841
    kernel_params = {}
1842

    
1843
    Device=InputCL['Device']
1844
    GpuStyle=InputCL['GpuStyle']
1845
    VariableType=InputCL['VariableType']
1846
    Size=InputCL['Size']
1847
    Mass=InputCL['Mass']
1848
    InternalRadius=InputCL['InternalRadius']
1849
    ExternalRadius=InputCL['ExternalRadius']
1850
    Angle=InputCL['Angle']
1851
    Method=InputCL['Method']
1852
    TrackPoints=InputCL['TrackPoints']
1853
    Physics=InputCL['Physics']
1854
    Threads=InputCL['Threads']
1855

    
1856
    PhysicsList=DictionariesAPI()
1857
    
1858
    if InputCL['BlackBody']:
1859
        # Spectrum is Black Body one
1860
        Line=0
1861
    else:
1862
        # Spectrum is Monochromatic Line one
1863
        Line=1
1864

    
1865
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1866
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1867

    
1868
    try:
1869
        # For PyCUDA import
1870
        import pycuda.driver as cuda
1871
        from pycuda.compiler import SourceModule
1872
        
1873
        cuda.init()
1874
        for Id in range(cuda.Device.count()):
1875
            if Id==Device:
1876
                XPU=cuda.Device(Id)
1877
                print("GPU selected %s" % XPU.name())
1878
        print
1879

    
1880
    except ImportError:
1881
        print("Platform does not seem to support CUDA")
1882

    
1883
    Context=XPU.make_context()
1884

    
1885
    try:
1886
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1887
        print("Compilation seems to be OK")
1888
    except:
1889
        print("Compilation seems to break")
1890

    
1891
    EachPixelCU=mod.get_function("EachPixel")
1892
    OriginalCU=mod.get_function("Original")
1893
    EachCircleCU=mod.get_function("EachCircle")
1894
    TrajectoryCU=mod.get_function("Trajectory")
1895
    PixelCU=mod.get_function("Pixel")
1896
    CircleCU=mod.get_function("Circle")
1897
    
1898
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1899
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1900
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1901
    cuda.memcpy_htod(zImageCU, zImage)
1902
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1903
    cuda.memcpy_htod(zImageCU, fImage)
1904
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1905
    cuda.memcpy_htod(IdLastCU, IdLast)
1906

    
1907
    start_time=time.time()
1908

    
1909
    if Method=='EachPixel':
1910
        EachPixelCU(zImageCU,fImageCU,
1911
                    numpy.float32(Mass),
1912
                    numpy.float32(InternalRadius),
1913
                    numpy.float32(ExternalRadius),
1914
                    numpy.float32(Angle),
1915
                    numpy.int32(Line),
1916
                    grid=(int(zImage.shape[0]/Threads),
1917
                          int(zImage.shape[1]/Threads)),
1918
                    block=(Threads,Threads,1))
1919
    elif Method=='EachCircle':
1920
        EachCircleCU(zImageCU,fImageCU,
1921
                   numpy.float32(Mass),
1922
                   numpy.float32(InternalRadius),
1923
                   numpy.float32(ExternalRadius),
1924
                   numpy.float32(Angle),
1925
                   numpy.int32(Line),
1926
                   grid=(int(zImage.shape[0]/Threads/2),1),
1927
                   block=(Threads,1,1))
1928
    elif Method=='Original':
1929
        OriginalCU(zImageCU,fImageCU,
1930
                   numpy.uint32(zImage.shape[0]/2),
1931
                   numpy.float32(Mass),
1932
                   numpy.float32(InternalRadius),
1933
                   numpy.float32(ExternalRadius),
1934
                   numpy.float32(Angle),
1935
                   numpy.int32(Line),
1936
                   grid=(1,1),
1937
                   block=(1,1,1))
1938
    elif Method=='TrajectoCircle':
1939
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1940
                     numpy.float32(Mass),
1941
                     numpy.float32(InternalRadius),
1942
                     numpy.float32(ExternalRadius),
1943
                     numpy.float32(Angle),
1944
                     numpy.int32(Line),
1945
                     grid=(int(Trajectories.shape[0]/Threads),1),
1946
                     block=(Threads,1,1))
1947
        
1948
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1949
                 numpy.float32(Mass),
1950
                 numpy.float32(InternalRadius),
1951
                 numpy.float32(ExternalRadius),
1952
                 numpy.float32(Angle),
1953
                 numpy.int32(Line),
1954
                 grid=(int(Trajectories.shape[0]/Threads),
1955
                       int(zImage.shape[0]*4/Threads)),
1956
                 block=(Threads,Threads,1))
1957
    else:
1958
        # Default method: TrajectoPixel
1959
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1960
                     numpy.float32(Mass),
1961
                     numpy.float32(InternalRadius),
1962
                     numpy.float32(ExternalRadius),
1963
                     numpy.float32(Angle),
1964
                     numpy.int32(Line),
1965
                     grid=(int(Trajectories.shape[0]/Threads),1),
1966
                     block=(Threads,1,1))
1967
        
1968
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1969
                numpy.uint32(Trajectories.shape[0]),
1970
                numpy.float32(Mass),
1971
                numpy.float32(InternalRadius),
1972
                numpy.float32(ExternalRadius),
1973
                numpy.float32(Angle),
1974
                numpy.int32(Line),
1975
                grid=(int(zImage.shape[0]/Threads),
1976
                      int(zImage.shape[1]/Threads),1),
1977
                block=(Threads,Threads,1))
1978

    
1979
    Context.synchronize()
1980

    
1981
    compute = time.time()-start_time
1982

    
1983
    cuda.memcpy_dtoh(zImage,zImageCU)
1984
    cuda.memcpy_dtoh(fImage,fImageCU)
1985
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1986
        cuda.memcpy_dtoh(Trajectories,TrajectoriesCU)
1987
    elapsed = time.time()-start_time
1988
    print("\nCompute Time : %f" % compute)
1989
    print("Elapsed Time : %f\n" % elapsed)
1990

    
1991
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1992
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1993
    print("Z max @(%f,%f) : %f" % ((1.*zMaxPosition[1][0]/zImage.shape[1]-0.5,1.*zMaxPosition[0][0]/zImage.shape[0]-0.5,zImage.max())))
1994
    print("Flux max @(%f,%f) : %f" % ((1.*fMaxPosition[1][0]/fImage.shape[1]-0.5,1.*fMaxPosition[0][0]/fImage.shape[0]-0.5,fImage.max())))
1995

    
1996
    
1997
    Context.pop()
1998

    
1999
    Context.detach()
2000

    
2001
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
2002
        if not NoImage:
2003
            AngleStep=4*numpy.pi/TrackPoints
2004
            Angles=numpy.arange(0.,4*numpy.pi,AngleStep)
2005
            Angles.shape=(1,TrackPoints)
2006
            Hostname=gethostname()
2007
            Date=time.strftime("%Y%m%d_%H%M%S")
2008
            ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2009
            
2010
            # numpy.savetxt("TrouNoirTrajectories_%s.csv" % ImageInfo,
2011
            #               numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
2012
            #               delimiter=' ', fmt='%.2e')
2013
            numpy.savetxt("TrouNoirTrajectories.csv",
2014
                          numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
2015
                          delimiter=' ', fmt='%.2e')
2016

    
2017
    
2018
    return(elapsed)
2019

    
2020
if __name__=='__main__':
2021

    
2022
    # Default device: first one!
2023
    Device=0
2024
    # Default implementation: OpenCL, most versatile!
2025
    GpuStyle = 'OpenCL'
2026
    Mass = 1.
2027
    # Internal Radius 3 times de Schwarzschild Radius
2028
    InternalRadius=6.*Mass
2029
    #
2030
    ExternalRadius=12.
2031
    #
2032
    # Angle with normal to disc 10 degrees
2033
    Angle = numpy.pi/180.*(90.-10.)
2034
    # Radiation of disc : BlackBody or Monochromatic
2035
    BlackBody = False
2036
    # Size of image
2037
    Size=1024
2038
    # Variable Type
2039
    VariableType='FP32'
2040
    # ?
2041
    q=-2
2042
    # Method of resolution
2043
    Method='TrajectoPixel'
2044
    # Colors for output image
2045
    Colors='Greyscale'
2046
    # Physics
2047
    Physics='Einstein'
2048
    # No output as image
2049
    NoImage = False
2050
    # Threads in CUDA
2051
    Threads = 32
2052
    # Trackpoints of trajectories
2053
    TrackPoints=2048
2054
    # Tracksave of trajectories
2055
    TrackSave=False
2056
    
2057
    HowToUse='%s -h [Help] -b [BlackBodyEmission] -j [TrackSave] -n [NoImage] -p <Einstein/Newton> -s <SizeInPixels> -m <Mass> -i <DiscInternalRadius> -e <DiscExternalRadius> -a <AngleAboveDisc> -d <DeviceId> -c <Greyscale/Red2Yellow> -g <CUDA/OpenCL> -o <EachPixel/TrajectoCircle/TrajectoPixel/EachCircle/Original> -t <ThreadsInCuda> -v <FP32/FP64> -k <TrackPoints>'
2058

    
2059
    try:
2060
        opts, args = getopt.getopt(sys.argv[1:],"hbnjs:m:i:e:a:d:g:v:o:t:c:p:k:",["tracksave","blackbody","noimage","camera","size=","mass=","internal=","external=","angle=","device=","gpustyle=","variabletype=","method=","threads=","colors=","physics=","trackpoints="])
2061
    except getopt.GetoptError:
2062
        print(HowToUse % sys.argv[0])
2063
        sys.exit(2)
2064

    
2065
    # List of Devices
2066
    Devices=[]
2067
    Alu={}
2068
        
2069
    for opt, arg in opts:
2070
        if opt == '-h':
2071
            print(HowToUse % sys.argv[0])
2072
            
2073
            print("\nInformations about devices detected under OpenCL API:")
2074
            # For PyOpenCL import
2075
            try:
2076
                import pyopencl as cl
2077
                Id=0
2078
                for platform in cl.get_platforms():
2079
                    for device in platform.get_devices():
2080
                        #deviceType=cl.device_type.to_string(device.type)
2081
                        deviceType="xPU"
2082
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
2083
                        Id=Id+1
2084
                        
2085
            except:
2086
                print("Your platform does not seem to support OpenCL")
2087
                
2088
            print("\nInformations about devices detected under CUDA API:")
2089
                # For PyCUDA import
2090
            try:
2091
                import pycuda.driver as cuda
2092
                cuda.init()
2093
                for Id in range(cuda.Device.count()):
2094
                    device=cuda.Device(Id)
2095
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
2096
                print
2097
            except:
2098
                print("Your platform does not seem to support CUDA")
2099
        
2100
            sys.exit()
2101
        
2102
        elif opt in ("-d", "--device"):
2103
#            Devices.append(int(arg))
2104
            Device=int(arg)
2105
        elif opt in ("-g", "--gpustyle"):
2106
            GpuStyle = arg
2107
        elif opt in ("-v", "--variabletype"):
2108
            VariableType = arg
2109
        elif opt in ("-s", "--size"):
2110
            Size = int(arg)
2111
        elif opt in ("-k", "--trackpoints"):
2112
            TrackPoints = int(arg)
2113
        elif opt in ("-m", "--mass"):
2114
            Mass = float(arg)
2115
        elif opt in ("-i", "--internal"):
2116
            InternalRadius = float(arg)
2117
        elif opt in ("-e", "--external"):
2118
            ExternalRadius = float(arg)
2119
        elif opt in ("-a", "--angle"):
2120
            Angle = numpy.pi/180.*(90.-float(arg))
2121
        elif opt in ("-b", "--blackbody"):
2122
            BlackBody = True
2123
        elif opt in ("-j", "--tracksave"):
2124
            TrackSave = True
2125
        elif opt in ("-n", "--noimage"):
2126
            NoImage = True
2127
        elif opt in ("-o", "--method"):
2128
            Method = arg
2129
        elif opt in ("-t", "--threads"):
2130
            Threads = int(arg)
2131
        elif opt in ("-c", "--colors"):
2132
            Colors = arg
2133
        elif opt in ("-p", "--physics"):
2134
            Physics = arg
2135

    
2136
    print("Device Identification selected : %s" % Device)
2137
    print("GpuStyle used : %s" % GpuStyle)
2138
    print("VariableType : %s" % VariableType)
2139
    print("Size : %i" % Size)
2140
    print("Mass : %f" % Mass)
2141
    print("Internal Radius : %f" % InternalRadius)
2142
    print("External Radius : %f" % ExternalRadius)
2143
    print("Angle with normal of (in radians) : %f" % Angle)
2144
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
2145
    print("Method of resolution : %s" % Method)
2146
    print("Colors for output images : %s" % Colors)
2147
    print("Physics used for Trajectories : %s" % Physics)
2148
    print("Trackpoints of Trajectories : %i" % TrackPoints)
2149
    print("Tracksave of Trajectories : %i" % TrackSave)
2150

    
2151
    if GpuStyle=='CUDA':
2152
        print("\nSelection of CUDA device") 
2153
        try:
2154
            # For PyCUDA import
2155
            import pycuda.driver as cuda
2156
            
2157
            cuda.init()
2158
            for Id in range(cuda.Device.count()):
2159
                device=cuda.Device(Id)
2160
                print("Device #%i of type GPU : %s" % (Id,device.name()))
2161
                if Id in Devices:
2162
                    Alu[Id]='GPU'
2163
            
2164
        except ImportError:
2165
            print("Platform does not seem to support CUDA")
2166

    
2167
    if GpuStyle=='OpenCL':
2168
        print("\nSelection of OpenCL device") 
2169
        try:
2170
            # For PyOpenCL import
2171
            import pyopencl as cl
2172
            Id=0
2173
            for platform in cl.get_platforms():
2174
                for device in platform.get_devices():
2175
                    #deviceType=cl.device_type.to_string(device.type)
2176
                    deviceType="xPU"
2177
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
2178

    
2179
                    if Id in Devices:
2180
                    # Set the Alu as detected Device Type
2181
                        Alu[Id]=deviceType
2182
                    Id=Id+1
2183
        except ImportError:
2184
            print("Platform does not seem to support OpenCL")
2185

    
2186

    
2187
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2188
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2189

    
2190
    InputCL={}
2191
    InputCL['Device']=Device
2192
    InputCL['GpuStyle']=GpuStyle
2193
    InputCL['VariableType']=VariableType
2194
    InputCL['Size']=Size
2195
    InputCL['Mass']=Mass
2196
    InputCL['InternalRadius']=InternalRadius
2197
    InputCL['ExternalRadius']=ExternalRadius
2198
    InputCL['Angle']=Angle
2199
    InputCL['BlackBody']=BlackBody
2200
    InputCL['Method']=Method
2201
    InputCL['TrackPoints']=TrackPoints
2202
    InputCL['Physics']=Physics
2203
    InputCL['Threads']=Threads
2204
    InputCL['NoImage']=NoImage
2205
    InputCL['TrackSave']=TrackSave
2206

    
2207
    if GpuStyle=='OpenCL':
2208
        duration=BlackHoleCL(zImage,fImage,InputCL)
2209
    else:
2210
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
2211
        
2212
    Hostname=gethostname()
2213
    Date=time.strftime("%Y%m%d_%H%M%S")
2214
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2215
    
2216
    
2217
    if not NoImage:
2218
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
2219
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)