Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 219

Historique | Voir | Annoter | Télécharger (39,34 ko)

1 217 equemene
#!/usr/bin/env python3
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 204 equemene
#
18 204 equemene
# Thanks to :
19 204 equemene
#
20 204 equemene
# - Herve Aussel for his part of code of black body spectrum
21 204 equemene
# - Didier Pelat for his help to perform this work
22 204 equemene
# - Jean-Pierre Luminet for his article published in 1979
23 204 equemene
# - Numerical Recipies for Runge Kutta recipies
24 204 equemene
# - Luc Blanchet for his disponibility about my questions in General Relativity
25 204 equemene
# - Pierre Lena for his passion about science and vulgarisation
26 199 equemene
27 199 equemene
import pyopencl as cl
28 199 equemene
import numpy
29 199 equemene
import time,string
30 199 equemene
from numpy.random import randint as nprnd
31 199 equemene
import sys
32 199 equemene
import getopt
33 199 equemene
import matplotlib.pyplot as plt
34 211 equemene
from socket import gethostname
35 199 equemene
36 211 equemene
def DictionariesAPI():
37 211 equemene
    PhysicsList={'Einstein':0,'Newton':1}
38 211 equemene
    return(PhysicsList)
39 204 equemene
40 211 equemene
BlobOpenCL= """
41 204 equemene

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

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

48 217 equemene
#define TRACKPOINTS 2048
49 217 equemene

50 199 equemene
float atanp(float x,float y)
51 199 equemene
{
52 199 equemene
  float angle;
53 199 equemene

54 199 equemene
  angle=atan2(y,x);
55 199 equemene

56 204 equemene
  if (angle<0.e0f)
57 199 equemene
    {
58 199 equemene
      angle+=(float)2.e0f*PI;
59 199 equemene
    }
60 199 equemene

61 199 equemene
  return angle;
62 199 equemene
}
63 199 equemene

64 199 equemene
float f(float v)
65 199 equemene
{
66 199 equemene
  return v;
67 199 equemene
}
68 199 equemene

69 211 equemene
#if PHYSICS == NEWTON
70 199 equemene
float g(float u,float m,float b)
71 199 equemene
{
72 211 equemene
  return (-u);
73 211 equemene
}
74 211 equemene
#else
75 211 equemene
float g(float u,float m,float b)
76 211 equemene
{
77 204 equemene
  return (3.e0f*m/b*pow(u,2)-u);
78 199 equemene
}
79 211 equemene
#endif
80 199 equemene

81 199 equemene
void calcul(float *us,float *vs,float up,float vp,
82 199 equemene
            float h,float m,float b)
83 199 equemene
{
84 199 equemene
  float c0,c1,c2,c3,d0,d1,d2,d3;
85 199 equemene

86 199 equemene
  c0=h*f(vp);
87 199 equemene
  c1=h*f(vp+c0/2.);
88 199 equemene
  c2=h*f(vp+c1/2.);
89 199 equemene
  c3=h*f(vp+c2);
90 199 equemene
  d0=h*g(up,m,b);
91 199 equemene
  d1=h*g(up+d0/2.,m,b);
92 199 equemene
  d2=h*g(up+d1/2.,m,b);
93 199 equemene
  d3=h*g(up+d2,m,b);
94 199 equemene

95 199 equemene
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
96 199 equemene
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
97 199 equemene
}
98 199 equemene

99 199 equemene
void rungekutta(float *ps,float *us,float *vs,
100 199 equemene
                float pp,float up,float vp,
101 199 equemene
                float h,float m,float b)
102 199 equemene
{
103 199 equemene
  calcul(us,vs,up,vp,h,m,b);
104 199 equemene
  *ps=pp+h;
105 199 equemene
}
106 199 equemene

107 199 equemene
float decalage_spectral(float r,float b,float phi,
108 199 equemene
                         float tho,float m)
109 199 equemene
{
110 199 equemene
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
111 199 equemene
}
112 199 equemene

113 199 equemene
float spectre(float rf,int q,float b,float db,
114 199 equemene
             float h,float r,float m,float bss)
115 199 equemene
{
116 199 equemene
  float flx;
117 199 equemene

118 209 equemene
  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
119 199 equemene
  return(flx);
120 199 equemene
}
121 199 equemene

122 209 equemene
float spectre_cn(float rf32,float b32,float db32,
123 209 equemene
                 float h32,float r32,float m32,float bss32)
124 199 equemene
{
125 209 equemene

126 213 equemene
#define MYFLOAT float
127 209 equemene

128 209 equemene
  MYFLOAT rf=(MYFLOAT)(rf32);
129 209 equemene
  MYFLOAT b=(MYFLOAT)(b32);
130 209 equemene
  MYFLOAT db=(MYFLOAT)(db32);
131 209 equemene
  MYFLOAT h=(MYFLOAT)(h32);
132 209 equemene
  MYFLOAT r=(MYFLOAT)(r32);
133 209 equemene
  MYFLOAT m=(MYFLOAT)(m32);
134 209 equemene
  MYFLOAT bss=(MYFLOAT)(bss32);
135 209 equemene

136 209 equemene
  MYFLOAT flx;
137 209 equemene
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
138 199 equemene
  int fi,posfreq;
139 199 equemene

140 209 equemene
#define planck 6.62e-34
141 209 equemene
#define k 1.38e-23
142 209 equemene
#define c2 9.e16
143 209 equemene
#define temp 3.e7
144 209 equemene
#define m_point 1.
145 199 equemene

146 209 equemene
#define lplanck (log(6.62)-34.*log(10.))
147 209 equemene
#define lk (log(1.38)-23.*log(10.))
148 209 equemene
#define lc2 (log(9.)+16.*log(10.))
149 199 equemene

150 209 equemene
  MYFLOAT v=1.-3./r;
151 199 equemene

152 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 ));
153 199 equemene

154 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)));
155 209 equemene

156 209 equemene
  flux_int=0.;
157 209 equemene
  flx=0.;
158 209 equemene

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

174 199 equemene
          flx+=flux_int;
175 199 equemene
        }
176 199 equemene
    }
177 209 equemene

178 209 equemene
  return((float)(flx));
179 199 equemene
}
180 199 equemene

181 199 equemene
void impact(float phi,float r,float b,float tho,float m,
182 199 equemene
            float *zp,float *fp,
183 199 equemene
            int q,float db,
184 204 equemene
            float h,int raie)
185 199 equemene
{
186 204 equemene
  float flx,rf,bss;
187 199 equemene

188 199 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
189 199 equemene

190 199 equemene
  if (raie==0)
191 199 equemene
    {
192 209 equemene
      bss=1.e19;
193 209 equemene
      flx=spectre_cn(rf,b,db,h,r,m,bss);
194 209 equemene
    }
195 209 equemene
  else
196 209 equemene
    {
197 204 equemene
      bss=2.;
198 199 equemene
      flx=spectre(rf,q,b,db,h,r,m,bss);
199 199 equemene
    }
200 199 equemene

201 199 equemene
  *zp=1./rf;
202 199 equemene
  *fp=flx;
203 199 equemene

204 199 equemene
}
205 199 equemene

206 204 equemene
__kernel void EachPixel(__global float *zImage,__global float *fImage,
207 204 equemene
                        float Mass,float InternalRadius,
208 204 equemene
                        float ExternalRadius,float Angle,
209 209 equemene
                        int Line)
210 199 equemene
{
211 199 equemene
   uint xi=(uint)get_global_id(0);
212 199 equemene
   uint yi=(uint)get_global_id(1);
213 199 equemene
   uint sizex=(uint)get_global_size(0);
214 199 equemene
   uint sizey=(uint)get_global_size(1);
215 199 equemene

216 204 equemene
   // Perform trajectory for each pixel, exit on hit
217 199 equemene

218 217 equemene
  private float m,rs,ri,re,tho;
219 217 equemene
  private int q,raie;
220 199 equemene

221 204 equemene
  m=Mass;
222 204 equemene
  rs=2.*m;
223 204 equemene
  ri=InternalRadius;
224 204 equemene
  re=ExternalRadius;
225 204 equemene
  tho=Angle;
226 204 equemene
  q=-2;
227 209 equemene
  raie=Line;
228 204 equemene

229 217 equemene
  private float d,bmx,db,b,h;
230 218 equemene
  private float rp0,rpp,rps;
231 217 equemene
  private float phi,thi,phd,php,nr,r;
232 217 equemene
  private int nh;
233 217 equemene
  private float zp,fp;
234 199 equemene

235 199 equemene
  // Autosize for image
236 199 equemene
  bmx=1.25*re;
237 199 equemene
  b=0.;
238 199 equemene

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

241 199 equemene
  // set origin as center of image
242 199 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
243 201 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
244 199 equemene
  // angle extracted from cylindric symmetry
245 199 equemene
  phi=atanp(x,y);
246 199 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
247 199 equemene

248 204 equemene
  float up,vp,pp,us,vs,ps;
249 204 equemene

250 204 equemene
  // impact parameter
251 204 equemene
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
252 204 equemene
  // step of impact parameter;
253 209 equemene
  db=bmx/(float)(sizex);
254 204 equemene

255 209 equemene
  up=0.;
256 209 equemene
  vp=1.;
257 199 equemene

258 199 equemene
  pp=0.;
259 199 equemene
  nh=0;
260 199 equemene

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

263 218 equemene
  rps=fabs(b/us);
264 218 equemene
  rp0=rps;
265 199 equemene

266 204 equemene
  int ExitOnImpact=0;
267 199 equemene

268 199 equemene
  do
269 199 equemene
  {
270 199 equemene
     nh++;
271 199 equemene
     pp=ps;
272 199 equemene
     up=us;
273 199 equemene
     vp=vs;
274 218 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
275 218 equemene
     rpp=rps;
276 218 equemene
     rps=fabs(b/us);
277 218 equemene
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;
278 199 equemene

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

281 199 equemene
  if (ExitOnImpact==1) {
282 218 equemene
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
283 199 equemene
  }
284 199 equemene
  else
285 199 equemene
  {
286 199 equemene
     zp=0.;
287 201 equemene
     fp=0.;
288 199 equemene
  }
289 199 equemene

290 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
291 204 equemene

292 201 equemene
  zImage[yi+sizex*xi]=(float)zp;
293 204 equemene
  fImage[yi+sizex*xi]=(float)fp;
294 204 equemene
}
295 199 equemene

296 204 equemene
__kernel void Pixel(__global float *zImage,__global float *fImage,
297 204 equemene
                    __global float *Trajectories,__global int *IdLast,
298 204 equemene
                    uint ImpactParameter,uint TrackPoints,
299 204 equemene
                    float Mass,float InternalRadius,
300 204 equemene
                    float ExternalRadius,float Angle,
301 209 equemene
                    int Line)
302 204 equemene
{
303 204 equemene
   uint xi=(uint)get_global_id(0);
304 204 equemene
   uint yi=(uint)get_global_id(1);
305 204 equemene
   uint sizex=(uint)get_global_size(0);
306 204 equemene
   uint sizey=(uint)get_global_size(1);
307 204 equemene

308 204 equemene
   // Perform trajectory for each pixel
309 204 equemene

310 209 equemene
  float m,rs,ri,re,tho;
311 209 equemene
  int q,raie;
312 204 equemene

313 204 equemene
  m=Mass;
314 204 equemene
  rs=2.*m;
315 204 equemene
  ri=InternalRadius;
316 204 equemene
  re=ExternalRadius;
317 204 equemene
  tho=Angle;
318 204 equemene
  q=-2;
319 209 equemene
  raie=Line;
320 204 equemene

321 204 equemene
  float d,bmx,db,b,h;
322 204 equemene
  float phi,thi,phd,php,nr,r;
323 204 equemene
  int nh;
324 204 equemene
  float zp=0,fp=0;
325 204 equemene

326 209 equemene
  // Autosize for image, 25% greater than external radius
327 204 equemene
  bmx=1.25*re;
328 204 equemene

329 204 equemene
  // Angular step of integration
330 204 equemene
  h=4.e0f*PI/(float)TrackPoints;
331 204 equemene

332 204 equemene
  // Step of Impact Parameter
333 209 equemene
  db=bmx/(2.e0*(float)ImpactParameter);
334 204 equemene

335 204 equemene
  // set origin as center of image
336 204 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
337 204 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
338 204 equemene

339 204 equemene
  // angle extracted from cylindric symmetry
340 204 equemene
  phi=atanp(x,y);
341 204 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
342 204 equemene

343 204 equemene
  // Real Impact Parameter
344 204 equemene
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
345 204 equemene

346 204 equemene
  // Integer Impact Parameter
347 204 equemene
  uint bi=(uint)sqrt(x*x+y*y);
348 204 equemene

349 204 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
350 204 equemene

351 204 equemene
  if (bi<ImpactParameter)
352 204 equemene
  {
353 204 equemene
    do
354 204 equemene
    {
355 204 equemene
      php=phd+(float)HalfLap*PI;
356 204 equemene
      nr=php/h;
357 204 equemene
      ni=(int)nr;
358 204 equemene

359 204 equemene
      if (ni<IdLast[bi])
360 204 equemene
      {
361 204 equemene
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
362 204 equemene
      }
363 204 equemene
      else
364 204 equemene
      {
365 204 equemene
        r=Trajectories[bi*TrackPoints+ni];
366 204 equemene
      }
367 204 equemene

368 204 equemene
      if ((r<=re)&&(r>=ri))
369 204 equemene
      {
370 204 equemene
        ExitOnImpact=1;
371 204 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
372 204 equemene
      }
373 204 equemene

374 204 equemene
      HalfLap++;
375 204 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
376 204 equemene

377 204 equemene
  }
378 204 equemene

379 201 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
380 204 equemene

381 204 equemene
  zImage[yi+sizex*xi]=zp;
382 204 equemene
  fImage[yi+sizex*xi]=fp;
383 204 equemene
}
384 204 equemene

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

401 204 equemene
   // Perform trajectory for each pixel
402 204 equemene

403 209 equemene
  float m,rs,ri,re,tho;
404 209 equemene
  int q,raie;
405 204 equemene

406 204 equemene
  m=Mass;
407 204 equemene
  rs=2.*m;
408 204 equemene
  ri=InternalRadius;
409 204 equemene
  re=ExternalRadius;
410 204 equemene
  tho=Angle;
411 209 equemene
  raie=Line;
412 204 equemene

413 204 equemene
  float bmx,db,b,h;
414 204 equemene
  float phi,thi,phd;
415 204 equemene
  int nh;
416 204 equemene
  float zp=0,fp=0;
417 204 equemene

418 204 equemene
  // Autosize for image
419 204 equemene
  bmx=1.25*re;
420 204 equemene

421 204 equemene
  // Angular step of integration
422 204 equemene
  h=4.e0f*PI/(float)TrackPoints;
423 204 equemene

424 204 equemene
  // impact parameter
425 204 equemene
  b=(float)bi/(float)bmaxi*bmx;
426 209 equemene
  db=bmx/(2.e0*(float)bmaxi);
427 204 equemene

428 204 equemene
  phi=2.*PI/(float)imx*(float)i;
429 204 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
430 204 equemene
  int yi=(int)((float)bi*sin(phi))+bmaxi;
431 204 equemene
  int xi=(int)((float)bi*cos(phi))+bmaxi;
432 204 equemene

433 204 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
434 204 equemene
  float php,nr,r;
435 204 equemene

436 204 equemene
  do
437 204 equemene
  {
438 204 equemene
     php=phd+(float)HalfLap*PI;
439 204 equemene
     nr=php/h;
440 204 equemene
     ni=(int)nr;
441 204 equemene

442 204 equemene
     if (ni<IdLast[bi])
443 204 equemene
     {
444 204 equemene
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
445 204 equemene
     }
446 204 equemene
     else
447 204 equemene
     {
448 204 equemene
        r=Trajectories[bi*TrackPoints+ni];
449 204 equemene
     }
450 204 equemene

451 204 equemene
     if ((r<=re)&&(r>=ri))
452 204 equemene
     {
453 204 equemene
        ExitOnImpact=1;
454 204 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
455 204 equemene
     }
456 204 equemene

457 204 equemene
     HalfLap++;
458 204 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
459 204 equemene

460 204 equemene
  zImage[yi+2*bmaxi*xi]=zp;
461 204 equemene
  fImage[yi+2*bmaxi*xi]=fp;
462 204 equemene

463 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
464 204 equemene

465 204 equemene
}
466 204 equemene

467 204 equemene
__kernel void Trajectory(__global float *Trajectories,
468 204 equemene
                         __global int *IdLast,int TrackPoints,
469 204 equemene
                         float Mass,float InternalRadius,
470 204 equemene
                         float ExternalRadius,float Angle,
471 209 equemene
                         int Line)
472 204 equemene
{
473 204 equemene
  // Integer Impact Parameter ID
474 204 equemene
  int bi=get_global_id(0);
475 204 equemene
  // Integer Impact Parameter Size (half of image)
476 204 equemene
  int bmaxi=get_global_size(0);
477 204 equemene

478 204 equemene
  // Perform trajectory for each pixel
479 204 equemene

480 209 equemene
  float m,rs,ri,re,tho;
481 209 equemene
  int raie,q;
482 204 equemene

483 204 equemene
  m=Mass;
484 204 equemene
  rs=2.*m;
485 204 equemene
  ri=InternalRadius;
486 204 equemene
  re=ExternalRadius;
487 204 equemene
  tho=Angle;
488 204 equemene
  q=-2;
489 209 equemene
  raie=Line;
490 204 equemene

491 204 equemene
  float d,bmx,db,b,h;
492 204 equemene
  float phi,thi,phd,php,nr,r;
493 204 equemene
  int nh;
494 204 equemene
  float zp,fp;
495 204 equemene

496 204 equemene
  // Autosize for image
497 204 equemene
  bmx=1.25*re;
498 204 equemene

499 204 equemene
  // Angular step of integration
500 204 equemene
  h=4.e0f*PI/(float)TrackPoints;
501 204 equemene

502 204 equemene
  // impact parameter
503 204 equemene
  b=(float)bi/(float)bmaxi*bmx;
504 204 equemene

505 204 equemene
  float up,vp,pp,us,vs,ps;
506 204 equemene

507 209 equemene
  up=0.;
508 209 equemene
  vp=1.;
509 204 equemene

510 204 equemene
  pp=0.;
511 204 equemene
  nh=0;
512 204 equemene

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

515 204 equemene
  // b versus us
516 204 equemene
  float bvus=fabs(b/us);
517 204 equemene
  float bvus0=bvus;
518 204 equemene
  Trajectories[bi*TrackPoints+nh]=bvus;
519 204 equemene

520 204 equemene
  do
521 204 equemene
  {
522 204 equemene
     nh++;
523 204 equemene
     pp=ps;
524 204 equemene
     up=us;
525 204 equemene
     vp=vs;
526 204 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
527 204 equemene
     bvus=fabs(b/us);
528 204 equemene
     Trajectories[bi*TrackPoints+nh]=bvus;
529 204 equemene

530 204 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
531 204 equemene

532 204 equemene
  IdLast[bi]=nh;
533 204 equemene

534 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
535 199 equemene

536 199 equemene
}
537 211 equemene
"""
538 199 equemene
539 217 equemene
def KernelCodeCuda():
540 217 equemene
    BlobCUDA= """
541 217 equemene

542 217 equemene
#define PI (float)3.14159265359
543 217 equemene
#define nbr 256
544 217 equemene

545 217 equemene
#define EINSTEIN 0
546 217 equemene
#define NEWTON 1
547 217 equemene

548 217 equemene
#define TRACKPOINTS 2048
549 217 equemene

550 217 equemene
__device__ float nothing(float x)
551 217 equemene
{
552 217 equemene
  return(x);
553 217 equemene
}
554 217 equemene

555 217 equemene
__device__ float atanp(float x,float y)
556 217 equemene
{
557 217 equemene
  float angle;
558 217 equemene

559 217 equemene
  angle=atan2(y,x);
560 217 equemene

561 217 equemene
  if (angle<0.e0f)
562 217 equemene
    {
563 217 equemene
      angle+=(float)2.e0f*PI;
564 217 equemene
    }
565 217 equemene

566 217 equemene
  return(angle);
567 217 equemene
}
568 217 equemene

569 217 equemene
__device__ float f(float v)
570 217 equemene
{
571 217 equemene
  return(v);
572 217 equemene
}
573 217 equemene

574 217 equemene
#if PHYSICS == NEWTON
575 217 equemene
__device__ float g(float u,float m,float b)
576 217 equemene
{
577 217 equemene
  return (-u);
578 217 equemene
}
579 217 equemene
#else
580 217 equemene
__device__ float g(float u,float m,float b)
581 217 equemene
{
582 217 equemene
  return (3.e0f*m/b*pow(u,2)-u);
583 217 equemene
}
584 217 equemene
#endif
585 217 equemene

586 217 equemene
__device__ void calcul(float *us,float *vs,float up,float vp,
587 217 equemene
                       float h,float m,float b)
588 217 equemene
{
589 217 equemene
  float c0,c1,c2,c3,d0,d1,d2,d3;
590 217 equemene

591 217 equemene
  c0=h*f(vp);
592 217 equemene
  c1=h*f(vp+c0/2.);
593 217 equemene
  c2=h*f(vp+c1/2.);
594 217 equemene
  c3=h*f(vp+c2);
595 217 equemene
  d0=h*g(up,m,b);
596 217 equemene
  d1=h*g(up+d0/2.,m,b);
597 217 equemene
  d2=h*g(up+d1/2.,m,b);
598 217 equemene
  d3=h*g(up+d2,m,b);
599 217 equemene

600 217 equemene
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
601 217 equemene
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
602 217 equemene
}
603 217 equemene

604 217 equemene
__device__ void rungekutta(float *ps,float *us,float *vs,
605 217 equemene
                           float pp,float up,float vp,
606 217 equemene
                           float h,float m,float b)
607 217 equemene
{
608 217 equemene
  calcul(us,vs,up,vp,h,m,b);
609 217 equemene
  *ps=pp+h;
610 217 equemene
}
611 217 equemene

612 217 equemene
__device__ float decalage_spectral(float r,float b,float phi,
613 217 equemene
                                   float tho,float m)
614 217 equemene
{
615 217 equemene
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
616 217 equemene
}
617 217 equemene

618 217 equemene
__device__ float spectre(float rf,int q,float b,float db,
619 217 equemene
                         float h,float r,float m,float bss)
620 217 equemene
{
621 217 equemene
  float flx;
622 217 equemene

623 217 equemene
  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
624 217 equemene
  return(flx);
625 217 equemene
}
626 217 equemene

627 217 equemene
__device__ float spectre_cn(float rf32,float b32,float db32,
628 217 equemene
                            float h32,float r32,float m32,float bss32)
629 217 equemene
{
630 217 equemene

631 217 equemene
#define MYFLOAT float
632 217 equemene

633 217 equemene
  MYFLOAT rf=(MYFLOAT)(rf32);
634 217 equemene
  MYFLOAT b=(MYFLOAT)(b32);
635 217 equemene
  MYFLOAT db=(MYFLOAT)(db32);
636 217 equemene
  MYFLOAT h=(MYFLOAT)(h32);
637 217 equemene
  MYFLOAT r=(MYFLOAT)(r32);
638 217 equemene
  MYFLOAT m=(MYFLOAT)(m32);
639 217 equemene
  MYFLOAT bss=(MYFLOAT)(bss32);
640 217 equemene

641 217 equemene
  MYFLOAT flx;
642 217 equemene
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
643 217 equemene
  int fi,posfreq;
644 217 equemene

645 217 equemene
#define planck 6.62e-34
646 217 equemene
#define k 1.38e-23
647 217 equemene
#define c2 9.e16
648 217 equemene
#define temp 3.e7
649 217 equemene
#define m_point 1.
650 217 equemene

651 217 equemene
#define lplanck (log(6.62)-34.*log(10.))
652 217 equemene
#define lk (log(1.38)-23.*log(10.))
653 217 equemene
#define lc2 (log(9.)+16.*log(10.))
654 217 equemene

655 217 equemene
  MYFLOAT v=1.-3./r;
656 217 equemene

657 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 ));
658 217 equemene

659 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)));
660 217 equemene

661 217 equemene
  flux_int=0.;
662 217 equemene
  flx=0.;
663 217 equemene

664 217 equemene
  for (fi=0;fi<nbr;fi++)
665 217 equemene
    {
666 217 equemene
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
667 217 equemene
      nu_rec=nu_em*rf;
668 217 equemene
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
669 217 equemene
      if ((posfreq>0)&&(posfreq<nbr))
670 217 equemene
        {
671 217 equemene
          // Initial version
672 217 equemene
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
673 217 equemene
          // Version with log used
674 217 equemene
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
675 217 equemene
          // flux_int*=pow(rf,3)*b*db*h;
676 217 equemene
          //flux_int*=exp(3.*log(rf))*b*db*h;
677 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;
678 217 equemene

679 217 equemene
          flx+=flux_int;
680 217 equemene
        }
681 217 equemene
    }
682 217 equemene

683 217 equemene
  return((float)(flx));
684 217 equemene
}
685 217 equemene

686 217 equemene
__device__ void impact(float phi,float r,float b,float tho,float m,
687 217 equemene
                       float *zp,float *fp,
688 217 equemene
                       int q,float db,
689 217 equemene
                       float h,int raie)
690 217 equemene
{
691 217 equemene
  float flx,rf,bss;
692 217 equemene

693 217 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
694 217 equemene

695 217 equemene
  if (raie==0)
696 217 equemene
    {
697 217 equemene
      bss=1.e19;
698 217 equemene
      flx=spectre_cn(rf,b,db,h,r,m,bss);
699 217 equemene
    }
700 217 equemene
  else
701 217 equemene
    {
702 217 equemene
      bss=2.;
703 217 equemene
      flx=spectre(rf,q,b,db,h,r,m,bss);
704 217 equemene
    }
705 217 equemene

706 217 equemene
  *zp=1./rf;
707 217 equemene
  *fp=flx;
708 217 equemene

709 217 equemene
}
710 217 equemene

711 217 equemene
__global__ void EachPixel(float *zImage,float *fImage,
712 217 equemene
                          float Mass,float InternalRadius,
713 217 equemene
                          float ExternalRadius,float Angle,
714 217 equemene
                          int Line)
715 217 equemene
{
716 218 equemene
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
717 218 equemene
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
718 217 equemene
   uint sizex=(uint)gridDim.x*blockDim.x;
719 217 equemene
   uint sizey=(uint)gridDim.y*blockDim.y;
720 217 equemene

721 217 equemene
   // Perform trajectory for each pixel, exit on hit
722 217 equemene

723 217 equemene
  float m,rs,ri,re,tho;
724 217 equemene
  int q,raie;
725 217 equemene

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

734 217 equemene
  float d,bmx,db,b,h;
735 218 equemene
  float rp0,rpp,rps;
736 217 equemene
  float phi,thi,phd,php,nr,r;
737 217 equemene
  int nh;
738 217 equemene
  float zp,fp;
739 217 equemene

740 217 equemene
  // Autosize for image
741 217 equemene
  bmx=1.25*re;
742 217 equemene
  b=0.;
743 217 equemene

744 217 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
745 217 equemene

746 217 equemene
  // set origin as center of image
747 217 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
748 217 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
749 217 equemene
  // angle extracted from cylindric symmetry
750 217 equemene
  phi=atanp(x,y);
751 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
752 217 equemene

753 217 equemene
  float up,vp,pp,us,vs,ps;
754 217 equemene

755 217 equemene
  // impact parameter
756 217 equemene
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
757 217 equemene
  // step of impact parameter;
758 217 equemene
//  db=bmx/(float)(sizex/2);
759 217 equemene
  db=bmx/(float)(sizex);
760 217 equemene

761 217 equemene
  up=0.;
762 217 equemene
  vp=1.;
763 217 equemene

764 217 equemene
  pp=0.;
765 217 equemene
  nh=0;
766 217 equemene

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

769 218 equemene
  rps=fabs(b/us);
770 218 equemene
  rp0=rps;
771 217 equemene

772 217 equemene
  int ExitOnImpact=0;
773 217 equemene

774 217 equemene
  do
775 217 equemene
  {
776 217 equemene
     nh++;
777 217 equemene
     pp=ps;
778 217 equemene
     up=us;
779 217 equemene
     vp=vs;
780 218 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
781 218 equemene
     rpp=rps;
782 218 equemene
     rps=fabs(b/us);
783 218 equemene
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;
784 217 equemene

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

787 217 equemene
  if (ExitOnImpact==1) {
788 218 equemene
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
789 217 equemene
  }
790 217 equemene
  else
791 217 equemene
  {
792 217 equemene
     zp=0.;
793 217 equemene
     fp=0.;
794 217 equemene
  }
795 217 equemene

796 218 equemene
  __syncthreads();
797 218 equemene

798 217 equemene
  zImage[yi+sizex*xi]=(float)zp;
799 217 equemene
  fImage[yi+sizex*xi]=(float)fp;
800 217 equemene
}
801 217 equemene

802 217 equemene
__global__ void Pixel(float *zImage,float *fImage,
803 217 equemene
                      float *Trajectories,int *IdLast,
804 217 equemene
                      uint ImpactParameter,uint TrackPoints,
805 217 equemene
                      float Mass,float InternalRadius,
806 217 equemene
                      float ExternalRadius,float Angle,
807 217 equemene
                      int Line)
808 217 equemene
{
809 219 equemene
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
810 219 equemene
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
811 219 equemene
   uint sizex=(uint)gridDim.x*blockDim.x;
812 219 equemene
   uint sizey=(uint)gridDim.y*blockDim.y;
813 219 equemene
//  uint xi=(uint)blockIdx.x;
814 219 equemene
//  uint yi=(uint)blockIdx.y;
815 219 equemene
//  uint sizex=(uint)gridDim.x*blockDim.x;
816 219 equemene
//  uint sizey=(uint)gridDim.y*blockDim.y;
817 217 equemene

818 217 equemene
  // Perform trajectory for each pixel
819 217 equemene

820 217 equemene
  float m,rs,ri,re,tho;
821 217 equemene
  int q,raie;
822 217 equemene

823 217 equemene
  m=Mass;
824 217 equemene
  rs=2.*m;
825 217 equemene
  ri=InternalRadius;
826 217 equemene
  re=ExternalRadius;
827 217 equemene
  tho=Angle;
828 217 equemene
  q=-2;
829 217 equemene
  raie=Line;
830 217 equemene

831 217 equemene
  float d,bmx,db,b,h;
832 217 equemene
  float phi,thi,phd,php,nr,r;
833 217 equemene
  int nh;
834 217 equemene
  float zp=0,fp=0;
835 217 equemene
  // Autosize for image, 25% greater than external radius
836 217 equemene
  bmx=1.25*re;
837 217 equemene

838 217 equemene
  // Angular step of integration
839 217 equemene
  h=4.e0f*PI/(float)TrackPoints;
840 217 equemene

841 217 equemene
  // Step of Impact Parameter
842 217 equemene
  db=bmx/(2.e0*(float)ImpactParameter);
843 217 equemene

844 217 equemene
  // set origin as center of image
845 217 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
846 217 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
847 217 equemene
  // angle extracted from cylindric symmetry
848 217 equemene
  phi=atanp(x,y);
849 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
850 217 equemene

851 217 equemene
  // Real Impact Parameter
852 217 equemene
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
853 217 equemene

854 217 equemene
  // Integer Impact Parameter
855 217 equemene
  uint bi=(uint)sqrt(x*x+y*y);
856 217 equemene

857 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
858 217 equemene

859 217 equemene
  if (bi<ImpactParameter)
860 217 equemene
  {
861 217 equemene
    do
862 217 equemene
    {
863 217 equemene
      php=phd+(float)HalfLap*PI;
864 217 equemene
      nr=php/h;
865 217 equemene
      ni=(int)nr;
866 217 equemene

867 217 equemene
      if (ni<IdLast[bi])
868 217 equemene
      {
869 217 equemene
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
870 217 equemene
      }
871 217 equemene
      else
872 217 equemene
      {
873 217 equemene
        r=Trajectories[bi*TrackPoints+ni];
874 217 equemene
      }
875 217 equemene

876 217 equemene
      if ((r<=re)&&(r>=ri))
877 217 equemene
      {
878 217 equemene
        ExitOnImpact=1;
879 217 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
880 217 equemene
      }
881 217 equemene

882 217 equemene
      HalfLap++;
883 217 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
884 217 equemene

885 217 equemene
  }
886 217 equemene

887 217 equemene
  zImage[yi+sizex*xi]=zp;
888 217 equemene
  fImage[yi+sizex*xi]=fp;
889 217 equemene
}
890 217 equemene

891 217 equemene
__global__ void Circle(float *Trajectories,int *IdLast,
892 217 equemene
                       float *zImage,float *fImage,
893 217 equemene
                       int TrackPoints,
894 217 equemene
                       float Mass,float InternalRadius,
895 217 equemene
                       float ExternalRadius,float Angle,
896 217 equemene
                       int Line)
897 217 equemene
{
898 217 equemene
   // Integer Impact Parameter ID
899 219 equemene
   int bi=blockIdx.x*blockDim.x+threadIdx.x;
900 217 equemene
   // Integer points on circle
901 219 equemene
   int i=blockIdx.y*blockDim.y+threadIdx.y;
902 217 equemene
   // Integer Impact Parameter Size (half of image)
903 217 equemene
   int bmaxi=gridDim.x*blockDim.x;
904 217 equemene
   // Integer Points on circle
905 217 equemene
   int imx=gridDim.y*blockDim.y;
906 217 equemene

907 219 equemene

908 219 equemene
   // Integer Impact Parameter ID
909 219 equemene
//   int bi=blockIdx.x;
910 219 equemene
   // Integer points on circle
911 219 equemene
//   int i=blockIdx.y;
912 219 equemene
   // Integer Impact Parameter Size (half of image)
913 219 equemene
//   int bmaxi=gridDim.x*blockDim.x;
914 219 equemene
   // Integer Points on circle
915 219 equemene
//   int imx=gridDim.y*blockDim.y;
916 219 equemene

917 217 equemene
   // Perform trajectory for each pixel
918 217 equemene

919 217 equemene
  float m,rs,ri,re,tho;
920 217 equemene
  int q,raie;
921 217 equemene

922 217 equemene
  m=Mass;
923 217 equemene
  rs=2.*m;
924 217 equemene
  ri=InternalRadius;
925 217 equemene
  re=ExternalRadius;
926 217 equemene
  tho=Angle;
927 217 equemene
  raie=Line;
928 217 equemene

929 217 equemene
  float bmx,db,b,h;
930 217 equemene
  float phi,thi,phd;
931 217 equemene
  int nh;
932 217 equemene
  float zp=0,fp=0;
933 217 equemene

934 217 equemene
  // Autosize for image
935 217 equemene
  bmx=1.25*re;
936 217 equemene

937 217 equemene
  // Angular step of integration
938 217 equemene
  h=4.e0f*PI/(float)TrackPoints;
939 217 equemene

940 217 equemene
  // impact parameter
941 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
942 217 equemene
  db=bmx/(2.e0*(float)bmaxi);
943 217 equemene

944 217 equemene
  phi=2.*PI/(float)imx*(float)i;
945 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
946 217 equemene
  int yi=(int)((float)bi*sin(phi))+bmaxi;
947 217 equemene
  int xi=(int)((float)bi*cos(phi))+bmaxi;
948 217 equemene

949 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
950 217 equemene
  float php,nr,r;
951 217 equemene

952 217 equemene
  do
953 217 equemene
  {
954 217 equemene
     php=phd+(float)HalfLap*PI;
955 217 equemene
     nr=php/h;
956 217 equemene
     ni=(int)nr;
957 217 equemene

958 217 equemene
     if (ni<IdLast[bi])
959 217 equemene
     {
960 217 equemene
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
961 217 equemene
     }
962 217 equemene
     else
963 217 equemene
     {
964 217 equemene
        r=Trajectories[bi*TrackPoints+ni];
965 217 equemene
     }
966 217 equemene

967 217 equemene
     if ((r<=re)&&(r>=ri))
968 217 equemene
     {
969 217 equemene
        ExitOnImpact=1;
970 217 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
971 217 equemene
     }
972 217 equemene

973 217 equemene
     HalfLap++;
974 217 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
975 217 equemene

976 217 equemene
  zImage[yi+2*bmaxi*xi]=zp;
977 217 equemene
  fImage[yi+2*bmaxi*xi]=fp;
978 217 equemene

979 217 equemene
}
980 217 equemene

981 217 equemene
__global__ void Trajectory(float *Trajectories,
982 217 equemene
                           int *IdLast,int TrackPoints,
983 217 equemene
                           float Mass,float InternalRadius,
984 217 equemene
                           float ExternalRadius,float Angle,
985 217 equemene
                           int Line)
986 217 equemene
{
987 217 equemene
  // Integer Impact Parameter ID
988 219 equemene
  //  int bi=blockIdx.x;
989 219 equemene
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
990 217 equemene
  // Integer Impact Parameter Size (half of image)
991 217 equemene
  int bmaxi=gridDim.x*blockDim.x;
992 217 equemene

993 217 equemene
  // Perform trajectory for each pixel
994 217 equemene

995 217 equemene
  float m,rs,ri,re,tho;
996 217 equemene
  int raie,q;
997 217 equemene

998 217 equemene
  m=Mass;
999 217 equemene
  rs=2.*m;
1000 217 equemene
  ri=InternalRadius;
1001 217 equemene
  re=ExternalRadius;
1002 217 equemene
  tho=Angle;
1003 217 equemene
  q=-2;
1004 217 equemene
  raie=Line;
1005 217 equemene

1006 217 equemene
  float d,bmx,db,b,h;
1007 217 equemene
  float phi,thi,phd,php,nr,r;
1008 217 equemene
  int nh;
1009 217 equemene
  float zp,fp;
1010 217 equemene

1011 217 equemene
  // Autosize for image
1012 217 equemene
  bmx=1.25*re;
1013 217 equemene

1014 217 equemene
  // Angular step of integration
1015 217 equemene
  h=4.e0f*PI/(float)TrackPoints;
1016 217 equemene

1017 217 equemene
  // impact parameter
1018 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
1019 217 equemene

1020 217 equemene
  float up,vp,pp,us,vs,ps;
1021 217 equemene

1022 217 equemene
  up=0.;
1023 217 equemene
  vp=1.;
1024 217 equemene

1025 217 equemene
  pp=0.;
1026 217 equemene
  nh=0;
1027 217 equemene

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

1030 217 equemene
  // b versus us
1031 217 equemene
  float bvus=fabs(b/us);
1032 217 equemene
  float bvus0=bvus;
1033 217 equemene
  Trajectories[bi*TrackPoints+nh]=bvus;
1034 217 equemene

1035 217 equemene
  do
1036 217 equemene
  {
1037 217 equemene
     nh++;
1038 217 equemene
     pp=ps;
1039 217 equemene
     up=us;
1040 217 equemene
     vp=vs;
1041 217 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1042 217 equemene
     bvus=fabs(b/us);
1043 217 equemene
     Trajectories[bi*TrackPoints+nh]=bvus;
1044 217 equemene

1045 217 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
1046 217 equemene

1047 217 equemene
  IdLast[bi]=nh;
1048 217 equemene

1049 217 equemene
}
1050 217 equemene
"""
1051 217 equemene
    return(BlobCUDA)
