Statistiques
| Révision :

root / TrouNoir / TrouNoir.py @ 203

Historique | Voir | Annoter | Télécharger (23,64 ko)

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

45 203 equemene
#define PI (float)3.14159265359
46 203 equemene
#define nbr 200
47 203 equemene

48 203 equemene
float atanp(float x,float y)
49 203 equemene
{
50 203 equemene
  float angle;
51 203 equemene

52 203 equemene
  angle=atan2(y,x);
53 203 equemene

54 203 equemene
  if (angle<0.e0f)
55 203 equemene
    {
56 203 equemene
      angle+=(float)2.e0f*PI;
57 203 equemene
    }
58 203 equemene

59 203 equemene
  return angle;
60 203 equemene
}
61 203 equemene

62 203 equemene
float f(float v)
63 203 equemene
{
64 203 equemene
  return v;
65 203 equemene
}
66 203 equemene

67 203 equemene
float g(float u,float m,float b)
68 203 equemene
{
69 203 equemene
  return (3.e0f*m/b*pow(u,2)-u);
70 203 equemene
}
71 203 equemene

72 203 equemene

73 203 equemene
void calcul(float *us,float *vs,float up,float vp,
74 203 equemene
            float h,float m,float b)
75 203 equemene
{
76 203 equemene
  float c0,c1,c2,c3,d0,d1,d2,d3;
77 203 equemene

78 203 equemene
  c0=h*f(vp);
79 203 equemene
  c1=h*f(vp+c0/2.);
80 203 equemene
  c2=h*f(vp+c1/2.);
81 203 equemene
  c3=h*f(vp+c2);
82 203 equemene
  d0=h*g(up,m,b);
83 203 equemene
  d1=h*g(up+d0/2.,m,b);
84 203 equemene
  d2=h*g(up+d1/2.,m,b);
85 203 equemene
  d3=h*g(up+d2,m,b);
86 203 equemene

87 203 equemene
  *us=up+(c0+2.*c1+2.*c2+c3)/6.;
88 203 equemene
  *vs=vp+(d0+2.*d1+2.*d2+d3)/6.;
89 203 equemene
}
90 203 equemene

91 203 equemene
void rungekutta(float *ps,float *us,float *vs,
92 203 equemene
                float pp,float up,float vp,
93 203 equemene
                float h,float m,float b)
94 203 equemene
{
95 203 equemene
  calcul(us,vs,up,vp,h,m,b);
96 203 equemene
  *ps=pp+h;
97 203 equemene
}
98 203 equemene

99 203 equemene
float decalage_spectral(float r,float b,float phi,
100 203 equemene
                         float tho,float m)
101 203 equemene
{
102 203 equemene
  return (sqrt(1-3*m/r)/(1+sqrt(m/pow(r,3))*b*sin(tho)*sin(phi)));
103 203 equemene
}
104 203 equemene

105 203 equemene
float spectre(float rf,int q,float b,float db,
106 203 equemene
             float h,float r,float m,float bss)
107 203 equemene
{
108 203 equemene
  float flx;
109 203 equemene

110 203 equemene
  flx=pow(r/m,q)*pow(rf,4)*b;
111 203 equemene
  return(flx);
112 203 equemene
}
113 203 equemene

114 203 equemene
float spectre_cn(float rf,float b,float db,float h,float r,float m,float bss)
115 203 equemene
{
116 203 equemene
  double flx;
117 203 equemene
  double nu_rec,nu_em,qu,v,temp,temp_em,flux_int,m_point,planck,c,k,psc2,psk;
118 203 equemene
  int fi,posfreq;
119 203 equemene

120 203 equemene
  planck=6.62e-34;
121 203 equemene
  k=1.38e-23;
122 203 equemene
  psk=5.38e-11;
123 203 equemene
  psc2=7.35e-51;
124 203 equemene
  temp=3.e7;
125 203 equemene
  m_point=1.e14;
126 203 equemene
  v=1.-3./r;
127 203 equemene

128 203 equemene
  qu=1./sqrt(1.-3./r)/sqrt(r)*(sqrt(r)-sqrt(6.)+sqrt(3.)/2.*log((sqrt(r)+sqrt(3.))/(sqrt(r)-sqrt(3.))*(sqrt(6.)-sqrt(3.))/(sqrt(6.)+sqrt(3.))));
129 203 equemene

130 203 equemene
  temp_em=temp*sqrt(m)*exp(0.25*log(m_point))*exp(-0.75*log(r))*exp(-0.125*log(v))*exp(0.25*log(qu));
131 203 equemene

132 203 equemene
  flux_int=0;
133 203 equemene
  flx=0;
134 203 equemene

135 203 equemene
  for (fi=0;fi<nbr;fi++)
136 203 equemene
    {
137 203 equemene
      nu_em=bss*(float)fi/(float)nbr;
138 203 equemene
      nu_rec=nu_em*rf;
139 203 equemene
      posfreq=(int)(1./bss*nu_rec*(float)nbr);
140 203 equemene
      if ((posfreq>0)&&(posfreq<nbr))
141 203 equemene
        {
142 203 equemene
          //flux_int=2*planck/9.e16f*pow(nu_em,3)/(exp(planck*nu_em/k/temp_em)-1.)*pow(rf,3)*b*db*h;
143 203 equemene
          //flux_int=pow(nu_em,3)/(exp(planck*nu_em/k/temp_em)-1.)*pow(rf,3)*b;
144 203 equemene
          flux_int=2.*psc2*pow(nu_em,3)/(exp(psk*nu_em/temp_em)-1.)*pow(rf,3)*b;
145 203 equemene
          flx+=flux_int;
146 203 equemene
        }
147 203 equemene
    }
148 203 equemene
  //return(pow(rf,3)*b*temp_em);
149 203 equemene
  return((float)flx);
150 203 equemene
  return(flx);
151 203 equemene
}
152 203 equemene

