Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 241

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

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

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

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

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

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

    
50

    
51

    
52

    
53

    
54

    
55

    
56

    
57

    
58

    
59

    
60

    
61

    
62

    
63

    
64

    
65

    
66

    
67

    
68

    
69

    
70

    
71

    
72

    
73

    
74

    
75

    
76

    
77

    
78

    
79

    
80

    
81

    
82

    
83

    
84

    
85

    
86

    
87

    
88

    
89

    
90

    
91

    
92

    
93

    
94

    
95

    
96

    
97

    
98

    
99

    
100

    
101

    
102

    
103
BlobOpenCL= """
104

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

108
#define EINSTEIN 0
109
#define NEWTON 1
110

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

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

121
  angle=atan2(y,x);
122

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

128
  return angle;
129
}
130

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

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

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

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

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

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

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

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

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

190
float spectre_cn(float rf32,float b32,float db32,
191
                 float h32,float r32,float m32,float bss32)
192
{
193

194
#define MYFLOAT float
195

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

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

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

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

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

220
  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 ));
221

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

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

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

242
          flx+=flux_int;
243
        }
244
    }
245

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

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

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

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

269
  *zp=1.e0f/rf;
270
  *fp=flx;
271

272
}
273

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

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

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

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

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

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

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

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

315

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

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

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

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

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

332
  int ExitOnImpact=0;
333

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

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

346

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

356
  barrier(CLK_GLOBAL_MEM_FENCE);
357

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

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

374
   // Perform trajectory for each pixel
375

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

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

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

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

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

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

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

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

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

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

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

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

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

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

438
      HalfLap++;
439
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
440

441
  }
442

443
  barrier(CLK_GLOBAL_MEM_FENCE);
444

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

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

464
   // Perform trajectory for each pixel
465

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

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

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

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

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

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

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

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

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

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

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

518
     HalfLap++;
519
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
520

521
  zImage[yi+2*bmaxi*xi]=zp;
522
  fImage[yi+2*bmaxi*xi]=fp;
523

524
  barrier(CLK_GLOBAL_MEM_FENCE);
525

526
}
527

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

538
  // Perform trajectory for each pixel
539

540
  float m,rs,re;
541

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

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

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

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

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

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

560
  up=0.e0f;
561
  vp=1.e0f;
562

563
  pp=0.e0f;
564
  nh=0;
565

566
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
567

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

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

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

585
  IdLast[bi]=nh;
586

587
  barrier(CLK_GLOBAL_MEM_FENCE);
588

589
}
590

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

601
  private float Trajectory[TRACKPOINTS];
602

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

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

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

617

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

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

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

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

630
  up=0.e0f;
631
  vp=1.e0f;
632

633
  pp=0.e0f;
634
  nh=0;
635

636
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
637

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

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

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

655

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

660

661
  uint imx=(uint)(16*bi);
662

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

671
     uint HalfLap=0,ExitOnImpact=0,ni;
672
     float php,nr,r;
673

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

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

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

695
        HalfLap++;
696

697
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
698

699
     zImage[yi+2*bmaxi*xi]=zp;
700
     fImage[yi+2*bmaxi*xi]=fp;
701

702
  }
703

704
  barrier(CLK_GLOBAL_MEM_FENCE);
705

706
}
707

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

716
   float Trajectory[TRACKPOINTS];
717

718
   // Perform trajectory for each pixel
719

720
   float m,rs,ri,re,tho;
721
   int raie,q;
722

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

731
   float bmx,db,b,h;
732
   uint nh;
733

734
   // Autosize for image
735
   bmx=1.25e0f*re;
736

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

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

747
      float up,vp,pp,us,vs,ps;
748

749
      up=0.e0f;
750
      vp=1.e0f;
751

752
      pp=0.e0f;
753
      nh=0;
754

755
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
756

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

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

772
      } while ((bvus>=rs)&&(bvus<=bvus0));
773

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

778
      int imx=(int)(16*bi);
779

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

788
         uint HalfLap=0,ExitOnImpact=0,ni;
789
         float php,nr,r;
790

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

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

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

812
            HalfLap++;
813

814
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
815

816
         zImage[yi+2*bmaxi*xi]=zp;
817
         fImage[yi+2*bmaxi*xi]=fp;
818

819
      }
820

821
   }
822

823
   barrier(CLK_GLOBAL_MEM_FENCE);
824

825
}
826
"""
827

    
828
def KernelCodeCuda():
829
    BlobCUDA= """
830

831
#define PI (float)3.14159265359
832
#define nbr 256
833

834
#define EINSTEIN 0
835
#define NEWTON 1
836

837
#ifdef SETTRACKPOINTS
838
#define TRACKPOINTS SETTRACKPOINTS
839
#else
840
#define TRACKPOINTS
841
#endif
842

843
__device__ float nothing(float x)
844
{
845
  return(x);
846
}
847

848
__device__ float atanp(float x,float y)
849
{
850
  float angle;
851

852
  angle=atan2(y,x);
853

854
  if (angle<0.e0f)
855
    {
856
      angle+=(float)2.e0f*PI;
857
    }
858

859
  return(angle);
860
}
861

862
__device__ float f(float v)
863
{
864
  return(v);
865
}
866

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

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

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

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

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

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

911
__device__ float spectre(float rf,int q,float b,float db,
912
                         float h,float r,float m,float bss)
913
{
914
  float flx;
915

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

921
__device__ float spectre_cn(float rf32,float b32,float db32,
922
                            float h32,float r32,float m32,float bss32)
923
{
924

925
#define MYFLOAT float
926

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

935
  MYFLOAT flx;
936
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
937
  int fi,posfreq;
938

939
#define planck 6.62e-34
940
#define k 1.38e-23
941
#define c2 9.e16
942
#define temp 3.e7
943
#define m_point 1.
944

945
#define lplanck (log(6.62)-34.*log(10.))
946
#define lk (log(1.38)-23.*log(10.))
947
#define lc2 (log(9.)+16.*log(10.))
948

949
  MYFLOAT v=1.-3./r;
950

951
  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 ));
952

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

955
  flux_int=0.;
956
  flx=0.;
957

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

973
          flx+=flux_int;
974
        }
975
    }
976

977
  return((float)(flx));
978
}
979

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

987
  rf=decalage_spectral(r,b,phi,tho,m);
988

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

1000
  *zp=1./rf;
1001
  *fp=flx;
1002

1003
}
1004

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

1015
   // Perform trajectory for each pixel, exit on hit
1016

1017
  float m,rs,ri,re,tho;
1018
  int q,raie;
1019

1020
  m=Mass;
1021
  rs=2.*m;
1022
  ri=InternalRadius;
1023
  re=ExternalRadius;
1024
  tho=Angle;
1025
  q=-2;
1026
  raie=Line;
1027

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

1034
  // Autosize for image
1035
  bmx=1.25*re;
1036
  b=0.;
1037

1038
  h=4.e0f*PI/(float)TRACKPOINTS;
1039

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

1047
  float up,vp,pp,us,vs,ps;
1048

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

1055
  up=0.;
1056
  vp=1.;
1057
  pp=0.;
1058
  nh=0;
1059

1060
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1061

1062
  rps=fabs(b/us);
1063
  rp0=rps;
1064

1065
  int ExitOnImpact=0;
1066

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

1078
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
1079

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

1089
  __syncthreads();
1090

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

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

1107
  // Perform trajectory for each pixel
1108

1109
  float m,rs,ri,re,tho;
1110
  int q,raie;
1111

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

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

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

1130
  // Step of Impact Parameter
1131
  db=bmx/(2.e0f*(float)ImpactParameter);
1132

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

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

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

1146
  int HalfLap=0,ExitOnImpact=0,ni;
1147

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

1156
      if (ni<IdLast[bi])
1157
      {
1158
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.e0f)+Trajectories[bi*TRACKPOINTS+ni];
1159
      }
1160
      else
1161
      {
1162
        r=Trajectories[bi*TRACKPOINTS+ni];
1163
      }
1164

1165
      if ((r<=re)&&(r>=ri))
1166
      {
1167
        ExitOnImpact=1;
1168
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1169
      }
1170

1171
      HalfLap++;
1172
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1173

1174
  }
1175

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

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

1195
   // Perform trajectory for each pixel
1196

1197
  float m,rs,ri,re,tho;
1198
  int q,raie;
1199

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

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

1212
  // Autosize for image
1213
  bmx=1.25e0f*re;
1214

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

1218
  // impact parameter
1219
  b=(float)bi/(float)bmaxi*bmx;
1220
  db=bmx/(2.e0f*(float)bmaxi);
1221

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

1227
  int HalfLap=0,ExitOnImpact=0,ni;
1228
  float php,nr,r;
1229

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

1236
     if (ni<IdLast[bi])
1237
     {
1238
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.e0f)+Trajectories[bi*TRACKPOINTS+ni];
1239
     }
1240
     else
1241
     {
1242
        r=Trajectories[bi*TRACKPOINTS+ni];
1243
     }
1244

1245
     if ((r<=re)&&(r>=ri))
1246
     {
1247
        ExitOnImpact=1;
1248
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1249
     }
1250

1251
     HalfLap++;
1252
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1253

1254
  zImage[yi+2*bmaxi*xi]=zp;
1255
  fImage[yi+2*bmaxi*xi]=fp;
1256

1257
}
1258

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

1269
  // Perform trajectory for each pixel
1270

1271
  float m,rs,ri,re,tho;
1272
  int raie,q;
1273

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

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

1287
  // Autosize for image
1288
  bmx=1.25e0f*re;
1289

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

1293
  // impact parameter
1294
  b=(float)bi/(float)bmaxi*bmx;
1295

1296
  float up,vp,pp,us,vs,ps;
1297

1298
  up=0.e0f;
1299
  vp=1.e0f;
1300
  pp=0.e0f;
1301
  nh=0;
1302

1303
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1304

1305
  // b versus us
1306
  float bvus=fabs(b/us);
1307
  float bvus0=bvus;
1308
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
1309

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

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

1322
  IdLast[bi]=nh;
1323

1324
}
1325

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

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

1337
  float Trajectory[2048];
1338

1339
  // Perform trajectory for each pixel
1340

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

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

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

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

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

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

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

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

1372
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1373

1374
  // b versus us
1375
  float bvus=fabs(b/us);
1376
  float bvus0=bvus;
1377
  Trajectory[nh]=bvus;
1378

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

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

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

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

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

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

1410
        if (ni<nh)
1411
        {
1412
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1413
        }
1414
        else
1415
        {
1416
           r=Trajectory[ni];
1417
        }
1418

1419
        if ((r<=re)&&(r>=ri))
1420
        {
1421
           ExitOnImpact=1;
1422
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1423
        }
1424

1425
        HalfLap++;
1426

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

1429
   __syncthreads();
1430

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

1434
  }
1435

1436
}
1437

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

1446
   float Trajectory[TRACKPOINTS];
1447

1448
   // Perform trajectory for each pixel
1449

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

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

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

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

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

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

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

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

1484
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1485

1486
      // b versus us
1487
      float bvus=fabs(b/us);
1488
      float bvus0=bvus;
1489
      Trajectory[nh]=bvus;
1490

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

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

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

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

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

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

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

1526
            if (ni<nh)
1527
            {
1528
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1529
            }
1530
            else
1531
            {
1532
               r=Trajectory[ni];
1533
            }
1534

1535
            if ((r<=re)&&(r>=ri))
1536
            {
1537
               ExitOnImpact=1;
1538
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1539
            }
1540

1541
            HalfLap++;
1542

1543
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
1544

1545
         zImage[yi+2*bmaxi*xi]=zp;
1546
         fImage[yi+2*bmaxi*xi]=fp;
1547

1548
      }
1549

1550
   }
1551

1552
}
1553
"""
1554
    return(BlobCUDA)
