Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 258

Historique | Voir | Annoter | Télécharger (53,06 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
# Thanks to Andreas Klockner for PyOpenCL and PyCUDA:
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

    
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
def KernelCodeCuda():
900
    BlobCUDA= """
901

902
#define PI (float)3.14159265359
903
#define nbr 256
904

905
#define EINSTEIN 0
906
#define NEWTON 1
907

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

918
__device__ float atanp(float x,float y)
919
{
920
  float angle;
921

922
  angle=atan2(y,x);
923

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

929
  return(angle);
930
}
931

932
__device__ float f(float v)
933
{
934
  return(v);
935
}
936

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

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

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

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

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

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

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

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

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

995
#define MYFLOAT float
996

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

1005
  MYFLOAT flx;
1006
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
1007
  int fi,posfreq;
1008

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

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

1019
  MYFLOAT v=1.-3./r;
1020

1021
  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 ));
1022

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

1025
  flux_int=0.;
1026
  flx=0.;
1027

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

1043
          flx+=flux_int;
1044
        }
1045
    }
1046

1047
  return((float)(flx));
1048
}
1049

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

1057
  rf=decalage_spectral(r,b,phi,tho,m);
1058

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

1070
  *zp=1./rf;
1071
  *fp=flx;
1072

1073
}
1074

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

1085

1086
   // Perform trajectory for each pixel, exit on hit
1087

1088
  float m,rs,ri,re,tho;
1089
  int q,raie;
1090

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

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

1105
  // Autosize for image
1106
  bmx=1.25*re;
1107
  b=0.;
1108

1109
  h=4.e0f*PI/(float)TRACKPOINTS;
1110

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

1118
  float up,vp,pp,us,vs,ps;
1119

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

1126
  up=0.;
1127
  vp=1.;
1128
  pp=0.;
1129
  nh=0;
1130

1131
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1132

1133
  rps=fabs(b/us);
1134
  rp0=rps;
1135

1136
  int ExitOnImpact=0;
1137

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

1149
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
1150

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

1160
  __syncthreads();
1161

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

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

1178
  // Perform trajectory for each pixel
1179

1180
  float m,ri,re,tho;
1181
  int q,raie;
1182

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

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

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

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

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

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

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

1215
  int HalfLap=0,ExitOnImpact=0,ni;
1216

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

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

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

1240
      HalfLap++;
1241
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1242

1243
  }
1244

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

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

1264
   // Perform trajectory for each pixel
1265

1266
  float m,ri,re,tho;
1267
  int q,raie;
1268

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

1275
  float bmx,db,b,h;
1276
  float phi,phd;
1277
  float zp=0,fp=0;
1278

1279
  // Autosize for image
1280
  bmx=1.25e0f*re;
1281

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

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

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

1294
  int HalfLap=0,ExitOnImpact=0,ni;
1295
  float php,nr,r;
1296

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

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

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

1318
     HalfLap++;
1319
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1320

1321
  zImage[yi+2*bmaxi*xi]=zp;
1322
  fImage[yi+2*bmaxi*xi]=fp;
1323

1324
}
1325

1326
__global__ void Trajectory(float *Trajectories,int *IdLast,
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
  // Integer Impact Parameter Size (half of image)
1334
  int bmaxi=gridDim.x*blockDim.x;
1335

1336
  // Perform trajectory for each pixel
1337

1338
  float m,rs,re;
1339

1340
  m=Mass;
1341
  rs=2.e0f*m;
1342
  re=ExternalRadius;
1343

1344
  float bmx,b,h;
1345
  int nh;
1346

1347
  // Autosize for image
1348
  bmx=1.25e0f*re;
1349

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

1353
  // impact parameter
1354
  b=(float)bi/(float)bmaxi*bmx;
1355

1356
  float up,vp,pp,us,vs,ps;
1357

1358
  up=0.e0f;
1359
  vp=1.e0f;
1360
  pp=0.e0f;
1361
  nh=0;
1362

1363
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1364

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

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

1380
  } while ((bvus>=rs)&&(bvus<=bvus0));
1381

1382
  IdLast[bi]=nh;
1383

1384
}
1385

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

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

1397
  float Trajectory[2048];
1398

1399
  // Perform trajectory for each pixel
1400

1401
  float m,rs,ri,re,tho;
1402
  int raie,q;
1403

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

1412
  float bmx,db,b,h;
1413
  int nh;
1414

1415
  // Autosize for image
1416
  bmx=1.25e0f*re;
1417

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

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

1425
  float up,vp,pp,us,vs,ps;
1426

1427
  up=0.;
1428
  vp=1.;
1429
  pp=0.;
1430
  nh=0;
1431

1432
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1433

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

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

1449
  } while ((bvus>=rs)&&(bvus<=bvus0));
1450

1451
  int imx=(int)(16*bi);
1452

1453
  for (int i=0;i<imx;i++)
1454
  {
1455
     float zp=0,fp=0;
1456
     float phi=2.*PI/(float)imx*(float)i;
1457
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
1458
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1459
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1460

1461
     int HalfLap=0,ExitOnImpact=0,ni;
1462
     float php,nr,r;
1463

1464
     do
1465
     {
1466
        php=phd+(float)HalfLap*PI;
1467
        nr=php/h;
1468
        ni=(int)nr;
1469

1470
        if (ni<nh)
1471
        {
1472
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1473
        }
1474
        else
1475
        {
1476
           r=Trajectory[ni];
1477
        }
1478

1479
        if ((r<=re)&&(r>=ri))
1480
        {
1481
           ExitOnImpact=1;
1482
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1483
        }
1484

1485
        HalfLap++;
1486

1487
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
1488

1489
   __syncthreads();
1490

1491
   zImage[yi+2*bmaxi*xi]=zp;
1492
   fImage[yi+2*bmaxi*xi]=fp;
1493

1494
  }
1495

1496
}
1497

1498
__global__ void Original(float *zImage,float *fImage,
1499
                         uint Size,float Mass,float InternalRadius,
1500
                         float ExternalRadius,float Angle,
1501
                         int Line)
1502
{
1503
   // Integer Impact Parameter Size (half of image)
1504
   uint bmaxi=(uint)Size;
1505

1506
   float Trajectory[TRACKPOINTS];
1507

1508
   // Perform trajectory for each pixel
1509

1510
   float m,rs,ri,re,tho;
1511
   int raie,q;
1512

1513
   m=Mass;
1514
   rs=2.e0f*m;
1515
   ri=InternalRadius;
1516
   re=ExternalRadius;
1517
   tho=Angle;
1518
   q=-2;
1519
   raie=Line;
1520

1521
   float bmx,db,b,h;
1522
   int nh;
1523

1524
   // Autosize for image
1525
   bmx=1.25e0f*re;
1526

1527
   // Angular step of integration
1528
   h=4.e0f*PI/(float)TRACKPOINTS;
1529

1530
   // Integer Impact Parameter ID
1531
   for (int bi=0;bi<bmaxi;bi++)
1532
   {
1533
      // impact parameter
1534
      b=(float)bi/(float)bmaxi*bmx;
1535
      db=bmx/(2.e0f*(float)bmaxi);
1536

1537
      float up,vp,pp,us,vs,ps;
1538

1539
      up=0.;
1540
      vp=1.;
1541
      pp=0.;
1542
      nh=0;
1543

1544
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1545

1546
      // b versus us
1547
      float bvus=fabs(b/us);
1548
      float bvus0=bvus;
1549
      Trajectory[nh]=bvus;
1550

1551
      do
1552
      {
1553
         nh++;
1554
         pp=ps;
1555
         up=us;
1556
         vp=vs;
1557
         rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1558
         bvus=fabs(b/us);
1559
         Trajectory[nh]=bvus;
1560

1561
      } while ((bvus>=rs)&&(bvus<=bvus0));
1562

1563
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
1564
         Trajectory[i]=0.e0f;
1565
      }
1566

1567
      int imx=(int)(16*bi);
1568

1569
      for (int i=0;i<imx;i++)
1570
      {
1571
         float zp=0,fp=0;
1572
         float phi=2.e0f*PI/(float)imx*(float)i;
1573
         float phd=atanp(cos(phi)*sin(tho),cos(tho));
1574
         uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1575
         uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1576

1577
         int HalfLap=0,ExitOnImpact=0,ni;
1578
         float php,nr,r;
1579

1580
         do
1581
         {
1582
            php=phd+(float)HalfLap*PI;
1583
            nr=php/h;
1584
            ni=(int)nr;
1585

1586
            if (ni<nh)
1587
            {
1588
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1589
            }
1590
            else
1591
            {
1592
               r=Trajectory[ni];
1593
            }
1594

1595
            if ((r<=re)&&(r>=ri))
1596
            {
1597
               ExitOnImpact=1;
1598
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1599
            }
1600

1601
            HalfLap++;
1602

1603
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
1604

1605
         zImage[yi+2*bmaxi*xi]=zp;
1606
         fImage[yi+2*bmaxi*xi]=fp;
1607

1608
      }
1609

1610
   }
1611

1612
}
1613
"""
1614
    return(BlobCUDA)
1615

    
1616
# def ImageOutput(sigma,prefix):
1617
#     from PIL import Image
1618
#     Max=sigma.max()
1619
#     Min=sigma.min()
1620
#     # Normalize value as 8bits Integer
1621
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1622
#     image = Image.fromarray(SigmaInt)
1623
#     image.save("%s.jpg" % prefix)
1624

    
1625
def ImageOutput(sigma,prefix,Colors):
1626
    import matplotlib.pyplot as plt
1627
    start_time=time.time()
1628
    if Colors == 'Red2Yellow':
1629
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1630
    else:
1631
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1632
    save_time = time.time()-start_time
1633
    print("Save image as %s.png file" % prefix)
1634
    print("Save Time : %f" % save_time)
1635

    
1636
def BlackHoleCL(zImage,fImage,InputCL):
1637

    
1638
    kernel_params = {}
1639

    
1640
    Device=InputCL['Device']
1641
    GpuStyle=InputCL['GpuStyle']
1642
    VariableType=InputCL['VariableType']
1643
    Size=InputCL['Size']
1644
    Mass=InputCL['Mass']
1645
    InternalRadius=InputCL['InternalRadius']
1646
    ExternalRadius=InputCL['ExternalRadius']
1647
    Angle=InputCL['Angle']
1648
    Method=InputCL['Method']
1649
    TrackPoints=InputCL['TrackPoints']
1650
    Physics=InputCL['Physics']
1651
    NoImage=InputCL['NoImage']
1652
    TrackSave=InputCL['TrackSave']
1653

    
1654
    PhysicsList=DictionariesAPI()
1655
    
1656
    if InputCL['BlackBody']:
1657
        # Spectrum is Black Body one
1658
        Line=0
1659
    else:
1660
        # Spectrum is Monochromatic Line one
1661
        Line=1
1662

    
1663
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1664
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1665

    
1666
    # Je detecte un peripherique GPU dans la liste des peripheriques
1667
    Id=0
1668
    HasXPU=False
1669
    for platform in cl.get_platforms():
1670
        for device in platform.get_devices():
1671
            if Id==Device:
1672
                PF4XPU=platform.name
1673
                XPU=device
1674
                print("CPU/GPU selected: ",device.name.lstrip())
1675
                HasXPU=True
1676
            Id+=1
1677

    
1678
    if HasXPU==False:
1679
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1680
        sys.exit()
1681

    
1682
    ctx = cl.Context([XPU])
1683
    queue = cl.CommandQueue(ctx,
1684
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1685

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

    
1688
    print('My Platform is ',PF4XPU)
1689
    
1690
    if 'Intel' in PF4XPU or 'Experimental' in PF4XPU or 'Clover' in PF4XPU or 'Portable' in PF4XPU :
1691
        print('No extra options for Intel and Clover!')
1692
    else:
1693
        BuildOptions = BuildOptions+" -cl-mad-enable"
1694

    
1695
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1696
    
1697
    # Je recupere les flag possibles pour les buffers
1698
    mf = cl.mem_flags
1699

    
1700
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1701
        TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories) 
1702
        IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1703

    
1704
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1705
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1706
    
1707
    start_time=time.time()
1708

    
1709
    if Method=='EachPixel':
1710
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1711
                                       None,zImageCL,fImageCL,
1712
                                       numpy.float32(Mass),
1713
                                       numpy.float32(InternalRadius),
1714
                                       numpy.float32(ExternalRadius),
1715
                                       numpy.float32(Angle),
1716
                                       numpy.int32(Line))
