Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 225

Historique | Voir | Annoter | Télécharger (49,36 ko)

1 221 equemene
#!/usr/bin/env python
2 199 equemene
#
3 199 equemene
# TrouNoir model using PyOpenCL
4 199 equemene
#
5 204 equemene
# CC BY-NC-SA 2019 : <emmanuel.quemener@ens-lyon.fr>
6 199 equemene
#
7 199 equemene
# Thanks to Andreas Klockner for PyOpenCL:
8 199 equemene
# http://mathema.tician.de/software/pyopencl
9 199 equemene
#
10 204 equemene
# Original code programmed in Fortran 77 in mars 1994
11 204 equemene
# for Practical Work of Numerical Simulation
12 204 equemene
# DEA (old Master2) in astrophysics and spatial techniques in Meudon
13 204 equemene
# by Herve Aussel & Emmanuel Quemener
14 204 equemene
#
15 204 equemene
# Conversion in C done by Emmanuel Quemener in august 1997
16 204 equemene
# GPUfication in OpenCL under Python in july 2019
17 221 equemene
# GPUfication in CUDA under Python in august 2019
18 204 equemene
#
19 204 equemene
# Thanks to :
20 204 equemene
#
21 204 equemene
# - Herve Aussel for his part of code of black body spectrum
22 204 equemene
# - Didier Pelat for his help to perform this work
23 204 equemene
# - Jean-Pierre Luminet for his article published in 1979
24 204 equemene
# - Numerical Recipies for Runge Kutta recipies
25 204 equemene
# - Luc Blanchet for his disponibility about my questions in General Relativity
26 204 equemene
# - Pierre Lena for his passion about science and vulgarisation
27 199 equemene
28 199 equemene
import pyopencl as cl
29 199 equemene
import numpy
30 199 equemene
import time,string
31 199 equemene
from numpy.random import randint as nprnd
32 199 equemene
import sys
33 199 equemene
import getopt
34 199 equemene
import matplotlib.pyplot as plt
35 211 equemene
from socket import gethostname
36 199 equemene
37 211 equemene
def DictionariesAPI():
38 211 equemene
    PhysicsList={'Einstein':0,'Newton':1}
39 211 equemene
    return(PhysicsList)