153 203 equemene
void impact(float phi,float r,float b,float tho,float m,
154 203 equemene
            float *zp,float *fp,
155 203 equemene
            int q,float db,
156 203 equemene
            float h,int raie)
157 203 equemene
{
158 203 equemene
  float flx,rf,bss;
159 203 equemene

160 203 equemene
  rf=decalage_spectral(r,b,phi,tho,m);
161 203 equemene

162 203 equemene
  if (raie==0)
163 203 equemene
    {
164 203 equemene
      bss=2.;
165 203 equemene
      flx=spectre(rf,q,b,db,h,r,m,bss);
166 203 equemene
    }
167 203 equemene
  else
168 203 equemene
    {
169 203 equemene
      bss=1.e6;
170 203 equemene
      flx=spectre_cn(rf,b,db,h,r,m,bss);
171 203 equemene
    }
172 203 equemene

173 203 equemene
  *zp=1./rf;
174 203 equemene
  *fp=flx;
175 203 equemene

176 203 equemene
}
177 203 equemene

178 203 equemene
__kernel void EachPixel(__global float *zImage,__global float *fImage,
179 203 equemene
                        float Mass,float InternalRadius,
180 203 equemene
                        float ExternalRadius,float Angle,
181 203 equemene
                        float ObserverDistance,
182 203 equemene
                        int BlackBody,int AngularCamera)
183 203 equemene
{
184 203 equemene
   uint xi=(uint)get_global_id(0);
185 203 equemene
   uint yi=(uint)get_global_id(1);
186 203 equemene
   uint sizex=(uint)get_global_size(0);
187 203 equemene
   uint sizey=(uint)get_global_size(1);
188 203 equemene

189 203 equemene
   // Perform trajectory for each pixel, exit on hit
190 203 equemene

191 203 equemene
  float m,rs,ri,re,tho,ro;
192 203 equemene
  int q,dist,raie;
193 203 equemene

194 203 equemene
  m=Mass;
195 203 equemene
  rs=2.*m;
196 203 equemene
  ri=InternalRadius;
197 203 equemene
  re=ExternalRadius;
198 203 equemene
  ro=ObserverDistance;
199 203 equemene
  tho=Angle;
200 203 equemene
  q=-2;
201 203 equemene
  dist=AngularCamera;
202 203 equemene
  raie=BlackBody;
203 203 equemene

204 203 equemene
  float d,bmx,db,b,h;
205 203 equemene
  float rp[2048];
206 203 equemene
  float phi,thi,phd,php,nr,r;
207 203 equemene
  int nh;
208 203 equemene
  float zp,fp;
209 203 equemene

210 203 equemene
  // Autosize for image
211 203 equemene
  bmx=1.25*re;
212 203 equemene
  b=0.;
213 203 equemene

214 203 equemene
  h=4.e0f*PI/(float)2048;
215 203 equemene

216 203 equemene
  // set origin as center of image
217 203 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
218 203 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
219 203 equemene
  // angle extracted from cylindric symmetry
220 203 equemene
  phi=atanp(x,y);
221 203 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
222 203 equemene

223 203 equemene
  float up,vp,pp,us,vs,ps;
224 203 equemene

225 203 equemene
  // impact parameter
226 203 equemene
  b=sqrt(x*x+y*y)*(float)2.e0f/(float)sizex*bmx;
227 203 equemene
  // step of impact parameter;
228 203 equemene
  db=bmx/(float)(sizex/2);
229 203 equemene

230 203 equemene
  if (dist==1)
231 203 equemene
  {
232 203 equemene
     up=0.;
233 203 equemene
     vp=1.;
234 203 equemene
  }
235 203 equemene
  else
236 203 equemene
  {
237 203 equemene
     up=sin(thi);
238 203 equemene
     vp=cos(thi);
239 203 equemene
  }
240 203 equemene

241 203 equemene
  pp=0.;
242 203 equemene
  nh=0;
243 203 equemene

244 203 equemene
  rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
245 203 equemene

246 203 equemene
  rp[nh]=fabs(b/us);
247 203 equemene

248 203 equemene
  int ExitOnImpact=0;
249 203 equemene

250 203 equemene
  do
251 203 equemene
  {
252 203 equemene
     nh++;
253 203 equemene
     pp=ps;
254 203 equemene
     up=us;
255 203 equemene
     vp=vs;
256 203 equemene
     rungekutta(&ps,&us,&vs,pp,up,vp,h,m,b);
257 203 equemene
     rp[nh]=fabs(b/us);
258 203 equemene
     ExitOnImpact = ((fmod(pp,PI)<fmod(phd,PI))&&(fmod(ps,PI)>fmod(phd,PI)))&&(rp[nh]>ri)&&(rp[nh]<re)?1:0;
259 203 equemene

260 203 equemene
  } while ((rp[nh]>=rs)&&(rp[nh]<=rp[0])&&(ExitOnImpact==0));
261 203 equemene

262 203 equemene
  if (ExitOnImpact==1) {
263 203 equemene
     impact(phi,rp[nh-1],b,tho,m,&zp,&fp,q,db,h,raie);
264 203 equemene
  }
265 203 equemene
  else
266 203 equemene
  {
267 203 equemene
     zp=0.;
268 203 equemene
     fp=0.;
269 203 equemene
  }
270 203 equemene

271 203 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
272 203 equemene

273 203 equemene
  zImage[yi+sizex*xi]=(float)zp;
274 203 equemene
  fImage[yi+sizex*xi]=(float)fp;
275 203 equemene
}
276 203 equemene