1555

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

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

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

    
1578
    kernel_params = {}
1579

    
1580
    print(InputCL)
1581

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

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

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

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

    
1619
    if HasXPU==False:
1620
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1621
        sys.exit()
1622

    
1623
    ctx = cl.Context([XPU])
1624
    queue = cl.CommandQueue(ctx,
1625
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1626

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

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

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

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

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

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

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

    
1732
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1733
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1734
    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())))
1735
    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())))
1736
    zImageCL.release()
1737
    fImageCL.release()
1738

    
1739
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1740
        if not NoImage:
1741
            AngleStep=4*numpy.pi/TrackPoints
1742
            Angles=numpy.arange(0.,4*numpy.pi,AngleStep)
1743
            Angles.shape=(1,TrackPoints)
1744
            Hostname=gethostname()
1745
            Date=time.strftime("%Y%m%d_%H%M%S")
1746
            ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1747
            
1748
            # numpy.savetxt("TrouNoirTrajectories_%s.csv" % ImageInfo,
1749
            #               numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1750
            #               delimiter=' ', fmt='%.2e')
1751
            numpy.savetxt("TrouNoirTrajectories.csv",
1752
                          numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1753
                          delimiter=' ', fmt='%.2e')
1754
        
1755
        TrajectoriesCL.release()