40 204 equemene
41 211 equemene
BlobOpenCL= """
42 204 equemene

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

46 211 equemene
#define EINSTEIN 0
47 211 equemene
#define NEWTON 1
48 211 equemene

49 224 equemene
#ifdef SETTRACKPOINTS
50 224 equemene
#define TRACKPOINTS SETTRACKPOINTS
51 224 equemene
#else
52 217 equemene
#define TRACKPOINTS 2048
53 224 equemene
#endif
54 217 equemene

55 199 equemene
float atanp(float x,float y)
56 199 equemene
{
57 199 equemene
  float angle;
58 199 equemene

59 199 equemene
  angle=atan2(y,x);
60 199 equemene

61 204 equemene
  if (angle<0.e0f)
62 199 equemene
    {
63 199 equemene
      angle+=(float)2.e0f*PI;
64 199 equemene
    }
65 199 equemene

66 199 equemene
  return angle;
67 199 equemene
}
68 199 equemene

69 199 equemene
float f(float v)
70 199 equemene
{
71 199 equemene
  return v;
72 199 equemene
}
73 199 equemene

74 211 equemene
#if PHYSICS == NEWTON
75 199 equemene
float g(float u,float m,float b)
76 199 equemene
{
77 211 equemene
  return (-u);
78 211 equemene
}
79 211 equemene
#else
80 211 equemene
float g(float u,float m,float b)
81 211 equemene
{
82 204 equemene
  return (3.e0f*m/b*pow(u,2)-u);
83 199 equemene
}
84 211 equemene
#endif
85 199 equemene

86 199 equemene
void calcul(float *us,float *vs,float up,float vp,
87 199 equemene
            float h,float m,float b)
88 199 equemene
{
89 199 equemene
  float c0,c1,c2,c3,d0,d1,d2,d3;
90 199 equemene

91 199 equemene
  c0=h*f(vp);
92 199 equemene
  c1=h*f(vp+c0/2.);
93 199 equemene
  c2=h*f(vp+c1/2.);
94 199 equemene
  c3=h*f(vp+c2);
95 199 equemene
  d0=h*g(up,m,b);
96 199 equemene
  d1=h*g(up+d0/2.,m,b);
97 199 equemene
  d2=h*g(up+d1/2.,m,b);
98 199 equemene
  d3=h*g(up+d2,m,b);
99 199 equemene

100 199 equemene
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
101 199 equemene
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
102 199 equemene
}
103 199 equemene

104 199 equemene
void rungekutta(float *ps,float *us,float *vs,
105 199 equemene
                float pp,float up,float vp,
106 199 equemene
                float h,float m,float b)
107 199 equemene
{
108 199 equemene
  calcul(us,vs,up,vp,h,m,b);
109 199 equemene
  *ps=pp+h;
110 199 equemene
}
111 199 equemene

112 199 equemene
float decalage_spectral(float r,float b,float phi,
113 199 equemene
                         float tho,float m)
114 199 equemene
{
115 199 equemene
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
116 199 equemene
}
117 199 equemene

118 199 equemene
float spectre(float rf,int q,float b,float db,
119 199 equemene
             float h,float r,float m,float bss)
120 199 equemene
{
121 199 equemene
  float flx;
122 199 equemene

123 221 equemene
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
124 221 equemene
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
125 199 equemene
  return(flx);
126 199 equemene
}
127 199 equemene

128 209 equemene
float spectre_cn(float rf32,float b32,float db32,
129 209 equemene
                 float h32,float r32,float m32,float bss32)
130 199 equemene
{
131 209 equemene

132 213 equemene
#define MYFLOAT float
133 209 equemene

134 209 equemene
  MYFLOAT rf=(MYFLOAT)(rf32);
135 209 equemene
  MYFLOAT b=(MYFLOAT)(b32);
136 209 equemene
  MYFLOAT db=(MYFLOAT)(db32);
137 209 equemene
  MYFLOAT h=(MYFLOAT)(h32);
138 209 equemene
  MYFLOAT r=(MYFLOAT)(r32);
139 209 equemene
  MYFLOAT m=(MYFLOAT)(m32);
140 209 equemene
  MYFLOAT bss=(MYFLOAT)(bss32);
141 209 equemene

142 209 equemene
  MYFLOAT flx;
143 209 equemene
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
144 199 equemene
  int fi,posfreq;
145 199 equemene

146 209 equemene
#define planck 6.62e-34
147 209 equemene
#define k 1.38e-23
148 209 equemene
#define c2 9.e16
149 209 equemene
#define temp 3.e7
150 209 equemene
#define m_point 1.
151 199 equemene

152 209 equemene
#define lplanck (log(6.62)-34.*log(10.))
153 209 equemene
#define lk (log(1.38)-23.*log(10.))
154 209 equemene
#define lc2 (log(9.)+16.*log(10.))
155 199 equemene

156 209 equemene
  MYFLOAT v=1.-3./r;
157 199 equemene

158 209 equemene
  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 ));
159 199 equemene

160 209 equemene
  temp_em=temp*sqrt(m)*exp(0.25*log(m_point)-0.75*log(r)-0.125*log(v)+0.25*log(fabs(qu)));
161 209 equemene

162 209 equemene
  flux_int=0.;
163 209 equemene
  flx=0.;
164 209 equemene

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

180 199 equemene
          flx+=flux_int;
181 199 equemene
        }
182 199 equemene
    }
183 209 equemene

184 209 equemene
  return((float)(flx));
185 199 equemene
}
186 199 equemene

187 199 equemene
void impact(float phi,float r,float b,float tho,float m,
188 199 equemene
            float *zp,float *fp,
189 199 equemene
            int q,float db,
190 204 equemene
            float h,int raie)
191 199 equemene
{
192 204 equemene
  float flx,rf,bss;
193 199 equemene

194 199 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
195 199 equemene

196 199 equemene
  if (raie==0)
197 199 equemene
    {
198 209 equemene
      bss=1.e19;
199 209 equemene
      flx=spectre_cn(rf,b,db,h,r,m,bss);
200 209 equemene
    }
201 209 equemene
  else
202 209 equemene
    {
203 204 equemene
      bss=2.;
204 199 equemene
      flx=spectre(rf,q,b,db,h,r,m,bss);
205 199 equemene
    }
206 199 equemene

207 199 equemene
  *zp=1./rf;
208 199 equemene
  *fp=flx;
209 199 equemene

210 199 equemene
}
211 199 equemene

212 204 equemene
__kernel void EachPixel(__global float *zImage,__global float *fImage,
213 204 equemene
                        float Mass,float InternalRadius,
214 204 equemene
                        float ExternalRadius,float Angle,
215 209 equemene
                        int Line)
216 199 equemene
{
217 199 equemene
   uint xi=(uint)get_global_id(0);
218 199 equemene
   uint yi=(uint)get_global_id(1);
219 199 equemene
   uint sizex=(uint)get_global_size(0);
220 199 equemene
   uint sizey=(uint)get_global_size(1);
221 199 equemene

222 204 equemene
   // Perform trajectory for each pixel, exit on hit
223 199 equemene

224 217 equemene
  private float m,rs,ri,re,tho;
225 217 equemene
  private int q,raie;
226 199 equemene

227 204 equemene
  m=Mass;
228 204 equemene
  rs=2.*m;
229 204 equemene
  ri=InternalRadius;
230 204 equemene
  re=ExternalRadius;
231 204 equemene
  tho=Angle;
232 204 equemene
  q=-2;
233 209 equemene
  raie=Line;
234 204 equemene

235 217 equemene
  private float d,bmx,db,b,h;
236 218 equemene
  private float rp0,rpp,rps;
237 217 equemene
  private float phi,thi,phd,php,nr,r;
238 217 equemene
  private int nh;
239 217 equemene
  private float zp,fp;
240 199 equemene

241 199 equemene
  // Autosize for image
242 199 equemene
  bmx=1.25*re;
243 199 equemene
  b=0.;
244 199 equemene

245 217 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
246 201 equemene

247 199 equemene
  // set origin as center of image
248 199 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
249 201 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
250 199 equemene
  // angle extracted from cylindric symmetry
251 199 equemene
  phi=atanp(x,y);
252 199 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
253 199 equemene

254 204 equemene
  float up,vp,pp,us,vs,ps;
255 204 equemene

256 204 equemene
  // impact parameter
257 204 equemene
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
258 204 equemene
  // step of impact parameter;
259 209 equemene
  db=bmx/(float)(sizex);
260 204 equemene

261 209 equemene
  up=0.;
262 209 equemene
  vp=1.;
263 199 equemene

264 199 equemene
  pp=0.;
265 199 equemene
  nh=0;
266 199 equemene

267 199 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
268 199 equemene

269 218 equemene
  rps=fabs(b/us);
270 218 equemene
  rp0=rps;
271 199 equemene

272 204 equemene
  int ExitOnImpact=0;
273 199 equemene

274 199 equemene
  do
275 199 equemene
  {
276 199 equemene
     nh++;
277 199 equemene
     pp=ps;
278 199 equemene
     up=us;
279 199 equemene
     vp=vs;
280 218 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
281 218 equemene
     rpp=rps;
282 218 equemene
     rps=fabs(b/us);
283 218 equemene
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;
284 199 equemene

285 218 equemene
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
286 199 equemene

287 199 equemene
  if (ExitOnImpact==1) {
288 218 equemene
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
289 199 equemene
  }
290 199 equemene
  else
291 199 equemene
  {
292 199 equemene
     zp=0.;
293 201 equemene
     fp=0.;
294 199 equemene
  }
295 199 equemene

296 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
297 204 equemene

298 201 equemene
  zImage[yi+sizex*xi]=(float)zp;
299 204 equemene
  fImage[yi+sizex*xi]=(float)fp;
300 204 equemene
}
301 199 equemene

302 204 equemene
__kernel void Pixel(__global float *zImage,__global float *fImage,
303 204 equemene
                    __global float *Trajectories,__global int *IdLast,
304 224 equemene
                    uint ImpactParameter,
305 204 equemene
                    float Mass,float InternalRadius,
306 204 equemene
                    float ExternalRadius,float Angle,
307 209 equemene
                    int Line)
308 204 equemene
{
309 204 equemene
   uint xi=(uint)get_global_id(0);
310 204 equemene
   uint yi=(uint)get_global_id(1);
311 204 equemene
   uint sizex=(uint)get_global_size(0);
312 204 equemene
   uint sizey=(uint)get_global_size(1);
313 204 equemene

314 204 equemene
   // Perform trajectory for each pixel
315 204 equemene

316 209 equemene
  float m,rs,ri,re,tho;
317 209 equemene
  int q,raie;
318 204 equemene

319 204 equemene
  m=Mass;
320 204 equemene
  rs=2.*m;
321 204 equemene
  ri=InternalRadius;
322 204 equemene
  re=ExternalRadius;
323 204 equemene
  tho=Angle;
324 204 equemene
  q=-2;
325 209 equemene
  raie=Line;
326 204 equemene

327 204 equemene
  float d,bmx,db,b,h;
328 204 equemene
  float phi,thi,phd,php,nr,r;
329 204 equemene
  int nh;
330 204 equemene
  float zp=0,fp=0;
331 204 equemene

332 209 equemene
  // Autosize for image, 25% greater than external radius
333 204 equemene
  bmx=1.25*re;
334 204 equemene

335 204 equemene
  // Angular step of integration
336 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
337 204 equemene

338 204 equemene
  // Step of Impact Parameter
339 209 equemene
  db=bmx/(2.e0*(float)ImpactParameter);
340 204 equemene

341 204 equemene
  // set origin as center of image
342 204 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
343 204 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
344 204 equemene

345 204 equemene
  // angle extracted from cylindric symmetry
346 204 equemene
  phi=atanp(x,y);
347 204 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
348 204 equemene

349 204 equemene
  // Real Impact Parameter
350 204 equemene
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
351 204 equemene

352 204 equemene
  // Integer Impact Parameter
353 204 equemene
  uint bi=(uint)sqrt(x*x+y*y);
354 204 equemene

355 204 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
356 204 equemene

357 204 equemene
  if (bi<ImpactParameter)
358 204 equemene
  {
359 204 equemene
    do
360 204 equemene
    {
361 204 equemene
      php=phd+(float)HalfLap*PI;
362 204 equemene
      nr=php/h;
363 204 equemene
      ni=(int)nr;
364 204 equemene

365 204 equemene
      if (ni<IdLast[bi])
366 204 equemene
      {
367 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
368 204 equemene
      }
369 204 equemene
      else
370 204 equemene
      {
371 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
372 204 equemene
      }
373 204 equemene

374 204 equemene
      if ((r<=re)&&(r>=ri))
375 204 equemene
      {
376 204 equemene
        ExitOnImpact=1;
377 204 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
378 204 equemene
      }
379 204 equemene

380 204 equemene
      HalfLap++;
381 204 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
382 204 equemene

383 204 equemene
  }
384 204 equemene

385 201 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
386 204 equemene

387 204 equemene
  zImage[yi+sizex*xi]=zp;
388 204 equemene
  fImage[yi+sizex*xi]=fp;
389 204 equemene
}
390 204 equemene

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

406 204 equemene
   // Perform trajectory for each pixel
407 204 equemene

408 209 equemene
  float m,rs,ri,re,tho;
409 209 equemene
  int q,raie;
410 204 equemene

411 204 equemene
  m=Mass;
412 204 equemene
  rs=2.*m;
413 204 equemene
  ri=InternalRadius;
414 204 equemene
  re=ExternalRadius;
415 204 equemene
  tho=Angle;
416 209 equemene
  raie=Line;
417 204 equemene

418 204 equemene
  float bmx,db,b,h;
419 204 equemene
  float phi,thi,phd;
420 204 equemene
  int nh;
421 204 equemene
  float zp=0,fp=0;
422 204 equemene

423 204 equemene
  // Autosize for image
424 204 equemene
  bmx=1.25*re;
425 204 equemene

426 204 equemene
  // Angular step of integration
427 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
428 204 equemene

429 204 equemene
  // impact parameter
430 204 equemene
  b=(float)bi/(float)bmaxi*bmx;
431 224 equemene
  db=bmx/(2.e0f*(float)bmaxi);
432 204 equemene

433 204 equemene
  phi=2.*PI/(float)imx*(float)i;
434 204 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
435 204 equemene
  int yi=(int)((float)bi*sin(phi))+bmaxi;
436 204 equemene
  int xi=(int)((float)bi*cos(phi))+bmaxi;
437 204 equemene

438 204 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
439 204 equemene
  float php,nr,r;
440 204 equemene

441 204 equemene
  do
442 204 equemene
  {
443 204 equemene
     php=phd+(float)HalfLap*PI;
444 204 equemene
     nr=php/h;
445 204 equemene
     ni=(int)nr;
446 204 equemene

447 204 equemene
     if (ni<IdLast[bi])
448 204 equemene
     {
449 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
450 204 equemene
     }
451 204 equemene
     else
452 204 equemene
     {
453 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
454 204 equemene
     }
455 204 equemene

456 204 equemene
     if ((r<=re)&&(r>=ri))
457 204 equemene
     {
458 204 equemene
        ExitOnImpact=1;
459 204 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
460 204 equemene
     }
461 204 equemene

462 204 equemene
     HalfLap++;
463 204 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
464 204 equemene

465 204 equemene
  zImage[yi+2*bmaxi*xi]=zp;
466 204 equemene
  fImage[yi+2*bmaxi*xi]=fp;
467 204 equemene

468 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
469 204 equemene

470 204 equemene
}
471 204 equemene

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

482 204 equemene
  // Perform trajectory for each pixel
483 204 equemene

484 224 equemene
  float m,rs,re;
485 224 equemene

486 224 equemene
  m=Mass;
487 224 equemene
  rs=2.*m;
488 224 equemene
  re=ExternalRadius;
489 224 equemene

490 224 equemene
  float bmx,b,h;
491 224 equemene
  int nh;
492 224 equemene

493 224 equemene
  // Autosize for image
494 224 equemene
  bmx=1.25*re;
495 224 equemene

496 224 equemene
  // Angular step of integration
497 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
498 224 equemene

499 224 equemene
  // impact parameter
500 224 equemene
  b=(float)bi/(float)bmaxi*bmx;
501 224 equemene

502 224 equemene
  float up,vp,pp,us,vs,ps;
503 224 equemene

504 224 equemene
  up=0.;
505 224 equemene
  vp=1.;
506 224 equemene

507 224 equemene
  pp=0.;
508 224 equemene
  nh=0;
509 224 equemene

510 224 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
511 224 equemene

512 224 equemene
  // b versus us
513 224 equemene
  float bvus=fabs(b/us);
514 224 equemene
  float bvus0=bvus;
515 224 equemene
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
516 224 equemene

517 224 equemene
  do
518 224 equemene
  {
519 224 equemene
     nh++;
520 224 equemene
     pp=ps;
521 224 equemene
     up=us;
522 224 equemene
     vp=vs;
523 224 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
524 224 equemene
     bvus=fabs(b/us);
525 224 equemene
     Trajectories[bi*TRACKPOINTS+nh]=bvus;
526 224 equemene

527 224 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
528 224 equemene

529 224 equemene
  IdLast[bi]=nh;
530 224 equemene

531 224 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
532 224 equemene

533 224 equemene
}
534 224 equemene

535 225 equemene
__kernel void EachCircle(__global float *zImage,__global float *fImage,
536 225 equemene
                         float Mass,float InternalRadius,
537 225 equemene
                         float ExternalRadius,float Angle,
538 225 equemene
                         int Line)
539 224 equemene
{
540 224 equemene
  // Integer Impact Parameter ID
541 224 equemene
  int bi=get_global_id(0);
542 224 equemene
  // Integer Impact Parameter Size (half of image)
543 224 equemene
  int bmaxi=get_global_size(0);
544 224 equemene

545 225 equemene
  float Trajectory[TRACKPOINTS];
546 224 equemene

547 209 equemene
  float m,rs,ri,re,tho;
548 209 equemene
  int raie,q;
549 204 equemene

550 204 equemene
  m=Mass;
551 204 equemene
  rs=2.*m;
552 204 equemene
  ri=InternalRadius;
553 204 equemene
  re=ExternalRadius;
554 204 equemene
  tho=Angle;
555 204 equemene
  q=-2;
556 209 equemene
  raie=Line;
557 204 equemene

558 224 equemene
  float bmx,db,b,h;
559 204 equemene
  int nh;
560 204 equemene

561 204 equemene
  // Autosize for image
562 224 equemene
  bmx=1.25e0f*re;
563 204 equemene

564 204 equemene
  // Angular step of integration
565 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
566 204 equemene

567 204 equemene
  // impact parameter
568 204 equemene
  b=(float)bi/(float)bmaxi*bmx;
569 224 equemene
  db=bmx/(2.e0f*(float)bmaxi);
570 204 equemene

571 204 equemene
  float up,vp,pp,us,vs,ps;
572 204 equemene

573 209 equemene
  up=0.;
574 209 equemene
  vp=1.;
575 204 equemene

576 204 equemene
  pp=0.;
577 204 equemene
  nh=0;
578 204 equemene

579 204 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
580 204 equemene

581 204 equemene
  // b versus us
582 204 equemene
  float bvus=fabs(b/us);
583 204 equemene
  float bvus0=bvus;
584 225 equemene
  Trajectory[nh]=bvus;
585 204 equemene

586 204 equemene
  do
587 204 equemene
  {
588 204 equemene
     nh++;
589 204 equemene
     pp=ps;
590 204 equemene
     up=us;
591 204 equemene
     vp=vs;
592 204 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
593 204 equemene
     bvus=fabs(b/us);
594 225 equemene
     Trajectory[nh]=bvus;
595 204 equemene

596 204 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
597 204 equemene

598 225 equemene
  for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
599 225 equemene
     Trajectory[i]=0.e0f;
600 225 equemene
  }
601 225 equemene

602 224 equemene
  int imx=(int)(16*bi);
603 204 equemene

604 224 equemene
  for (int i=0;i<imx;i++)
605 224 equemene
  {
606 225 equemene
     float zp=0.,fp=0.;
607 224 equemene
     float phi=2.*PI/(float)imx*(float)i;
608 224 equemene
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
609 224 equemene
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
610 224 equemene
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
611 224 equemene

612 224 equemene
     int HalfLap=0,ExitOnImpact=0,ni;
613 224 equemene
     float php,nr,r;
614 224 equemene

615 224 equemene
     do
616 224 equemene
     {
617 224 equemene
        php=phd+(float)HalfLap*PI;
618 224 equemene
        nr=php/h;
619 224 equemene
        ni=(int)nr;
620 224 equemene

621 224 equemene
        if (ni<nh)
622 224 equemene
        {
623 225 equemene
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
624 224 equemene
        }
625 224 equemene
        else
626 224 equemene
        {
627 225 equemene
           r=Trajectory[ni];
628 224 equemene
        }
629 224 equemene

630 224 equemene
        if ((r<=re)&&(r>=ri))
631 224 equemene
        {
632 224 equemene
           ExitOnImpact=1;
633 224 equemene
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
634 224 equemene
        }
635 224 equemene

636 224 equemene
        HalfLap++;
637 224 equemene

638 224 equemene
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
639 224 equemene

640 224 equemene
     zImage[yi+2*bmaxi*xi]=zp;
641 224 equemene
     fImage[yi+2*bmaxi*xi]=fp;
642 224 equemene

643 224 equemene
  }
644 224 equemene

645 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
646 199 equemene

647 199 equemene
}
648 225 equemene

649 225 equemene
__kernel void Original(__global float *zImage,__global float *fImage,
650 225 equemene
                       uint Size,float Mass,float InternalRadius,
651 225 equemene
                       float ExternalRadius,float Angle,
652 225 equemene
                       int Line)
653 225 equemene
{
654 225 equemene
   // Integer Impact Parameter Size (half of image)
655 225 equemene
   uint bmaxi=(uint)Size;
656 225 equemene

657 225 equemene
   float Trajectory[TRACKPOINTS];
658 225 equemene

659 225 equemene
   // Perform trajectory for each pixel
660 225 equemene

661 225 equemene
   float m,rs,ri,re,tho;
662 225 equemene
   int raie,q;
663 225 equemene

664 225 equemene
   m=Mass;
665 225 equemene
   rs=2.*m;
666 225 equemene
   ri=InternalRadius;
667 225 equemene
   re=ExternalRadius;
668 225 equemene
   tho=Angle;
669 225 equemene
   q=-2;
670 225 equemene
   raie=Line;
671 225 equemene

672 225 equemene
   float bmx,db,b,h;
673 225 equemene
   int nh;
674 225 equemene

675 225 equemene
   // Autosize for image
676 225 equemene
   bmx=1.25e0f*re;
677 225 equemene

678 225 equemene
   // Angular step of integration
679 225 equemene
   h=4.e0f*PI/(float)TRACKPOINTS;
680 225 equemene

681 225 equemene
   // Integer Impact Parameter ID
682 225 equemene
   for (int bi=0;bi<bmaxi;bi++)
683 225 equemene
   {
684 225 equemene
      // impact parameter
685 225 equemene
      b=(float)bi/(float)bmaxi*bmx;
686 225 equemene
      db=bmx/(2.e0f*(float)bmaxi);
687 225 equemene

688 225 equemene
      float up,vp,pp,us,vs,ps;
689 225 equemene

690 225 equemene
      up=0.;
691 225 equemene
      vp=1.;
692 225 equemene

693 225 equemene
      pp=0.;
694 225 equemene
      nh=0;
695 225 equemene

696 225 equemene
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
697 225 equemene

698 225 equemene
      // b versus us
699 225 equemene
      float bvus=fabs(b/us);
700 225 equemene
      float bvus0=bvus;
701 225 equemene
      Trajectory[nh]=bvus;
702 225 equemene

703 225 equemene
      do
704 225 equemene
      {
705 225 equemene
         nh++;
706 225 equemene
         pp=ps;
707 225 equemene
         up=us;
708 225 equemene
         vp=vs;
709 225 equemene
         rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
710 225 equemene
         bvus=fabs(b/us);
711 225 equemene
         Trajectory[nh]=bvus;
712 225 equemene

713 225 equemene
      } while ((bvus>=rs)&&(bvus<=bvus0));
714 225 equemene

715 225 equemene
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
716 225 equemene
         Trajectory[i]=0.e0f;
717 225 equemene
      }
718 225 equemene

719 225 equemene
      int imx=(int)(16*bi);
720 225 equemene

721 225 equemene
      for (int i=0;i<imx;i++)
722 225 equemene
      {
723 225 equemene
         float zp=0,fp=0;
724 225 equemene
         float phi=2.e0f*PI/(float)imx*(float)i;
725 225 equemene
         float phd=atanp(cos(phi)*sin(tho),cos(tho));
726 225 equemene
         uint yi=(uint)((float)bi*sin(phi)+bmaxi);
727 225 equemene
         uint xi=(uint)((float)bi*cos(phi)+bmaxi);
728 225 equemene

729 225 equemene
         int HalfLap=0,ExitOnImpact=0,ni;
730 225 equemene
         float php,nr,r;
731 225 equemene

732 225 equemene
         do
733 225 equemene
         {
734 225 equemene
            php=phd+(float)HalfLap*PI;
735 225 equemene
            nr=php/h;
736 225 equemene
            ni=(int)nr;
737 225 equemene

738 225 equemene
            if (ni<nh)
739 225 equemene
            {
740 225 equemene
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
741 225 equemene
            }
742 225 equemene
            else
743 225 equemene
            {
744 225 equemene
               r=Trajectory[ni];
745 225 equemene
            }
746 225 equemene

747 225 equemene
            if ((r<=re)&&(r>=ri))
748 225 equemene
            {
749 225 equemene
               ExitOnImpact=1;
750 225 equemene
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
751 225 equemene
            }
752 225 equemene

753 225 equemene
            HalfLap++;
754 225 equemene

755 225 equemene
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
756 225 equemene

757 225 equemene
         zImage[yi+2*bmaxi*xi]=zp;
758 225 equemene
         fImage[yi+2*bmaxi*xi]=fp;
759 225 equemene

760 225 equemene
      }
761 225 equemene

762 225 equemene
   }
763 225 equemene

764 225 equemene
   barrier(CLK_GLOBAL_MEM_FENCE);
765 225 equemene

766 225 equemene
}
767 211 equemene
"""
768 199 equemene
769 217 equemene
def KernelCodeCuda():
770 217 equemene
    BlobCUDA= """
771 217 equemene

772 217 equemene
#define PI (float)3.14159265359
773 217 equemene
#define nbr 256
774 217 equemene

775 217 equemene
#define EINSTEIN 0
776 217 equemene
#define NEWTON 1
777 217 equemene

778 224 equemene
#ifdef SETTRACKPOINTS
779 224 equemene
#define TRACKPOINTS SETTRACKPOINTS
780 224 equemene
#else
781 224 equemene
#define TRACKPOINTS
782 224 equemene
#endif
783 217 equemene

784 217 equemene
__device__ float nothing(float x)
785 217 equemene
{
786 217 equemene
  return(x);
787 217 equemene
}
788 217 equemene

789 217 equemene
__device__ float atanp(float x,float y)
790 217 equemene
{
791 217 equemene
  float angle;
792 217 equemene

793 217 equemene
  angle=atan2(y,x);
794 217 equemene

795 217 equemene
  if (angle<0.e0f)
796 217 equemene
    {
797 217 equemene
      angle+=(float)2.e0f*PI;
798 217 equemene
    }
799 217 equemene

800 217 equemene
  return(angle);
801 217 equemene
}
802 217 equemene

803 217 equemene
__device__ float f(float v)
804 217 equemene
{
805 217 equemene
  return(v);
806 217 equemene
}
807 217 equemene

808 217 equemene
#if PHYSICS == NEWTON
809 217 equemene
__device__ float g(float u,float m,float b)
810 217 equemene
{
811 217 equemene
  return (-u);
812 217 equemene
}
813 217 equemene
#else
814 217 equemene
__device__ float g(float u,float m,float b)
815 217 equemene
{
816 217 equemene
  return (3.e0f*m/b*pow(u,2)-u);
817 217 equemene
}
818 217 equemene
#endif
819 217 equemene

820 217 equemene
__device__ void calcul(float *us,float *vs,float up,float vp,
821 217 equemene
                       float h,float m,float b)
822 217 equemene
{
823 217 equemene
  float c0,c1,c2,c3,d0,d1,d2,d3;
824 217 equemene

825 217 equemene
  c0=h*f(vp);
826 217 equemene
  c1=h*f(vp+c0/2.);
827 217 equemene
  c2=h*f(vp+c1/2.);
828 217 equemene
  c3=h*f(vp+c2);
829 217 equemene
  d0=h*g(up,m,b);
830 217 equemene
  d1=h*g(up+d0/2.,m,b);
831 217 equemene
  d2=h*g(up+d1/2.,m,b);
832 217 equemene
  d3=h*g(up+d2,m,b);
833 217 equemene

834 217 equemene
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
835 217 equemene
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
836 217 equemene
}
837 217 equemene

838 217 equemene
__device__ void rungekutta(float *ps,float *us,float *vs,
839 217 equemene
                           float pp,float up,float vp,
840 217 equemene
                           float h,float m,float b)
841 217 equemene
{
842 217 equemene
  calcul(us,vs,up,vp,h,m,b);
843 217 equemene
  *ps=pp+h;
844 217 equemene
}
845 217 equemene

846 217 equemene
__device__ float decalage_spectral(float r,float b,float phi,
847 217 equemene
                                   float tho,float m)
848 217 equemene
{
849 217 equemene
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
850 217 equemene
}
851 217 equemene

852 217 equemene
__device__ float spectre(float rf,int q,float b,float db,
853 217 equemene
                         float h,float r,float m,float bss)
854 217 equemene
{
855 217 equemene
  float flx;
856 217 equemene

857 221 equemene
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
858 221 equemene
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
859 217 equemene
  return(flx);
860 217 equemene
}
861 217 equemene

862 217 equemene
__device__ float spectre_cn(float rf32,float b32,float db32,
863 217 equemene
                            float h32,float r32,float m32,float bss32)
864 217 equemene
{
865 217 equemene

866 217 equemene
#define MYFLOAT float
867 217 equemene

868 217 equemene
  MYFLOAT rf=(MYFLOAT)(rf32);
869 217 equemene
  MYFLOAT b=(MYFLOAT)(b32);
870 217 equemene
  MYFLOAT db=(MYFLOAT)(db32);
871 217 equemene
  MYFLOAT h=(MYFLOAT)(h32);
872 217 equemene
  MYFLOAT r=(MYFLOAT)(r32);
873 217 equemene
  MYFLOAT m=(MYFLOAT)(m32);
874 217 equemene
  MYFLOAT bss=(MYFLOAT)(bss32);
875 217 equemene

876 217 equemene
  MYFLOAT flx;
877 217 equemene
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
878 217 equemene
  int fi,posfreq;
879 217 equemene

880 217 equemene
#define planck 6.62e-34
881 217 equemene
#define k 1.38e-23
882 217 equemene
#define c2 9.e16
883 217 equemene
#define temp 3.e7
884 217 equemene
#define m_point 1.
885 217 equemene

886 217 equemene
#define lplanck (log(6.62)-34.*log(10.))
887 217 equemene
#define lk (log(1.38)-23.*log(10.))
888 217 equemene
#define lc2 (log(9.)+16.*log(10.))
889 217 equemene

890 217 equemene
  MYFLOAT v=1.-3./r;
891 217 equemene

892 217 equemene
  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 ));
893 217 equemene

894 217 equemene
  temp_em=temp*sqrt(m)*exp(0.25*log(m_point)-0.75*log(r)-0.125*log(v)+0.25*log(fabs(qu)));
895 217 equemene

896 217 equemene
  flux_int=0.;
897 217 equemene
  flx=0.;
898 217 equemene

899 217 equemene
  for (fi=0;fi<nbr;fi++)
900 217 equemene
    {
901 217 equemene
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
902 217 equemene
      nu_rec=nu_em*rf;
903 217 equemene
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
904 217 equemene
      if ((posfreq>0)&&(posfreq<nbr))
905 217 equemene
        {
906 217 equemene
          // Initial version
907 217 equemene
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
908 217 equemene
          // Version with log used
909 217 equemene
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
910 217 equemene
          // flux_int*=pow(rf,3)*b*db*h;
911 217 equemene
          //flux_int*=exp(3.*log(rf))*b*db*h;
912 217 equemene
          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;
913 217 equemene

914 217 equemene
          flx+=flux_int;
915 217 equemene
        }
916 217 equemene
    }
917 217 equemene

918 217 equemene
  return((float)(flx));
919 217 equemene
}
920 217 equemene

921 217 equemene
__device__ void impact(float phi,float r,float b,float tho,float m,
922 217 equemene
                       float *zp,float *fp,
923 217 equemene
                       int q,float db,
924 217 equemene
                       float h,int raie)
925 217 equemene
{
926 217 equemene
  float flx,rf,bss;
927 217 equemene

928 217 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
929 217 equemene

930 217 equemene
  if (raie==0)
931 217 equemene
    {
932 217 equemene
      bss=1.e19;
933 217 equemene
      flx=spectre_cn(rf,b,db,h,r,m,bss);
934 217 equemene
    }
935 217 equemene
  else
936 217 equemene
    {
937 217 equemene
      bss=2.;
938 217 equemene
      flx=spectre(rf,q,b,db,h,r,m,bss);
939 217 equemene
    }
940 217 equemene

941 217 equemene
  *zp=1./rf;
942 217 equemene
  *fp=flx;
943 217 equemene

944 217 equemene
}
945 217 equemene

946 217 equemene
__global__ void EachPixel(float *zImage,float *fImage,
947 217 equemene
                          float Mass,float InternalRadius,
948 217 equemene
                          float ExternalRadius,float Angle,
949 217 equemene
                          int Line)
950 217 equemene
{
951 218 equemene
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
952 218 equemene
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
953 217 equemene
   uint sizex=(uint)gridDim.x*blockDim.x;
954 217 equemene
   uint sizey=(uint)gridDim.y*blockDim.y;
955 217 equemene

956 217 equemene
   // Perform trajectory for each pixel, exit on hit
957 217 equemene

958 217 equemene
  float m,rs,ri,re,tho;
959 217 equemene
  int q,raie;
960 217 equemene

961 217 equemene
  m=Mass;
962 217 equemene
  rs=2.*m;
963 217 equemene
  ri=InternalRadius;
964 217 equemene
  re=ExternalRadius;
965 217 equemene
  tho=Angle;
966 217 equemene
  q=-2;
967 217 equemene
  raie=Line;
968 217 equemene

969 217 equemene
  float d,bmx,db,b,h;
970 218 equemene
  float rp0,rpp,rps;
971 217 equemene
  float phi,thi,phd,php,nr,r;
972 217 equemene
  int nh;
973 217 equemene
  float zp,fp;
974 217 equemene

975 217 equemene
  // Autosize for image
976 217 equemene
  bmx=1.25*re;
977 217 equemene
  b=0.;
978 217 equemene

979 217 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
980 217 equemene

981 217 equemene
  // set origin as center of image
982 217 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
983 217 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
984 217 equemene
  // angle extracted from cylindric symmetry
985 217 equemene
  phi=atanp(x,y);
986 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
987 217 equemene

988 217 equemene
  float up,vp,pp,us,vs,ps;
989 217 equemene

990 217 equemene
  // impact parameter
991 217 equemene
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
992 217 equemene
  // step of impact parameter;
993 217 equemene
//  db=bmx/(float)(sizex/2);
994 217 equemene
  db=bmx/(float)(sizex);
995 217 equemene

996 217 equemene
  up=0.;
997 217 equemene
  vp=1.;
998 217 equemene

999 217 equemene
  pp=0.;
1000 217 equemene
  nh=0;
1001 217 equemene

1002 217 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1003 217 equemene

1004 218 equemene
  rps=fabs(b/us);
1005 218 equemene
  rp0=rps;
1006 217 equemene

1007 217 equemene
  int ExitOnImpact=0;
1008 217 equemene

1009 217 equemene
  do
1010 217 equemene
  {
1011 217 equemene
     nh++;
1012 217 equemene
     pp=ps;
1013 217 equemene
     up=us;
1014 217 equemene
     vp=vs;
1015 218 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1016 218 equemene
     rpp=rps;
1017 218 equemene
     rps=fabs(b/us);
1018 218 equemene
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;
1019 217 equemene

1020 218 equemene
  } while ((rps>=rs)&&(rps<=rp0)&&(ExitOnImpact==0));
1021 217 equemene

1022 217 equemene
  if (ExitOnImpact==1) {
1023 218 equemene
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
1024 217 equemene
  }
1025 217 equemene
  else
1026 217 equemene
  {
1027 217 equemene
     zp=0.;
1028 217 equemene
     fp=0.;
1029 217 equemene
  }
1030 217 equemene

1031 218 equemene
  __syncthreads();
1032 218 equemene

1033 217 equemene
  zImage[yi+sizex*xi]=(float)zp;
1034 217 equemene
  fImage[yi+sizex*xi]=(float)fp;
1035 217 equemene
}
1036 217 equemene

1037 217 equemene
__global__ void Pixel(float *zImage,float *fImage,
1038 217 equemene
                      float *Trajectories,int *IdLast,
1039 224 equemene
                      uint ImpactParameter,
1040 217 equemene
                      float Mass,float InternalRadius,
1041 217 equemene
                      float ExternalRadius,float Angle,
1042 217 equemene
                      int Line)
1043 217 equemene
{
1044 219 equemene
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
1045 219 equemene
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
1046 219 equemene
   uint sizex=(uint)gridDim.x*blockDim.x;
1047 219 equemene
   uint sizey=(uint)gridDim.y*blockDim.y;
1048 217 equemene

1049 217 equemene
  // Perform trajectory for each pixel
1050 217 equemene

1051 217 equemene
  float m,rs,ri,re,tho;
1052 217 equemene
  int q,raie;
1053 217 equemene

1054 217 equemene
  m=Mass;
1055 217 equemene
  rs=2.*m;
1056 217 equemene
  ri=InternalRadius;
1057 217 equemene
  re=ExternalRadius;
1058 217 equemene
  tho=Angle;
1059 217 equemene
  q=-2;
1060 217 equemene
  raie=Line;
1061 217 equemene

1062 217 equemene
  float d,bmx,db,b,h;
1063 217 equemene
  float phi,thi,phd,php,nr,r;
1064 217 equemene
  int nh;
1065 217 equemene
  float zp=0,fp=0;
1066 217 equemene
  // Autosize for image, 25% greater than external radius
1067 217 equemene
  bmx=1.25*re;
1068 217 equemene

1069 217 equemene
  // Angular step of integration
1070 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1071 217 equemene

1072 217 equemene
  // Step of Impact Parameter
1073 217 equemene
  db=bmx/(2.e0*(float)ImpactParameter);
1074 217 equemene

1075 217 equemene
  // set origin as center of image
1076 217 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
1077 217 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
1078 217 equemene
  // angle extracted from cylindric symmetry
1079 217 equemene
  phi=atanp(x,y);
1080 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1081 217 equemene

1082 217 equemene
  // Real Impact Parameter
1083 217 equemene
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
1084 217 equemene

1085 217 equemene
  // Integer Impact Parameter
1086 217 equemene
  uint bi=(uint)sqrt(x*x+y*y);
1087 217 equemene

1088 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
1089 217 equemene

1090 217 equemene
  if (bi<ImpactParameter)
1091 217 equemene
  {
1092 217 equemene
    do
1093 217 equemene
    {
1094 217 equemene
      php=phd+(float)HalfLap*PI;
1095 217 equemene
      nr=php/h;
1096 217 equemene
      ni=(int)nr;
1097 217 equemene

1098 217 equemene
      if (ni<IdLast[bi])
1099 217 equemene
      {
1100 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1101 217 equemene
      }
1102 217 equemene
      else
1103 217 equemene
      {
1104 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
1105 217 equemene
      }
1106 217 equemene

1107 217 equemene
      if ((r<=re)&&(r>=ri))
1108 217 equemene
      {
1109 217 equemene
        ExitOnImpact=1;
1110 217 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1111 217 equemene
      }
1112 217 equemene

1113 217 equemene
      HalfLap++;
1114 217 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1115 217 equemene

1116 217 equemene
  }
1117 217 equemene

1118 217 equemene
  zImage[yi+sizex*xi]=zp;
1119 217 equemene
  fImage[yi+sizex*xi]=fp;
1120 217 equemene
}
1121 217 equemene

1122 217 equemene
__global__ void Circle(float *Trajectories,int *IdLast,
1123 217 equemene
                       float *zImage,float *fImage,
1124 217 equemene
                       float Mass,float InternalRadius,
1125 217 equemene
                       float ExternalRadius,float Angle,
1126 217 equemene
                       int Line)
1127 217 equemene
{
1128 217 equemene
   // Integer Impact Parameter ID
1129 219 equemene
   int bi=blockIdx.x*blockDim.x+threadIdx.x;
1130 217 equemene
   // Integer points on circle
1131 219 equemene
   int i=blockIdx.y*blockDim.y+threadIdx.y;
1132 217 equemene
   // Integer Impact Parameter Size (half of image)
1133 217 equemene
   int bmaxi=gridDim.x*blockDim.x;
1134 217 equemene
   // Integer Points on circle
1135 217 equemene
   int imx=gridDim.y*blockDim.y;
1136 217 equemene

1137 217 equemene
   // Perform trajectory for each pixel
1138 217 equemene

1139 217 equemene
  float m,rs,ri,re,tho;
1140 217 equemene
  int q,raie;
1141 217 equemene

1142 217 equemene
  m=Mass;
1143 217 equemene
  rs=2.*m;
1144 217 equemene
  ri=InternalRadius;
1145 217 equemene
  re=ExternalRadius;
1146 217 equemene
  tho=Angle;
1147 217 equemene
  raie=Line;
1148 217 equemene

1149 217 equemene
  float bmx,db,b,h;
1150 217 equemene
  float phi,thi,phd;
1151 217 equemene
  int nh;
1152 217 equemene
  float zp=0,fp=0;
1153 217 equemene

1154 217 equemene
  // Autosize for image
1155 217 equemene
  bmx=1.25*re;
1156 217 equemene

1157 217 equemene
  // Angular step of integration
1158 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1159 217 equemene

1160 217 equemene
  // impact parameter
1161 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
1162 217 equemene
  db=bmx/(2.e0*(float)bmaxi);
1163 217 equemene

1164 217 equemene
  phi=2.*PI/(float)imx*(float)i;
1165 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1166 217 equemene
  int yi=(int)((float)bi*sin(phi))+bmaxi;
1167 217 equemene
  int xi=(int)((float)bi*cos(phi))+bmaxi;
1168 217 equemene

1169 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
1170 217 equemene
  float php,nr,r;
1171 217 equemene

1172 217 equemene
  do
1173 217 equemene
  {
1174 217 equemene
     php=phd+(float)HalfLap*PI;
1175 217 equemene
     nr=php/h;
1176 217 equemene
     ni=(int)nr;
1177 217 equemene

1178 217 equemene
     if (ni<IdLast[bi])
1179 217 equemene
     {
1180 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1181 217 equemene
     }
1182 217 equemene
     else
1183 217 equemene
     {
1184 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
1185 217 equemene
     }
1186 217 equemene

1187 217 equemene
     if ((r<=re)&&(r>=ri))
1188 217 equemene
     {
1189 217 equemene
        ExitOnImpact=1;
1190 217 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1191 217 equemene
     }
1192 217 equemene

1193 217 equemene
     HalfLap++;
1194 217 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1195 217 equemene

1196 217 equemene
  zImage[yi+2*bmaxi*xi]=zp;
1197 217 equemene
  fImage[yi+2*bmaxi*xi]=fp;
1198 217 equemene

1199 217 equemene
}
1200 217 equemene

1201 224 equemene
__global__ void Trajectory(float *Trajectories,int *IdLast,
1202 217 equemene
                           float Mass,float InternalRadius,
1203 217 equemene
                           float ExternalRadius,float Angle,
1204 217 equemene
                           int Line)
1205 217 equemene
{
1206 217 equemene
  // Integer Impact Parameter ID
1207 219 equemene
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1208 217 equemene
  // Integer Impact Parameter Size (half of image)
1209 217 equemene
  int bmaxi=gridDim.x*blockDim.x;
1210 217 equemene

1211 217 equemene
  // Perform trajectory for each pixel
1212 217 equemene

1213 217 equemene
  float m,rs,ri,re,tho;
1214 217 equemene
  int raie,q;
1215 217 equemene

1216 217 equemene
  m=Mass;
1217 217 equemene
  rs=2.*m;
1218 217 equemene
  ri=InternalRadius;
1219 217 equemene
  re=ExternalRadius;
1220 217 equemene
  tho=Angle;
1221 217 equemene
  q=-2;
1222 217 equemene
  raie=Line;
1223 217 equemene

1224 217 equemene
  float d,bmx,db,b,h;
1225 217 equemene
  float phi,thi,phd,php,nr,r;
1226 217 equemene
  int nh;
1227 217 equemene
  float zp,fp;
1228 217 equemene

1229 217 equemene
  // Autosize for image
1230 217 equemene
  bmx=1.25*re;
1231 217 equemene

1232 217 equemene
  // Angular step of integration
1233 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1234 217 equemene

1235 217 equemene
  // impact parameter
1236 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
1237 217 equemene

1238 217 equemene
  float up,vp,pp,us,vs,ps;
1239 217 equemene

1240 217 equemene
  up=0.;
1241 217 equemene
  vp=1.;
1242 217 equemene

1243 217 equemene
  pp=0.;
1244 217 equemene
  nh=0;
1245 217 equemene

1246 217 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1247 217 equemene

1248 217 equemene
  // b versus us
1249 217 equemene
  float bvus=fabs(b/us);
1250 217 equemene
  float bvus0=bvus;
1251 224 equemene
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
1252 217 equemene

1253 217 equemene
  do
1254 217 equemene
  {
1255 217 equemene
     nh++;
1256 217 equemene
     pp=ps;
1257 217 equemene
     up=us;
1258 217 equemene
     vp=vs;
1259 217 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1260 217 equemene
     bvus=fabs(b/us);
1261 224 equemene
     Trajectories[bi*TRACKPOINTS+nh]=bvus;
1262 217 equemene

1263 217 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
1264 217 equemene

1265 217 equemene
  IdLast[bi]=nh;
1266 217 equemene

1267 217 equemene
}
1268 224 equemene

1269 225 equemene
__global__ void EachCircle(float *zImage,float *fImage,
1270 225 equemene
                           float Mass,float InternalRadius,
1271 225 equemene
                           float ExternalRadius,float Angle,
1272 225 equemene
                           int Line)
1273 224 equemene
{
1274 224 equemene
  // Integer Impact Parameter ID
1275 224 equemene
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1276 224 equemene

1277 224 equemene
  // Integer Impact Parameter Size (half of image)
1278 224 equemene
  int bmaxi=gridDim.x*blockDim.x;
1279 224 equemene

1280 225 equemene
  float Trajectory[2048];
1281 224 equemene

1282 224 equemene
  // Perform trajectory for each pixel
1283 224 equemene

1284 224 equemene
  float m,rs,ri,re,tho;
1285 224 equemene
  int raie,q;
1286 224 equemene

1287 224 equemene
  m=Mass;
1288 224 equemene
  rs=2.*m;
1289 224 equemene
  ri=InternalRadius;
1290 224 equemene
  re=ExternalRadius;
1291 224 equemene
  tho=Angle;
1292 224 equemene
  q=-2;
1293 224 equemene
  raie=Line;
1294 224 equemene

1295 224 equemene
  float bmx,db,b,h;
1296 224 equemene
  int nh;
1297 224 equemene

1298 224 equemene
  // Autosize for image
1299 224 equemene
  bmx=1.25e0f*re;
1300 224 equemene

1301 224 equemene
  // Angular step of integration
1302 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1303 224 equemene

1304 224 equemene
  // impact parameter
1305 224 equemene
  b=(float)bi/(float)bmaxi*bmx;
1306 224 equemene
  db=bmx/(2.e0f*(float)bmaxi);
1307 224 equemene

1308 224 equemene
  float up,vp,pp,us,vs,ps;
1309 224 equemene

1310 224 equemene
  up=0.;
1311 224 equemene
  vp=1.;
1312 224 equemene

1313 224 equemene
  pp=0.;
1314 224 equemene
  nh=0;
1315 224 equemene

1316 224 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1317 224 equemene

1318 224 equemene
  // b versus us
1319 224 equemene
  float bvus=fabs(b/us);
1320 224 equemene
  float bvus0=bvus;
1321 225 equemene
  Trajectory[nh]=bvus;
1322 224 equemene

1323 224 equemene
  do
1324 224 equemene
  {
1325 224 equemene
     nh++;
1326 224 equemene
     pp=ps;
1327 224 equemene
     up=us;
1328 224 equemene
     vp=vs;
1329 224 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1330 224 equemene
     bvus=fabs(b/us);
1331 225 equemene
     Trajectory[nh]=bvus;
1332 224 equemene

1333 224 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
1334 224 equemene

1335 224 equemene
  int imx=(int)(16*bi);
1336 224 equemene

1337 224 equemene
  for (int i=0;i<imx;i++)
1338 224 equemene
  {
1339 224 equemene
     float zp=0,fp=0;
1340 224 equemene
     float phi=2.*PI/(float)imx*(float)i;
1341 224 equemene
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
1342 224 equemene
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1343 224 equemene
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1344 224 equemene

1345 224 equemene
     int HalfLap=0,ExitOnImpact=0,ni;
1346 224 equemene
     float php,nr,r;
1347 224 equemene

1348 224 equemene
     do
1349 224 equemene
     {
1350 224 equemene
        php=phd+(float)HalfLap*PI;
1351 224 equemene
        nr=php/h;
1352 224 equemene
        ni=(int)nr;
1353 224 equemene

1354 224 equemene
        if (ni<nh)
1355 224 equemene
        {
1356 225 equemene
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
1357 224 equemene
        }
1358 224 equemene
        else
1359 224 equemene
        {
1360 225 equemene
           r=Trajectory[ni];
1361 224 equemene
        }
1362 224 equemene

1363 224 equemene
        if ((r<=re)&&(r>=ri))
1364 224 equemene
        {
1365 224 equemene
           ExitOnImpact=1;
1366 224 equemene
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1367 224 equemene
        }
1368 224 equemene

1369 224 equemene
        HalfLap++;
1370 224 equemene

1371 224 equemene
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
1372 224 equemene

1373 224 equemene
   __syncthreads();
1374 224 equemene

1375 224 equemene
   zImage[yi+2*bmaxi*xi]=zp;
1376 224 equemene
   fImage[yi+2*bmaxi*xi]=fp;
1377 224 equemene

1378 224 equemene
  }
1379 224 equemene

1380 224 equemene
}
1381 225 equemene

1382 225 equemene
__global__ void Original(float *zImage,float *fImage,
1383 225 equemene
                         uint Size,float Mass,float InternalRadius,
1384 225 equemene
                         float ExternalRadius,float Angle,
1385 225 equemene
                         int Line)
1386 225 equemene
{
1387 225 equemene
   // Integer Impact Parameter Size (half of image)
1388 225 equemene
   uint bmaxi=(uint)Size;
1389 225 equemene

1390 225 equemene
   float Trajectory[TRACKPOINTS];
1391 225 equemene

1392 225 equemene
   // Perform trajectory for each pixel
1393 225 equemene

1394 225 equemene
   float m,rs,ri,re,tho;
1395 225 equemene
   int raie,q;
1396 225 equemene

1397 225 equemene
   m=Mass;
1398 225 equemene
   rs=2.*m;
1399 225 equemene
   ri=InternalRadius;
1400 225 equemene
   re=ExternalRadius;
1401 225 equemene
   tho=Angle;
1402 225 equemene
   q=-2;
1403 225 equemene
   raie=Line;
1404 225 equemene

1405 225 equemene
   float bmx,db,b,h;
1406 225 equemene
   int nh;
1407 225 equemene

1408 225 equemene
   // Autosize for image
1409 225 equemene
   bmx=1.25e0f*re;
1410 225 equemene

1411 225 equemene
   // Angular step of integration
1412 225 equemene
   h=4.e0f*PI/(float)TRACKPOINTS;
1413 225 equemene

1414 225 equemene
   // Integer Impact Parameter ID
1415 225 equemene
   for (int bi=0;bi<bmaxi;bi++)
1416 225 equemene
   {
1417 225 equemene
      // impact parameter
1418 225 equemene
      b=(float)bi/(float)bmaxi*bmx;
1419 225 equemene
      db=bmx/(2.e0f*(float)bmaxi);
1420 225 equemene

1421 225 equemene
      float up,vp,pp,us,vs,ps;
1422 225 equemene

1423 225 equemene
      up=0.;
1424 225 equemene
      vp=1.;
1425 225 equemene

1426 225 equemene
      pp=0.;
1427 225 equemene
      nh=0;
1428 225 equemene

1429 225 equemene
      rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1430 225 equemene

1431 225 equemene
      // b versus us
1432 225 equemene
      float bvus=fabs(b/us);
1433 225 equemene
      float bvus0=bvus;
1434 225 equemene
      Trajectory[nh]=bvus;
1435 225 equemene

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

1446 225 equemene
      } while ((bvus>=rs)&&(bvus<=bvus0));
1447 225 equemene

1448 225 equemene
      for (uint i=(uint)nh+1;i<TRACKPOINTS;i++) {
1449 225 equemene
         Trajectory[i]=0.e0f;
1450 225 equemene
      }
1451 225 equemene

1452 225 equemene
      int imx=(int)(16*bi);
1453 225 equemene

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

1462 225 equemene
         int HalfLap=0,ExitOnImpact=0,ni;
1463 225 equemene
         float php,nr,r;
1464 225 equemene

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

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

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

1486 225 equemene
            HalfLap++;
1487 225 equemene

1488 225 equemene
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
1489 225 equemene

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

1493 225 equemene
      }
1494 225 equemene

1495 225 equemene
   }
1496 225 equemene

1497 225 equemene
}
1498 217 equemene
"""
1499 217 equemene
    return(BlobCUDA)