277 203 equemene
__kernel void Pixel(__global float *zImage,__global float *fImage,
278 203 equemene
                    __global float *Trajectories,__global int *IdLast,
279 203 equemene
                    uint ImpactParameter,uint TrackPoints,
280 203 equemene
                    float Mass,float InternalRadius,
281 203 equemene
                    float ExternalRadius,float Angle,
282 203 equemene
                    float ObserverDistance,
283 203 equemene
                    int BlackBody,int AngularCamera)
284 203 equemene
{
285 203 equemene
   uint xi=(uint)get_global_id(0);
286 203 equemene
   uint yi=(uint)get_global_id(1);
287 203 equemene
   uint sizex=(uint)get_global_size(0);
288 203 equemene
   uint sizey=(uint)get_global_size(1);
289 203 equemene

290 203 equemene
   // Perform trajectory for each pixel
291 203 equemene

292 203 equemene
  float m,rs,ri,re,tho,ro;
293 203 equemene
  int q,dist,raie;
294 203 equemene

295 203 equemene
  m=Mass;
296 203 equemene
  rs=2.*m;
297 203 equemene
  ri=InternalRadius;
298 203 equemene
  re=ExternalRadius;
299 203 equemene
  ro=ObserverDistance;
300 203 equemene
  tho=Angle;
301 203 equemene
  q=-2;
302 203 equemene
  dist=AngularCamera;
303 203 equemene
  raie=BlackBody;
304 203 equemene

305 203 equemene
  float d,bmx,db,b,h;
306 203 equemene
  float phi,thi,phd,php,nr,r;
307 203 equemene
  int nh;
308 203 equemene
  float zp=0,fp=0;
309 203 equemene

310 203 equemene
  // Autosize for image, 25% greater than externa radius
311 203 equemene
  bmx=1.25*re;
312 203 equemene

313 203 equemene
  // Angular step of integration
314 203 equemene
  h=4.e0f*PI/(float)TrackPoints;
315 203 equemene

316 203 equemene
  // Step of Impact Parameter
317 203 equemene
  db=bmx/(float)ImpactParameter;
318 203 equemene

319 203 equemene
  // set origin as center of image
320 203 equemene
  float x=(float)xi-(float)(sizex/2)+(float)5e-1f;
321 203 equemene
  float y=(float)yi-(float)(sizey/2)+(float)5e-1f;
322 203 equemene

323 203 equemene
  // angle extracted from cylindric symmetry
324 203 equemene
  phi=atanp(x,y);
325 203 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
326 203 equemene

327 203 equemene
  // Real Impact Parameter
328 203 equemene
  b=sqrt(x*x+y*y)*bmx/(float)ImpactParameter;
329 203 equemene

330 203 equemene
  // Integer Impact Parameter
331 203 equemene
  uint bi=(uint)sqrt(x*x+y*y);
332 203 equemene

333 203 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
334 203 equemene

335 203 equemene
  if (bi<ImpactParameter)
336 203 equemene
  {
337 203 equemene
    do
338 203 equemene
    {
339 203 equemene
      php=phd+(float)HalfLap*PI;
340 203 equemene
      nr=php/h;
341 203 equemene
      ni=(int)nr;
342 203 equemene

343 203 equemene
      if (ni<IdLast[bi])
344 203 equemene
      {
345 203 equemene
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
346 203 equemene
      }
347 203 equemene
      else
348 203 equemene
      {
349 203 equemene
        r=Trajectories[bi*TrackPoints+ni];
350 203 equemene
      }
351 203 equemene

352 203 equemene
      if ((r<=re)&&(r>=ri))
353 203 equemene
      {
354 203 equemene
        ExitOnImpact=1;
355 203 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
356 203 equemene
      }
357 203 equemene

358 203 equemene
      HalfLap++;
359 203 equemene
    } while ((HalfLap<=2)&&(ExitOnImpact==0));
360 203 equemene

361 203 equemene
  }
362 203 equemene

363 203 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
364 203 equemene

365 203 equemene
  zImage[yi+sizex*xi]=zp;
366 203 equemene
  fImage[yi+sizex*xi]=fp;
367 203 equemene
}
368 203 equemene