1717
        CLLaunch.wait()
1718
    elif Method=='Original':
1719
        CLLaunch=BlackHoleCL.Original(queue,(1,),
1720
                                      None,zImageCL,fImageCL,
1721
                                      numpy.uint32(zImage.shape[0]/2),
1722
                                      numpy.float32(Mass),
1723
                                      numpy.float32(InternalRadius),
1724
                                      numpy.float32(ExternalRadius),
1725
                                      numpy.float32(Angle),
1726
                                      numpy.int32(Line))
1727
        CLLaunch.wait()
1728
    elif Method=='EachCircle':
1729
        CLLaunch=BlackHoleCL.EachCircle(queue,(int(zImage.shape[0]/2),),
1730
                                        None,zImageCL,fImageCL,
1731
                                        numpy.float32(Mass),
1732
                                        numpy.float32(InternalRadius),
1733
                                        numpy.float32(ExternalRadius),
1734
                                        numpy.float32(Angle),
1735
                                        numpy.int32(Line))
1736
        CLLaunch.wait()
1737
    elif Method=='TrajectoCircle':
1738
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1739
                                        None,TrajectoriesCL,IdLastCL,
1740
                                        numpy.float32(Mass),
1741
                                        numpy.float32(InternalRadius),
1742
                                        numpy.float32(ExternalRadius),