1756
        IdLastCL.release()
1757
        
1758
    return(elapsed)
1759

    
1760
def BlackHoleCUDA(zImage,fImage,InputCL):
1761

    
1762
    kernel_params = {}
1763

    
1764
    print(InputCL)
1765

    
1766
    Device=InputCL['Device']
1767
    GpuStyle=InputCL['GpuStyle']
1768
    VariableType=InputCL['VariableType']
1769
    Size=InputCL['Size']
1770
    Mass=InputCL['Mass']
1771
    InternalRadius=InputCL['InternalRadius']
1772
    ExternalRadius=InputCL['ExternalRadius']
1773
    Angle=InputCL['Angle']
1774
    Method=InputCL['Method']
1775
    TrackPoints=InputCL['TrackPoints']
1776
    Physics=InputCL['Physics']
1777
    Threads=InputCL['Threads']
1778

    
1779
    PhysicsList=DictionariesAPI()
1780
    
1781
    if InputCL['BlackBody']:
1782
        # Spectrum is Black Body one
1783
        Line=0
1784
    else:
1785
        # Spectrum is Monochromatic Line one
1786
        Line=1
1787

    
1788
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1789
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1790

    
1791
    try:
1792
        # For PyCUDA import
1793
        import pycuda.driver as cuda
1794
        from pycuda.compiler import SourceModule
