Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 227

Historique | Voir | Annoter | Télécharger (49,45 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, disable VECTORIZER with
29
# export CL_CONFIG_USE_VECTORIZER=0
30

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

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

    
44
#
45
# Blank space below to simplify debugging on OpenCL code
46
#
47

    
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
BlobOpenCL= """
102

103
#define PI (float)3.14159265359
104
#define nbr 256
105

106
#define EINSTEIN 0
107
#define NEWTON 1
108

109
#ifdef SETTRACKPOINTS
110
#define TRACKPOINTS SETTRACKPOINTS
111
#else
112
#define TRACKPOINTS 2048
113
#endif
114

115
float atanp(float x,float y)
116
{
117
  float angle;
118

119
  angle=atan2(y,x);
120

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

126
  return angle;
127
}
128

129
float f(float v)
130
{
131
  return v;
132
}
133

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

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

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

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

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

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

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

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

188
float spectre_cn(float rf32,float b32,float db32,
189
                 float h32,float r32,float m32,float bss32)
190
{
191
  
192
#define MYFLOAT float
193

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

202
  MYFLOAT flx;
203
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
204
  int fi,posfreq;
205

206
#define planck 6.62e-34
207
#define k 1.38e-23
208
#define c2 9.e16
209
#define temp 3.e7
210
#define m_point 1.
211

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

216
  MYFLOAT v=1.-3./r;
217

218
  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 ));
219

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

222
  flux_int=0.;
223
  flx=0.;
224

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

240
          flx+=flux_int;
241
        }
242
    }
243

244
  return((float)(flx));
245
}
246

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

254
  rf=decalage_spectral(r,b,phi,tho,m);
255

256
  if (raie==0)
257
    {
258
      bss=1.e19;
259
      flx=spectre_cn(rf,b,db,h,r,m,bss);
260
    }
261
  else
262
    { 
263
      bss=2.;
264
      flx=spectre(rf,q,b,db,h,r,m,bss);
265
    }
266
  
267
  *zp=1./rf;
268
  *fp=flx;
269

270
}
271

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

282
   // Perform trajectory for each pixel, exit on hit
283

284
  float m,rs,ri,re,tho;
285
  int q,raie;
286

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

295
  float bmx,db,b,h;
296
  float rp0,rpp,rps;
297
  float phi,phd;
298
  int nh;
299
  float zp=0.,fp=0.;
300

301
  // Autosize for image
302
  bmx=1.25*re;
303

304
  h=4.e0f*PI/(float)TRACKPOINTS;
305

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

313
  float up,vp,pp,us,vs,ps;
314

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

320
  up=0.;
321
  vp=1.;
322
      
323
  pp=0.;
324
  nh=0;
325

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

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

331
  int ExitOnImpact=0;
332

333
  do
334
  {
335
     nh++;
336
     pp=ps;
337
     up=us;
338
     vp=vs;
339
     rpp=rps;
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));
345

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

355
  barrier(CLK_GLOBAL_MEM_FENCE);
356

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

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

373
   // Perform trajectory for each pixel
374

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

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

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

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

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

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

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

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

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

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

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

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

422
      if (ni<IdLast[bi])
423
      {
424
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
425
      }
426
      else
427
      {
428
        r=Trajectories[bi*TRACKPOINTS+ni];
429
      }
430
           
431
      if ((r<=re)&&(r>=ri))
432
      {
433
        ExitOnImpact=1;
434
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
435
      }
436
              
437
      HalfLap++;
438
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
439

440
  }
441

442
  barrier(CLK_GLOBAL_MEM_FENCE);
443

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

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

463
   // Perform trajectory for each pixel
464

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

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

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

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

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

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

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

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

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

502
     if (ni<IdLast[bi])
503
     {
504
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
505
     }
506
     else
507
     {
508
        r=Trajectories[bi*TRACKPOINTS+ni];
509
     }
510
           
511
     if ((r<=re)&&(r>=ri))
512
     {
513
        ExitOnImpact=1;
514
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
515
     }
516
              
517
     HalfLap++;
518
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
519
  
520
  zImage[yi+2*bmaxi*xi]=zp;
521
  fImage[yi+2*bmaxi*xi]=fp;  
522

523
  barrier(CLK_GLOBAL_MEM_FENCE);
524

525
}
526

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

537
  // Perform trajectory for each pixel
538

539
  float m,rs,re;
540

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

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

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

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

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

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

559
  up=0.;
560
  vp=1.;
561
      
562
  pp=0.;
563
  nh=0;
564

565
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
566
  
567
  // b versus us
568
  float bvus=fabs(b/us);
569
  float bvus0=bvus;
570
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
571

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

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

584
  IdLast[bi]=nh;
585

586
  barrier(CLK_GLOBAL_MEM_FENCE);
587
 
588
}
589

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

600
  float Trajectory[TRACKPOINTS];
601

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

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

613
  float bmx,db,b,h;
614
  int nh;
615

616
  // Autosize for image
617
  bmx=1.25e0f*re;
618

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

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

626
  float up,vp,pp,us,vs,ps;
627

628
  up=0.;
629
  vp=1.;
630
      
631
  pp=0.;
632
  nh=0;
633

634
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
635
  
636
  // b versus us
637
  float bvus=fabs(b/us);
638
  float bvus0=bvus;
639
  Trajectory[nh]=bvus;
640

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

651
  } while ((bvus>=rs)&&(bvus<=bvus0));
652

653
  for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
654
     Trajectory[i]=0.e0f;
655
  }
656

657
  int imx=(int)(16*bi);
658

659
  for (int i=0;i<imx;i++)
660
  {
661
     float zp=0.,fp=0.;
662
     float phi=2.*PI/(float)imx*(float)i;
663
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
664
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
665
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
666

667
     int HalfLap=0,ExitOnImpact=0,ni;
668
     float php,nr,r;
669

670
     do
671
     {
672
        php=phd+(float)HalfLap*PI;
673
        nr=php/h;
674
        ni=(int)nr;
675

676
        if (ni<nh)
677
        {
678
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
679
        }
680
        else
681
        {
682
           r=Trajectory[ni];
683
        }
684
           
685
        if ((r<=re)&&(r>=ri))
686
        {
687
           ExitOnImpact=1;
688
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
689
        }
690
              
691
        HalfLap++;
692

693
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
694
 
695
     zImage[yi+2*bmaxi*xi]=zp;
696
     fImage[yi+2*bmaxi*xi]=fp;
697

698
  }
699

700
  barrier(CLK_GLOBAL_MEM_FENCE);
701
 
702
}
703

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

712
   float Trajectory[TRACKPOINTS];
713

714
   // Perform trajectory for each pixel
715

716
   float m,rs,ri,re,tho;
717
   int raie,q;
718

719
   m=Mass;
720
   rs=2.*m;
721
   ri=InternalRadius;
722
   re=ExternalRadius;
723
   tho=Angle;
724
   q=-2;
725
   raie=Line;
726

727
   float bmx,db,b,h;
728
   int nh;
729

730
   // Autosize for image
731
   bmx=1.25e0f*re;
732

733
   // Angular step of integration
734
   h=4.e0f*PI/(float)TRACKPOINTS;
735

736
   // Integer Impact Parameter ID 
737
   for (int bi=0;bi<bmaxi;bi++)
738
   {
739
      // impact parameter
740
      b=(float)bi/(float)bmaxi*bmx;
741
      db=bmx/(2.e0f*(float)bmaxi);
742

743
      float up,vp,pp,us,vs,ps;
744

745
      up=0.;
746
      vp=1.;
747
      
748
      pp=0.;
749
      nh=0;
750

751
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
752
  
753
      // b versus us
754
      float bvus=fabs(b/us);
755
      float bvus0=bvus;
756
      Trajectory[nh]=bvus;
757

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

768
      } while ((bvus>=rs)&&(bvus<=bvus0));
769

770
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
771
         Trajectory[i]=0.e0f;
772
      }
773

774
      int imx=(int)(16*bi);
775

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

784
         int HalfLap=0,ExitOnImpact=0,ni;
785
         float php,nr,r;
786

787
         do
788
         {
789
            php=phd+(float)HalfLap*PI;
790
            nr=php/h;
791
            ni=(int)nr;
792

793
            if (ni<nh)
794
            {
795
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
796
            }
797
            else
798
            {
799
               r=Trajectory[ni];
800
            }
801
           
802
            if ((r<=re)&&(r>=ri))
803
            {
804
               ExitOnImpact=1;
805
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
806
            }
807
               
808
            HalfLap++;
809
 
810
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
811
  
812
         zImage[yi+2*bmaxi*xi]=zp;
813
         fImage[yi+2*bmaxi*xi]=fp;
814
  
815
      }
816

817
   } 
818
   
819
   barrier(CLK_GLOBAL_MEM_FENCE);
820
 
821
}
822
"""
823

    
824
def KernelCodeCuda():
825
    BlobCUDA= """
826

827
#define PI (float)3.14159265359
828
#define nbr 256
829

830
#define EINSTEIN 0
831
#define NEWTON 1
832

833
#ifdef SETTRACKPOINTS
834
#define TRACKPOINTS SETTRACKPOINTS
835
#else
836
#define TRACKPOINTS
837
#endif
838

839
__device__ float nothing(float x)
840
{
841
  return(x);
842
}
843

844
__device__ float atanp(float x,float y)
845
{
846
  float angle;
847

848
  angle=atan2(y,x);
849

850
  if (angle<0.e0f)
851
    {
852
      angle+=(float)2.e0f*PI;
853
    }
854

855
  return(angle);
856
}
857

858
__device__ float f(float v)
859
{
860
  return(v);
861
}
862

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

875
__device__ void calcul(float *us,float *vs,float up,float vp,
876
                       float h,float m,float b)
877
{
878
  float c0,c1,c2,c3,d0,d1,d2,d3;
879

880
  c0=h*f(vp);
881
  c1=h*f(vp+c0/2.);
882
  c2=h*f(vp+c1/2.);
883
  c3=h*f(vp+c2);
884
  d0=h*g(up,m,b);
885
  d1=h*g(up+d0/2.,m,b);
886
  d2=h*g(up+d1/2.,m,b);
887
  d3=h*g(up+d2,m,b);
888

889
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
890
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
891
}
892

893
__device__ void rungekutta(float *ps,float *us,float *vs,
894
                           float pp,float up,float vp,
895
                           float h,float m,float b)
896
{
897
  calcul(us,vs,up,vp,h,m,b);
898
  *ps=pp+h;
899
}
900

901
__device__ float decalage_spectral(float r,float b,float phi,
902
                                   float tho,float m)
903
{
904
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
905
}
906

907
__device__ float spectre(float rf,int q,float b,float db,
908
                         float h,float r,float m,float bss)
909
{
910
  float flx;
911

912
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
913
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
914
  return(flx);
915
}
916

917
__device__ float spectre_cn(float rf32,float b32,float db32,
918
                            float h32,float r32,float m32,float bss32)
919
{
920
  
921
#define MYFLOAT float
922

923
  MYFLOAT rf=(MYFLOAT)(rf32);
924
  MYFLOAT b=(MYFLOAT)(b32);
925
  MYFLOAT db=(MYFLOAT)(db32);
926
  MYFLOAT h=(MYFLOAT)(h32);
927
  MYFLOAT r=(MYFLOAT)(r32);
928
  MYFLOAT m=(MYFLOAT)(m32);
929
  MYFLOAT bss=(MYFLOAT)(bss32);
930

931
  MYFLOAT flx;
932
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
933
  int fi,posfreq;
934

935
#define planck 6.62e-34
936
#define k 1.38e-23
937
#define c2 9.e16
938
#define temp 3.e7
939
#define m_point 1.
940

941
#define lplanck (log(6.62)-34.*log(10.))
942
#define lk (log(1.38)-23.*log(10.))
943
#define lc2 (log(9.)+16.*log(10.))
944

945
  MYFLOAT v=1.-3./r;
946

947
  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 ));
948

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

951
  flux_int=0.;
952
  flx=0.;
953

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

969
          flx+=flux_int;
970
        }
971
    }
972

973
  return((float)(flx));
974
}
975

976
__device__ void impact(float phi,float r,float b,float tho,float m,
977
                       float *zp,float *fp,
978
                       int q,float db,
979
                       float h,int raie)
980
{
981
  float flx,rf,bss;
982

983
  rf=decalage_spectral(r,b,phi,tho,m);
984

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

999
}
1000

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

1011
   // Perform trajectory for each pixel, exit on hit
1012

1013
  float m,rs,ri,re,tho;
1014
  int q,raie;
1015

1016
  m=Mass;
1017
  rs=2.*m;
1018
  ri=InternalRadius;
1019
  re=ExternalRadius;
1020
  tho=Angle;
1021
  q=-2;
1022
  raie=Line;
1023

1024
  float d,bmx,db,b,h;
1025
  float rp0,rpp,rps;
1026
  float phi,thi,phd,php,nr,r;
1027
  int nh;
1028
  float zp,fp;
1029

1030
  // Autosize for image
1031
  bmx=1.25*re;
1032
  b=0.;
1033

1034
  h=4.e0f*PI/(float)TRACKPOINTS;
1035

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

1043
  float up,vp,pp,us,vs,ps;
1044

1045
  // impact parameter
1046
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
1047
  // step of impact parameter;
1048
//  db=bmx/(float)(sizex/2);
1049
  db=bmx/(float)(sizex);
1050

1051
  up=0.;
1052
  vp=1.;
1053
      
1054
  pp=0.;
1055
  nh=0;
1056

1057
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1058
    
1059
  rps=fabs(b/us);
1060
  rp0=rps;
1061

1062
  int ExitOnImpact=0;
1063

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

1075
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
1076

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

1086
  __syncthreads();
1087

1088
  zImage[yi+sizex*xi]=(float)zp;
1089
  fImage[yi+sizex*xi]=(float)fp;
1090
}
1091

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

1104
  // Perform trajectory for each pixel
1105

1106
  float m,rs,ri,re,tho;
1107
  int q,raie;
1108

1109
  m=Mass;
1110
  rs=2.*m;
1111
  ri=InternalRadius;
1112
  re=ExternalRadius;
1113
  tho=Angle;
1114
  q=-2;
1115
  raie=Line;
1116

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

1124
  // Angular step of integration
1125
  h=4.e0f*PI/(float)TRACKPOINTS;
1126

1127
  // Step of Impact Parameter
1128
  db=bmx/(2.e0*(float)ImpactParameter);
1129

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

1137
  // Real Impact Parameter
1138
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
1139

1140
  // Integer Impact Parameter
1141
  uint bi=(uint)sqrt(x*x+y*y);
1142

1143
  int HalfLap=0,ExitOnImpact=0,ni;
1144

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

1153
      if (ni<IdLast[bi])
1154
      {
1155
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1156
      }
1157
      else
1158
      {
1159
        r=Trajectories[bi*TRACKPOINTS+ni];
1160
      }
1161
           
1162
      if ((r<=re)&&(r>=ri))
1163
      {
1164
        ExitOnImpact=1;
1165
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1166
      }
1167
              
1168
      HalfLap++;
1169
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1170

1171
  }
1172

1173
  zImage[yi+sizex*xi]=zp;
1174
  fImage[yi+sizex*xi]=fp;
1175
}
1176

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

1192
   // Perform trajectory for each pixel
1193

1194
  float m,rs,ri,re,tho;
1195
  int q,raie;
1196

1197
  m=Mass;
1198
  rs=2.*m;
1199
  ri=InternalRadius;
1200
  re=ExternalRadius;
1201
  tho=Angle;
1202
  raie=Line;
1203

1204
  float bmx,db,b,h;
1205
  float phi,thi,phd;
1206
  int nh;
1207
  float zp=0,fp=0;
1208

1209
  // Autosize for image
1210
  bmx=1.25*re;
1211

1212
  // Angular step of integration
1213
  h=4.e0f*PI/(float)TRACKPOINTS;
1214

1215
  // impact parameter
1216
  b=(float)bi/(float)bmaxi*bmx;
1217
  db=bmx/(2.e0*(float)bmaxi);
1218

1219
  phi=2.*PI/(float)imx*(float)i;
1220
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1221
  int yi=(int)((float)bi*sin(phi))+bmaxi;
1222
  int xi=(int)((float)bi*cos(phi))+bmaxi;
1223

1224
  int HalfLap=0,ExitOnImpact=0,ni;
1225
  float php,nr,r;
1226

1227
  do
1228
  {
1229
     php=phd+(float)HalfLap*PI;
1230
     nr=php/h;
1231
     ni=(int)nr;
1232

1233
     if (ni<IdLast[bi])
1234
     {
1235
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1236
     }
1237
     else
1238
     {
1239
        r=Trajectories[bi*TRACKPOINTS+ni];
1240
     }
1241
           
1242
     if ((r<=re)&&(r>=ri))
1243
     {
1244
        ExitOnImpact=1;
1245
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1246
     }
1247
              
1248
     HalfLap++;
1249
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1250
  
1251
  zImage[yi+2*bmaxi*xi]=zp;
1252
  fImage[yi+2*bmaxi*xi]=fp;  
1253

1254
}
1255

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

1266
  // Perform trajectory for each pixel
1267

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

1271
  m=Mass;
1272
  rs=2.*m;
1273
  ri=InternalRadius;
1274
  re=ExternalRadius;
1275
  tho=Angle;
1276
  q=-2;
1277
  raie=Line;
1278

1279
  float d,bmx,db,b,h;
1280
  float phi,thi,phd,php,nr,r;
1281
  int nh;
1282
  float zp,fp;
1283

1284
  // Autosize for image
1285
  bmx=1.25*re;
1286

1287
  // Angular step of integration
1288
  h=4.e0f*PI/(float)TRACKPOINTS;
1289

1290
  // impact parameter
1291
  b=(float)bi/(float)bmaxi*bmx;
1292

1293
  float up,vp,pp,us,vs,ps;
1294

1295
  up=0.;
1296
  vp=1.;
1297
      
1298
  pp=0.;
1299
  nh=0;
1300

1301
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1302
  
1303
  // b versus us
1304
  float bvus=fabs(b/us);
1305
  float bvus0=bvus;
1306
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
1307

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

1318
  } while ((bvus>=rs)&&(bvus<=bvus0));
1319

1320
  IdLast[bi]=nh;
1321

1322
}
1323

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

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

1335
  float Trajectory[2048];
1336

1337
  // Perform trajectory for each pixel
1338

1339
  float m,rs,ri,re,tho;
1340
  int raie,q;
1341

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

1350
  float bmx,db,b,h;
1351
  int nh;
1352

1353
  // Autosize for image
1354
  bmx=1.25e0f*re;
1355

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

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

1363
  float up,vp,pp,us,vs,ps;
1364

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

1371
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1372
  
1373
  // b versus us
1374
  float bvus=fabs(b/us);
1375
  float bvus0=bvus;
1376
  Trajectory[nh]=bvus;
1377

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

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

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

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

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

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

1409
        if (ni<nh)
1410
        {
1411
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1412
        }
1413
        else
1414
        {
1415
           r=Trajectory[ni];
1416
        }
1417
           
1418
        if ((r<=re)&&(r>=ri))
1419
        {
1420
           ExitOnImpact=1;
1421
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1422
        }
1423
              
1424
        HalfLap++;
1425

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

1428
   __syncthreads();
1429

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

1433
  }
1434
 
1435
}
1436

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

1445
   float Trajectory[TRACKPOINTS];
1446

1447
   // Perform trajectory for each pixel
1448

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

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

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

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

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

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

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

1478
      up=0.;
1479
      vp=1.;
1480
      
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
    
1561
#     # Normalize value as 8bits Integer
1562
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1563
#     image = Image.fromarray(SigmaInt)
1564
#     image.save("%s.jpg" % prefix)
1565

    
1566
def ImageOutput(sigma,prefix,Colors):
1567
    import matplotlib.pyplot as plt
1568
    start_time=time.time() 
1569
    if Colors == 'Red2Yellow':
1570
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1571
    else:
1572
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1573
    save_time = time.time()-start_time
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

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

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

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

    
1617
    if HasXPU==False:
1618
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1619
        sys.exit()                
1620
        
1621
    ctx = cl.Context([XPU])
1622
    queue = cl.CommandQueue(ctx,
1623
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1624

    
1625
    
1626
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1627

    
1628
    BuildOptions="-cl-mad-enable -DPHYSICS=%i -DSETTRACKPOINTS=%i " % (PhysicsList[Physics],InputCL['TrackPoints'])
1629

    
1630
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1631
    
1632
    # Je recupere les flag possibles pour les buffers
1633
    mf = cl.mem_flags
1634

    
1635
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
1636
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1637
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1638
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1639
    
1640
    start_time=time.time() 
1641

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

    
1719
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1720
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1721
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1722
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1723
    zImageCL.release()
1724
    fImageCL.release()
1725

    
1726
    TrajectoriesCL.release()
1727
    IdLastCL.release()
1728

    
1729
    return(elapsed)
1730

    
1731
def BlackHoleCUDA(zImage,fImage,InputCL):
1732

    
1733
    kernel_params = {}
1734

    
1735
    print(InputCL)
1736

    
1737
    Device=InputCL['Device']
1738
    GpuStyle=InputCL['GpuStyle']
1739
    VariableType=InputCL['VariableType']
1740
    Size=InputCL['Size']
1741
    Mass=InputCL['Mass']
1742
    InternalRadius=InputCL['InternalRadius']
1743
    ExternalRadius=InputCL['ExternalRadius']
1744
    Angle=InputCL['Angle']
1745
    Method=InputCL['Method']
1746
    TrackPoints=InputCL['TrackPoints']
1747
    Physics=InputCL['Physics']
1748

    
1749
    PhysicsList=DictionariesAPI()
1750
    
1751
    if InputCL['BlackBody']:
1752
        # Spectrum is Black Body one
1753
        Line=0
1754
    else:
1755
        # Spectrum is Monochromatic Line one
1756
        Line=1
1757

    
1758
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1759
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1760

    
1761
    try:
1762
        # For PyCUDA import
1763
        import pycuda.driver as cuda
1764
        from pycuda.compiler import SourceModule
1765
        
1766
        cuda.init()
1767
        for Id in range(cuda.Device.count()):
1768
            if Id==Device:
1769
                XPU=cuda.Device(Id)
1770
                print("GPU selected %s" % XPU.name())
1771
        print
1772

    
1773
    except ImportError:
1774
        print("Platform does not seem to support CUDA")
1775

    
1776
    Context=XPU.make_context()
1777

    
1778
    try:
1779
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1780
        print("Compilation seems to be OK")
1781
    except:
1782
        print("Compilation seems to break")
1783

    
1784
    EachPixelCU=mod.get_function("EachPixel")
1785
    OriginalCU=mod.get_function("Original")
1786
    EachCircleCU=mod.get_function("EachCircle")
1787
    TrajectoryCU=mod.get_function("Trajectory")
1788
    PixelCU=mod.get_function("Pixel")
1789
    CircleCU=mod.get_function("Circle")
1790
    
1791
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1792
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1793
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1794
    cuda.memcpy_htod(zImageCU, zImage)
1795
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1796
    cuda.memcpy_htod(zImageCU, fImage)
1797
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1798
    cuda.memcpy_htod(IdLastCU, IdLast)
1799

    
1800
    start_time=time.time() 
1801

    
1802
    if Method=='EachPixel':
1803
        EachPixelCU(zImageCU,fImageCU,
1804
                    numpy.float32(Mass),
1805
                    numpy.float32(InternalRadius),
1806
                    numpy.float32(ExternalRadius),
1807
                    numpy.float32(Angle),
1808
                    numpy.int32(Line),
1809
                    grid=(zImage.shape[0]/32,zImage.shape[1]/32),
1810
                    block=(32,32,1))
1811
    elif Method=='EachCircle':
1812
        EachCircleCU(zImageCU,fImageCU,
1813
                   numpy.float32(Mass),
1814
                   numpy.float32(InternalRadius),
1815
                   numpy.float32(ExternalRadius),
1816
                   numpy.float32(Angle),
1817
                   numpy.int32(Line),
1818
                   grid=(zImage.shape[0]/32/2,1),
1819
                   block=(32,1,1))
1820
    elif Method=='Original':
1821
        OriginalCU(zImageCU,fImageCU,
1822
                   numpy.uint32(zImage.shape[0]/2),
1823
                   numpy.float32(Mass),
1824
                   numpy.float32(InternalRadius),
1825
                   numpy.float32(ExternalRadius),
1826
                   numpy.float32(Angle),
1827
                   numpy.int32(Line),
1828
                   grid=(1,1),
1829
                   block=(1,1,1))
1830
    elif Method=='TrajectoCircle':
1831
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1832
                     numpy.float32(Mass),
1833
                     numpy.float32(InternalRadius),
1834
                     numpy.float32(ExternalRadius),
1835
                     numpy.float32(Angle),
1836
                     numpy.int32(Line),
1837
                     grid=(Trajectories.shape[0]/32,1),
1838
                     block=(32,1,1))
1839
        
1840
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1841
                 numpy.float32(Mass),
1842
                 numpy.float32(InternalRadius),
1843
                 numpy.float32(ExternalRadius),
1844
                 numpy.float32(Angle),
1845
                 numpy.int32(Line),
1846
                 grid=(Trajectories.shape[0]/32,zImage.shape[0]*4/32),
1847
                 block=(32,32,1))
1848
    else:
1849
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1850
                     numpy.float32(Mass),
1851
                     numpy.float32(InternalRadius),
1852
                     numpy.float32(ExternalRadius),
1853
                     numpy.float32(Angle),
1854
                     numpy.int32(Line),
1855
                     grid=(Trajectories.shape[0]/32,1),
1856
                     block=(32,1,1))
1857
        
1858
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1859
                numpy.uint32(Trajectories.shape[0]),
1860
                numpy.float32(Mass),
1861
                numpy.float32(InternalRadius),
1862
                numpy.float32(ExternalRadius),
1863
                numpy.float32(Angle),
1864
                numpy.int32(Line),
1865
                grid=(zImage.shape[0]/32,zImage.shape[1]/32,1),
1866
                block=(32,32,1))
1867

    
1868
    Context.synchronize()
1869

    
1870
    compute = time.time()-start_time
1871

    
1872
    cuda.memcpy_dtoh(zImage,zImageCU)
1873
    cuda.memcpy_dtoh(fImage,fImageCU)
1874
    elapsed = time.time()-start_time
1875
    print("\nCompute Time : %f" % compute)
1876
    print("Elapsed Time : %f\n" % elapsed)
1877

    
1878
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1879
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1880
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1881
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1882

    
1883
    
1884
    Context.pop()
1885

    
1886
    Context.detach()
1887

    
1888
    return(elapsed)
1889

    
1890
if __name__=='__main__':
1891
        
1892
    GpuStyle = 'OpenCL'
1893
    Mass = 1.
1894
    # Internal Radius 3 times de Schwarzschild Radius
1895
    InternalRadius=6.*Mass
1896
    #
1897
    ExternalRadius=12.
1898
    #
1899
    # Angle with normal to disc 10 degrees
1900
    Angle = numpy.pi/180.*(90.-10.)
1901
    # Radiation of disc : BlackBody or Monochromatic
1902
    BlackBody = False
1903
    # Size of image
1904
    Size=256
1905
    # Variable Type
1906
    VariableType='FP32'
1907
    # ?
1908
    q=-2
1909
    # Method of resolution
1910
    Method='TrajectoPixel'
1911
    # Colors for output image
1912
    Colors='Greyscale'
1913
    # Physics
1914
    Physics='Einstein'
1915
    # No output as image
1916
    NoImage = False
1917
    
1918
    HowToUse='%s -h [Help] -b [BlackBodyEmission] -n [NoImage] -p <Einstein/Newton> -s <SizeInPixels> -m <Mass> -i <DiscInternalRadius> -x <DiscExternalRadius> -a <AngleAboveDisc> -d <DeviceId> -c <Greyscale/Red2Yellow> -g <CUDA/OpenCL> -t <EachPixel/TrajectoCircle/TrajectoPixel/EachCircle/Original> -v <FP32/FP64>'
1919

    
1920
    try:
1921
        opts, args = getopt.getopt(sys.argv[1:],"hbns:m:i:x:a:d:g:v:t:c:p:",["blackbody","noimage","camera","size=","mass=","internal=","external=","angle=","device=","gpustyle=","variabletype=","method=","colors=","physics="])
1922
    except getopt.GetoptError:
1923
        print(HowToUse % sys.argv[0])
1924
        sys.exit(2)
1925

    
1926
    # List of Devices
1927
    Devices=[]
1928
    Alu={}
1929
        
1930
    for opt, arg in opts:
1931
        if opt == '-h':
1932
            print(HowToUse % sys.argv[0])
1933
            
1934
            print("\nInformations about devices detected under OpenCL API:")
1935
            # For PyOpenCL import
1936
            try:
1937
                import pyopencl as cl
1938
                Id=0
1939
                for platform in cl.get_platforms():
1940
                    for device in platform.get_devices():
1941
                        #deviceType=cl.device_type.to_string(device.type)
1942
                        deviceType="xPU"
1943
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
1944
                        Id=Id+1
1945
                        
1946
            except:
1947
                print("Your platform does not seem to support OpenCL")
1948
                
1949
            print("\nInformations about devices detected under CUDA API:")
1950
                # For PyCUDA import
1951
            try:
1952
                import pycuda.driver as cuda
1953
                cuda.init()
1954
                for Id in range(cuda.Device.count()):
1955
                    device=cuda.Device(Id)
1956
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
1957
                print
1958
            except:
1959
                print("Your platform does not seem to support CUDA")
1960
        
1961
            sys.exit()
1962
        
1963
        elif opt in ("-d", "--device"):
1964
#            Devices.append(int(arg))
1965
            Device=int(arg)
1966
        elif opt in ("-g", "--gpustyle"):
1967
            GpuStyle = arg
1968
        elif opt in ("-v", "--variabletype"):
1969
            VariableType = arg
1970
        elif opt in ("-s", "--size"):
1971
            Size = int(arg)
1972
        elif opt in ("-m", "--mass"):
1973
            Mass = float(arg)
1974
        elif opt in ("-i", "--internal"):
1975
            InternalRadius = float(arg)
1976
        elif opt in ("-e", "--external"):
1977
            ExternalRadius = float(arg)
1978
        elif opt in ("-a", "--angle"):
1979
            Angle = numpy.pi/180.*(90.-float(arg))
1980
        elif opt in ("-b", "--blackbody"):
1981
            BlackBody = True
1982
        elif opt in ("-n", "--noimage"):
1983
            NoImage = True
1984
        elif opt in ("-t", "--method"):
1985
            Method = arg
1986
        elif opt in ("-c", "--colors"):
1987
            Colors = arg
1988
        elif opt in ("-p", "--physics"):
1989
            Physics = arg
1990

    
1991
    print("Device Identification selected : %s" % Device)
1992
    print("GpuStyle used : %s" % GpuStyle)
1993
    print("VariableType : %s" % VariableType)
1994
    print("Size : %i" % Size)
1995
    print("Mass : %f" % Mass)
1996
    print("Internal Radius : %f" % InternalRadius)
1997
    print("External Radius : %f" % ExternalRadius)
1998
    print("Angle with normal of (in radians) : %f" % Angle)
1999
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
2000
    print("Method of resolution : %s" % Method)
2001
    print("Colors for output images : %s" % Colors)
2002
    print("Physics used for Trajectories : %s" % Physics)
2003

    
2004
    if GpuStyle=='CUDA':
2005
        print("\nSelection of CUDA device") 
2006
        try:
2007
            # For PyCUDA import
2008
            import pycuda.driver as cuda
2009
            
2010
            cuda.init()
2011
            for Id in range(cuda.Device.count()):
2012
                device=cuda.Device(Id)
2013
                print("Device #%i of type GPU : %s" % (Id,device.name()))
2014
                if Id in Devices:
2015
                    Alu[Id]='GPU'
2016
            
2017
        except ImportError:
2018
            print("Platform does not seem to support CUDA")
2019

    
2020
    if GpuStyle=='OpenCL':
2021
        print("\nSelection of OpenCL device") 
2022
        try:
2023
            # For PyOpenCL import
2024
            import pyopencl as cl
2025
            Id=0
2026
            for platform in cl.get_platforms():
2027
                for device in platform.get_devices():
2028
                    #deviceType=cl.device_type.to_string(device.type)
2029
                    deviceType="xPU"
2030
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
2031

    
2032
                    if Id in Devices:
2033
                    # Set the Alu as detected Device Type
2034
                        Alu[Id]=deviceType
2035
                    Id=Id+1
2036
        except ImportError:
2037
            print("Platform does not seem to support OpenCL")
2038

    
2039
#    print(Devices,Alu)
2040

    
2041
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
2042
    TrackPoints=2048
2043
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2044
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
2045

    
2046
    InputCL={}
2047
    InputCL['Device']=Device
2048
    InputCL['GpuStyle']=GpuStyle
2049
    InputCL['VariableType']=VariableType
2050
    InputCL['Size']=Size
2051
    InputCL['Mass']=Mass
2052
    InputCL['InternalRadius']=InternalRadius
2053
    InputCL['ExternalRadius']=ExternalRadius
2054
    InputCL['Angle']=Angle
2055
    InputCL['BlackBody']=BlackBody
2056
    InputCL['Method']=Method
2057
    InputCL['TrackPoints']=TrackPoints
2058
    InputCL['Physics']=Physics
2059

    
2060
    if GpuStyle=='OpenCL':
2061
        duration=BlackHoleCL(zImage,fImage,InputCL)
2062
    else:
2063
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
2064
        
2065
    Hostname=gethostname()
2066
    Date=time.strftime("%Y%m%d_%H%M%S")
2067
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2068
    
2069
    
2070
    if not NoImage:
2071
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
2072
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)