1743
                                        numpy.float32(Angle),
1744
                                        numpy.int32(Line))
1745
        
1746
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1747
                                           int(zImage.shape[0]*4)),None,
1748
                                    TrajectoriesCL,IdLastCL,
1749
                                    zImageCL,fImageCL,
1750
                                    numpy.float32(Mass),
1751
                                    numpy.float32(InternalRadius),
1752
                                    numpy.float32(ExternalRadius),
1753
                                    numpy.float32(Angle),
1754
                                    numpy.int32(Line))
1755
        CLLaunch.wait()
1756
    else:
1757
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1758
                                        None,TrajectoriesCL,IdLastCL,
1759
                                        numpy.float32(Mass),
1760
                                        numpy.float32(InternalRadius),
1761
                                        numpy.float32(ExternalRadius),
1762
                                        numpy.float32(Angle),
1763
                                        numpy.int32(Line))
1764
        
1765
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],zImage.shape[1]),None,
1766
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1767
                                   numpy.uint32(Trajectories.shape[0]),
1768
                                   numpy.float32(Mass),
1769
                                   numpy.float32(InternalRadius),
1770
                                   numpy.float32(ExternalRadius),
1771
                                   numpy.float32(Angle),
1772
                                   numpy.int32(Line))
1773
        CLLaunch.wait()