1795
        
1796
        cuda.init()
1797
        for Id in range(cuda.Device.count()):
1798
            if Id==Device:
1799
                XPU=cuda.Device(Id)
1800
                print("GPU selected %s" % XPU.name())
1801
        print
1802

    
1803
    except ImportError:
1804
        print("Platform does not seem to support CUDA")
1805

    
1806
    Context=XPU.make_context()
1807

    
1808
    try:
1809
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1810
        print("Compilation seems to be OK")
1811
    except:
1812
        print("Compilation seems to break")
1813

    
1814
    EachPixelCU=mod.get_function("EachPixel")
1815
    OriginalCU=mod.get_function("Original")
1816
    EachCircleCU=mod.get_function("EachCircle")
1817
    TrajectoryCU=mod.get_function("Trajectory")
1818
    PixelCU=mod.get_function("Pixel")
1819
    CircleCU=mod.get_function("Circle")
1820
    
1821
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1822
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1823
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1824
    cuda.memcpy_htod(zImageCU, zImage)
1825
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1826
    cuda.memcpy_htod(zImageCU, fImage)
1827
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1828
    cuda.memcpy_htod(IdLastCU, IdLast)
1829

    
1830
    start_time=time.time()
1831

    
1832
    if Method=='EachPixel':
1833
        EachPixelCU(zImageCU,fImageCU,
1834
                    numpy.float32(Mass),
1835
                    numpy.float32(InternalRadius),
1836
                    numpy.float32(ExternalRadius),
1837
                    numpy.float32(Angle),
1838
                    numpy.int32(Line),
1839
                    grid=(zImage.shape[0]/Threads,zImage.shape[1]/Threads),
1840
                    block=(Threads,Threads,1))