1500 217 equemene
1501 211 equemene
# def ImageOutput(sigma,prefix):
1502 211 equemene
#     from PIL import Image
1503 211 equemene
#     Max=sigma.max()
1504 211 equemene
#     Min=sigma.min()
1505 199 equemene
1506 211 equemene
#     # Normalize value as 8bits Integer
1507 211 equemene
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1508 211 equemene
#     image = Image.fromarray(SigmaInt)
1509 211 equemene
#     image.save("%s.jpg" % prefix)
1510 211 equemene
1511 211 equemene
def ImageOutput(sigma,prefix,Colors):
1512 211 equemene
    import matplotlib.pyplot as plt
1513 222 equemene
    start_time=time.time()
1514 211 equemene
    if Colors == 'Red2Yellow':
1515 211 equemene
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1516 211 equemene
    else:
1517 211 equemene
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1518 222 equemene
    save_time = time.time()-start_time
1519 222 equemene
    print("Save Time : %f" % save_time)
1520 211 equemene
1521 204 equemene
def BlackHoleCL(zImage,fImage,InputCL):
1522 199 equemene
1523 199 equemene
    kernel_params = {}
1524 199 equemene
1525 204 equemene
    print(InputCL)
1526 204 equemene
1527 204 equemene
    Device=InputCL['Device']