1774
    
1775
    compute = time.time()-start_time
1776
    
1777
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1778
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1779
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1780
        cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1781
        cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1782
    elapsed = time.time()-start_time
1783
    print("\nCompute Time : %f" % compute)
1784
    print("Elapsed Time : %f\n" % elapsed)
1785

    
1786
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1787
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1788
    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())))
1789
    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())))
1790
    zImageCL.release()
1791
    fImageCL.release()
1792

    
1793
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1794
        if not NoImage:
1795
            AngleStep=4*numpy.pi/TrackPoints
1796
            Angles=numpy.arange(0.,4*numpy.pi,AngleStep)
1797
            Angles.shape=(1,TrackPoints)
1798
            Hostname=gethostname()
1799
            Date=time.strftime("%Y%m%d_%H%M%S")
1800
            ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1801
            
1802
            if TrackSave:
1803
                # numpy.savetxt("TrouNoirTrajectories_%s.csv" % ImageInfo,
1804
                #               numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1805
                #               delimiter=' ', fmt='%.2e')
1806
                numpy.savetxt("TrouNoirTrajectories.csv",
1807
                            numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),delimiter=' ', fmt='%.2e')
1808
        
1809
        TrajectoriesCL.release()
1810
        IdLastCL.release()
1811
        
1812
    return(elapsed)