1841
    elif Method=='EachCircle':
1842
        EachCircleCU(zImageCU,fImageCU,
1843
                   numpy.float32(Mass),
1844
                   numpy.float32(InternalRadius),
1845
                   numpy.float32(ExternalRadius),
1846
                   numpy.float32(Angle),
1847
                   numpy.int32(Line),
1848
                   grid=(int(zImage.shape[0]/Threads/2),1),
1849
                   block=(Threads,1,1))
1850
    elif Method=='Original':
1851
        OriginalCU(zImageCU,fImageCU,
1852
                   numpy.uint32(zImage.shape[0]/2),
1853
                   numpy.float32(Mass),
1854
                   numpy.float32(InternalRadius),
1855
                   numpy.float32(ExternalRadius),
1856
                   numpy.float32(Angle),
1857
                   numpy.int32(Line),
1858
                   grid=(1,1),
1859
                   block=(1,1,1))
1860
    elif Method=='TrajectoCircle':
1861
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1862
                     numpy.float32(Mass),
1863
                     numpy.float32(InternalRadius),
1864
                     numpy.float32(ExternalRadius),
1865
                     numpy.float32(Angle),
1866
                     numpy.int32(Line),
1867
                     grid=(int(Trajectories.shape[0]/Threads),1),
1868
                     block=(Threads,1,1))
1869
        
1870
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1871
                 numpy.float32(Mass),
1872
                 numpy.float32(InternalRadius),
1873
                 numpy.float32(ExternalRadius),
1874
                 numpy.float32(Angle),
1875
                 numpy.int32(Line),
1876
                 grid=(int(Trajectories.shape[0]/Threads),
1877
                       int(zImage.shape[0]*4/Threads)),
1878
                 block=(Threads,Threads,1))
1879
    else:
1880
        # Default method: TrajectoPixel
1881
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1882
                     numpy.float32(Mass),
1883
                     numpy.float32(InternalRadius),
1884
                     numpy.float32(ExternalRadius),
1885
                     numpy.float32(Angle),
1886
                     numpy.int32(Line),
1887
                     grid=(int(Trajectories.shape[0]/Threads),1),
1888
                     block=(Threads,1,1))
1889
        
1890
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1891
                numpy.uint32(Trajectories.shape[0]),
1892
                numpy.float32(Mass),
1893
                numpy.float32(InternalRadius),
1894
                numpy.float32(ExternalRadius),
1895
                numpy.float32(Angle),
1896
                numpy.int32(Line),
1897
                grid=(int(zImage.shape[0]/Threads),
1898
                      int(zImage.shape[1]/Threads),1),
1899
                block=(Threads,Threads,1))
1900

    
1901
    Context.synchronize()
1902

    
1903
    compute = time.time()-start_time
1904

    
1905
    cuda.memcpy_dtoh(zImage,zImageCU)
1906
    cuda.memcpy_dtoh(fImage,fImageCU)
1907
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1908
        cuda.memcpy_dtoh(Trajectories,TrajectoriesCU)
1909
    elapsed = time.time()-start_time
1910
    print("\nCompute Time : %f" % compute)
1911
    print("Elapsed Time : %f\n" % elapsed)
1912

    
1913
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1914
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1915
    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())))
1916
    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())))
1917

    
1918
    
1919
    Context.pop()
1920

    
1921
    Context.detach()
1922

    
1923
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1924
        if not NoImage:
1925
            AngleStep=4*numpy.pi/TrackPoints
1926
            Angles=numpy.arange(0.,4*numpy.pi,AngleStep)
1927
            Angles.shape=(1,TrackPoints)
1928
            Hostname=gethostname()
1929
            Date=time.strftime("%Y%m%d_%H%M%S")