1528 204 equemene
    GpuStyle=InputCL['GpuStyle']
1529 204 equemene
    VariableType=InputCL['VariableType']
1530 204 equemene
    Size=InputCL['Size']
1531 204 equemene
    Mass=InputCL['Mass']
1532 204 equemene
    InternalRadius=InputCL['InternalRadius']
1533 204 equemene
    ExternalRadius=InputCL['ExternalRadius']
1534 204 equemene
    Angle=InputCL['Angle']
1535 204 equemene
    Method=InputCL['Method']
1536 204 equemene
    TrackPoints=InputCL['TrackPoints']
1537 211 equemene
    Physics=InputCL['Physics']
1538 204 equemene
1539 211 equemene
    PhysicsList=DictionariesAPI()
1540 211 equemene
1541 204 equemene
    if InputCL['BlackBody']:
1542 209 equemene
        # Spectrum is Black Body one
1543 209 equemene
        Line=0
1544 204 equemene
    else:
1545 209 equemene
        # Spectrum is Monochromatic Line one
1546 209 equemene
        Line=1
1547 204 equemene
1548 217 equemene
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1549 217 equemene
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1550 204 equemene
1551 199 equemene
    # Je detecte un peripherique GPU dans la liste des peripheriques
1552 199 equemene
    Id=0