1813

    
1814
def BlackHoleCUDA(zImage,fImage,InputCL):
1815

    
1816
    kernel_params = {}
1817

    
1818
    Device=InputCL['Device']
1819
    GpuStyle=InputCL['GpuStyle']
1820
    VariableType=InputCL['VariableType']
1821
    Size=InputCL['Size']
1822
    Mass=InputCL['Mass']
1823
    InternalRadius=InputCL['InternalRadius']
1824
    ExternalRadius=InputCL['ExternalRadius']
1825
    Angle=InputCL['Angle']
1826
    Method=InputCL['Method']
1827
    TrackPoints=InputCL['TrackPoints']
1828
    Physics=InputCL['Physics']
1829
    Threads=InputCL['Threads']
1830

    
1831
    PhysicsList=DictionariesAPI()
1832
    
1833
    if InputCL['BlackBody']:
1834
        # Spectrum is Black Body one
1835
        Line=0
1836
    else:
1837
        # Spectrum is Monochromatic Line one
1838
        Line=1
1839

    
1840
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1841
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1842

    
1843
    try:
1844
        # For PyCUDA import
1845
        import pycuda.driver as cuda
1846
        from pycuda.compiler import SourceModule
1847
        
1848
        cuda.init()
1849
        for Id in range(cuda.Device.count()):
1850
            if Id==Device:
1851
                XPU=cuda.Device(Id)
1852
                print("GPU selected %s" % XPU.name())
1853
        print
1854

    
1855
    except ImportError:
1856
        print("Platform does not seem to support CUDA")
1857

    
1858
    Context=XPU.make_context()
1859

    
1860
    try:
1861
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1862
        print("Compilation seems to be OK")
1863
    except:
1864
        print("Compilation seems to break")
1865

    
1866
    EachPixelCU=mod.get_function("EachPixel")
1867
    OriginalCU=mod.get_function("Original")
1868
    EachCircleCU=mod.get_function("EachCircle")
1869
    TrajectoryCU=mod.get_function("Trajectory")
1870
    PixelCU=mod.get_function("Pixel")
1871
    CircleCU=mod.get_function("Circle")
1872
    
1873
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1874
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1875
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1876
    cuda.memcpy_htod(zImageCU, zImage)
1877
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1878
    cuda.memcpy_htod(zImageCU, fImage)
1879
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1880
    cuda.memcpy_htod(IdLastCU, IdLast)
1881

    
1882
    start_time=time.time()
1883

    
1884
    if Method=='EachPixel':
1885
        EachPixelCU(zImageCU,fImageCU,
1886
                    numpy.float32(Mass),
1887
                    numpy.float32(InternalRadius),
1888
                    numpy.float32(ExternalRadius),
1889
                    numpy.float32(Angle),
1890
                    numpy.int32(Line),
1891
                    grid=(int(zImage.shape[0]/Threads),
1892
                          int(zImage.shape[1]/Threads)),
1893
                    block=(Threads,Threads,1))
1894
    elif Method=='EachCircle':