1052 217 equemene
1053 211 equemene
# def ImageOutput(sigma,prefix):
1054 211 equemene
#     from PIL import Image
1055 211 equemene
#     Max=sigma.max()
1056 211 equemene
#     Min=sigma.min()
1057 199 equemene
1058 211 equemene
#     # Normalize value as 8bits Integer
1059 211 equemene
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1060 211 equemene
#     image = Image.fromarray(SigmaInt)
1061 211 equemene
#     image.save("%s.jpg" % prefix)
1062 211 equemene
1063 211 equemene
def ImageOutput(sigma,prefix,Colors):
1064 211 equemene
    import matplotlib.pyplot as plt
1065 211 equemene
    if Colors == 'Red2Yellow':
1066 211 equemene
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1067 211 equemene
    else:
1068 211 equemene
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1069 211 equemene
1070 204 equemene
def BlackHoleCL(zImage,fImage,InputCL):
1071 199 equemene
1072 199 equemene
    kernel_params = {}
1073 199 equemene
1074 204 equemene
    print(InputCL)
1075 204 equemene
1076 204 equemene
    Device=InputCL['Device']
1077 204 equemene
    GpuStyle=InputCL['GpuStyle']
1078 204 equemene
    VariableType=InputCL['VariableType']
1079 204 equemene
    Size=InputCL['Size']
