Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 224

Historique | Voir | Annoter | Télécharger (43,55 ko)

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

132 213 equemene
#define MYFLOAT float
133 209 equemene

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

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

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

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

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

158 209 equemene
  qu=1./sqrt((1.-3./r)*r)*(sqrt(r)-sqrt(6.)+sqrt(3.)/2.*log((sqrt(r)+sqrt(3.))/(sqrt(r)-sqrt(3.))* 0.17157287525380988 ));
159 199 equemene

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

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

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

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

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

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

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

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

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

210 199 equemene
}
211 199 equemene

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

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

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

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

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

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

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

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

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

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

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

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

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

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

272 204 equemene
  int ExitOnImpact=0;
273 199 equemene

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

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

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

296 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
297 204 equemene

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

383 204 equemene
  }
384 204 equemene

385 201 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
386 204 equemene

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

468 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
469 204 equemene

470 204 equemene
}
471 204 equemene

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

531 224 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
532 224 equemene

533 224 equemene
}
534 224 equemene

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

545 224 equemene
  float Trajectories[2048];
546 224 equemene

547 224 equemene
  // Perform trajectory for each pixel
548 224 equemene

549 209 equemene
  float m,rs,ri,re,tho;
550 209 equemene
  int raie,q;
551 204 equemene

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

560 224 equemene
  float bmx,db,b,h;
561 204 equemene
  int nh;
562 204 equemene

563 204 equemene
  // Autosize for image
564 224 equemene
  bmx=1.25e0f*re;
565 204 equemene

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

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

573 204 equemene
  float up,vp,pp,us,vs,ps;
574 204 equemene

575 209 equemene
  up=0.;
576 209 equemene
  vp=1.;
577 204 equemene

578 204 equemene
  pp=0.;
579 204 equemene
  nh=0;
580 204 equemene

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

583 204 equemene
  // b versus us
584 204 equemene
  float bvus=fabs(b/us);
585 204 equemene
  float bvus0=bvus;
586 224 equemene
  Trajectories[nh]=bvus;
587 204 equemene

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

598 204 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
599 204 equemene

600 224 equemene
  int imx=(int)(16*bi);
601 204 equemene

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

610 224 equemene
     int HalfLap=0,ExitOnImpact=0,ni;
611 224 equemene
     float php,nr,r;
612 224 equemene

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

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

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

634 224 equemene
        HalfLap++;
635 224 equemene

636 224 equemene
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
637 224 equemene

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

641 224 equemene
  }
642 224 equemene

643 204 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
644 199 equemene

645 199 equemene
}
646 211 equemene
"""
647 199 equemene
648 217 equemene
def KernelCodeCuda():
649 217 equemene
    BlobCUDA= """
650 217 equemene

651 217 equemene
#define PI (float)3.14159265359
652 217 equemene
#define nbr 256
653 217 equemene

654 217 equemene
#define EINSTEIN 0
655 217 equemene
#define NEWTON 1
656 217 equemene

657 224 equemene
#ifdef SETTRACKPOINTS
658 224 equemene
#define TRACKPOINTS SETTRACKPOINTS
659 224 equemene
#else
660 224 equemene
#define TRACKPOINTS
661 224 equemene
#endif
662 217 equemene

663 217 equemene
__device__ float nothing(float x)
664 217 equemene
{
665 217 equemene
  return(x);
666 217 equemene
}
667 217 equemene

668 217 equemene
__device__ float atanp(float x,float y)
669 217 equemene
{
670 217 equemene
  float angle;
671 217 equemene

672 217 equemene
  angle=atan2(y,x);
673 217 equemene

674 217 equemene
  if (angle<0.e0f)
675 217 equemene
    {
676 217 equemene
      angle+=(float)2.e0f*PI;
677 217 equemene
    }
678 217 equemene

679 217 equemene
  return(angle);
680 217 equemene
}
681 217 equemene

682 217 equemene
__device__ float f(float v)
683 217 equemene
{
684 217 equemene
  return(v);
685 217 equemene
}
686 217 equemene

687 217 equemene
#if PHYSICS == NEWTON
688 217 equemene
__device__ float g(float u,float m,float b)
689 217 equemene
{
690 217 equemene
  return (-u);
691 217 equemene
}
692 217 equemene
#else
693 217 equemene
__device__ float g(float u,float m,float b)
694 217 equemene
{
695 217 equemene
  return (3.e0f*m/b*pow(u,2)-u);
696 217 equemene
}
697 217 equemene
#endif
698 217 equemene

699 217 equemene
__device__ void calcul(float *us,float *vs,float up,float vp,
700 217 equemene
                       float h,float m,float b)
701 217 equemene
{
702 217 equemene
  float c0,c1,c2,c3,d0,d1,d2,d3;
703 217 equemene

704 217 equemene
  c0=h*f(vp);
705 217 equemene
  c1=h*f(vp+c0/2.);
706 217 equemene
  c2=h*f(vp+c1/2.);
707 217 equemene
  c3=h*f(vp+c2);
708 217 equemene
  d0=h*g(up,m,b);
709 217 equemene
  d1=h*g(up+d0/2.,m,b);
710 217 equemene
  d2=h*g(up+d1/2.,m,b);
711 217 equemene
  d3=h*g(up+d2,m,b);
712 217 equemene

713 217 equemene
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
714 217 equemene
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
715 217 equemene
}
716 217 equemene

717 217 equemene
__device__ void rungekutta(float *ps,float *us,float *vs,
718 217 equemene
                           float pp,float up,float vp,
719 217 equemene
                           float h,float m,float b)
720 217 equemene
{
721 217 equemene
  calcul(us,vs,up,vp,h,m,b);
722 217 equemene
  *ps=pp+h;
723 217 equemene
}
724 217 equemene

725 217 equemene
__device__ float decalage_spectral(float r,float b,float phi,
726 217 equemene
                                   float tho,float m)
727 217 equemene
{
728 217 equemene
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
729 217 equemene
}
730 217 equemene

731 217 equemene
__device__ float spectre(float rf,int q,float b,float db,
732 217 equemene
                         float h,float r,float m,float bss)
733 217 equemene
{
734 217 equemene
  float flx;
735 217 equemene

736 221 equemene
//  flx=exp(q*log(r/m))*pow(rf,4)*b*db*h;
737 221 equemene
  flx=exp(q*log(r/m)+4.*log(rf))*b*db*h;
738 217 equemene
  return(flx);
739 217 equemene
}
740 217 equemene

741 217 equemene
__device__ float spectre_cn(float rf32,float b32,float db32,
742 217 equemene
                            float h32,float r32,float m32,float bss32)