1895
        EachCircleCU(zImageCU,fImageCU,
1896
                   numpy.float32(Mass),
1897
                   numpy.float32(InternalRadius),
1898
                   numpy.float32(ExternalRadius),
1899
                   numpy.float32(Angle),
1900
                   numpy.int32(Line),
1901
                   grid=(int(zImage.shape[0]/Threads/2),1),
1902
                   block=(Threads,1,1))
1903
    elif Method=='Original':
1904
        OriginalCU(zImageCU,fImageCU,
1905
                   numpy.uint32(zImage.shape[0]/2),
1906
                   numpy.float32(Mass),
1907
                   numpy.float32(InternalRadius),
1908
                   numpy.float32(ExternalRadius),
1909
                   numpy.float32(Angle),
1910
                   numpy.int32(Line),
1911
                   grid=(1,1),
1912
                   block=(1,1,1))
1913
    elif Method=='TrajectoCircle':
1914
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1915
                     numpy.float32(Mass),
1916
                     numpy.float32(InternalRadius),
1917
                     numpy.float32(ExternalRadius),
1918
                     numpy.float32(Angle),
1919
                     numpy.int32(Line),
1920
                     grid=(int(Trajectories.shape[0]/Threads),1),
1921
                     block=(Threads,1,1))
1922
        
1923
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1924
                 numpy.float32(Mass),
1925
                 numpy.float32(InternalRadius),
1926
                 numpy.float32(ExternalRadius),
1927
                 numpy.float32(Angle),
1928
                 numpy.int32(Line),
1929
                 grid=(int(Trajectories.shape[0]/Threads),
1930
                       int(zImage.shape[0]*4/Threads)),
1931
                 block=(Threads,Threads,1))
1932
    else:
1933
        # Default method: TrajectoPixel
1934
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1935
                     numpy.float32(Mass),
1936
                     numpy.float32(InternalRadius),
1937
                     numpy.float32(ExternalRadius),
1938
                     numpy.float32(Angle),
1939
                     numpy.int32(Line),
1940
                     grid=(int(Trajectories.shape[0]/Threads),1),
1941
                     block=(Threads,1,1))
1942
        
1943
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1944
                numpy.uint32(Trajectories.shape[0]),
1945
                numpy.float32(Mass),
1946
                numpy.float32(InternalRadius),
1947
                numpy.float32(ExternalRadius),
1948
                numpy.float32(Angle),
1949
                numpy.int32(Line),
1950
                grid=(int(zImage.shape[0]/Threads),
1951
                      int(zImage.shape[1]/Threads),1),
1952
                block=(Threads,Threads,1))
1953

    
1954
    Context.synchronize()
1955

    
1956
    compute = time.time()-start_time
1957

    
1958
    cuda.memcpy_dtoh(zImage,zImageCU)
1959
    cuda.memcpy_dtoh(fImage,fImageCU)
1960
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1961
        cuda.memcpy_dtoh(Trajectories,TrajectoriesCU)
1962
    elapsed = time.time()-start_time
1963
    print("\nCompute Time : %f" % compute)
1964
    print("Elapsed Time : %f\n" % elapsed)
1965

    
1966
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1967
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1968
    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())))
1969
    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())))
1970

    
1971
    
1972
    Context.pop()
1973

    
1974
    Context.detach()
1975

    
1976
    if Method=='TrajectoPixel' or Method=='TrajectoCircle':
1977
        if not NoImage:
1978
            AngleStep=4*numpy.pi/TrackPoints
1979
            Angles=numpy.arange(0.,4*numpy.pi,AngleStep)
1980
            Angles.shape=(1,TrackPoints)
1981
            Hostname=gethostname()
1982
            Date=time.strftime("%Y%m%d_%H%M%S")
1983
            ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1984
            
1985
            # numpy.savetxt("TrouNoirTrajectories_%s.csv" % ImageInfo,
1986
            #               numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1987
            #               delimiter=' ', fmt='%.2e')