1930
            ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1931
            
1932
            # numpy.savetxt("TrouNoirTrajectories_%s.csv" % ImageInfo,
1933
            #               numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1934
            #               delimiter=' ', fmt='%.2e')
1935
            numpy.savetxt("TrouNoirTrajectories.csv",
1936
                          numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1937
                          delimiter=' ', fmt='%.2e')
1938

    
1939
    
1940
    return(elapsed)
1941

    
1942
if __name__=='__main__':
1943
        
1944
    GpuStyle = 'OpenCL'
1945
    Mass = 1.
1946
    # Internal Radius 3 times de Schwarzschild Radius
1947
    InternalRadius=6.*Mass
1948
    #
1949
    ExternalRadius=12.
1950
    #
1951
    # Angle with normal to disc 10 degrees
1952
    Angle = numpy.pi/180.*(90.-10.)
1953
    # Radiation of disc : BlackBody or Monochromatic
1954
    BlackBody = False
1955
    # Size of image
1956
    Size=256
1957
    # Variable Type
1958
    VariableType='FP32'
1959
    # ?
1960
    q=-2
1961
    # Method of resolution
1962
    Method='TrajectoPixel'
1963
    # Colors for output image
1964
    Colors='Greyscale'
1965
    # Physics
1966
    Physics='Einstein'
1967
    # No output as image
1968
    NoImage = False
1969
    # Threads in CUDA
1970
    Threads = 32
1971
    # Trackpoints of trajectories
1972
    TrackPoints=2048
1973
    
1974
    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> -o <EachPixel/TrajectoCircle/TrajectoPixel/EachCircle/Original> -t <ThreadsInCuda> -v <FP32/FP64> -k <TrackPoints>'
1975

    
1976
    try:
1977
        opts, args = getopt.getopt(sys.argv[1:],"hbns:m:i:x:a:d:g:v:o:t:c:p:k:",["blackbody","noimage","camera","size=","mass=","internal=","external=","angle=","device=","gpustyle=","variabletype=","method=","threads=","colors=","physics=","trackpoints="])
1978
    except getopt.GetoptError:
1979
        print(HowToUse % sys.argv[0])
1980
        sys.exit(2)
1981

    
1982
    # List of Devices
1983
    Devices=[]
1984
    Alu={}
1985
        
1986
    for opt, arg in opts:
1987
        if opt == '-h':
1988
            print(HowToUse % sys.argv[0])
1989
            
1990
            print("\nInformations about devices detected under OpenCL API:")
1991
            # For PyOpenCL import
1992
            try:
1993
                import pyopencl as cl
1994
                Id=0
1995
                for platform in cl.get_platforms():
1996
                    for device in platform.get_devices():
1997
                        #deviceType=cl.device_type.to_string(device.type)
1998
                        deviceType="xPU"
1999
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
2000
                        Id=Id+1
2001
                        
2002
            except:
2003
                print("Your platform does not seem to support OpenCL")
2004
                
2005
            print("\nInformations about devices detected under CUDA API:")
2006
                # For PyCUDA import
2007
            try:
2008
                import pycuda.driver as cuda
2009
                cuda.init()
2010
                for Id in range(cuda.Device.count()):
2011
                    device=cuda.Device(Id)
2012
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
2013
                print
2014
            except:
2015
                print("Your platform does not seem to support CUDA")
2016
        
2017
            sys.exit()
2018
        
2019
        elif opt in ("-d", "--device"):
2020
#            Devices.append(int(arg))
2021
            Device=int(arg)
2022
        elif opt in ("-g", "--gpustyle"):
2023
            GpuStyle = arg
2024
        elif opt in ("-v", "--variabletype"):
2025
            VariableType = arg
2026
        elif opt in ("-s", "--size"):
2027
            Size = int(arg)
2028
        elif opt in ("-k", "--trackpoints"):
2029
            TrackPoints = int(arg)
2030
        elif opt in ("-m", "--mass"):
2031
            Mass = float(arg)
2032
        elif opt in ("-i", "--internal"):