369 203 equemene
__kernel void Circle(__global float *Trajectories,__global int *IdLast,
370 203 equemene
                     __global float *zImage,__global float *fImage,
371 203 equemene
                     int TrackPoints,
372 203 equemene
                     float Mass,float InternalRadius,
373 203 equemene
                     float ExternalRadius,float Angle,
374 203 equemene
                     float ObserverDistance,
375 203 equemene
                     int BlackBody,int AngularCamera)
376 203 equemene
{
377 203 equemene
   // Integer Impact Parameter ID
378 203 equemene
   int bi=get_global_id(0);
379 203 equemene
   // Integer points on circle
380 203 equemene
   int i=get_global_id(1);
381 203 equemene
   // Integer Impact Parameter Size (half of image)
382 203 equemene
   int bmaxi=get_global_size(0);
383 203 equemene
   // Integer Points on circle
384 203 equemene
   int imx=get_global_size(1);
385 203 equemene

386 203 equemene
   // Perform trajectory for each pixel
387 203 equemene

388 203 equemene
  float m,rs,ri,re,tho,ro;
389 203 equemene
  int q,dist,raie;
390 203 equemene

391 203 equemene
  m=Mass;
392 203 equemene
  rs=2.*m;
393 203 equemene
  ri=InternalRadius;
394 203 equemene
  re=ExternalRadius;
395 203 equemene
  ro=ObserverDistance;
396 203 equemene
  tho=Angle;
397 203 equemene
  q=-2;
398 203 equemene
  dist=AngularCamera;
399 203 equemene
  raie=BlackBody;
400 203 equemene

401 203 equemene
  float bmx,db,b,h;
402 203 equemene
  float phi,thi,phd;
403 203 equemene
  int nh;
404 203 equemene
  float zp=0,fp=0;
405 203 equemene

406 203 equemene
  // Autosize for image
407 203 equemene
  bmx=1.25*re;
408 203 equemene

409 203 equemene
  // Angular step of integration
410 203 equemene
  h=4.e0f*PI/(float)TrackPoints;
411 203 equemene

412 203 equemene
  // impact parameter
413 203 equemene
  b=(float)bi/(float)bmaxi*bmx;
414 203 equemene
  db=bmx/(float)bmaxi;
415 203 equemene

416 203 equemene
  phi=2.*PI/(float)imx*(float)i;
417 203 equemene
  phd=atanp(cos(phi)*sin(tho),cos(tho));
418 203 equemene
  int yi=(int)((float)bi*sin(phi))+bmaxi;
419 203 equemene
  int xi=(int)((float)bi*cos(phi))+bmaxi;
420 203 equemene

421 203 equemene
  int HalfLap=0,ExitOnImpact=0,ni;
422 203 equemene
  float php,nr,r;
423 203 equemene

424 203 equemene
  do
425 203 equemene
  {
426 203 equemene
     php=phd+(float)HalfLap*PI;
427 203 equemene
     nr=php/h;
428 203 equemene
     ni=(int)nr;
429 203 equemene

430 203 equemene
     if (ni<IdLast[bi])
431 203 equemene
     {
432 203 equemene
        r=(Trajectories[bi*TrackPoints+ni+1]-Trajectories[bi*TrackPoints+ni])*(nr-ni*1.)+Trajectories[bi*TrackPoints+ni];
433 203 equemene
     }
434 203 equemene
     else
435 203 equemene
     {
436 203 equemene
        r=Trajectories[bi*TrackPoints+ni];
437 203 equemene
     }
438 203 equemene

439 203 equemene
     if ((r<=re)&&(r>=ri))
440 203 equemene
     {
441 203 equemene
        ExitOnImpact=1;
442 203 equemene
        impact(phi,r,b,tho,m,&zp,&fp,q,db,h,raie);
443 203 equemene
     }
444 203 equemene

445 203 equemene
     HalfLap++;
446 203 equemene
  } while ((HalfLap<=2)&&(ExitOnImpact==0));
447 203 equemene

448 203 equemene
  zImage[yi+2*bmaxi*xi]=zp;
449 203 equemene
  fImage[yi+2*bmaxi*xi]=fp;
450 203 equemene

451 203 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
452 203 equemene

453 203 equemene
}
454 203 equemene

455 203 equemene
__kernel void Trajectory(__global float *Trajectories,
456 203 equemene
                         __global int *IdLast,int TrackPoints,
457 203 equemene
                         float Mass,float InternalRadius,
458 203 equemene
                         float ExternalRadius,float Angle,
459 203 equemene
                         float ObserverDistance,
460 203 equemene
                         int BlackBody,int AngularCamera)