1080 204 equemene
    Mass=InputCL['Mass']
1081 204 equemene
    InternalRadius=InputCL['InternalRadius']
1082 204 equemene
    ExternalRadius=InputCL['ExternalRadius']
1083 204 equemene
    Angle=InputCL['Angle']
1084 204 equemene
    Method=InputCL['Method']
1085 204 equemene
    TrackPoints=InputCL['TrackPoints']
1086 211 equemene
    Physics=InputCL['Physics']
1087 204 equemene
1088 211 equemene
    PhysicsList=DictionariesAPI()
1089 211 equemene
1090 204 equemene
    if InputCL['BlackBody']:
1091 209 equemene
        # Spectrum is Black Body one
1092 209 equemene
        Line=0
1093 204 equemene
    else:
1094 209 equemene
        # Spectrum is Monochromatic Line one
1095 209 equemene
        Line=1
1096 204 equemene
1097 217 equemene
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1098 217 equemene
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1099 204 equemene
1100 199 equemene
    # Je detecte un peripherique GPU dans la liste des peripheriques
1101 199 equemene
    Id=0
1102 199 equemene
    HasXPU=False
1103 199 equemene
    for platform in cl.get_platforms():
1104 199 equemene
        for device in platform.get_devices():
1105 199 equemene
            if Id==Device:
1106 199 equemene
                XPU=device
1107 217 equemene
                print("CPU/GPU selected: ",device.name.lstrip())
1108 199 equemene
                HasXPU=True
1109 199 equemene
            Id+=1
1110 199 equemene
1111 199 equemene
    if HasXPU==False:
1112 217 equemene
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1113 199 equemene
        sys.exit()
1114 199 equemene
1115 199 equemene
    ctx = cl.Context([XPU])
1116 199 equemene
    queue = cl.CommandQueue(ctx,
1117 199 equemene
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1118 199 equemene
1119 199 equemene
1120 211 equemene
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1121 211 equemene
1122 211 equemene
    BuildOptions="-cl-mad-enable -DPHYSICS=%i " % (PhysicsList[Physics])
1123 211 equemene
1124 211 equemene
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1125 211 equemene
1126 199 equemene
    # Je recupere les flag possibles pour les buffers
1127 199 equemene
    mf = cl.mem_flags
1128 199 equemene
1129 204 equemene
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
1130 201 equemene
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1131 201 equemene
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1132 204 equemene
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1133 199 equemene
1134 199 equemene
    start_time=time.time()
1135 199 equemene
1136 204 equemene
    if Method=='EachPixel':
1137 204 equemene
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1138 204 equemene
                                       None,zImageCL,fImageCL,
1139 204 equemene
                                       numpy.float32(Mass),
1140 204 equemene
                                       numpy.float32(InternalRadius),
1141 204 equemene
                                       numpy.float32(ExternalRadius),
1142 204 equemene
                                       numpy.float32(Angle),
1143 209 equemene
                                       numpy.int32(Line))