743 217 equemene
{
744 217 equemene

745 217 equemene
#define MYFLOAT float
746 217 equemene

747 217 equemene
  MYFLOAT rf=(MYFLOAT)(rf32);
748 217 equemene
  MYFLOAT b=(MYFLOAT)(b32);
749 217 equemene
  MYFLOAT db=(MYFLOAT)(db32);
750 217 equemene
  MYFLOAT h=(MYFLOAT)(h32);
751 217 equemene
  MYFLOAT r=(MYFLOAT)(r32);
752 217 equemene
  MYFLOAT m=(MYFLOAT)(m32);
753 217 equemene
  MYFLOAT bss=(MYFLOAT)(bss32);
754 217 equemene

755 217 equemene
  MYFLOAT flx;
756 217 equemene
  MYFLOAT nu_rec,nu_em,qu,temp_em,flux_int;
757 217 equemene
  int fi,posfreq;
758 217 equemene

759 217 equemene
#define planck 6.62e-34
760 217 equemene
#define k 1.38e-23
761 217 equemene
#define c2 9.e16
762 217 equemene
#define temp 3.e7
763 217 equemene
#define m_point 1.
764 217 equemene

765 217 equemene
#define lplanck (log(6.62)-34.*log(10.))
766 217 equemene
#define lk (log(1.38)-23.*log(10.))
767 217 equemene
#define lc2 (log(9.)+16.*log(10.))
768 217 equemene

769 217 equemene
  MYFLOAT v=1.-3./r;
770 217 equemene

771 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 ));
772 217 equemene

773 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)));
774 217 equemene

775 217 equemene
  flux_int=0.;
776 217 equemene
  flx=0.;
777 217 equemene

778 217 equemene
  for (fi=0;fi<nbr;fi++)
779 217 equemene
    {
780 217 equemene
      nu_em=bss*(MYFLOAT)fi/(MYFLOAT)nbr;
781 217 equemene
      nu_rec=nu_em*rf;
782 217 equemene
      posfreq=(int)(nu_rec*(MYFLOAT)nbr/bss);
783 217 equemene
      if ((posfreq>0)&&(posfreq<nbr))
784 217 equemene
        {
785 217 equemene
          // Initial version
786 217 equemene
          // flux_int=2.*planck/c2*pow(nu_em,3)/(exp(planck*nu_em/(k*temp_em))-1.);
787 217 equemene
          // Version with log used
788 217 equemene
          //flux_int=2.*exp(lplanck-lc2+3.*log(nu_em))/(exp(exp(lplanck-lk+log(nu_em/temp_em)))-1.);
789 217 equemene
          // flux_int*=pow(rf,3)*b*db*h;
790 217 equemene
          //flux_int*=exp(3.*log(rf))*b*db*h;
791 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;
792 217 equemene

793 217 equemene
          flx+=flux_int;
794 217 equemene
        }
795 217 equemene
    }
796 217 equemene

797 217 equemene
  return((float)(flx));
798 217 equemene
}
799 217 equemene

800 217 equemene
__device__ void impact(float phi,float r,float b,float tho,float m,
801 217 equemene
                       float *zp,float *fp,
802 217 equemene
                       int q,float db,
803 217 equemene
                       float h,int raie)
804 217 equemene
{
805 217 equemene
  float flx,rf,bss;
806 217 equemene

807 217 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
808 217 equemene

809 217 equemene
  if (raie==0)
810 217 equemene
    {
811 217 equemene
      bss=1.e19;
812 217 equemene
      flx=spectre_cn(rf,b,db,h,r,m,bss);
813 217 equemene
    }
814 217 equemene
  else
815 217 equemene
    {
816 217 equemene
      bss=2.;
817 217 equemene
      flx=spectre(rf,q,b,db,h,r,m,bss);
818 217 equemene
    }
819 217 equemene

820 217 equemene
  *zp=1./rf;
821 217 equemene
  *fp=flx;
822 217 equemene

823 217 equemene
}
824 217 equemene

825 217 equemene
__global__ void EachPixel(float *zImage,float *fImage,
826 217 equemene
                          float Mass,float InternalRadius,
827 217 equemene
                          float ExternalRadius,float Angle,
828 217 equemene
                          int Line)
829 217 equemene
{
830 218 equemene
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
831 218 equemene
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
832 217 equemene
   uint sizex=(uint)gridDim.x*blockDim.x;
833 217 equemene
   uint sizey=(uint)gridDim.y*blockDim.y;
834 217 equemene

835 217 equemene
   // Perform trajectory for each pixel, exit on hit
836 217 equemene

837 217 equemene
  float m,rs,ri,re,tho;
838 217 equemene
  int q,raie;
839 217 equemene

840 217 equemene
  m=Mass;
841 217 equemene
  rs=2.*m;
842 217 equemene
  ri=InternalRadius;
843 217 equemene
  re=ExternalRadius;
844 217 equemene
  tho=Angle;
845 217 equemene
  q=-2;
846 217 equemene
  raie=Line;
847 217 equemene

848 217 equemene
  float d,bmx,db,b,h;
849 218 equemene
  float rp0,rpp,rps;
850 217 equemene
  float phi,thi,phd,php,nr,r;
851 217 equemene
  int nh;
852 217 equemene
  float zp,fp;
853 217 equemene

854 217 equemene
  // Autosize for image
855 217 equemene
  bmx=1.25*re;
856 217 equemene
  b=0.;
857 217 equemene

858 217 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
859 217 equemene

860 217 equemene
  // set origin as center of image
861 217 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
862 217 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
863 217 equemene
  // angle extracted from cylindric symmetry
864 217 equemene
  phi=atanp(x,y);
865 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
866 217 equemene

867 217 equemene
  float up,vp,pp,us,vs,ps;
868 217 equemene

869 217 equemene
  // impact parameter
870 217 equemene
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
871 217 equemene
  // step of impact parameter;
872 217 equemene
//  db=bmx/(float)(sizex/2);
873 217 equemene
  db=bmx/(float)(sizex);
874 217 equemene

875 217 equemene
  up=0.;
876 217 equemene
  vp=1.;
877 217 equemene

878 217 equemene
  pp=0.;
879 217 equemene
  nh=0;
880 217 equemene

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

883 218 equemene
  rps=fabs(b/us);
884 218 equemene
  rp0=rps;
885 217 equemene

886 217 equemene
  int ExitOnImpact=0;
887 217 equemene

888 217 equemene
  do
889 217 equemene
  {
890 217 equemene
     nh++;
891 217 equemene
     pp=ps;
892 217 equemene
     up=us;
893 217 equemene
     vp=vs;
894 218 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
895 218 equemene
     rpp=rps;
896 218 equemene
     rps=fabs(b/us);
897 218 equemene
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rps>ri)&&(rps<re)?1:0;
898 217 equemene

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

901 217 equemene
  if (ExitOnImpact==1) {
902 218 equemene
     impact(phi,rpp,b,tho,m,&zp,&fp,q,db,h,raie);
903 217 equemene
  }
904 217 equemene
  else
905 217 equemene
  {
906 217 equemene
     zp=0.;
907 217 equemene
     fp=0.;
908 217 equemene
  }
909 217 equemene

910 218 equemene
  __syncthreads();
911 218 equemene

912 217 equemene
  zImage[yi+sizex*xi]=(float)zp;
913 217 equemene
  fImage[yi+sizex*xi]=(float)fp;
914 217 equemene
}
915 217 equemene