461 203 equemene
{
462 203 equemene
  // Integer Impact Parameter ID
463 203 equemene
  int bi=get_global_id(0);
464 203 equemene
  // Integer Impact Parameter Size (half of image)
465 203 equemene
  int bmaxi=get_global_size(0);
466 203 equemene

467 203 equemene
  // Perform trajectory for each pixel
468 203 equemene

469 203 equemene
  float m,rs,ri,re,tho,ro;
470 203 equemene
  int dist,raie,q;
471 203 equemene

472 203 equemene
  m=Mass;
473 203 equemene
  rs=2.*m;
474 203 equemene
  ri=InternalRadius;
475 203 equemene
  re=ExternalRadius;
476 203 equemene
  ro=ObserverDistance;
477 203 equemene
  tho=Angle;
478 203 equemene
  q=-2;
479 203 equemene
  dist=AngularCamera;
480 203 equemene
  raie=BlackBody;
481 203 equemene

482 203 equemene
  float d,bmx,db,b,h;
483 203 equemene
  float phi,thi,phd,php,nr,r;
484 203 equemene
  int nh;
485 203 equemene
  float zp,fp;
486 203 equemene

487 203 equemene
  // Autosize for image
488 203 equemene
  bmx=1.25*re;
489 203 equemene

490 203 equemene
  // Angular step of integration
491 203 equemene
  h=4.e0f*PI/(float)TrackPoints;
492 203 equemene

493 203 equemene
  // impact parameter
494 203 equemene
  b=(float)bi/(float)bmaxi*bmx;
495 203 equemene

496 203 equemene
  float up,vp,pp,us,vs,ps;
497 203 equemene

498 203 equemene
  if (dist==1)
499 203 equemene
  {
500 203 equemene
     up=0.;
501 203 equemene
     vp=1.;
502 203 equemene
  }
503 203 equemene
  else
504 203 equemene
  {
505 203 equemene
     thi=asin(b/ro);
506 203 equemene
     up=sin(thi);
507 203 equemene
     vp=cos(thi);
508 203 equemene
  }
509 203 equemene

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

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

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

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

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

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

534 203 equemene
  barrier(CLK_GLOBAL_MEM_FENCE);
535 203 equemene

536 203 equemene
}
537 203 equemene
""")
538 203 equemene
539 203 equemene
def ImageOutput(sigma,prefix):
540 203 equemene
    Max=sigma.max()
541 203 equemene
    Min=sigma.min()
542 203 equemene
543 203 equemene
    # Normalize value as 8bits Integer
544 203 equemene
    SigmaInt=(255*(sigma-Min)/(Max-Min)).astype('uint8')
545 203 equemene
    image = Image.fromarray(SigmaInt)
546 203 equemene
    image.save("%s.jpg" % prefix)
547 203 equemene
548 203 equemene
def BlackHoleCL(zImage,fImage,InputCL):
549 203 equemene
550 203 equemene
    kernel_params = {}
551 203 equemene
552 203 equemene
    print(InputCL)
553 203 equemene
554 203 equemene
    Device=InputCL['Device']
555 203 equemene
    GpuStyle=InputCL['GpuStyle']
556 203 equemene
    VariableType=InputCL['VariableType']
557 203 equemene
    Size=InputCL['Size']
558 203 equemene
    Mass=InputCL['Mass']
559 203 equemene
    InternalRadius=InputCL['InternalRadius']
560 203 equemene
    ExternalRadius=InputCL['ExternalRadius']
561 203 equemene
    Angle=InputCL['Angle']
562 203 equemene
    ObserverDistance=InputCL['ObserverDistance']
563 203 equemene
    Method=InputCL['Method']
564 203 equemene
    TrackPoints=InputCL['TrackPoints']
565 203 equemene
566 203 equemene
    if InputCL['BlackBody']:
567 203 equemene
        BlackBody=0
568 203 equemene
    else:
569 203 equemene
        BlackBody=1
570 203 equemene
571 203 equemene
    if InputCL['AngularCamera']:
572 203 equemene
        AngularCamera=1
573 203 equemene
    else:
574 203 equemene
        AngularCamera=0
575 203 equemene
576 203 equemene
    Trajectories=numpy.zeros((InputCL['Size']/2,InputCL['TrackPoints']),
577 203 equemene
                              dtype=numpy.float32)
578 203 equemene
    IdLast=numpy.zeros(InputCL['Size']/2,dtype=numpy.int32)
579 203 equemene
580 203 equemene
    # Je detecte un peripherique GPU dans la liste des peripheriques
581 203 equemene
    Id=0
582 203 equemene
    HasXPU=False
583 203 equemene
    for platform in cl.get_platforms():
584 203 equemene
        for device in platform.get_devices():
585 203 equemene
            if Id==Device:
586 203 equemene
                XPU=device
587 203 equemene
                print "CPU/GPU selected: ",device.name.lstrip()
588 203 equemene
                HasXPU=True
589 203 equemene
            Id+=1
590 203 equemene
591 203 equemene
    if HasXPU==False:
592 203 equemene
        print "No XPU #%i found in all of %i devices, sorry..." % (Device,Id-1)
593 203 equemene
        sys.exit()
594 203 equemene
595 203 equemene
    ctx = cl.Context([XPU])
596 203 equemene
    queue = cl.CommandQueue(ctx,
597 203 equemene
                            properties=cl.command_queue_properties.PROFILING_ENABLE)
598 203 equemene
599 203 equemene
    BlackHoleCL = cl.Program(ctx,KERNEL_CODE.substitute(kernel_params)).build()
600 203 equemene
601 203 equemene
    # Je recupere les flag possibles pour les buffers
602 203 equemene
    mf = cl.mem_flags
603 203 equemene
604 203 equemene
    print(zImage)
605 203 equemene
606 203 equemene
    TrajectoriesCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=Trajectories)
607 203 equemene
    zImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=zImage)
608 203 equemene
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
609 203 equemene
    fImageCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=fImage)
610 203 equemene
    IdLastCL = cl.Buffer(ctx, mf.WRITE_ONLY | mf.COPY_HOST_PTR, hostbuf=IdLast)
611 203 equemene
612 203 equemene
    start_time=time.time()
613 203 equemene
614 203 equemene
    print(Trajectories.shape)
615 203 equemene
    if Method=='EachPixel':
616 203 equemene
        CLLaunch=BlackHoleCL.EachPixel(queue,(zImage.shape[0],zImage.shape[1]),
617 203 equemene
                                       None,zImageCL,fImageCL,
618 203 equemene
                                       numpy.float32(Mass),
619 203 equemene
                                       numpy.float32(InternalRadius),
620 203 equemene
                                       numpy.float32(ExternalRadius),
621 203 equemene
                                       numpy.float32(Angle),
622 203 equemene
                                       numpy.float32(ObserverDistance),
623 203 equemene
                                       numpy.int32(BlackBody),
624 203 equemene
                                       numpy.int32(AngularCamera))
625 203 equemene
        CLLaunch.wait()
626 203 equemene
    elif Method=='TrajectoCircle':
627 203 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
628 203 equemene
                                        None,TrajectoriesCL,IdLastCL,
629 203 equemene
                                        numpy.uint32(Trajectories.shape[1]),
630 203 equemene
                                        numpy.float32(Mass),
631 203 equemene
                                        numpy.float32(InternalRadius),
632 203 equemene
                                        numpy.float32(ExternalRadius),
633 203 equemene
                                        numpy.float32(Angle),
634 203 equemene
                                        numpy.float32(ObserverDistance),
635 203 equemene
                                        numpy.int32(BlackBody),
636 203 equemene
                                        numpy.int32(AngularCamera))
637 203 equemene
638 203 equemene
        CLLaunch=BlackHoleCL.Circle(queue,(Trajectories.shape[0],
639 203 equemene
                                           zImage.shape[0]*4),None,
640 203 equemene
                                    TrajectoriesCL,IdLastCL,
641 203 equemene
                                    zImageCL,fImageCL,
642 203 equemene
                                    numpy.uint32(Trajectories.shape[1]),
643 203 equemene
                                    numpy.float32(Mass),
644 203 equemene
                                    numpy.float32(InternalRadius),
645 203 equemene
                                    numpy.float32(ExternalRadius),
646 203 equemene
                                    numpy.float32(Angle),
647 203 equemene
                                    numpy.float32(ObserverDistance),
648 203 equemene
                                    numpy.int32(BlackBody),
649 203 equemene
                                    numpy.int32(AngularCamera))
650 203 equemene
        CLLaunch.wait()
651 203 equemene
    else:
652 203 equemene
        CLLaunch=BlackHoleCL.Trajectory(queue,(Trajectories.shape[0],),
653 203 equemene
                                        None,TrajectoriesCL,IdLastCL,
654 203 equemene
                                        numpy.uint32(Trajectories.shape[1]),
655 203 equemene
                                        numpy.float32(Mass),
656 203 equemene
                                        numpy.float32(InternalRadius),
657 203 equemene
                                        numpy.float32(ExternalRadius),
658 203 equemene
                                        numpy.float32(Angle),
659 203 equemene
                                        numpy.float32(ObserverDistance),
660 203 equemene
                                        numpy.int32(BlackBody),
661 203 equemene
                                        numpy.int32(AngularCamera))
662 203 equemene
663 203 equemene
        CLLaunch=BlackHoleCL.Pixel(queue,(zImage.shape[0],
664 203 equemene
                                                  zImage.shape[1]),None,
665 203 equemene
                                   zImageCL,fImageCL,TrajectoriesCL,IdLastCL,
666 203 equemene
                                   numpy.uint32(Trajectories.shape[0]),
667 203 equemene
                                   numpy.uint32(Trajectories.shape[1]),
668 203 equemene
                                        numpy.float32(Mass),
669 203 equemene
                                        numpy.float32(InternalRadius),
670 203 equemene
                                        numpy.float32(ExternalRadius),
671 203 equemene
                                        numpy.float32(Angle),
672 203 equemene
                                        numpy.float32(ObserverDistance),
673 203 equemene
                                        numpy.int32(BlackBody),
674 203 equemene
                                        numpy.int32(AngularCamera))
675 203 equemene
        CLLaunch.wait()
676 203 equemene
677 203 equemene
    elapsed = time.time()-start_time
678 203 equemene
    print("Elapsed %f: " % elapsed)
679 203 equemene
680 203 equemene
    cl.enqueue_copy(queue,zImage,zImageCL).wait()
681 203 equemene
    cl.enqueue_copy(queue,fImage,fImageCL).wait()
682 203 equemene
    cl.enqueue_copy(queue,Trajectories,TrajectoriesCL).wait()
683 203 equemene
    cl.enqueue_copy(queue,IdLast,IdLastCL).wait()
684 203 equemene
685 203 equemene
    print(zImage.max())
686 203 equemene
    print(fImage.max())
687 203 equemene
    zImageCL.release()
688 203 equemene
    fImageCL.release()
689 203 equemene
690 203 equemene
    TrajectoriesCL.release()
691 203 equemene
    IdLastCL.release()
692 203 equemene
693 203 equemene
    return(elapsed)
694 203 equemene
695 203 equemene
if __name__=='__main__':
696 203 equemene
697 203 equemene
    GpuStyle = 'OpenCL'
698 203 equemene
    Mass = 1.
699 203 equemene
    # Internal Radius 3 times de Schwarzschild Radius
700 203 equemene
    InternalRadius=6.*Mass
701 203 equemene
    #
702 203 equemene
    ExternalRadius=12.
703 203 equemene
    #
704 203 equemene
    ObserverDistance=100.
705 203 equemene
    # Angle with normal to disc 10 degrees
706 203 equemene
    Angle = numpy.pi/180.*(90.-10.)
707 203 equemene
    # Radiation of disc : BlackBody or Monochromatic
708 203 equemene
    BlackBody = True
709 203 equemene
    # Camera : Angular Camera or plate with the size of camera
710 203 equemene
    AngularCamera = True
711 203 equemene
    # Size of image
712 203 equemene
    Size=256
713 203 equemene
    # Variable Type
714 203 equemene
    VariableType='FP32'
715 203 equemene
    # ?
716 203 equemene
    q=-2
717 203 equemene
    # Method of resolution
718 203 equemene
    Method='Trajecto'
719 203 equemene
720 203 equemene
    HowToUse='%s -h [Help] -b [BlackBodyEmission] -c [AngularCamera] -s <SizeInPixels> -m <Mass> -i <DiscInternalRadius> -x <DiscExternalRadius> -o <ObservatorDistance> -a <AngleAboveDisc> -d <DeviceId> -g <CUDA/OpenCL> -t <EachPixel/TrajecoCircle/TrajectoPixel> -v <FP32/FP64>'
721 203 equemene
722 203 equemene
    try:
723 203 equemene
        opts, args = getopt.getopt(sys.argv[1:],"hbcs:m:i:x:o:a:d:g:v:t:",["blackbody","camera","size=","mass=","internal=","external=","observer=","angle=","device=","gpustyle=","variabletype=","method="])
724 203 equemene
    except getopt.GetoptError:
725 203 equemene
        print(HowToUse % sys.argv[0])
726 203 equemene
        sys.exit(2)
727 203 equemene
728 203 equemene
    # List of Devices
729 203 equemene
    Devices=[]
730 203 equemene
    Alu={}
731 203 equemene
732 203 equemene
    for opt, arg in opts:
733 203 equemene
        if opt == '-h':
734 203 equemene
            print(HowToUse % sys.argv[0])
735 203 equemene
736 203 equemene
            print("\nInformations about devices detected under OpenCL API:")
737 203 equemene
            # For PyOpenCL import
738 203 equemene
            try:
739 203 equemene
                import pyopencl as cl
740 203 equemene
                Id=0
741 203 equemene
                for platform in cl.get_platforms():
742 203 equemene
                    for device in platform.get_devices():
743 203 equemene
                        #deviceType=cl.device_type.to_string(device.type)
744 203 equemene
                        deviceType="xPU"
745 203 equemene
                        print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip(),deviceType,device.name.lstrip()))
746 203 equemene
                        Id=Id+1
747 203 equemene
748 203 equemene
            except:
749 203 equemene
                print("Your platform does not seem to support OpenCL")
750 203 equemene
751 203 equemene
            print("\nInformations about devices detected under CUDA API:")
752 203 equemene
                # For PyCUDA import
753 203 equemene
            try:
754 203 equemene
                import pycuda.driver as cuda
755 203 equemene
                cuda.init()
756 203 equemene
                for Id in range(cuda.Device.count()):
757 203 equemene
                    device=cuda.Device(Id)
758 203 equemene
                    print("Device #%i of type GPU : %s" % (Id,device.name()))
759 203 equemene
                print
760 203 equemene
            except:
761 203 equemene
                print("Your platform does not seem to support CUDA")
762 203 equemene
763 203 equemene
            sys.exit()
764 203 equemene
765 203 equemene
        elif opt in ("-d", "--device"):
766 203 equemene
#            Devices.append(int(arg))
767 203 equemene
            Device=int(arg)
768 203 equemene
        elif opt in ("-g", "--gpustyle"):
769 203 equemene
            GpuStyle = arg
770 203 equemene
        elif opt in ("-v", "--variabletype"):
771 203 equemene
            VariableType = arg
772 203 equemene
        elif opt in ("-s", "--size"):
773 203 equemene
            Size = int(arg)
774 203 equemene
        elif opt in ("-m", "--mass"):
775 203 equemene
            Mass = float(arg)
776 203 equemene
        elif opt in ("-i", "--internal"):
777 203 equemene
            InternalRadius = float(arg)
778 203 equemene
        elif opt in ("-e", "--external"):
779 203 equemene
            ExternalRadius = float(arg)
780 203 equemene
        elif opt in ("-o", "--observer"):
781 203 equemene
            ObserverDistance = float(arg)
782 203 equemene
        elif opt in ("-a", "--angle"):
783 203 equemene
            Angle = numpy.pi/180.*(90.-float(arg))
784 203 equemene
        elif opt in ("-b", "--blackbody"):
785 203 equemene
            BlackBody = True
786 203 equemene
        elif opt in ("-c", "--camera"):
787 203 equemene
            AngularCamera = True
788 203 equemene
        elif opt in ("-t", "--method"):
789 203 equemene
            Method = arg
790 203 equemene
791 203 equemene
    print("Device Identification selected : %s" % Device)
792 203 equemene
    print("GpuStyle used : %s" % GpuStyle)
793 203 equemene
    print("VariableType : %s" % VariableType)
794 203 equemene
    print("Size : %i" % Size)
795 203 equemene
    print("Mass : %f" % Mass)
796 203 equemene
    print("Internal Radius : %f" % InternalRadius)
797 203 equemene
    print("External Radius : %f" % ExternalRadius)
798 203 equemene
    print("Observer Distance : %f" % ObserverDistance)
799 203 equemene
    print("Angle with normal of (in radians) : %f" % Angle)
800 203 equemene
    print("Black Body Disc Emission (monochromatic instead) : %s" % BlackBody)
801 203 equemene
    print("Angular Camera (dimension of object instead) : %s" % AngularCamera)
802 203 equemene
    print("Method of resolution : %s" % Method)
803 203 equemene
804 203 equemene
    if GpuStyle=='CUDA':
805 203 equemene
        print("\nSelection of CUDA device")
806 203 equemene
        try:
807 203 equemene
            # For PyCUDA import
808 203 equemene
            import pycuda.driver as cuda
809 203 equemene
810 203 equemene
            cuda.init()
811 203 equemene
            for Id in range(cuda.Device.count()):
812 203 equemene
                device=cuda.Device(Id)
813 203 equemene
                print("Device #%i of type GPU : %s" % (Id,device.name()))
814 203 equemene
                if Id in Devices:
815 203 equemene
                    Alu[Id]='GPU'
816 203 equemene
817 203 equemene
        except ImportError:
818 203 equemene
            print("Platform does not seem to support CUDA")
819 203 equemene
820 203 equemene
    if GpuStyle=='OpenCL':
821 203 equemene
        print("\nSelection of OpenCL device")
822 203 equemene
        try:
823 203 equemene
            # For PyOpenCL import
824 203 equemene
            import pyopencl as cl
825 203 equemene
            Id=0
826 203 equemene
            for platform in cl.get_platforms():
827 203 equemene
                for device in platform.get_devices():
828 203 equemene
                    #deviceType=cl.device_type.to_string(device.type)
829 203 equemene
                    deviceType="xPU"
830 203 equemene
                    print("Device #%i from %s of type %s : %s" % (Id,platform.vendor.lstrip().rstrip(),deviceType,device.name.lstrip().rstrip()))
831 203 equemene
832 203 equemene
                    if Id in Devices:
833 203 equemene
                    # Set the Alu as detected Device Type
834 203 equemene
                        Alu[Id]=deviceType
835 203 equemene
                    Id=Id+1
836 203 equemene
        except ImportError:
837 203 equemene
            print("Platform does not seem to support OpenCL")
838 203 equemene
839 203 equemene
#    print(Devices,Alu)
840 203 equemene
841 203 equemene
#    MyImage=numpy.where(numpy.random.zeros(Size,Size)>0,1,-1).astype(numpy.float32)
842 203 equemene
    TrackPoints=2048
843 203 equemene
    zImage=numpy.zeros((Size,Size),dtype=numpy.float32)
844 203 equemene
    fImage=numpy.zeros((Size,Size),dtype=numpy.float32)
845 203 equemene
846 203 equemene
    InputCL={}
847 203 equemene
    InputCL['Device']=Device
848 203 equemene
    InputCL['GpuStyle']=GpuStyle
849 203 equemene
    InputCL['VariableType']=VariableType
850 203 equemene
    InputCL['Size']=Size
851 203 equemene
    InputCL['Mass']=Mass
852 203 equemene
    InputCL['InternalRadius']=InternalRadius
853 203 equemene
    InputCL['ExternalRadius']=ExternalRadius
854 203 equemene
    InputCL['ObserverDistance']=ObserverDistance
855 203 equemene
    InputCL['Angle']=Angle
856 203 equemene
    InputCL['BlackBody']=BlackBody
857 203 equemene
    InputCL['AngularCamera']=AngularCamera
858 203 equemene
    InputCL['Method']=Method
859 203 equemene
    InputCL['TrackPoints']=TrackPoints
860 203 equemene
861 203 equemene
    duration=BlackHoleCL(zImage,fImage,InputCL)
862 203 equemene
863 203 equemene
    ImageOutput(zImage,"TrouNoirZ_%s" % Method)
864 203 equemene
    ImageOutput(fImage,"TrouNoirF_%s" % Method)