1553 199 equemene
    HasXPU=False
1554 199 equemene
    for platform in cl.get_platforms():
1555 199 equemene
        for device in platform.get_devices():
1556 199 equemene
            if Id==Device:
1557 199 equemene
                XPU=device
1558 217 equemene
                print("CPU/GPU selected: ",device.name.lstrip())
1559 199 equemene
                HasXPU=True
1560 199 equemene
            Id+=1
1561 199 equemene
1562 199 equemene
    if HasXPU==False:
1563 217 equemene
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1564 199 equemene
        sys.exit()
1565 199 equemene
1566 199 equemene
    ctx = cl.Context([XPU])
1567 199 equemene
    queue = cl.CommandQueue(ctx,
1568 199 equemene
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1569 199 equemene
1570 199 equemene
1571 211 equemene
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1572 211 equemene
1573 224 equemene
    BuildOptions="-cl-mad-enable -DPHYSICS=%i -DSETTRACKPOINTS=%i " % (PhysicsList[Physics],InputCL['TrackPoints'])
1574 211 equemene
1575 211 equemene
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1576 211 equemene
1577 199 equemene
    # Je recupere les flag possibles pour les buffers
1578 199 equemene
    mf = cl.mem_flags
1579 199 equemene
1580 204 equemene
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
1581 201 equemene
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1582 201 equemene
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1583 204 equemene
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1584 199 equemene
1585 199 equemene
    start_time=time.time()
1586 199 equemene
1587 204 equemene
    if Method=='EachPixel':
1588 204 equemene
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1589 204 equemene
                                       None,zImageCL,fImageCL,
1590 204 equemene
                                       numpy.float32(Mass),
1591 204 equemene
                                       numpy.float32(InternalRadius),
1592 204 equemene
                                       numpy.float32(ExternalRadius),
1593 204 equemene
                                       numpy.float32(Angle),
1594 209 equemene
                                       numpy.int32(Line))