1144 204 equemene
        CLLaunch.wait()
1145 204 equemene
    elif Method=='TrajectoCircle':
1146 204 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1147 204 equemene
                                        None,TrajectoriesCL,IdLastCL,
1148 204 equemene
                                        numpy.uint32(Trajectories.shape[1]),
1149 204 equemene
                                        numpy.float32(Mass),
1150 204 equemene
                                        numpy.float32(InternalRadius),
1151 204 equemene
                                        numpy.float32(ExternalRadius),
1152 204 equemene
                                        numpy.float32(Angle),
1153 209 equemene
                                        numpy.int32(Line))
1154 204 equemene
1155 204 equemene
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1156 204 equemene
                                           zImage.shape[0]*4),None,
1157 204 equemene
                                    TrajectoriesCL,IdLastCL,
1158 204 equemene
                                    zImageCL,fImageCL,
1159 204 equemene
                                    numpy.uint32(Trajectories.shape[1]),
1160 204 equemene
                                    numpy.float32(Mass),
1161 204 equemene
                                    numpy.float32(InternalRadius),
1162 204 equemene
                                    numpy.float32(ExternalRadius),
1163 204 equemene
                                    numpy.float32(Angle),
1164 209 equemene
                                    numpy.int32(Line))
1165 204 equemene
        CLLaunch.wait()