2033
            InternalRadius = float(arg)
2034
        elif opt in ("-e", "--external"):
2035
            ExternalRadius = float(arg)
2036
        elif opt in ("-a", "--angle"):
2037
            Angle = numpy.pi/180.*(90.-float(arg))
2038
        elif opt in ("-b", "--blackbody"):
2039
            BlackBody = True
2040
        elif opt in ("-n", "--noimage"):
2041
            NoImage = True
2042
        elif opt in ("-o", "--method"):
2043
            Method = arg
2044
        elif opt in ("-t", "--threads"):
2045
            Threads = int(arg)
2046
        elif opt in ("-c", "--colors"):
2047
            Colors = arg
2048
        elif opt in ("-p", "--physics"):
2049
            Physics = arg
2050

    
2051
    print("Device Identification selected : %s" % Device)
2052
    print("GpuStyle used : %s" % GpuStyle)
2053
    print("VariableType : %s" % VariableType)
2054
    print("Size : %i" % Size)
2055
    print("Mass : %f" % Mass)
2056
    print("Internal Radius : %f" % InternalRadius)
2057
    print("External Radius : %f" % ExternalRadius)
2058
    print("Angle with normal of (in radians) : %f" % Angle)
2059
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
2060
    print("Method of resolution : %s" % Method)
2061
    print("Colors for output images : %s" % Colors)
2062
    print("Physics used for Trajectories : %s" % Physics)
2063
    print("Trackpoints of Trajectories : %i" % TrackPoints)
2064

    
2065
    if GpuStyle=='CUDA':
2066
        print("\nSelection of CUDA device") 
2067
        try:
2068
            # For PyCUDA import
2069
            import pycuda.driver as cuda
2070
            
2071
            cuda.init()
2072
            for Id in range(cuda.Device.count()):
2073
                device=cuda.Device(Id)
2074
                print("Device #%i of type GPU : %s" % (Id,device.name()))
2075
                if Id in Devices:
2076
                    Alu[Id]='GPU'
2077
            
2078
        except ImportError:
2079
            print("Platform does not seem to support CUDA")
2080

    
2081
    if GpuStyle=='OpenCL':
2082
        print("\nSelection of OpenCL device") 
2083
        try:
2084
            # For PyOpenCL import
2085
            import pyopencl as cl
2086
            Id=0
2087
            for platform in cl.get_platforms():
2088
                for device in platform.get_devices():
2089
                    #deviceType=cl.device_type.to_string(device.type)
2090
                    deviceType="xPU"
2091
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
2092

    
2093
                    if Id in Devices:
2094
                    # Set the Alu as detected Device Type
2095
                        Alu[Id]=deviceType
2096
                    Id=Id+1
2097
        except ImportError:
2098
            print("Platform does not seem to support OpenCL")
2099

    
2100
#    print(Devices,Alu)
2101

    
2102
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
2103
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2104
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2105

    
2106
    InputCL={}
2107
    InputCL['Device']=Device
2108
    InputCL['GpuStyle']=GpuStyle
2109
    InputCL['VariableType']=VariableType
2110
    InputCL['Size']=Size
2111
    InputCL['Mass']=Mass
2112
    InputCL['InternalRadius']=InternalRadius
2113
    InputCL['ExternalRadius']=ExternalRadius
2114
    InputCL['Angle']=Angle
2115
    InputCL['BlackBody']=BlackBody
2116
    InputCL['Method']=Method
2117
    InputCL['TrackPoints']=TrackPoints
2118
    InputCL['Physics']=Physics
2119
    InputCL['Threads']=Threads
2120
    InputCL['NoImage']=NoImage
2121

    
2122
    if GpuStyle=='OpenCL':
2123
        duration=BlackHoleCL(zImage,fImage,InputCL)
2124
    else:
2125
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
2126
        
2127
    Hostname=gethostname()
2128
    Date=time.strftime("%Y%m%d_%H%M%S")
2129
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2130
    
2131
    
2132
    if not NoImage:
2133
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
2134
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)