1595 204 equemene
        CLLaunch.wait()
1596 224 equemene
    elif Method=='Original':
1597 225 equemene
        CLLaunch=BlackHoleCL.Original(queue,(1,),
1598 224 equemene
                                      None,zImageCL,fImageCL,
1599 225 equemene
                                      numpy.uint32(zImage.shape[0]/2),
1600 224 equemene
                                      numpy.float32(Mass),
1601 224 equemene
                                      numpy.float32(InternalRadius),
1602 224 equemene
                                      numpy.float32(ExternalRadius),
1603 224 equemene
                                      numpy.float32(Angle),
1604 224 equemene
                                      numpy.int32(Line))
1605 224 equemene
        CLLaunch.wait()
1606 225 equemene
    elif Method=='EachCircle':
1607 225 equemene
        CLLaunch=BlackHoleCL.EachCircle(queue,(zImage.shape[0]/2,),
1608 225 equemene
                                        None,zImageCL,fImageCL,
1609 225 equemene
                                        numpy.float32(Mass),
1610 225 equemene
                                        numpy.float32(InternalRadius),
1611 225 equemene
                                        numpy.float32(ExternalRadius),
1612 225 equemene
                                        numpy.float32(Angle),
1613 225 equemene
                                        numpy.int32(Line))
1614 225 equemene
        CLLaunch.wait()