1166 204 equemene
    else:
1167 204 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1168 204 equemene
                                        None,TrajectoriesCL,IdLastCL,
1169 204 equemene
                                        numpy.uint32(Trajectories.shape[1]),
1170 204 equemene
                                        numpy.float32(Mass),
1171 204 equemene
                                        numpy.float32(InternalRadius),
1172 204 equemene
                                        numpy.float32(ExternalRadius),
1173 204 equemene
                                        numpy.float32(Angle),
1174 209 equemene
                                        numpy.int32(Line))
1175 204 equemene
1176 204 equemene
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],
1177 217 equemene
                                          zImage.shape[1]),None,
1178 204 equemene
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1179 204 equemene
                                   numpy.uint32(Trajectories.shape[0]),
1180 204 equemene
                                   numpy.uint32(Trajectories.shape[1]),
1181 204 equemene
                                        numpy.float32(Mass),
1182 204 equemene
                                        numpy.float32(InternalRadius),
1183 204 equemene
                                        numpy.float32(ExternalRadius),
1184 204 equemene
                                        numpy.float32(Angle),
1185 209 equemene
                                        numpy.int32(Line))
1186 204 equemene
        CLLaunch.wait()
1187 204 equemene
1188 218 equemene
    compute = time.time()-start_time
1189 199 equemene
1190 201 equemene
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1191 201 equemene
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1192 204 equemene
    cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1193 204 equemene
    cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1194 218 equemene
    elapsed = time.time()-start_time
1195 218 equemene
    print("\nCompute Time : %f" % compute)
1196 218 equemene
    print("Elapsed Time : %f\n" % elapsed)
1197 211 equemene
1198 211 equemene
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1199 211 equemene
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1200 219 equemene
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1201 219 equemene
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1202 201 equemene
    zImageCL.release()
1203 201 equemene
    fImageCL.release()
1204 204 equemene
1205 204 equemene
    TrajectoriesCL.release()
1206 204 equemene
    IdLastCL.release()
1207 204 equemene
1208 199 equemene
    return(elapsed)
1209 199 equemene
1210 217 equemene
def BlackHoleCUDA(zImage,fImage,InputCL):
1211 217 equemene
1212 217 equemene
    kernel_params = {}
1213 217 equemene
1214 217 equemene
    print(InputCL)