1988
            numpy.savetxt("TrouNoirTrajectories.csv",
1989
                          numpy.transpose(numpy.concatenate((Angles,Trajectories),axis=0)),
1990
                          delimiter=' ', fmt='%.2e')
1991

    
1992
    
1993
    return(elapsed)
1994

    
1995
if __name__=='__main__':
1996

    
1997
    # Default device: first one!
1998
    Device=0
1999
    # Default implementation: OpenCL, most versatile!
2000
    GpuStyle = 'OpenCL'
2001
    Mass = 1.
2002
    # Internal Radius 3 times de Schwarzschild Radius
2003
    InternalRadius=6.*Mass
2004
    #
2005
    ExternalRadius=12.
2006
    #
2007
    # Angle with normal to disc 10 degrees
2008
    Angle = numpy.pi/180.*(90.-10.)
2009
    # Radiation of disc : BlackBody or Monochromatic
2010
    BlackBody = False
2011
    # Size of image
2012
    Size=1024
2013
    # Variable Type
2014
    VariableType='FP32'
2015
    # ?
2016
    q=-2
2017
    # Method of resolution
2018
    Method='TrajectoPixel'
2019
    # Colors for output image
2020
    Colors='Greyscale'
2021
    # Physics
2022
    Physics='Einstein'
2023
    # No output as image
2024
    NoImage = False
2025
    # Threads in CUDA
2026
    Threads = 32
2027
    # Trackpoints of trajectories
2028
    TrackPoints=2048
2029
    # Tracksave of trajectories
2030
    TrackSave=False
2031
    
2032
    HowToUse='%s -h [Help] -b [BlackBodyEmission] -j [TrackSave] -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>'
2033

    
2034
    try:
2035
        opts, args = getopt.getopt(sys.argv[1:],"hbnjs:m:i:x: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="])
2036
    except getopt.GetoptError:
2037
        print(HowToUse % sys.argv[0])
2038
        sys.exit(2)
2039

    
2040
    # List of Devices
2041
    Devices=[]
2042
    Alu={}
2043
        
2044
    for opt, arg in opts:
2045
        if opt == '-h':
2046
            print(HowToUse % sys.argv[0])
2047
            
2048
            print("\nInformations about devices detected under OpenCL API:")
2049
            # For PyOpenCL import
2050
            try:
2051
                import pyopencl as cl
2052
                Id=0
2053
                for platform in cl.get_platforms():
2054
                    for device in platform.get_devices():
2055
                        #deviceType=cl.device_type.to_string(device.type)
2056
                        deviceType="xPU"
2057
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
2058
                        Id=Id+1
2059
                        
2060
            except:
2061
                print("Your platform does not seem to support OpenCL")
2062
                
2063
            print("\nInformations about devices detected under CUDA API:")
2064
                # For PyCUDA import
2065
            try:
2066
                import pycuda.driver as cuda
2067
                cuda.init()
2068
                for Id in range(cuda.Device.count()):
2069
                    device=cuda.Device(Id)
2070
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
2071
                print
2072
            except:
2073
                print("Your platform does not seem to support CUDA")
2074
        
2075
            sys.exit()
2076
        
2077
        elif opt in ("-d", "--device"):
2078
#            Devices.append(int(arg))
2079
            Device=int(arg)
2080
        elif opt in ("-g", "--gpustyle"):
2081
            GpuStyle = arg
2082
        elif opt in ("-v", "--variabletype"):
2083
            VariableType = arg
2084
        elif opt in ("-s", "--size"):
2085
            Size = int(arg)
2086
        elif opt in ("-k", "--trackpoints"):
2087
            TrackPoints = int(arg)
2088
        elif opt in ("-m", "--mass"):
2089
            Mass = float(arg)
2090
        elif opt in ("-i", "--internal"):
2091
            InternalRadius = float(arg)
2092
        elif opt in ("-e", "--external"):
2093
            ExternalRadius = float(arg)
2094
        elif opt in ("-a", "--angle"):