1615 204 equemene
    elif Method=='TrajectoCircle':
1616 204 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1617 204 equemene
                                        None,TrajectoriesCL,IdLastCL,
1618 204 equemene
                                        numpy.float32(Mass),
1619 204 equemene
                                        numpy.float32(InternalRadius),
1620 204 equemene
                                        numpy.float32(ExternalRadius),
1621 204 equemene
                                        numpy.float32(Angle),
1622 209 equemene
                                        numpy.int32(Line))
1623 204 equemene
1624 204 equemene
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1625 204 equemene
                                           zImage.shape[0]*4),None,
1626 204 equemene
                                    TrajectoriesCL,IdLastCL,
1627 204 equemene
                                    zImageCL,fImageCL,
1628 204 equemene
                                    numpy.float32(Mass),
1629 204 equemene
                                    numpy.float32(InternalRadius),
1630 204 equemene
                                    numpy.float32(ExternalRadius),
1631 204 equemene
                                    numpy.float32(Angle),
1632 209 equemene
                                    numpy.int32(Line))
1633 204 equemene
        CLLaunch.wait()
1634 204 equemene
    else:
1635 204 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1636 204 equemene
                                        None,TrajectoriesCL,IdLastCL,
1637 204 equemene
                                        numpy.float32(Mass),
1638 204 equemene
                                        numpy.float32(InternalRadius),
1639 204 equemene
                                        numpy.float32(ExternalRadius),
1640 204 equemene
                                        numpy.float32(Angle),
1641 209 equemene
                                        numpy.int32(Line))
1642 204 equemene
1643 204 equemene
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],
1644 217 equemene
                                          zImage.shape[1]),None,
1645 204 equemene
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1646 204 equemene
                                   numpy.uint32(Trajectories.shape[0]),
1647 204 equemene
                                        numpy.float32(Mass),
1648 204 equemene
                                        numpy.float32(InternalRadius),
1649 204 equemene
                                        numpy.float32(ExternalRadius),
1650 204 equemene
                                        numpy.float32(Angle),
1651 209 equemene
                                        numpy.int32(Line))
1652 204 equemene
        CLLaunch.wait()
1653 204 equemene
1654 218 equemene
    compute = time.time()-start_time
1655 199 equemene
1656 201 equemene
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1657 201 equemene
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1658 204 equemene
    cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1659 204 equemene
    cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1660 218 equemene
    elapsed = time.time()-start_time
1661 218 equemene
    print("\nCompute Time : %f" % compute)
1662 218 equemene
    print("Elapsed Time : %f\n" % elapsed)
1663 211 equemene
1664 211 equemene
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1665 211 equemene
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1666 219 equemene
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1667 219 equemene
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1668 201 equemene
    zImageCL.release()
1669 201 equemene
    fImageCL.release()
1670 204 equemene
1671 204 equemene
    TrajectoriesCL.release()
1672 204 equemene
    IdLastCL.release()
1673 204 equemene
1674 199 equemene
    return(elapsed)
1675 199 equemene
1676 217 equemene
def BlackHoleCUDA(zImage,fImage,InputCL):
1677 217 equemene
1678 217 equemene
    kernel_params = {}
1679 217 equemene
1680 217 equemene
    print(InputCL)
1681 217 equemene
1682 217 equemene
    Device=InputCL['Device']
1683 217 equemene
    GpuStyle=InputCL['GpuStyle']
1684 217 equemene
    VariableType=InputCL['VariableType']
1685 217 equemene
    Size=InputCL['Size']
1686 217 equemene
    Mass=InputCL['Mass']
1687 217 equemene
    InternalRadius=InputCL['InternalRadius']
1688 217 equemene
    ExternalRadius=InputCL['ExternalRadius']
1689 217 equemene
    Angle=InputCL['Angle']
1690 217 equemene
    Method=InputCL['Method']
1691 217 equemene
    TrackPoints=InputCL['TrackPoints']
1692 217 equemene
    Physics=InputCL['Physics']
1693 217 equemene
1694 217 equemene
    PhysicsList=DictionariesAPI()
1695 217 equemene
1696 217 equemene
    if InputCL['BlackBody']:
1697 217 equemene
        # Spectrum is Black Body one
1698 217 equemene
        Line=0
1699 217 equemene
    else:
1700 217 equemene
        # Spectrum is Monochromatic Line one
1701 217 equemene
        Line=1
1702 217 equemene
1703 217 equemene
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1704 217 equemene
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1705 217 equemene
1706 217 equemene
    try:
1707 217 equemene
        # For PyCUDA import
1708 217 equemene
        import pycuda.driver as cuda
1709 217 equemene
        from pycuda.compiler import SourceModule
1710 217 equemene
1711 217 equemene
        cuda.init()
1712 217 equemene
        for Id in range(cuda.Device.count()):
1713 217 equemene
            if Id==Device:
1714 217 equemene
                XPU=cuda.Device(Id)
1715 217 equemene
                print("GPU selected %s" % XPU.name())
1716 217 equemene
        print
1717 217 equemene
1718 217 equemene
    except ImportError:
1719 217 equemene
        print("Platform does not seem to support CUDA")
1720 217 equemene
1721 217 equemene
    Context=XPU.make_context()
1722 217 equemene
1723 217 equemene
    try:
1724 224 equemene
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1725 217 equemene
        print("Compilation seems to be OK")
1726 217 equemene
    except:
1727 217 equemene
        print("Compilation seems to break")
1728 217 equemene
1729 217 equemene
    EachPixelCU=mod.get_function("EachPixel")
1730 224 equemene
    OriginalCU=mod.get_function("Original")
1731 225 equemene
    EachCircleCU=mod.get_function("EachCircle")
1732 217 equemene
    TrajectoryCU=mod.get_function("Trajectory")
1733 217 equemene
    PixelCU=mod.get_function("Pixel")
1734 217 equemene
    CircleCU=mod.get_function("Circle")
1735 217 equemene
1736 217 equemene
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1737 217 equemene
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1738 217 equemene
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1739 217 equemene
    cuda.memcpy_htod(zImageCU, zImage)
1740 217 equemene
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1741 217 equemene
    cuda.memcpy_htod(zImageCU, fImage)
1742 217 equemene
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1743 217 equemene
    cuda.memcpy_htod(IdLastCU, IdLast)
1744 217 equemene
1745 217 equemene
    start_time=time.time()
1746 217 equemene
1747 217 equemene
    if Method=='EachPixel':
1748 217 equemene
        EachPixelCU(zImageCU,fImageCU,
1749 217 equemene
                    numpy.float32(Mass),
1750 217 equemene
                    numpy.float32(InternalRadius),
1751 217 equemene
                    numpy.float32(ExternalRadius),
1752 217 equemene
                    numpy.float32(Angle),
1753 217 equemene
                    numpy.int32(Line),
1754 219 equemene
                    grid=(zImage.shape[0]/32,zImage.shape[1]/32),
1755 219 equemene
                    block=(32,32,1))
1756 225 equemene
    elif Method=='EachCircle':
1757 225 equemene
        EachCircleCU(zImageCU,fImageCU,
1758 225 equemene
                   numpy.float32(Mass),
1759 225 equemene
                   numpy.float32(InternalRadius),
1760 225 equemene
                   numpy.float32(ExternalRadius),
1761 225 equemene
                   numpy.float32(Angle),
1762 225 equemene
                   numpy.int32(Line),
1763 225 equemene
                   grid=(zImage.shape[0]/32/2,1),
1764 225 equemene
                   block=(32,1,1))
1765 224 equemene
    elif Method=='Original':
1766 224 equemene
        OriginalCU(zImageCU,fImageCU,
1767 225 equemene
                   numpy.uint32(zImage.shape[0]/2),
1768 224 equemene
                   numpy.float32(Mass),
1769 224 equemene
                   numpy.float32(InternalRadius),
1770 224 equemene
                   numpy.float32(ExternalRadius),
1771 224 equemene
                   numpy.float32(Angle),
1772 224 equemene
                   numpy.int32(Line),
1773 225 equemene
                   grid=(1,1),
1774 225 equemene
                   block=(1,1,1))
1775 217 equemene
    elif Method=='TrajectoCircle':
1776 217 equemene
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1777 217 equemene
                     numpy.float32(Mass),
1778 217 equemene
                     numpy.float32(InternalRadius),
1779 217 equemene
                     numpy.float32(ExternalRadius),
1780 217 equemene
                     numpy.float32(Angle),
1781 217 equemene
                     numpy.int32(Line),
1782 219 equemene
                     grid=(Trajectories.shape[0]/32,1),
1783 219 equemene
                     block=(32,1,1))
1784 217 equemene
1785 217 equemene
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1786 217 equemene
                 numpy.float32(Mass),
1787 217 equemene
                 numpy.float32(InternalRadius),
1788 217 equemene
                 numpy.float32(ExternalRadius),
1789 217 equemene
                 numpy.float32(Angle),
1790 217 equemene
                 numpy.int32(Line),
1791 219 equemene
                 grid=(Trajectories.shape[0]/32,zImage.shape[0]*4/32),
1792 219 equemene
                 block=(32,32,1))
1793 217 equemene
    else:
1794 217 equemene
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1795 217 equemene
                     numpy.float32(Mass),
1796 217 equemene
                     numpy.float32(InternalRadius),
1797 217 equemene
                     numpy.float32(ExternalRadius),
1798 217 equemene
                     numpy.float32(Angle),
1799 217 equemene
                     numpy.int32(Line),
1800 219 equemene
                     grid=(Trajectories.shape[0]/32,1),
1801 219 equemene
                     block=(32,1,1))
1802 217 equemene
1803 217 equemene
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1804 217 equemene
                numpy.uint32(Trajectories.shape[0]),
1805 217 equemene
                numpy.float32(Mass),