1215 217 equemene
1216 217 equemene
    Device=InputCL['Device']
1217 217 equemene
    GpuStyle=InputCL['GpuStyle']
1218 217 equemene
    VariableType=InputCL['VariableType']
1219 217 equemene
    Size=InputCL['Size']
1220 217 equemene
    Mass=InputCL['Mass']
1221 217 equemene
    InternalRadius=InputCL['InternalRadius']
1222 217 equemene
    ExternalRadius=InputCL['ExternalRadius']
1223 217 equemene
    Angle=InputCL['Angle']
1224 217 equemene
    Method=InputCL['Method']
1225 217 equemene
    TrackPoints=InputCL['TrackPoints']
1226 217 equemene
    Physics=InputCL['Physics']
1227 217 equemene
1228 217 equemene
    PhysicsList=DictionariesAPI()
1229 217 equemene
1230 217 equemene
    if InputCL['BlackBody']:
1231 217 equemene
        # Spectrum is Black Body one
1232 217 equemene
        Line=0
1233 217 equemene
    else:
1234 217 equemene
        # Spectrum is Monochromatic Line one
1235 217 equemene
        Line=1
1236 217 equemene
1237 217 equemene
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1238 217 equemene
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1239 217 equemene
1240 217 equemene
    try:
1241 217 equemene
        # For PyCUDA import
1242 217 equemene
        import pycuda.driver as cuda
1243 217 equemene
        from pycuda.compiler import SourceModule
1244 217 equemene
1245 217 equemene
        cuda.init()
1246 217 equemene
        for Id in range(cuda.Device.count()):
1247 217 equemene
            if Id==Device:
1248 217 equemene
                XPU=cuda.Device(Id)
1249 217 equemene
                print("GPU selected %s" % XPU.name())
1250 217 equemene
        print
1251 217 equemene
1252 217 equemene
    except ImportError:
1253 217 equemene
        print("Platform does not seem to support CUDA")
1254 217 equemene
1255 217 equemene
    Context=XPU.make_context()
1256 217 equemene
1257 217 equemene
    try:
1258 217 equemene
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i' % (PhysicsList[Physics])])
1259 217 equemene
        print("Compilation seems to be OK")
1260 217 equemene
    except:
1261 217 equemene
        print("Compilation seems to break")
1262 217 equemene
1263 217 equemene
    EachPixelCU=mod.get_function("EachPixel")
1264 217 equemene
    TrajectoryCU=mod.get_function("Trajectory")
1265 217 equemene
    PixelCU=mod.get_function("Pixel")
1266 217 equemene
    CircleCU=mod.get_function("Circle")
1267 217 equemene
1268 217 equemene
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1269 217 equemene
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1270 217 equemene
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1271 217 equemene
    cuda.memcpy_htod(zImageCU, zImage)
1272 217 equemene
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1273 217 equemene
    cuda.memcpy_htod(zImageCU, fImage)
1274 217 equemene
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1275 217 equemene
    cuda.memcpy_htod(IdLastCU, IdLast)
1276 217 equemene
1277 217 equemene
    start_time=time.time()
1278 217 equemene
1279 217 equemene
    if Method=='EachPixel':
1280 217 equemene
        EachPixelCU(zImageCU,fImageCU,
1281 217 equemene
                    numpy.float32(Mass),
1282 217 equemene
                    numpy.float32(InternalRadius),
1283 217 equemene
                    numpy.float32(ExternalRadius),
1284 217 equemene
                    numpy.float32(Angle),
1285 217 equemene
                    numpy.int32(Line),
1286 219 equemene
                    grid=(zImage.shape[0]/32,zImage.shape[1]/32),
1287 219 equemene
                    block=(32,32,1))
1288 217 equemene
    elif Method=='TrajectoCircle':
1289 217 equemene
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1290 217 equemene
                     numpy.uint32(Trajectories.shape[1]),
1291 217 equemene
                     numpy.float32(Mass),
1292 217 equemene
                     numpy.float32(InternalRadius),
1293 217 equemene
                     numpy.float32(ExternalRadius),
1294 217 equemene
                     numpy.float32(Angle),
1295 217 equemene
                     numpy.int32(Line),
1296 219 equemene
                     grid=(Trajectories.shape[0]/32,1),
1297 219 equemene
                     block=(32,1,1))
1298 217 equemene
1299 217 equemene
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1300 217 equemene
                 numpy.uint32(Trajectories.shape[1]),
1301 217 equemene
                 numpy.float32(Mass),
1302 217 equemene
                 numpy.float32(InternalRadius),
1303 217 equemene
                 numpy.float32(ExternalRadius),
1304 217 equemene
                 numpy.float32(Angle),
1305 217 equemene
                 numpy.int32(Line),
1306 219 equemene
                 grid=(Trajectories.shape[0]/32,zImage.shape[0]*4/32),
1307 219 equemene
                 block=(32,32,1))
1308 217 equemene
    else:
1309 217 equemene
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1310 217 equemene
                     numpy.uint32(Trajectories.shape[1]),
1311 217 equemene
                     numpy.float32(Mass),
1312 217 equemene
                     numpy.float32(InternalRadius),
1313 217 equemene
                     numpy.float32(ExternalRadius),
1314 217 equemene
                     numpy.float32(Angle),
1315 217 equemene
                     numpy.int32(Line),
1316 219 equemene
                     grid=(Trajectories.shape[0]/32,1),
1317 219 equemene
                     block=(32,1,1))
1318 217 equemene
1319 217 equemene
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1320 217 equemene
                numpy.uint32(Trajectories.shape[0]),
1321 217 equemene
                numpy.uint32(Trajectories.shape[1]),
1322 217 equemene
                numpy.float32(Mass),
1323 217 equemene
                numpy.float32(InternalRadius),
1324 217 equemene
                numpy.float32(ExternalRadius),
1325 217 equemene
                numpy.float32(Angle),
1326 217 equemene
                numpy.int32(Line),
1327 219 equemene
                grid=(zImage.shape[0]/32,zImage.shape[1]/32,1),
1328 219 equemene
                block=(32,32,1))
1329 217 equemene
1330 217 equemene
1331 218 equemene
    compute = time.time()-start_time
1332 217 equemene
1333 217 equemene
    cuda.memcpy_dtoh(zImage,zImageCU)
1334 217 equemene
    cuda.memcpy_dtoh(fImage,fImageCU)
1335 218 equemene
    elapsed = time.time()-start_time
1336 218 equemene
    print("\nCompute Time : %f" % compute)
1337 218 equemene
    print("Elapsed Time : %f\n" % elapsed)
1338 217 equemene
1339 217 equemene
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1340 217 equemene
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1341 219 equemene
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1342 219 equemene
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1343 217 equemene
1344 218 equemene
    Context.pop()
1345 218 equemene
1346 217 equemene
    Context.detach()
1347 217 equemene
1348 217 equemene
    return(elapsed)
1349 217 equemene
1350 199 equemene
if __name__=='__main__':
1351 199 equemene
1352 199 equemene
    GpuStyle = 'OpenCL'
1353 199 equemene
    Mass = 1.
1354 199 equemene
    # Internal Radius 3 times de Schwarzschild Radius
1355 199 equemene
    InternalRadius=6.*Mass
1356 199 equemene
    #
1357 199 equemene
    ExternalRadius=12.
1358 199 equemene
    #
1359 199 equemene
    # Angle with normal to disc 10 degrees
1360 199 equemene
    Angle = numpy.pi/180.*(90.-10.)
1361 199 equemene
    # Radiation of disc : BlackBody or Monochromatic
1362 209 equemene
    BlackBody = False
1363 199 equemene
    # Size of image
1364 199 equemene
    Size=256
1365 199 equemene
    # Variable Type
1366 199 equemene
    VariableType='FP32'
1367 199 equemene
    # ?
1368 199 equemene
    q=-2
1369 204 equemene
    # Method of resolution
1370 209 equemene
    Method='TrajectoPixel'
1371 211 equemene
    # Colors for output image
1372 211 equemene
    Colors='Greyscale'