916 217 equemene
__global__ void Pixel(float *zImage,float *fImage,
917 217 equemene
                      float *Trajectories,int *IdLast,
918 224 equemene
                      uint ImpactParameter,
919 217 equemene
                      float Mass,float InternalRadius,
920 217 equemene
                      float ExternalRadius,float Angle,
921 217 equemene
                      int Line)
922 217 equemene
{
923 219 equemene
   uint xi=(uint)(blockIdx.x*blockDim.x+threadIdx.x);
924 219 equemene
   uint yi=(uint)(blockIdx.y*blockDim.y+threadIdx.y);
925 219 equemene
   uint sizex=(uint)gridDim.x*blockDim.x;
926 219 equemene
   uint sizey=(uint)gridDim.y*blockDim.y;
927 217 equemene

928 217 equemene
  // Perform trajectory for each pixel
929 217 equemene

930 217 equemene
  float m,rs,ri,re,tho;
931 217 equemene
  int q,raie;
932 217 equemene

933 217 equemene
  m=Mass;
934 217 equemene
  rs=2.*m;
935 217 equemene
  ri=InternalRadius;
936 217 equemene
  re=ExternalRadius;
937 217 equemene
  tho=Angle;
938 217 equemene
  q=-2;
939 217 equemene
  raie=Line;
940 217 equemene

941 217 equemene
  float d,bmx,db,b,h;
942 217 equemene
  float phi,thi,phd,php,nr,r;
943 217 equemene
  int nh;
944 217 equemene
  float zp=0,fp=0;
945 217 equemene
  // Autosize for image, 25% greater than external radius
946 217 equemene
  bmx=1.25*re;
947 217 equemene

948 217 equemene
  // Angular step of integration
949 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
950 217 equemene

951 217 equemene
  // Step of Impact Parameter
952 217 equemene
  db=bmx/(2.e0*(float)ImpactParameter);
953 217 equemene

954 217 equemene
  // set origin as center of image
955 217 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
956 217 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
957 217 equemene
  // angle extracted from cylindric symmetry
958 217 equemene
  phi=atanp(x,y);
959 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
960 217 equemene

961 217 equemene
  // Real Impact Parameter
962 217 equemene
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
963 217 equemene

964 217 equemene
  // Integer Impact Parameter
965 217 equemene
  uint bi=(uint)sqrt(x*x+y*y);
966 217 equemene

967 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
968 217 equemene

969 217 equemene
  if (bi<ImpactParameter)
970 217 equemene
  {
971 217 equemene
    do
972 217 equemene
    {
973 217 equemene
      php=phd+(float)HalfLap*PI;
974 217 equemene
      nr=php/h;
975 217 equemene
      ni=(int)nr;
976 217 equemene

977 217 equemene
      if (ni<IdLast[bi])
978 217 equemene
      {
979 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
980 217 equemene
      }
981 217 equemene
      else
982 217 equemene
      {
983 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
984 217 equemene
      }
985 217 equemene

986 217 equemene
      if ((r<=re)&&(r>=ri))
987 217 equemene
      {
988 217 equemene
        ExitOnImpact=1;
989 217 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
990 217 equemene
      }
991 217 equemene

992 217 equemene
      HalfLap++;
993 217 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
994 217 equemene

995 217 equemene
  }
996 217 equemene

997 217 equemene
  zImage[yi+sizex*xi]=zp;
998 217 equemene
  fImage[yi+sizex*xi]=fp;
999 217 equemene
}
1000 217 equemene

1001 217 equemene
__global__ void Circle(float *Trajectories,int *IdLast,
1002 217 equemene
                       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 217 equemene
   // Integer Impact Parameter ID
1008 219 equemene
   int bi=blockIdx.x*blockDim.x+threadIdx.x;
1009 217 equemene
   // Integer points on circle
1010 219 equemene
   int i=blockIdx.y*blockDim.y+threadIdx.y;
1011 217 equemene
   // Integer Impact Parameter Size (half of image)
1012 217 equemene
   int bmaxi=gridDim.x*blockDim.x;
1013 217 equemene
   // Integer Points on circle
1014 217 equemene
   int imx=gridDim.y*blockDim.y;
1015 217 equemene

1016 217 equemene
   // Perform trajectory for each pixel
1017 217 equemene

1018 217 equemene
  float m,rs,ri,re,tho;
1019 217 equemene
  int q,raie;
1020 217 equemene

1021 217 equemene
  m=Mass;
1022 217 equemene
  rs=2.*m;
1023 217 equemene
  ri=InternalRadius;
1024 217 equemene
  re=ExternalRadius;
1025 217 equemene
  tho=Angle;
1026 217 equemene
  raie=Line;
1027 217 equemene

1028 217 equemene
  float bmx,db,b,h;
1029 217 equemene
  float phi,thi,phd;
1030 217 equemene
  int nh;
1031 217 equemene
  float zp=0,fp=0;
1032 217 equemene

1033 217 equemene
  // Autosize for image
1034 217 equemene
  bmx=1.25*re;
1035 217 equemene

1036 217 equemene
  // Angular step of integration
1037 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1038 217 equemene

1039 217 equemene
  // impact parameter
1040 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
1041 217 equemene
  db=bmx/(2.e0*(float)bmaxi);
1042 217 equemene

1043 217 equemene
  phi=2.*PI/(float)imx*(float)i;
1044 217 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
1045 217 equemene
  int yi=(int)((float)bi*sin(phi))+bmaxi;
1046 217 equemene
  int xi=(int)((float)bi*cos(phi))+bmaxi;
1047 217 equemene

1048 217 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
1049 217 equemene
  float php,nr,r;
1050 217 equemene

1051 217 equemene
  do
1052 217 equemene
  {
1053 217 equemene
     php=phd+(float)HalfLap*PI;
1054 217 equemene
     nr=php/h;
1055 217 equemene
     ni=(int)nr;
1056 217 equemene

1057 217 equemene
     if (ni<IdLast[bi])
1058 217 equemene
     {
1059 224 equemene
        r=(Trajectories[bi*TRACKPOINTS+ni+1]-Trajectories[bi*TRACKPOINTS+ni])*(nr-ni*1.)+Trajectories[bi*TRACKPOINTS+ni];
1060 217 equemene
     }
1061 217 equemene
     else
1062 217 equemene
     {
1063 224 equemene
        r=Trajectories[bi*TRACKPOINTS+ni];
1064 217 equemene
     }
1065 217 equemene

1066 217 equemene
     if ((r<=re)&&(r>=ri))
1067 217 equemene
     {
1068 217 equemene
        ExitOnImpact=1;
1069 217 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
1070 217 equemene
     }
1071 217 equemene

1072 217 equemene
     HalfLap++;
1073 217 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
1074 217 equemene

1075 217 equemene
  zImage[yi+2*bmaxi*xi]=zp;
1076 217 equemene
  fImage[yi+2*bmaxi*xi]=fp;
1077 217 equemene

1078 217 equemene
}
1079 217 equemene

1080 224 equemene
__global__ void Trajectory(float *Trajectories,int *IdLast,
1081 217 equemene
                           float Mass,float InternalRadius,
1082 217 equemene
                           float ExternalRadius,float Angle,
1083 217 equemene
                           int Line)
1084 217 equemene
{
1085 217 equemene
  // Integer Impact Parameter ID
1086 219 equemene
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1087 217 equemene
  // Integer Impact Parameter Size (half of image)
1088 217 equemene
  int bmaxi=gridDim.x*blockDim.x;
1089 217 equemene

1090 217 equemene
  // Perform trajectory for each pixel
1091 217 equemene

1092 217 equemene
  float m,rs,ri,re,tho;
1093 217 equemene
  int raie,q;
1094 217 equemene

1095 217 equemene
  m=Mass;
1096 217 equemene
  rs=2.*m;
1097 217 equemene
  ri=InternalRadius;
1098 217 equemene
  re=ExternalRadius;
1099 217 equemene
  tho=Angle;
1100 217 equemene
  q=-2;
1101 217 equemene
  raie=Line;
1102 217 equemene

1103 217 equemene
  float d,bmx,db,b,h;
1104 217 equemene
  float phi,thi,phd,php,nr,r;
1105 217 equemene
  int nh;
1106 217 equemene
  float zp,fp;
1107 217 equemene

1108 217 equemene
  // Autosize for image
1109 217 equemene
  bmx=1.25*re;
1110 217 equemene

1111 217 equemene
  // Angular step of integration
1112 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1113 217 equemene

1114 217 equemene
  // impact parameter
1115 217 equemene
  b=(float)bi/(float)bmaxi*bmx;
1116 217 equemene

1117 217 equemene
  float up,vp,pp,us,vs,ps;
1118 217 equemene

1119 217 equemene
  up=0.;
1120 217 equemene
  vp=1.;
1121 217 equemene

1122 217 equemene
  pp=0.;
1123 217 equemene
  nh=0;
1124 217 equemene

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

1127 217 equemene
  // b versus us
1128 217 equemene
  float bvus=fabs(b/us);
1129 217 equemene
  float bvus0=bvus;
1130 224 equemene
  Trajectories[bi*TRACKPOINTS+nh]=bvus;
1131 217 equemene

1132 217 equemene
  do
1133 217 equemene
  {
1134 217 equemene
     nh++;
1135 217 equemene
     pp=ps;
1136 217 equemene
     up=us;
1137 217 equemene
     vp=vs;
1138 217 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1139 217 equemene
     bvus=fabs(b/us);
1140 224 equemene
     Trajectories[bi*TRACKPOINTS+nh]=bvus;
1141 217 equemene

1142 217 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
1143 217 equemene

1144 217 equemene
  IdLast[bi]=nh;
1145 217 equemene

1146 217 equemene
}
1147 224 equemene

1148 224 equemene
__global__ void Original(float *zImage,float *fImage,
1149 224 equemene
                         float Mass,float InternalRadius,
1150 224 equemene
                         float ExternalRadius,float Angle,
1151 224 equemene
                         int Line)
1152 224 equemene
{
1153 224 equemene
  // Integer Impact Parameter ID
1154 224 equemene
  int bi=blockIdx.x*blockDim.x+threadIdx.x;
1155 224 equemene

1156 224 equemene
  // Integer Impact Parameter Size (half of image)
1157 224 equemene
  int bmaxi=gridDim.x*blockDim.x;
1158 224 equemene

1159 224 equemene
  float Trajectories[2048];
1160 224 equemene

1161 224 equemene
  // Perform trajectory for each pixel
1162 224 equemene

1163 224 equemene
  float m,rs,ri,re,tho;
1164 224 equemene
  int raie,q;
1165 224 equemene

1166 224 equemene
  m=Mass;
1167 224 equemene
  rs=2.*m;
1168 224 equemene
  ri=InternalRadius;
1169 224 equemene
  re=ExternalRadius;
1170 224 equemene
  tho=Angle;
1171 224 equemene
  q=-2;
1172 224 equemene
  raie=Line;
1173 224 equemene

1174 224 equemene
  float bmx,db,b,h;
1175 224 equemene
  int nh;
1176 224 equemene

1177 224 equemene
  // Autosize for image
1178 224 equemene
  bmx=1.25e0f*re;
1179 224 equemene

1180 224 equemene
  // Angular step of integration
1181 224 equemene
  h=4.e0f*PI/(float)TRACKPOINTS;
1182 224 equemene

1183 224 equemene
  // impact parameter
1184 224 equemene
  b=(float)bi/(float)bmaxi*bmx;
1185 224 equemene
  db=bmx/(2.e0f*(float)bmaxi);
1186 224 equemene

1187 224 equemene
  float up,vp,pp,us,vs,ps;
1188 224 equemene

1189 224 equemene
  up=0.;
1190 224 equemene
  vp=1.;
1191 224 equemene

1192 224 equemene
  pp=0.;
1193 224 equemene
  nh=0;
1194 224 equemene

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

1197 224 equemene
  // b versus us
1198 224 equemene
  float bvus=fabs(b/us);
1199 224 equemene
  float bvus0=bvus;
1200 224 equemene
  Trajectories[nh]=bvus;
1201 224 equemene

1202 224 equemene
  do
1203 224 equemene
  {
1204 224 equemene
     nh++;
1205 224 equemene
     pp=ps;
1206 224 equemene
     up=us;
1207 224 equemene
     vp=vs;
1208 224 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
1209 224 equemene
     bvus=fabs(b/us);
1210 224 equemene
     Trajectories[nh]=bvus;
1211 224 equemene

1212 224 equemene
  } while ((bvus>=rs)&&(bvus<=bvus0));
1213 224 equemene

1214 224 equemene
  int imx=(int)(16*bi);
1215 224 equemene

1216 224 equemene
  for (int i=0;i<imx;i++)
1217 224 equemene
  {
1218 224 equemene
     float zp=0,fp=0;
1219 224 equemene
     float phi=2.*PI/(float)imx*(float)i;
1220 224 equemene
     float phd=atanp(cos(phi)*sin(tho),cos(tho));
1221 224 equemene
     uint yi=(uint)((float)bi*sin(phi)+bmaxi);
1222 224 equemene
     uint xi=(uint)((float)bi*cos(phi)+bmaxi);
1223 224 equemene

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

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

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

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

1248 224 equemene
        HalfLap++;
1249 224 equemene

1250 224 equemene
     } while ((HalfLap<=2)&&(ExitOnImpact==0));
1251 224 equemene

1252 224 equemene
   __syncthreads();
1253 224 equemene

1254 224 equemene
   zImage[yi+2*bmaxi*xi]=zp;
1255 224 equemene
   fImage[yi+2*bmaxi*xi]=fp;
1256 224 equemene

1257 224 equemene
  }
1258 224 equemene

1259 224 equemene
}
1260 217 equemene
"""
1261 217 equemene
    return(BlobCUDA)
1262 217 equemene
1263 211 equemene
# def ImageOutput(sigma,prefix):
1264 211 equemene
#     from PIL import Image
1265 211 equemene
#     Max=sigma.max()
1266 211 equemene
#     Min=sigma.min()
1267 199 equemene
1268 211 equemene
#     # Normalize value as 8bits Integer
1269 211 equemene
#     SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
1270 211 equemene
#     image = Image.fromarray(SigmaInt)
1271 211 equemene
#     image.save("%s.jpg" % prefix)
1272 211 equemene
1273 211 equemene
def ImageOutput(sigma,prefix,Colors):
1274 211 equemene
    import matplotlib.pyplot as plt
1275 222 equemene
    start_time=time.time()
1276 211 equemene
    if Colors == 'Red2Yellow':
1277 211 equemene
        plt.imsave("%s.png" % prefix, sigma, cmap='afmhot')
1278 211 equemene
    else:
1279 211 equemene
        plt.imsave("%s.png" % prefix, sigma, cmap='Greys_r')
1280 222 equemene
    save_time = time.time()-start_time
1281 222 equemene
    print("Save Time : %f" % save_time)
1282 211 equemene
1283 204 equemene
def BlackHoleCL(zImage,fImage,InputCL):
1284 199 equemene
1285 199 equemene
    kernel_params = {}
1286 199 equemene
1287 204 equemene
    print(InputCL)
1288 204 equemene
1289 204 equemene
    Device=InputCL['Device']
1290 204 equemene
    GpuStyle=InputCL['GpuStyle']
1291 204 equemene
    VariableType=InputCL['VariableType']
1292 204 equemene
    Size=InputCL['Size']
1293 204 equemene
    Mass=InputCL['Mass']
1294 204 equemene
    InternalRadius=InputCL['InternalRadius']
1295 204 equemene
    ExternalRadius=InputCL['ExternalRadius']
1296 204 equemene
    Angle=InputCL['Angle']
1297 204 equemene
    Method=InputCL['Method']
1298 204 equemene
    TrackPoints=InputCL['TrackPoints']
1299 211 equemene
    Physics=InputCL['Physics']
1300 204 equemene
1301 211 equemene
    PhysicsList=DictionariesAPI()
1302 211 equemene
1303 204 equemene
    if InputCL['BlackBody']:
1304 209 equemene
        # Spectrum is Black Body one
1305 209 equemene
        Line=0
1306 204 equemene
    else:
1307 209 equemene
        # Spectrum is Monochromatic Line one
1308 209 equemene
        Line=1
1309 204 equemene
1310 217 equemene
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1311 217 equemene
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1312 204 equemene
1313 199 equemene
    # Je detecte un peripherique GPU dans la liste des peripheriques
1314 199 equemene
    Id=0
1315 199 equemene
    HasXPU=False
1316 199 equemene
    for platform in cl.get_platforms():
1317 199 equemene
        for device in platform.get_devices():
1318 199 equemene
            if Id==Device:
1319 199 equemene
                XPU=device
1320 217 equemene
                print("CPU/GPU selected: ",device.name.lstrip())
1321 199 equemene
                HasXPU=True
1322 199 equemene
            Id+=1
1323 199 equemene
1324 199 equemene
    if HasXPU==False:
1325 217 equemene
        print("No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1))
1326 199 equemene
        sys.exit()
1327 199 equemene
1328 199 equemene
    ctx = cl.Context([XPU])
1329 199 equemene
    queue = cl.CommandQueue(ctx,
1330 199 equemene
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
1331 199 equemene
1332 199 equemene
1333 211 equemene
    #    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
1334 211 equemene
1335 224 equemene
    BuildOptions="-cl-mad-enable -DPHYSICS=%i -DSETTRACKPOINTS=%i " % (PhysicsList[Physics],InputCL['TrackPoints'])
1336 211 equemene
1337 211 equemene
    BlackHoleCL = cl.Program(ctx,BlobOpenCL).build(options = BuildOptions)
1338 211 equemene
1339 199 equemene
    # Je recupere les flag possibles pour les buffers
1340 199 equemene
    mf = cl.mem_flags
1341 199 equemene
1342 204 equemene
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
1343 201 equemene
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
1344 201 equemene
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
1345 204 equemene
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
1346 199 equemene
1347 199 equemene
    start_time=time.time()
1348 199 equemene
1349 204 equemene
    if Method=='EachPixel':
1350 204 equemene
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
1351 204 equemene
                                       None,zImageCL,fImageCL,
1352 204 equemene
                                       numpy.float32(Mass),
1353 204 equemene
                                       numpy.float32(InternalRadius),
1354 204 equemene
                                       numpy.float32(ExternalRadius),
1355 204 equemene
                                       numpy.float32(Angle),
1356 209 equemene
                                       numpy.int32(Line))
1357 204 equemene
        CLLaunch.wait()
1358 224 equemene
    elif Method=='Original':
1359 224 equemene
        CLLaunch=BlackHoleCL.Original(queue,(zImage.shape[0]/2,),
1360 224 equemene
                                      None,zImageCL,fImageCL,
1361 224 equemene
                                      numpy.float32(Mass),
1362 224 equemene
                                      numpy.float32(InternalRadius),
1363 224 equemene
                                      numpy.float32(ExternalRadius),
1364 224 equemene
                                      numpy.float32(Angle),
1365 224 equemene
                                      numpy.int32(Line))
1366 224 equemene
        CLLaunch.wait()
1367 204 equemene
    elif Method=='TrajectoCircle':
1368 204 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1369 204 equemene
                                        None,TrajectoriesCL,IdLastCL,
1370 204 equemene
                                        numpy.float32(Mass),
1371 204 equemene
                                        numpy.float32(InternalRadius),
1372 204 equemene
                                        numpy.float32(ExternalRadius),
1373 204 equemene
                                        numpy.float32(Angle),
1374 209 equemene
                                        numpy.int32(Line))
1375 204 equemene
1376 204 equemene
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
1377 204 equemene
                                           zImage.shape[0]*4),None,
1378 204 equemene
                                    TrajectoriesCL,IdLastCL,
1379 204 equemene
                                    zImageCL,fImageCL,
1380 204 equemene
                                    numpy.float32(Mass),
1381 204 equemene
                                    numpy.float32(InternalRadius),
1382 204 equemene
                                    numpy.float32(ExternalRadius),
1383 204 equemene
                                    numpy.float32(Angle),
1384 209 equemene
                                    numpy.int32(Line))
1385 204 equemene
        CLLaunch.wait()
1386 204 equemene
    else:
1387 204 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
1388 204 equemene
                                        None,TrajectoriesCL,IdLastCL,
1389 204 equemene
                                        numpy.float32(Mass),
1390 204 equemene
                                        numpy.float32(InternalRadius),
1391 204 equemene
                                        numpy.float32(ExternalRadius),
1392 204 equemene
                                        numpy.float32(Angle),
1393 209 equemene
                                        numpy.int32(Line))
1394 204 equemene
1395 204 equemene
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],
1396 217 equemene
                                          zImage.shape[1]),None,
1397 204 equemene
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
1398 204 equemene
                                   numpy.uint32(Trajectories.shape[0]),
1399 204 equemene
                                        numpy.float32(Mass),
1400 204 equemene
                                        numpy.float32(InternalRadius),
1401 204 equemene
                                        numpy.float32(ExternalRadius),
1402 204 equemene
                                        numpy.float32(Angle),
1403 209 equemene
                                        numpy.int32(Line))
1404 204 equemene
        CLLaunch.wait()
1405 204 equemene
1406 218 equemene
    compute = time.time()-start_time
1407 199 equemene
1408 201 equemene
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
1409 201 equemene
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
1410 204 equemene
    cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
1411 204 equemene
    cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
1412 218 equemene
    elapsed = time.time()-start_time
1413 218 equemene
    print("\nCompute Time : %f" % compute)
1414 218 equemene
    print("Elapsed Time : %f\n" % elapsed)
1415 211 equemene
1416 211 equemene
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1417 211 equemene
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1418 219 equemene
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1419 219 equemene
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1420 201 equemene
    zImageCL.release()
1421 201 equemene
    fImageCL.release()
1422 204 equemene
1423 204 equemene
    TrajectoriesCL.release()
1424 204 equemene
    IdLastCL.release()
1425 204 equemene
1426 199 equemene
    return(elapsed)
1427 199 equemene
1428 217 equemene
def BlackHoleCUDA(zImage,fImage,InputCL):
1429 217 equemene
1430 217 equemene
    kernel_params = {}
1431 217 equemene
1432 217 equemene
    print(InputCL)
1433 217 equemene
1434 217 equemene
    Device=InputCL['Device']
1435 217 equemene
    GpuStyle=InputCL['GpuStyle']
1436 217 equemene
    VariableType=InputCL['VariableType']
1437 217 equemene
    Size=InputCL['Size']
1438 217 equemene
    Mass=InputCL['Mass']
1439 217 equemene
    InternalRadius=InputCL['InternalRadius']
1440 217 equemene
    ExternalRadius=InputCL['ExternalRadius']
1441 217 equemene
    Angle=InputCL['Angle']
1442 217 equemene
    Method=InputCL['Method']
1443 217 equemene
    TrackPoints=InputCL['TrackPoints']
1444 217 equemene
    Physics=InputCL['Physics']
1445 217 equemene
1446 217 equemene
    PhysicsList=DictionariesAPI()
1447 217 equemene
1448 217 equemene
    if InputCL['BlackBody']:
1449 217 equemene
        # Spectrum is Black Body one
1450 217 equemene
        Line=0
1451 217 equemene
    else:
1452 217 equemene
        # Spectrum is Monochromatic Line one
1453 217 equemene
        Line=1
1454 217 equemene
1455 217 equemene
    Trajectories=numpy.zeros((int(InputCL['Size']/2),InputCL['TrackPoints']),dtype=numpy.float32)
1456 217 equemene
    IdLast=numpy.zeros(int(InputCL['Size']/2),dtype=numpy.int32)
1457 217 equemene
1458 217 equemene
    try:
1459 217 equemene
        # For PyCUDA import
1460 217 equemene
        import pycuda.driver as cuda
1461 217 equemene
        from pycuda.compiler import SourceModule
1462 217 equemene
1463 217 equemene
        cuda.init()
1464 217 equemene
        for Id in range(cuda.Device.count()):
1465 217 equemene
            if Id==Device:
1466 217 equemene
                XPU=cuda.Device(Id)
1467 217 equemene
                print("GPU selected %s" % XPU.name())
1468 217 equemene
        print
1469 217 equemene
1470 217 equemene
    except ImportError:
1471 217 equemene
        print("Platform does not seem to support CUDA")
1472 217 equemene
1473 217 equemene
    Context=XPU.make_context()
1474 217 equemene
1475 217 equemene
    try:
1476 224 equemene
        mod = SourceModule(KernelCodeCuda(),options=['--compiler-options','-DPHYSICS=%i -DSETTRACKPOINTS=%i' % (PhysicsList[Physics],TrackPoints)])
1477 217 equemene
        print("Compilation seems to be OK")
1478 217 equemene
    except:
1479 217 equemene
        print("Compilation seems to break")
1480 217 equemene
1481 217 equemene
    EachPixelCU=mod.get_function("EachPixel")
1482 224 equemene
    OriginalCU=mod.get_function("Original")
1483 217 equemene
    TrajectoryCU=mod.get_function("Trajectory")
1484 217 equemene
    PixelCU=mod.get_function("Pixel")
1485 217 equemene
    CircleCU=mod.get_function("Circle")
1486 217 equemene
1487 217 equemene
    TrajectoriesCU = cuda.mem_alloc(Trajectories.size*Trajectories.dtype.itemsize)
1488 217 equemene
    cuda.memcpy_htod(TrajectoriesCU, Trajectories)
1489 217 equemene
    zImageCU = cuda.mem_alloc(zImage.size*zImage.dtype.itemsize)
1490 217 equemene
    cuda.memcpy_htod(zImageCU, zImage)
1491 217 equemene
    fImageCU = cuda.mem_alloc(fImage.size*fImage.dtype.itemsize)
1492 217 equemene
    cuda.memcpy_htod(zImageCU, fImage)
1493 217 equemene
    IdLastCU = cuda.mem_alloc(IdLast.size*IdLast.dtype.itemsize)
1494 217 equemene
    cuda.memcpy_htod(IdLastCU, IdLast)
1495 217 equemene
1496 217 equemene
    start_time=time.time()
1497 217 equemene
1498 217 equemene
    if Method=='EachPixel':
1499 217 equemene
        EachPixelCU(zImageCU,fImageCU,
1500 217 equemene
                    numpy.float32(Mass),
1501 217 equemene
                    numpy.float32(InternalRadius),
1502 217 equemene
                    numpy.float32(ExternalRadius),
1503 217 equemene
                    numpy.float32(Angle),
1504 217 equemene
                    numpy.int32(Line),
1505 219 equemene
                    grid=(zImage.shape[0]/32,zImage.shape[1]/32),
1506 219 equemene
                    block=(32,32,1))
1507 224 equemene
    elif Method=='Original':
1508 224 equemene
        OriginalCU(zImageCU,fImageCU,
1509 224 equemene
                   numpy.float32(Mass),
1510 224 equemene
                   numpy.float32(InternalRadius),
1511 224 equemene
                   numpy.float32(ExternalRadius),
1512 224 equemene
                   numpy.float32(Angle),
1513 224 equemene
                   numpy.int32(Line),
1514 224 equemene
                   grid=(zImage.shape[0]/32/2,1),
1515 224 equemene
                   block=(32,1,1))
1516 217 equemene
    elif Method=='TrajectoCircle':
1517 217 equemene
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1518 217 equemene
                     numpy.float32(Mass),
1519 217 equemene
                     numpy.float32(InternalRadius),
1520 217 equemene
                     numpy.float32(ExternalRadius),
1521 217 equemene
                     numpy.float32(Angle),
1522 217 equemene
                     numpy.int32(Line),
1523 219 equemene
                     grid=(Trajectories.shape[0]/32,1),
1524 219 equemene
                     block=(32,1,1))
1525 217 equemene
1526 217 equemene
        CircleCU(TrajectoriesCU,IdLastCU,zImageCU,fImageCU,
1527 217 equemene
                 numpy.float32(Mass),
1528 217 equemene
                 numpy.float32(InternalRadius),
1529 217 equemene
                 numpy.float32(ExternalRadius),
1530 217 equemene
                 numpy.float32(Angle),
1531 217 equemene
                 numpy.int32(Line),
1532 219 equemene
                 grid=(Trajectories.shape[0]/32,zImage.shape[0]*4/32),
1533 219 equemene
                 block=(32,32,1))
1534 217 equemene
    else:
1535 217 equemene
        TrajectoryCU(TrajectoriesCU,IdLastCU,
1536 217 equemene
                     numpy.float32(Mass),
1537 217 equemene
                     numpy.float32(InternalRadius),
1538 217 equemene
                     numpy.float32(ExternalRadius),
1539 217 equemene
                     numpy.float32(Angle),
1540 217 equemene
                     numpy.int32(Line),
1541 219 equemene
                     grid=(Trajectories.shape[0]/32,1),
1542 219 equemene
                     block=(32,1,1))
1543 217 equemene
1544 217 equemene
        PixelCU(zImageCU,fImageCU,TrajectoriesCU,IdLastCU,
1545 217 equemene
                numpy.uint32(Trajectories.shape[0]),
1546 217 equemene
                numpy.float32(Mass),
1547 217 equemene
                numpy.float32(InternalRadius),
1548 217 equemene
                numpy.float32(ExternalRadius),
1549 217 equemene
                numpy.float32(Angle),
1550 217 equemene
                numpy.int32(Line),
1551 219 equemene
                grid=(zImage.shape[0]/32,zImage.shape[1]/32,1),
1552 219 equemene
                block=(32,32,1))
1553 217 equemene
1554 220 equemene
    Context.synchronize()
1555 220 equemene
1556 218 equemene
    compute = time.time()-start_time
1557 217 equemene
1558 217 equemene
    cuda.memcpy_dtoh(zImage,zImageCU)
1559 217 equemene
    cuda.memcpy_dtoh(fImage,fImageCU)
1560 218 equemene
    elapsed = time.time()-start_time
1561 218 equemene
    print("\nCompute Time : %f" % compute)
1562 218 equemene
    print("Elapsed Time : %f\n" % elapsed)
1563 217 equemene
1564 217 equemene
    zMaxPosition=numpy.where(zImage[:,:]==zImage.max())
1565 217 equemene
    fMaxPosition=numpy.where(fImage[:,:]==fImage.max())
1566 219 equemene
    print("Z max @(%i,%i) : %f" % (zMaxPosition[1][0],zMaxPosition[0][0],zImage.max()))
1567 219 equemene
    print("Flux max @(%i,%i) : %f\n" % (fMaxPosition[1][0],fMaxPosition[0][0],fImage.max()))
1568 217 equemene
1569 220 equemene
1570 218 equemene
    Context.pop()
1571 220 equemene
1572 217 equemene
    Context.detach()
1573 217 equemene
1574 217 equemene
    return(elapsed)
1575 217 equemene
1576 199 equemene
if __name__=='__main__':
1577 199 equemene
1578 199 equemene
    GpuStyle = 'OpenCL'
1579 199 equemene
    Mass = 1.
1580 199 equemene
    # Internal Radius 3 times de Schwarzschild Radius
1581 199 equemene
    InternalRadius=6.*Mass
1582 199 equemene
    #
1583 199 equemene
    ExternalRadius=12.
1584 199 equemene
    #
1585 199 equemene
    # Angle with normal to disc 10 degrees
1586 199 equemene
    Angle = numpy.pi/180.*(90.-10.)
1587 199 equemene
    # Radiation of disc : BlackBody or Monochromatic
1588 209 equemene
    BlackBody = False
1589 199 equemene
    # Size of image
1590 199 equemene
    Size=256
1591 199 equemene
    # Variable Type
1592 199 equemene
    VariableType='FP32'
1593 199 equemene
    # ?
1594 199 equemene
    q=-2
1595 204 equemene
    # Method of resolution
1596 209 equemene
    Method='TrajectoPixel'
1597 211 equemene
    # Colors for output image
1598 211 equemene
    Colors='Greyscale'
1599 211 equemene
    # Physics
1600 211 equemene
    Physics='Einstein'
1601 211 equemene
    # No output as image
1602 211 equemene
    NoImage = False
1603 211 equemene
1604 224 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/Original> -v <FP32/FP64>'
1605 199 equemene
1606 199 equemene
    try:
1607 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="])
1608 199 equemene
    except getopt.GetoptError:
1609 199 equemene
        print(HowToUse % sys.argv[0])
1610 199 equemene
        sys.exit(2)
1611 199 equemene
1612 199 equemene
    # List of Devices
1613 199 equemene
    Devices=[]
1614 199 equemene
    Alu={}
1615 199 equemene
1616 199 equemene
    for opt, arg in opts:
1617 199 equemene
        if opt == '-h':
1618 199 equemene
            print(HowToUse % sys.argv[0])
1619 199 equemene
1620 199 equemene
            print("\nInformations about devices detected under OpenCL API:")
1621 199 equemene
            # For PyOpenCL import
1622 199 equemene
            try:
1623 199 equemene
                import pyopencl as cl
1624 199 equemene
                Id=0
1625 199 equemene
                for platform in cl.get_platforms():
1626 199 equemene
                    for device in platform.get_devices():
1627 199 equemene
                        #deviceType=cl.device_type.to_string(device.type)
1628 199 equemene
                        deviceType="xPU"
1629 199 equemene
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
1630 199 equemene
                        Id=Id+1
1631 199 equemene
1632 199 equemene
            except:
1633 199 equemene
                print("Your platform does not seem to support OpenCL")
1634 199 equemene
1635 199 equemene
            print("\nInformations about devices detected under CUDA API:")
1636 199 equemene
                # For PyCUDA import
1637 199 equemene
            try:
1638 199 equemene
                import pycuda.driver as cuda
1639 199 equemene
                cuda.init()
1640 199 equemene
                for Id in range(cuda.Device.count()):
1641 199 equemene
                    device=cuda.Device(Id)
1642 199 equemene
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
1643 199 equemene
                print
1644 199 equemene
            except:
1645 199 equemene
                print("Your platform does not seem to support CUDA")
1646 199 equemene
1647 199 equemene
            sys.exit()
1648 199 equemene
1649 199 equemene
        elif opt in ("-d", "--device"):
1650 199 equemene
#            Devices.append(int(arg))
1651 199 equemene
            Device=int(arg)
1652 199 equemene
        elif opt in ("-g", "--gpustyle"):
1653 199 equemene
            GpuStyle = arg
1654 204 equemene
        elif opt in ("-v", "--variabletype"):
1655 199 equemene
            VariableType = arg
1656 199 equemene
        elif opt in ("-s", "--size"):
1657 199 equemene
            Size = int(arg)
1658 199 equemene
        elif opt in ("-m", "--mass"):
1659 199 equemene
            Mass = float(arg)
1660 199 equemene
        elif opt in ("-i", "--internal"):
1661 199 equemene
            InternalRadius = float(arg)
1662 199 equemene
        elif opt in ("-e", "--external"):
1663 199 equemene
            ExternalRadius = float(arg)
1664 199 equemene
        elif opt in ("-a", "--angle"):
1665 199 equemene
            Angle = numpy.pi/180.*(90.-float(arg))
1666 199 equemene
        elif opt in ("-b", "--blackbody"):
1667 199 equemene
            BlackBody = True
1668 211 equemene
        elif opt in ("-n", "--noimage"):
1669 211 equemene
            NoImage = True
1670 204 equemene
        elif opt in ("-t", "--method"):
1671 204 equemene
            Method = arg
1672 211 equemene
        elif opt in ("-c", "--colors"):
1673 211 equemene
            Colors = arg
1674 211 equemene
        elif opt in ("-p", "--physics"):
1675 211 equemene
            Physics = arg
1676 199 equemene
1677 199 equemene
    print("Device Identification selected : %s" % Device)
1678 199 equemene
    print("GpuStyle used : %s" % GpuStyle)
1679 199 equemene
    print("VariableType : %s" % VariableType)
1680 199 equemene
    print("Size : %i" % Size)
1681 199 equemene
    print("Mass : %f" % Mass)
1682 199 equemene
    print("Internal Radius : %f" % InternalRadius)
1683 199 equemene
    print("External Radius : %f" % ExternalRadius)
1684 199 equemene
    print("Angle with normal of (in radians) : %f" % Angle)
1685 199 equemene
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
1686 204 equemene
    print("Method of resolution : %s" % Method)
1687 211 equemene
    print("Colors for output images : %s" % Colors)
1688 211 equemene
    print("Physics used for Trajectories : %s" % Physics)
1689 199 equemene
1690 199 equemene
    if GpuStyle=='CUDA':
1691 199 equemene
        print("\nSelection of CUDA device")
1692 199 equemene
        try:
1693 199 equemene
            # For PyCUDA import
1694 199 equemene
            import pycuda.driver as cuda
1695 199 equemene
1696 199 equemene
            cuda.init()
1697 199 equemene
            for Id in range(cuda.Device.count()):
1698 199 equemene
                device=cuda.Device(Id)
1699 199 equemene
                print("Device #%i of type GPU : %s" % (Id,device.name()))
1700 199 equemene
                if Id in Devices:
1701 199 equemene
                    Alu[Id]='GPU'
1702 199 equemene
1703 199 equemene
        except ImportError:
1704 199 equemene
            print("Platform does not seem to support CUDA")
1705 199 equemene
1706 199 equemene
    if GpuStyle=='OpenCL':
1707 199 equemene
        print("\nSelection of OpenCL device")
1708 199 equemene
        try:
1709 199 equemene
            # For PyOpenCL import
1710 199 equemene
            import pyopencl as cl
1711 199 equemene
            Id=0
1712 199 equemene
            for platform in cl.get_platforms():
1713 199 equemene
                for device in platform.get_devices():
1714 199 equemene
                    #deviceType=cl.device_type.to_string(device.type)
1715 199 equemene
                    deviceType="xPU"
1716 199 equemene
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
1717 199 equemene
1718 199 equemene
                    if Id in Devices:
1719 199 equemene
                    # Set the Alu as detected Device Type
1720 199 equemene
                        Alu[Id]=deviceType
1721 199 equemene
                    Id=Id+1
1722 199 equemene
        except ImportError:
1723 199 equemene
            print("Platform does not seem to support OpenCL")
1724 199 equemene
1725 199 equemene
#    print(Devices,Alu)
1726 199 equemene
1727 199 equemene
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
1728 204 equemene
    TrackPoints=2048
1729 201 equemene
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1730 201 equemene
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
1731 199 equemene
1732 204 equemene
    InputCL={}
1733 204 equemene
    InputCL['Device']=Device
1734 204 equemene
    InputCL['GpuStyle']=GpuStyle
1735 204 equemene
    InputCL['VariableType']=VariableType
1736 204 equemene
    InputCL['Size']=Size
1737 204 equemene
    InputCL['Mass']=Mass
1738 204 equemene
    InputCL['InternalRadius']=InternalRadius
1739 204 equemene
    InputCL['ExternalRadius']=ExternalRadius
1740 204 equemene
    InputCL['Angle']=Angle
1741 204 equemene
    InputCL['BlackBody']=BlackBody
1742 204 equemene
    InputCL['Method']=Method
1743 204 equemene
    InputCL['TrackPoints']=TrackPoints
1744 211 equemene
    InputCL['Physics']=Physics
1745 199 equemene
1746 217 equemene
    if GpuStyle=='OpenCL':
1747 217 equemene
        duration=BlackHoleCL(zImage,fImage,InputCL)
1748 217 equemene
    else:
1749 217 equemene
        duration=BlackHoleCUDA(zImage,fImage,InputCL)
1750 217 equemene
1751 211 equemene
    Hostname=gethostname()
1752 211 equemene
    Date=time.strftime("%Y%m%d_%H%M%S")
1753 211 equemene
    ImageInfo="%s_Device%i_%s_%s" % (Method,Device,Hostname,Date)
1754 211 equemene
1755 211 equemene
1756 211 equemene
    if not NoImage:
1757 211 equemene
        ImageOutput(zImage,"TrouNoirZ_%s" % ImageInfo,Colors)
1758 211 equemene
        ImageOutput(fImage,"TrouNoirF_%s" % ImageInfo,Colors)