2095
            Angle = numpy.pi/180.*(90.-float(arg))
2096
        elif opt in ("-b", "--blackbody"):
2097
            BlackBody = True
2098
        elif opt in ("-j", "--tracksave"):
2099
            TrackSave = True
2100
        elif opt in ("-n", "--noimage"):
2101
            NoImage = True
2102
        elif opt in ("-o", "--method"):
2103
            Method = arg
2104
        elif opt in ("-t", "--threads"):
2105
            Threads = int(arg)
2106
        elif opt in ("-c", "--colors"):
2107
            Colors = arg
2108
        elif opt in ("-p", "--physics"):
2109
            Physics = arg
2110

    
2111
    print("Device Identification selected : %s" % Device)
2112
    print("GpuStyle used : %s" % GpuStyle)
2113
    print("VariableType : %s" % VariableType)
2114
    print("Size : %i" % Size)
2115
    print("Mass : %f" % Mass)
2116
    print("Internal Radius : %f" % InternalRadius)
2117
    print("External Radius : %f" % ExternalRadius)
2118
    print("Angle with normal of (in radians) : %f" % Angle)
2119
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
2120
    print("Method of resolution : %s" % Method)
2121
    print("Colors for output images : %s" % Colors)
2122
    print("Physics used for Trajectories : %s" % Physics)
2123
    print("Trackpoints of Trajectories : %i" % TrackPoints)
2124
    print("Tracksave of Trajectories : %i" % TrackSave)
2125

    
2126
    if GpuStyle=='CUDA':
2127
        print("\nSelection of CUDA device") 
2128
        try:
2129
            # For PyCUDA import
2130
            import pycuda.driver as cuda
2131
            
2132
            cuda.init()
2133
            for Id in range(cuda.Device.count()):
2134
                device=cuda.Device(Id)
2135
                print("Device #%i of type GPU : %s" % (Id,device.name()))
2136
                if Id in Devices:
2137
                    Alu[Id]='GPU'
2138
            
2139
        except ImportError:
2140
            print("Platform does not seem to support CUDA")
2141

    
2142
    if GpuStyle=='OpenCL':
2143
        print("\nSelection of OpenCL device") 
2144
        try:
2145
            # For PyOpenCL import
2146
            import pyopencl as cl
2147
            Id=0
2148
            for platform in cl.get_platforms():
2149
                for device in platform.get_devices():
2150
                    #deviceType=cl.device_type.to_string(device.type)
2151
                    deviceType="xPU"
2152
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
2153

    
2154
                    if Id in Devices:
2155
                    # Set the Alu as detected Device Type
2156
                        Alu[Id]=deviceType
2157
                    Id=Id+1
2158
        except ImportError:
2159
            print("Platform does not seem to support OpenCL")
2160

    
2161

    
2162
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2163
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2164

    
2165
    InputCL={}
2166
    InputCL['Device']=Device
2167
    InputCL['GpuStyle']=GpuStyle
2168
    InputCL['VariableType']=VariableType
2169
    InputCL['Size']=Size
2170
    InputCL['Mass']=Mass
2171
    InputCL['InternalRadius']=InternalRadius
2172
    InputCL['ExternalRadius']=ExternalRadius
2173
    InputCL['Angle']=Angle
2174
    InputCL['BlackBody']=BlackBody
2175
    InputCL['Method']=Method
2176
    InputCL['TrackPoints']=TrackPoints
2177
    InputCL['Physics']=Physics
2178
    InputCL['Threads']=Threads
2179
    InputCL['NoImage']=NoImage
2180
    InputCL['TrackSave']=TrackSave
2181

    
2182
    if GpuStyle=='OpenCL':
2183
        duration=BlackHoleCL(zImage,fImage,InputCL)
2184
    else:
2185
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
2186
        
2187
    Hostname=gethostname()
2188
    Date=time.strftime("%Y%m%d_%H%M%S")
2189
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2190
    
2191
    
2192
    if not NoImage:
2193
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
2194
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)