1373 211 equemene
    # Physics
1374 211 equemene
    Physics='Einstein'
1375 211 equemene
    # No output as image
1376 211 equemene
    NoImage = False
1377 211 equemene
1378 211 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> -v <FP32/FP64>'
1379 199 equemene
1380 199 equemene
    try:
1381 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="])
1382 199 equemene
    except getopt.GetoptError:
1383 199 equemene
        print(HowToUse % sys.argv[0])
1384 199 equemene
        sys.exit(2)
1385 199 equemene
1386 199 equemene
    # List of Devices
1387 199 equemene
    Devices=[]
1388 199 equemene
    Alu={}
1389 199 equemene
1390 199 equemene
    for opt, arg in opts:
1391 199 equemene
        if opt == '-h':
1392 199 equemene
            print(HowToUse % sys.argv[0])
1393 199 equemene
1394 199 equemene
            print("\nInformations about devices detected under OpenCL API:")
1395 199 equemene
            # For PyOpenCL import
1396 199 equemene
            try:
1397 199 equemene
                import pyopencl as cl
1398 199 equemene
                Id=0
1399 199 equemene
                for platform in cl.get_platforms():
1400 199 equemene
                    for device in platform.get_devices():
1401 199 equemene
                        #deviceType=cl.device_type.to_string(device.type)
1402 199 equemene
                        deviceType="xPU"
1403 199 equemene
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
1404 199 equemene
                        Id=Id+1
1405 199 equemene
1406 199 equemene
            except:
1407 199 equemene
                print("Your platform does not seem to support OpenCL")
1408 199 equemene
1409 199 equemene
            print("\nInformations about devices detected under CUDA API:")
1410 199 equemene
                # For PyCUDA import
1411 199 equemene
            try:
1412 199 equemene
                import pycuda.driver as cuda
1413 199 equemene
                cuda.init()
1414 199 equemene
                for Id in range(cuda.Device.count()):
1415 199 equemene
                    device=cuda.Device(Id)
1416 199 equemene
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
1417 199 equemene
                print
1418 199 equemene
            except:
1419 199 equemene
                print("Your platform does not seem to support CUDA")
1420 199 equemene
1421 199 equemene
            sys.exit()
1422 199 equemene
1423 199 equemene
        elif opt in ("-d", "--device"):
1424 199 equemene
#            Devices.append(int(arg))
1425 199 equemene
            Device=int(arg)
1426 199 equemene
        elif opt in ("-g", "--gpustyle"):
1427 199 equemene
            GpuStyle = arg
1428 204 equemene
        elif opt in ("-v", "--variabletype"):
1429 199 equemene
            VariableType = arg
1430 199 equemene
        elif opt in ("-s", "--size"):
1431 199 equemene
            Size = int(arg)
1432 199 equemene
        elif opt in ("-m", "--mass"):
1433 199 equemene
            Mass = float(arg)
1434 199 equemene
        elif opt in ("-i", "--internal"):
1435 199 equemene
            InternalRadius = float(arg)
1436 199 equemene
        elif opt in ("-e", "--external"):
1437 199 equemene
            ExternalRadius = float(arg)
1438 199 equemene
        elif opt in ("-a", "--angle"):
1439 199 equemene
            Angle = numpy.pi/180.*(90.-float(arg))
1440 199 equemene
        elif opt in ("-b", "--blackbody"):
1441 199 equemene
            BlackBody = True
1442 211 equemene
        elif opt in ("-n", "--noimage"):
1443 211 equemene
            NoImage = True
1444 204 equemene
        elif opt in ("-t", "--method"):
1445 204 equemene
            Method = arg
1446 211 equemene
        elif opt in ("-c", "--colors"):
1447 211 equemene
            Colors = arg
1448 211 equemene
        elif opt in ("-p", "--physics"):
1449 211 equemene
            Physics = arg
1450 199 equemene
1451 199 equemene
    print("Device Identification selected : %s" % Device)
1452 199 equemene
    print("GpuStyle used : %s" % GpuStyle)
1453 199 equemene
    print("VariableType : %s" % VariableType)
1454 199 equemene
    print("Size : %i" % Size)
1455 199 equemene
    print("Mass : %f" % Mass)
1456 199 equemene
    print("Internal Radius : %f" % InternalRadius)
1457 199 equemene
    print("External Radius : %f" % ExternalRadius)
1458 199 equemene
    print("Angle with normal of (in radians) : %f" % Angle)
1459 199 equemene
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
1460 204 equemene
    print("Method of resolution : %s" % Method)
1461 211 equemene
    print("Colors for output images : %s" % Colors)
1462 211 equemene
    print("Physics used for Trajectories : %s" % Physics)
1463 199 equemene
1464 199 equemene
    if GpuStyle=='CUDA':
1465 199 equemene
        print("\nSelection of CUDA device")
1466 199 equemene
        try:
1467 199 equemene
            # For PyCUDA import
1468 199 equemene
            import pycuda.driver as cuda
1469 199 equemene
1470 199 equemene
            cuda.init()
1471 199 equemene
            for Id in range(cuda.Device.count()):
1472 199 equemene
                device=cuda.Device(Id)
1473 199 equemene
                print("Device #%i of type GPU : %s" % (Id,device.name()))
1474 199 equemene
                if Id in Devices:
1475 199 equemene
                    Alu[Id]='GPU'
1476 199 equemene
1477 199 equemene
        except ImportError:
1478 199 equemene
            print("Platform does not seem to support CUDA")
1479 199 equemene
1480 199 equemene
    if GpuStyle=='OpenCL':
1481 199 equemene
        print("\nSelection of OpenCL device")
1482 199 equemene
        try:
1483 199 equemene
            # For PyOpenCL import
1484 199 equemene
            import pyopencl as cl
1485 199 equemene
            Id=0
1486 199 equemene
            for platform in cl.get_platforms():
1487 199 equemene
                for device in platform.get_devices():
1488 199 equemene
                    #deviceType=cl.device_type.to_string(device.type)
1489 199 equemene
                    deviceType="xPU"
1490 199 equemene
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
1491 199 equemene
1492 199 equemene
                    if Id in Devices:
1493 199 equemene
                    # Set the Alu as detected Device Type
1494 199 equemene
                        Alu[Id]=deviceType
1495 199 equemene
                    Id=Id+1
1496 199 equemene
        except ImportError:
1497 199 equemene
            print("Platform does not seem to support OpenCL")
1498 199 equemene
1499 199 equemene
#    print(Devices,Alu)
1500 199 equemene
1501 199 equemene
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
1502 204 equemene
    TrackPoints=2048
1503 201 equemene
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1504 201 equemene
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1505 199 equemene
1506 204 equemene
    InputCL={}
1507 204 equemene
    InputCL['Device']=Device
1508 204 equemene
    InputCL['GpuStyle']=GpuStyle
1509 204 equemene
    InputCL['VariableType']=VariableType
1510 204 equemene
    InputCL['Size']=Size
1511 204 equemene
    InputCL['Mass']=Mass
1512 204 equemene
    InputCL['InternalRadius']=InternalRadius
1513 204 equemene
    InputCL['ExternalRadius']=ExternalRadius
1514 204 equemene
    InputCL['Angle']=Angle
1515 204 equemene
    InputCL['BlackBody']=BlackBody
1516 204 equemene
    InputCL['Method']=Method
1517 204 equemene
    InputCL['TrackPoints']=TrackPoints
1518 211 equemene
    InputCL['Physics']=Physics
1519 199 equemene
1520 217 equemene
    if GpuStyle=='OpenCL':
1521 217 equemene
        duration=BlackHoleCL(zImage,fImage,InputCL)
1522 217 equemene
    else:
1523 217 equemene
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
1524 217 equemene
1525 211 equemene
    Hostname=gethostname()
1526 211 equemene
    Date=time.strftime("%Y%m%d_%H%M%S")
1527 211 equemene
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1528 211 equemene
1529 211 equemene
1530 211 equemene
    if not NoImage:
1531 211 equemene
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
1532 211 equemene
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)