Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 233

Historique | Voir | Annoter | Télécharger (50 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 229 equemene
# If crash on OpenCL Intel implementation, add following options and force
29 229 equemene
#export PYOPENCL_COMPILER_OUTPUT=1
30 229 equemene
#export CL_CONFIG_USE_VECTORIZER=True
31 229 equemene
#export CL_CONFIG_CPU_VECTORIZER_MODE=16
32 226 equemene
33 199 equemene
import pyopencl as cl
34 199 equemene
import numpy
35 199 equemene
import time,string
36 199 equemene
from numpy.random import randint as nprnd
37 199 equemene
import sys
38 199 equemene
import getopt
39 199 equemene
import matplotlib.pyplot as plt
40 211 equemene
from socket import gethostname
41 199 equemene
42 211 equemene
def DictionariesAPI():
43 211 equemene
    PhysicsList={'Einstein':0,'Newton':1}
44 211 equemene
    return(PhysicsList)
45 204 equemene
46 226 equemene
#
47 226 equemene
# Blank space below to simplify debugging on OpenCL code
48 226 equemene
#
49 226 equemene
50 226 equemene
51 226 equemene
52 226 equemene
53 226 equemene
54 226 equemene
55 226 equemene
56 226 equemene
57 226 equemene
58 226 equemene
59 226 equemene
60 226 equemene
61 226 equemene
62 226 equemene
63 226 equemene
64 226 equemene
65 226 equemene
66 226 equemene
67 226 equemene
68 226 equemene
69 226 equemene
70 226 equemene
71 226 equemene
72 226 equemene
73 226 equemene
74 226 equemene
75 226 equemene
76 226 equemene
77 226 equemene
78 226 equemene
79 226 equemene
80 226 equemene
81 226 equemene
82 226 equemene
83 226 equemene
84 226 equemene
85 226 equemene
86 226 equemene
87 226 equemene
88 226 equemene
89 226 equemene
90 226 equemene
91 226 equemene
92 226 equemene
93 226 equemene
94 226 equemene
95 226 equemene
96 226 equemene
97 226 equemene
98 226 equemene
99 226 equemene
100 226 equemene
101 226 equemene
102 226 equemene
103 211 equemene
BlobOpenCL= """
104 204 equemene

105 199 equemene
#define PI (float)3.14159265359
106 209 equemene
#define nbr 256
107 199 equemene

108 211 equemene
#define EINSTEIN 0
109 211 equemene
#define NEWTON 1
110 211 equemene

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

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

121 199 equemene
  angle=atan2(y,x);
122 199 equemene

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

128 199 equemene
  return angle;
129 199 equemene
}
130 199 equemene

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

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

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

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

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

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

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

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

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

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

194 213 equemene
#define MYFLOAT float
195 209 equemene

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

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

208 209 equemene
#define planck 6.62e-34
209 209 equemene
#define k 1.38e-23
210 209 equemene
#define c2 9.e16
211 209 equemene
#define temp 3.e7
212 209 equemene
#define m_point 1.
213 199 equemene

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

218 209 equemene
  MYFLOAT v=1.-3./r;
219 199 equemene

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

222 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)));
223 209 equemene

224 209 equemene
  flux_int=0.;
225 209 equemene
  flx=0.;
226 209 equemene

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

242 199 equemene
          flx+=flux_int;
243 199 equemene
        }
244 199 equemene
    }
245 209 equemene

246 209 equemene
  return((float)(flx));
247 199 equemene
}
248 199 equemene

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

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

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

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

272 199 equemene
}
273 199 equemene

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

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

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

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

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

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

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

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

315 228 equemene

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

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

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

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

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

332 204 equemene
  int ExitOnImpact=0;
333 199 equemene

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

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

346 228 equemene

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

356 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
357 204 equemene

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

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

374 204 equemene
   // Perform trajectory for each pixel
375 204 equemene

376 226 equemene
  float m,ri,re,tho;
377 209 equemene
  int q,raie;
378 204 equemene

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

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

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

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

396 204 equemene
  // Step of Impact Parameter
397 209 equemene
  db=bmx/(2.e0*(float)ImpactParameter);
398 204 equemene

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

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

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

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

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

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

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

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

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

441 204 equemene
  }
442 204 equemene

443 201 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
444 204 equemene

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

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

464 204 equemene
   // Perform trajectory for each pixel
465 204 equemene

466 226 equemene
  float m,ri,re,tho;
467 209 equemene
  int q,raie;
468 204 equemene

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

475 204 equemene
  float bmx,db,b,h;
476 226 equemene
  float phi,phd;
477 204 equemene
  float zp=0,fp=0;
478 204 equemene

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

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

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

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

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

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

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

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

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

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

524 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
525 204 equemene

526 204 equemene
}
527 204 equemene

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

538 204 equemene
  // Perform trajectory for each pixel
539 204 equemene

540 224 equemene
  float m,rs,re;
541 224 equemene

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

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

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

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

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

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

560 224 equemene
  up=0.;
561 224 equemene
  vp=1.;
562 224 equemene

563 224 equemene
  pp=0.;
564 224 equemene
  nh=0;
565 224 equemene

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

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

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

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

585 224 equemene
  IdLast[bi]=nh;
586 224 equemene

587 224 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
588 224 equemene

589 224 equemene
}
590 224 equemene

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

601 225 equemene
  float Trajectory[TRACKPOINTS];
602 224 equemene

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

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

614 224 equemene
  float bmx,db,b,h;
615 204 equemene
  int nh;
616 204 equemene

617 204 equemene
  // Autosize for image
618 224 equemene
  bmx=1.25e0f*re;
619 204 equemene

620 204 equemene
  // Angular step of integration
621 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
622 204 equemene

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

627 204 equemene
  float up,vp,pp,us,vs,ps;
628 204 equemene

629 209 equemene
  up=0.;
630 209 equemene
  vp=1.;
631 204 equemene

632 204 equemene
  pp=0.;
633 204 equemene
  nh=0;
634 204 equemene

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

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

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

652 204 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
653 204 equemene

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

658 224 equemene
  int imx=(int)(16*bi);
659 204 equemene

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

668 224 equemene
     int HalfLap=0,ExitOnImpact=0,ni;
669 224 equemene
     float php,nr,r;
670 224 equemene

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

677 224 equemene
        if (ni<nh)
678 224 equemene
        {
679 225 equemene
           r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
680 224 equemene
        }
681 224 equemene
        else
682 224 equemene
        {
683 225 equemene
           r=Trajectory[ni];
684 224 equemene
        }
685 224 equemene

686 224 equemene
        if ((r<=re)&&(r>=ri))
687 224 equemene
        {
688 224 equemene
           ExitOnImpact=1;
689 224 equemene
           impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
690 224 equemene
        }
691 224 equemene

692 224 equemene
        HalfLap++;
693 224 equemene

694 224 equemene
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
695 224 equemene

696 224 equemene
     zImage[yi+2*bmaxi*xi]=zp;
697 224 equemene
     fImage[yi+2*bmaxi*xi]=fp;
698 224 equemene

699 224 equemene
  }
700 224 equemene

701 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
702 199 equemene

703 199 equemene
}
704 225 equemene

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

713 225 equemene
   float Trajectory[TRACKPOINTS];
714 225 equemene

715 225 equemene
   // Perform trajectory for each pixel
716 225 equemene

717 225 equemene
   float m,rs,ri,re,tho;
718 225 equemene
   int raie,q;
719 225 equemene

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

728 225 equemene
   float bmx,db,b,h;
729 225 equemene
   int nh;
730 225 equemene

731 225 equemene
   // Autosize for image
732 225 equemene
   bmx=1.25e0f*re;
733 225 equemene

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

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

744 225 equemene
      float up,vp,pp,us,vs,ps;
745 225 equemene

746 225 equemene
      up=0.;
747 225 equemene
      vp=1.;
748 225 equemene

749 225 equemene
      pp=0.;
750 225 equemene
      nh=0;
751 225 equemene

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

754 225 equemene
      // b versus us
755 225 equemene
      float bvus=fabs(b/us);
756 225 equemene
      float bvus0=bvus;
757 225 equemene
      Trajectory[nh]=bvus;
758 225 equemene

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

769 225 equemene
      } while ((bvus>=rs)&&(bvus<=bvus0));
770 225 equemene

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

775 225 equemene
      int imx=(int)(16*bi);
776 225 equemene

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

785 225 equemene
         int HalfLap=0,ExitOnImpact=0,ni;
786 225 equemene
         float php,nr,r;
787 225 equemene

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

794 225 equemene
            if (ni<nh)
795 225 equemene
            {
796 225 equemene
               r=(Trajectory[ni+1]-Trajectory[ni])*(nr-ni*1.)+Trajectory[ni];
797 225 equemene
            }
798 225 equemene
            else
799 225 equemene
            {
800 225 equemene
               r=Trajectory[ni];
801 225 equemene
            }
802 225 equemene

803 225 equemene
            if ((r<=re)&&(r>=ri))
804 225 equemene
            {
805 225 equemene
               ExitOnImpact=1;
806 225 equemene
               impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
807 225 equemene
            }
808 225 equemene

809 225 equemene
            HalfLap++;
810 225 equemene

811 225 equemene
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
812 225 equemene

813 225 equemene
         zImage[yi+2*bmaxi*xi]=zp;
814 225 equemene
         fImage[yi+2*bmaxi*xi]=fp;
815 225 equemene

816 225 equemene
      }
817 225 equemene

818 225 equemene
   }
819 225 equemene

820 225 equemene
   barrier(CLK_GLOBAL_MEM_FENCE);
821 225 equemene

822 225 equemene
}
823 211 equemene
"""
824 199 equemene
825 217 equemene
def KernelCodeCuda():
826 217 equemene
    BlobCUDA= """
827 217 equemene

828 217 equemene
#define PI (float)3.14159265359
829 217 equemene
#define nbr 256
830 217 equemene

831 217 equemene
#define EINSTEIN 0
832 217 equemene
#define NEWTON 1
833 217 equemene

834 224 equemene
#ifdef SETTRACKPOINTS
835 224 equemene
#define TRACKPOINTS SETTRACKPOINTS
836 224 equemene
#else
837 224 equemene
#define TRACKPOINTS
838 224 equemene
#endif
839 217 equemene

840 217 equemene
__device__ float nothing(float x)
841 217 equemene
{
842 217 equemene
  return(x);
843 217 equemene
}
844 217 equemene

845 217 equemene
__device__ float atanp(float x,float y)
846 217 equemene
{
847 217 equemene
  float angle;
848 217 equemene

849 217 equemene
  angle=atan2(y,x);
850 217 equemene

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

856 217 equemene
  return(angle);
857 217 equemene
}
858 217 equemene

859 217 equemene
__device__ float f(float v)
860 217 equemene
{
861 217 equemene
  return(v);
862 217 equemene
}
863 217 equemene

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

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

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

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

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

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

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

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

918 217 equemene
__device__ float spectre_cn(float rf32,float b32,float db32,
919 217 equemene
                            float h32,float r32,float m32,float bss32)
920 217 equemene
{
921 217 equemene

922 217 equemene
#define MYFLOAT float
923 217 equemene

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

932 217 equemene
  MYFLOAT flx;
933 217 equemene
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
934 217 equemene
  int fi,posfreq;
935 217 equemene

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

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

946 217 equemene
  MYFLOAT v=1.-3./r;
947 217 equemene

948 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 ));
949 217 equemene

950 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)));
951 217 equemene

952 217 equemene
  flux_int=0.;
953 217 equemene
  flx=0.;
954 217 equemene

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

970 217 equemene
          flx+=flux_int;
971 217 equemene
        }
972 217 equemene
    }
973 217 equemene

974 217 equemene
  return((float)(flx));
975 217 equemene
}
976 217 equemene

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

984 217 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
985 217 equemene

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

997 217 equemene
  *zp=1./rf;
998 217 equemene
  *fp=flx;
999 217 equemene

1000 217 equemene
}
1001 217 equemene

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

1012 217 equemene
   // Perform trajectory for each pixel, exit on hit
1013 217 equemene

1014 217 equemene
  float m,rs,ri,re,tho;
1015 217 equemene
  int q,raie;
1016 217 equemene

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

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

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

1035 217 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1036 217 equemene

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

1044 217 equemene
  float up,vp,pp,us,vs,ps;
1045 217 equemene

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

1052 217 equemene
  up=0.;
1053 217 equemene
  vp=1.;
1054 217 equemene

1055 217 equemene
  pp=0.;
1056 217 equemene
  nh=0;
1057 217 equemene

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

1060 218 equemene
  rps=fabs(b/us);
1061 218 equemene
  rp0=rps;
1062 217 equemene

1063 217 equemene
  int ExitOnImpact=0;
1064 217 equemene

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

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

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

1087 218 equemene
  __syncthreads();
1088 218 equemene

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

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

1105 217 equemene
  // Perform trajectory for each pixel
1106 217 equemene

1107 217 equemene
  float m,rs,ri,re,tho;
1108 217 equemene
  int q,raie;
1109 217 equemene

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

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

1125 217 equemene
  // Angular step of integration
1126 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1127 217 equemene

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

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

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

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

1144 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
1145 217 equemene

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

1154 217 equemene
      if (ni<IdLast[bi])
1155 217 equemene
      {
1156 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1157 217 equemene
      }
1158 217 equemene
      else
1159 217 equemene
      {
1160 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
1161 217 equemene
      }
1162 217 equemene

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

1169 217 equemene
      HalfLap++;
1170 217 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
1171 217 equemene

1172 217 equemene
  }
1173 217 equemene

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

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

1193 217 equemene
   // Perform trajectory for each pixel
1194 217 equemene

1195 217 equemene
  float m,rs,ri,re,tho;
1196 217 equemene
  int q,raie;
1197 217 equemene

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

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

1210 217 equemene
  // Autosize for image
1211 217 equemene
  bmx=1.25*re;
1212 217 equemene

1213 217 equemene
  // Angular step of integration
1214 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1215 217 equemene

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

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

1225 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
1226 217 equemene
  float php,nr,r;
1227 217 equemene

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

1234 217 equemene
     if (ni<IdLast[bi])
1235 217 equemene
     {
1236 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1237 217 equemene
     }
1238 217 equemene
     else
1239 217 equemene
     {
1240 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
1241 217 equemene
     }
1242 217 equemene

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

1249 217 equemene
     HalfLap++;
1250 217 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1251 217 equemene

1252 217 equemene
  zImage[yi+2*bmaxi*xi]=zp;
1253 217 equemene
  fImage[yi+2*bmaxi*xi]=fp;
1254 217 equemene

1255 217 equemene
}
1256 217 equemene

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

1267 217 equemene
  // Perform trajectory for each pixel
1268 217 equemene

1269 217 equemene
  float m,rs,ri,re,tho;
1270 217 equemene
  int raie,q;
1271 217 equemene

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

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

1285 217 equemene
  // Autosize for image
1286 217 equemene
  bmx=1.25*re;
1287 217 equemene

1288 217 equemene
  // Angular step of integration
1289 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1290 217 equemene

1291 217 equemene
  // impact parameter
1292 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
1293 217 equemene

1294 217 equemene
  float up,vp,pp,us,vs,ps;
1295 217 equemene

1296 217 equemene
  up=0.;
1297 217 equemene
  vp=1.;
1298 217 equemene

1299 217 equemene
  pp=0.;
1300 217 equemene
  nh=0;
1301 217 equemene

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

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

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

1319 217 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
1320 217 equemene

1321 217 equemene
  IdLast[bi]=nh;
1322 217 equemene

1323 217 equemene
}
1324 224 equemene

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

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

1336 225 equemene
  float Trajectory[2048];
1337 224 equemene

1338 224 equemene
  // Perform trajectory for each pixel
1339 224 equemene

1340 224 equemene
  float m,rs,ri,re,tho;
1341 224 equemene
  int raie,q;
1342 224 equemene

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

1351 224 equemene
  float bmx,db,b,h;
1352 224 equemene
  int nh;
1353 224 equemene

1354 224 equemene
  // Autosize for image
1355 224 equemene
  bmx=1.25e0f*re;
1356 224 equemene

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

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

1364 224 equemene
  float up,vp,pp,us,vs,ps;
1365 224 equemene

1366 224 equemene
  up=0.;
1367 224 equemene
  vp=1.;
1368 224 equemene

1369 224 equemene
  pp=0.;
1370 224 equemene
  nh=0;
1371 224 equemene

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

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

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

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

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

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

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

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

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

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

1425 224 equemene
        HalfLap++;
1426 224 equemene

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

1429 224 equemene
   __syncthreads();
1430 224 equemene

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

1434 224 equemene
  }
1435 224 equemene

1436 224 equemene
}
1437 225 equemene

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

1446 225 equemene
   float Trajectory[TRACKPOINTS];
1447 225 equemene

1448 225 equemene
   // Perform trajectory for each pixel
1449 225 equemene

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

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

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

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

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

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

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

1479 225 equemene
      up=0.;
1480 225 equemene
      vp=1.;
1481 225 equemene

1482 225 equemene
      pp=0.;
1483 225 equemene
      nh=0;
1484 225 equemene

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

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

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

1502 225 equemene
      } while ((bvus>=rs)&&(bvus<=bvus0));
1503 225 equemene

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

1508 225 equemene
      int imx=(int)(16*bi);
1509 225 equemene

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

1518 225 equemene
         int HalfLap=0,ExitOnImpact=0,ni;
1519 225 equemene
         float php,nr,r;
1520 225 equemene

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

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

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

1542 225 equemene
            HalfLap++;
1543 225 equemene

1544 225 equemene
         } while ((HalfLap<=2)&&(ExitOnImpact==0));
1545 225 equemene

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

1549 225 equemene
      }
1550 225 equemene

1551 225 equemene
   }
1552 225 equemene

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