1806 217 equemene
                numpy.float32(InternalRadius),
1807 217 equemene
                numpy.float32(ExternalRadius),
1808 217 equemene
                numpy.float32(Angle),
1809 217 equemene
                numpy.int32(Line),
1810 219 equemene
                grid=(zImage.shape[0]/32,zImage.shape[1]/32,1),
1811 219 equemene
                block=(32,32,1))
1812 217 equemene
1813 220 equemene
    Context.synchronize()
1814 220 equemene
1815 218 equemene
    compute = time.time()-start_time
1816 217 equemene
1817 217 equemene
    cuda.memcpy_dtoh(zImage,zImageCU)
1818 217 equemene
    cuda.memcpy_dtoh(fImage,fImageCU)
1819 218 equemene
    elapsed = time.time()-start_time
1820 218 equemene
    print("\nCompute Time : %f" % compute)
1821 218 equemene
    print("Elapsed Time : %f\n" % elapsed)
1822 217 equemene
1823 217 equemene
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1824 217 equemene
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1825 219 equemene
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1826 219 equemene
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1827 217 equemene
1828 220 equemene
1829 218 equemene
    Context.pop()
1830 220 equemene
1831 217 equemene
    Context.detach()
1832 217 equemene
1833 217 equemene
    return(elapsed)
1834 217 equemene
1835 199 equemene
if __name__=='__main__':
1836 199 equemene
1837 199 equemene
    GpuStyle = 'OpenCL'
1838 199 equemene
    Mass = 1.
1839 199 equemene
    # Internal Radius 3 times de Schwarzschild Radius
1840 199 equemene
    InternalRadius=6.*Mass
1841 199 equemene
    #
1842 199 equemene
    ExternalRadius=12.
1843 199 equemene
    #
1844 199 equemene
    # Angle with normal to disc 10 degrees
1845 199 equemene
    Angle = numpy.pi/180.*(90.-10.)
1846 199 equemene
    # Radiation of disc : BlackBody or Monochromatic
1847 209 equemene
    BlackBody = False
1848 199 equemene
    # Size of image
1849 199 equemene
    Size=256
1850 199 equemene
    # Variable Type
1851 199 equemene
    VariableType='FP32'
1852 199 equemene
    # ?
1853 199 equemene
    q=-2
1854 204 equemene
    # Method of resolution
1855 209 equemene
    Method='TrajectoPixel'
1856 211 equemene
    # Colors for output image
1857 211 equemene
    Colors='Greyscale'
1858 211 equemene
    # Physics
1859 211 equemene
    Physics='Einstein'
1860 211 equemene
    # No output as image
1861 211 equemene
    NoImage = False
1862 211 equemene
1863 225 equemene
    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>'
1864 199 equemene
1865 199 equemene
    try:
1866 211 equemene
        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="])
1867 199 equemene
    except getopt.GetoptError:
1868 199 equemene
        print(HowToUse % sys.argv[0])
1869 199 equemene
        sys.exit(2)
1870 199 equemene
1871 199 equemene
    # List of Devices
1872 199 equemene
    Devices=[]
1873 199 equemene
    Alu={}
1874 199 equemene
1875 199 equemene
    for opt, arg in opts:
1876 199 equemene
        if opt == '-h':
1877 199 equemene
            print(HowToUse % sys.argv[0])
1878 199 equemene
1879 199 equemene
            print("\nInformations about devices detected under OpenCL API:")
1880 199 equemene
            # For PyOpenCL import
1881 199 equemene
            try:
1882 199 equemene
                import pyopencl as cl
1883 199 equemene
                Id=0
1884 199 equemene
                for platform in cl.get_platforms():
1885 199 equemene
                    for device in platform.get_devices():
1886 199 equemene
                        #deviceType=cl.device_type.to_string(device.type)
1887 199 equemene
                        deviceType="xPU"
1888 199 equemene
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
1889 199 equemene
                        Id=Id+1
1890 199 equemene
1891 199 equemene
            except:
1892 199 equemene
                print("Your platform does not seem to support OpenCL")
1893 199 equemene
1894 199 equemene
            print("\nInformations about devices detected under CUDA API:")
1895 199 equemene
                # For PyCUDA import
1896 199 equemene
            try:
1897 199 equemene
                import pycuda.driver as cuda
1898 199 equemene
                cuda.init()
1899 199 equemene
                for Id in range(cuda.Device.count()):
1900 199 equemene
                    device=cuda.Device(Id)
1901 199 equemene
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
1902 199 equemene
                print
1903 199 equemene
            except:
1904 199 equemene
                print("Your platform does not seem to support CUDA")
1905 199 equemene
1906 199 equemene
            sys.exit()
1907 199 equemene
1908 199 equemene
        elif opt in ("-d", "--device"):
1909 199 equemene
#            Devices.append(int(arg))
1910 199 equemene
            Device=int(arg)
1911 199 equemene
        elif opt in ("-g", "--gpustyle"):
1912 199 equemene
            GpuStyle = arg
1913 204 equemene
        elif opt in ("-v", "--variabletype"):
1914 199 equemene
            VariableType = arg
1915 199 equemene
        elif opt in ("-s", "--size"):
1916 199 equemene
            Size = int(arg)
1917 199 equemene
        elif opt in ("-m", "--mass"):
1918 199 equemene
            Mass = float(arg)
1919 199 equemene
        elif opt in ("-i", "--internal"):
1920 199 equemene
            InternalRadius = float(arg)
1921 199 equemene
        elif opt in ("-e", "--external"):
1922 199 equemene
            ExternalRadius = float(arg)
1923 199 equemene
        elif opt in ("-a", "--angle"):
1924 199 equemene
            Angle = numpy.pi/180.*(90.-float(arg))
1925 199 equemene
        elif opt in ("-b", "--blackbody"):
1926 199 equemene
            BlackBody = True
1927 211 equemene
        elif opt in ("-n", "--noimage"):
1928 211 equemene
            NoImage = True
1929 204 equemene
        elif opt in ("-t", "--method"):
1930 204 equemene
            Method = arg
1931 211 equemene
        elif opt in ("-c", "--colors"):
1932 211 equemene
            Colors = arg
1933 211 equemene
        elif opt in ("-p", "--physics"):
1934 211 equemene
            Physics = arg
1935 199 equemene
1936 199 equemene
    print("Device Identification selected : %s" % Device)
1937 199 equemene
    print("GpuStyle used : %s" % GpuStyle)
1938 199 equemene
    print("VariableType : %s" % VariableType)
1939 199 equemene
    print("Size : %i" % Size)
1940 199 equemene
    print("Mass : %f" % Mass)
1941 199 equemene
    print("Internal Radius : %f" % InternalRadius)
1942 199 equemene
    print("External Radius : %f" % ExternalRadius)
1943 199 equemene
    print("Angle with normal of (in radians) : %f" % Angle)
1944 199 equemene
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
1945 204 equemene
    print("Method of resolution : %s" % Method)
1946 211 equemene
    print("Colors for output images : %s" % Colors)
1947 211 equemene
    print("Physics used for Trajectories : %s" % Physics)
1948 199 equemene
1949 199 equemene
    if GpuStyle=='CUDA':
1950 199 equemene
        print("\nSelection of CUDA device")
1951 199 equemene
        try:
1952 199 equemene
            # For PyCUDA import
1953 199 equemene
            import pycuda.driver as cuda
1954 199 equemene
1955 199 equemene
            cuda.init()
1956 199 equemene
            for Id in range(cuda.Device.count()):
1957 199 equemene
                device=cuda.Device(Id)
1958 199 equemene
                print("Device #%i of type GPU : %s" % (Id,device.name()))
1959 199 equemene
                if Id in Devices:
1960 199 equemene
                    Alu[Id]='GPU'
1961 199 equemene
1962 199 equemene
        except ImportError:
1963 199 equemene
            print("Platform does not seem to support CUDA")
1964 199 equemene
1965 199 equemene
    if GpuStyle=='OpenCL':
1966 199 equemene
        print("\nSelection of OpenCL device")
1967 199 equemene
        try:
1968 199 equemene
            # For PyOpenCL import
1969 199 equemene
            import pyopencl as cl
1970 199 equemene
            Id=0
1971 199 equemene
            for platform in cl.get_platforms():
1972 199 equemene
                for device in platform.get_devices():
1973 199 equemene
                    #deviceType=cl.device_type.to_string(device.type)
1974 199 equemene
                    deviceType="xPU"
1975 199 equemene
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
1976 199 equemene
1977 199 equemene
                    if Id in Devices:
1978 199 equemene
                    # Set the Alu as detected Device Type
1979 199 equemene
                        Alu[Id]=deviceType
1980 199 equemene
                    Id=Id+1
1981 199 equemene
        except ImportError:
1982 199 equemene
            print("Platform does not seem to support OpenCL")
1983 199 equemene
1984 199 equemene
#    print(Devices,Alu)
1985 199 equemene
1986 199 equemene
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
1987 204 equemene
    TrackPoints=2048
1988 201 equemene
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1989 201 equemene
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1990 199 equemene
1991 204 equemene
    InputCL={}
1992 204 equemene
    InputCL['Device']=Device
1993 204 equemene
    InputCL['GpuStyle']=GpuStyle
1994 204 equemene
    InputCL['VariableType']=VariableType
1995 204 equemene
    InputCL['Size']=Size
1996 204 equemene
    InputCL['Mass']=Mass
1997 204 equemene
    InputCL['InternalRadius']=InternalRadius
1998 204 equemene
    InputCL['ExternalRadius']=ExternalRadius
1999 204 equemene
    InputCL['Angle']=Angle
2000 204 equemene
    InputCL['BlackBody']=BlackBody
2001 204 equemene
    InputCL['Method']=Method
2002 204 equemene
    InputCL['TrackPoints']=TrackPoints
2003 211 equemene
    InputCL['Physics']=Physics
2004 199 equemene
2005 217 equemene
    if GpuStyle=='OpenCL':
2006 217 equemene
        duration=BlackHoleCL(zImage,fImage,InputCL)
2007 217 equemene
    else:
2008 217 equemene
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
2009 217 equemene
2010 211 equemene
    Hostname=gethostname()
2011 211 equemene
    Date=time.strftime("%Y%m%d_%H%M%S")
2012 211 equemene
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
2013 211 equemene
2014 211 equemene
2015 211 equemene
    if not NoImage:
2016 211 equemene
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
2017 211 equemene
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)