source: Sophya/trunk/Cosmo/SimLSS/geneutils.cc@ 3325

Last change on this file since 3325 was 3325, checked in by cmv, 18 years ago

Mise en conformite / SOPHYA lib:

  1. on enleve #include "sopnamsp.h" dans les .cc de la librairie
  2. on encadre par "namespace SOPHYA { ... }" tout le code des .cc de la librairie y compris les fonctions
  3. on met les fcts des .h dans le "namespace SOPHYA { ... }"
  4. on met #include "sopnamsp.h" dans tous les cmv*.cc cad les main programs

cmv le mauvais eleve (sur les conseils de Reza) 13/09/2007

File size: 14.1 KB
Line 
1#include "machdefs.h"
2#include <iostream>
3#include <stdlib.h>
4#include <stdio.h>
5#include <string.h>
6#include <math.h>
7#include <vector>
8
9#include "pexceptions.h"
10
11#include "histos.h"
12#include "hisprof.h"
13#include "srandgen.h"
14
15#include "geneutils.h"
16
17namespace SOPHYA {
18
19//-------------------------------------------------------------------
20// Classe d'interpolation lineaire:
21// Le vecteur y a n elements y_i tels que y_i = f(x_i)
22// ou les x_i sont regulierement espaces
23// et x_0=xmin et x_{n-1}=xmax (xmax inclus!)
24InterpFunc::InterpFunc(double xmin,double xmax,vector<double>& y)
25 : _xmin(xmin), _xmax(xmax), _y(y)
26{
27 if(_xmin>=_xmax || _y.size()<=2) { // au moins 3 points!
28 cout<<"InterpFunc::InterpFunc : bad arguments values"<<endl;
29 throw ParmError("InterpFunc::InterpFunc : bad arguments values");
30 }
31 _nm1 = _y.size()-1;
32 _dx = (_xmax-_xmin)/(double)_nm1;
33}
34
35double InterpFunc::Linear(double x,unsigned short& ok)
36{
37 ok=0; if(x<_xmin) ok=1; else if(x>_xmax) ok=2;
38 x -= _xmin;
39 long i = long(x/_dx); // On prend le "i" juste en dessous
40 if(i<0) i=0; else if(i>=_nm1) i=_nm1-1;
41 return _y[i] + (_y[i+1]-_y[i])/_dx*(x-i*_dx);
42}
43
44double InterpFunc::Parab(double x,unsigned short& ok)
45{
46 ok=0; if(x<_xmin) ok=1; else if(x>_xmax) ok=2;
47 x -= _xmin;
48 long i = long(x/_dx+0.5); // On prend le "i" le + proche
49 if(i<1) i=1; else if(i>=_nm1-1) i=_nm1-2;
50 double a = (_y[i+1]-2.*_y[i]+_y[i-1])/(2.*_dx*_dx);
51 double b = (_y[i+1]-_y[i-1])/(2.*_dx);
52 return _y[i] + (x-i*_dx)*(a*(x-i*_dx)+b);
53}
54
55//-------------------------------------------------------------------
56// Classe d'inversion d'une fonction STRICTEMENT MONOTONE CROISSANTE
57// - Le vecteur y a "Nin" elements y_i tels que "y_i = f(x_i)"
58// - On a x(i) < x(i+1) et y(i) < y(i+1)
59// - La classe renvoie ymin=y(0) , ymax=y(Nin -1)
60// et le vecteur x = f^-1(y) de "Nout" elements
61// Les y_i sont regulierement espaces et ymin et ymax
62// La re-interpolation inverse est faite par lineaire
63InverseFunc::InverseFunc(vector<double>& x,vector<double>& y)
64 : _ymin(0.) , _ymax(0.) , _x(x) , _y(y)
65{
66 int_4 ns = _x.size();
67 if(ns<3 || _y.size()<=0 || ns!=_y.size())
68 throw ParmError("InverseFunc::InverseFunc_Error: bad array size");
69
70 // Check "x" strictement monotone croissant
71 for(int_4 i=0;i<ns-1;i++)
72 if(_x[i+1]<=_x[i]) {
73 cout<<"InverseFunc::InverseFunc_Error: _x array not stricly growing"<<endl;
74 throw ParmError("InverseFunc::InverseFunc_Error: _x array not stricly growing");
75 }
76
77 // Check "y" monotone croissant
78 for(int_4 i=0;i<ns-1;i++)
79 if(_y[i+1]<_y[i]) {
80 cout<<"InverseFunc::InverseFunc_Error: _y array not growing"<<endl;
81 throw ParmError("InverseFunc::InverseFunc_Error: _y array not growing");
82 }
83
84 // define limits
85 _ymin = _y[0];
86 _ymax = _y[ns-1];
87
88}
89
90InverseFunc::~InverseFunc(void)
91{
92}
93
94int InverseFunc::ComputeLinear(long n,vector<double>& xfcty)
95// Compute table "xfcty" by linear interpolation of "x" versus "y"
96// on "n" points from "ymin" to "ymax":
97// xfcty[i] = interpolation of function "x" for "ymin+i*(ymax-ymin)/(n-1.)"
98{
99 if(n<3) return -1;
100
101 xfcty.resize(n);
102
103 long i1,i2;
104 double x;
105 for(int_4 i=0;i<n;i++) {
106 double y = _ymin + i*(_ymax-_ymin)/(n-1.);
107 find_in_y(y,i1,i2);
108 double dy = _y[i2]-_y[i1];
109 if(dy==0.) {
110 x = (_x[i2]+_x[i1])/2.; // la fct a inverser est plate!
111 } else {
112 x = _x[i1] + (_x[i2]-_x[i1])/dy * (y-_y[i1]);
113 }
114 xfcty[i] = x;
115 }
116
117 return 0;
118}
119
120int InverseFunc::ComputeParab(long n,vector<double>& xfcty)
121{
122 if(n<3) return -1;
123
124 xfcty.resize(n);
125
126 long i1,i2,i3;
127 double x;
128 for(int_4 i=0;i<n;i++) {
129 double y = _ymin + i*(_ymax-_ymin)/(n-1.);
130 find_in_y(y,i1,i2);
131 // On cherche le 3ieme point selon la position de y / au 2 premiers
132 double my = (_y[i1]+_y[i2])/2.;
133 if(y<my) {i3=i2; i2=i1; i1--;} else {i3=i2+1;}
134 // Protection
135 if(i1<0) {i1++; i2++; i3++;}
136 if(i3==_y.size()) {i1--; i2--; i3--;}
137 // Interpolation parabolique
138 double dy = _y[i3]-_y[i1];
139 if(dy==0.) {
140 x = (_x[i3]+_x[i1])/2.; // la fct a inverser est plate!
141 } else {
142 double X1=_x[i1]-_x[i2], X3=_x[i3]-_x[i2];
143 double Y1=_y[i1]-_y[i2], Y3=_y[i3]-_y[i2];
144 double den = Y1*Y3*dy;
145 double a = (X3*Y1-X1*Y3)/den;
146 double b = (X1*Y3*Y3-X3*Y1*Y1)/den;
147 y -= _y[i2];
148 x = (a*y+b)*y + _x[i2];
149 }
150 xfcty[i] = x;
151 }
152
153 return 0;
154}
155
156//----------------------------------------------------
157double InterpTab(double x0,vector<double>& X,vector<double>& Y,unsigned short typint)
158// Interpole in x0 the table Y = f(X)
159// X doit etre ordonne par ordre croissant (strictement)
160// typint = 0 : nearest value
161// 1 : linear interpolation
162// 2 : parabolique interpolation
163{
164 long n = X.size();
165 if(Y.size()<n || n<2)
166 throw ParmError("InterpTab_Error : incompatible size between X and Y tables!");
167
168 if(x0<X[0] || x0>X[n-1]) return 0.;
169 if(typint>2) typint = 0;
170
171 // Recherche des indices encadrants par dichotomie
172 long k, klo=0, khi=n-1;
173 while (khi-klo > 1) {
174 k = (khi+klo) >> 1;
175 if (X[k] > x0) khi=k; else klo=k;
176 }
177
178 // Quel est le plus proche?
179 k = (x0-X[klo]<X[khi]-x0) ? klo: khi;
180
181 // On retourne le plus proche
182 if(typint==0) return Y[k];
183
184 // On retourne l'extrapolation lineaire
185 if(typint==1 || n<3)
186 return Y[klo] + (Y[khi]-Y[klo])/(X[khi]-X[klo])*(x0-X[klo]);
187
188 // On retourne l'extrapolation parabolique
189 if(k==0) k++; else if(k==n-1) k--;
190 klo = k-1; khi = k+1;
191 double x1 = X[klo]-X[k], x2 = X[khi]-X[k];
192 double y1 = Y[klo]-Y[k], y2 = Y[khi]-Y[k];
193 double den = x1*x2*(x1-x2);
194 double a = (y1*x2-y2*x1)/den;
195 double b = (y2*x1*x1-y1*x2*x2)/den;
196 x0 -= X[k];
197 return Y[k] + (a*x0+b)*x0;;
198
199}
200
201//-------------------------------------------------------------------
202int FuncToHisto(GenericFunc& func,Histo& h,bool logaxex)
203// Remplit l'histo 1D "h" avec la fonction "func"
204// INPUT:
205// logaxex = false : remplissage lineaire
206// les abscisses "x" des bins sont remplis avec f(x)
207// logaxex = true : remplissage logarithmique (base 10)
208// les abscisses "x" des bins sont remplis avec f(10^x)
209// RETURN:
210// 0 = OK
211// 1 = error
212{
213 if(h.NBins()<=0) return 1;
214
215 h.Zero();
216
217 for(int_4 i=0;i<h.NBins();i++) {
218 double x = h.BinCenter(i);
219 if(logaxex) x = pow(10.,x);
220 h.SetBin(i,func(x));
221 }
222
223 return 0;
224}
225
226int FuncToVec(GenericFunc& func,TVector<r_8>& v,double xmin,double xmax,bool logaxex)
227// Remplit le TVector avec la fonction "func"
228// INPUT:
229// logaxex = false : remplissage lineaire
230// les abscisses "x" des bins sont remplis avec f(x)
231// logaxex = true : remplissage logarithmique (base 10)
232// les abscisses "x" des bins sont remplis avec f(10^x)
233// RETURN:
234// 0 = OK
235// 1 = error
236// Remarque:
237// v(i) = f(xmin+i*dx) avec dx = (xmax-xmin)/v.NElts()
238{
239 if(v.NElts()<=0 || xmax<=xmin) return 1;
240
241 v = 0.;
242 double dx = (xmax-xmin)/v.NElts();
243
244 for(int_4 i=0;i<v.NElts();i++) {
245 double x = xmin + i * dx;;
246 if(logaxex) x = pow(10.,x);
247 v(i) = func(x);
248 }
249
250 return 0;
251}
252
253//-------------------------------------------------------------------
254double AngSol(double dtheta,double dphi,double theta0)
255// Retourne l'angle solide d'un "rectangle" et coordonnees spheriques
256// de DEMI-COTE "dtheta" x "dphi" et centre en "theta0"
257// Attention: Le "theta0" de l'equateur est Pi/2 (et non pas zero)
258// Les unites des angles sont en radians
259// theta0 in [0,Pi]
260// dtheta in [0,Pi]
261// dphi in [0,2Pi]
262// Return: l'angle solide en steradian
263{
264 double theta1 = theta0-dtheta, theta2 = theta0+dtheta;
265 if(theta1<0.) theta1=0.;
266 if(theta2>M_PI) theta2=M_PI;
267
268 return 2.*dphi * (cos(theta1)-cos(theta2));
269}
270
271double AngSol(double dtheta)
272// Retourne l'angle solide d'une calotte spherique de demi-ouverture "dtheta"
273// Attention: Les unites des angles sont en radians
274// dtheta in [0,Pi]
275// Return: l'angle solide en steradian
276// Approx pour theta petit: PI theta^2
277{
278 return 2.*M_PI * (1.-cos(dtheta));
279}
280
281//-------------------------------------------------------------------
282unsigned long PoissRandLimit(double mu,double mumax)
283{
284 double pp,ppi,x;
285 unsigned long n;
286
287 if(mu>=mumax) {
288 pp = sqrt(mu);
289 while( (x=pp*NorRand()) < -mu );
290 return (unsigned long)(mu+x+0.5);
291 }
292
293 ppi = pp = exp(-mu);
294 x = drand01();
295 n = 0;
296 while (x > ppi) {
297 n++;
298 pp = mu*pp/(double)n;
299 ppi += pp;
300 }
301 return n;
302}
303
304//-------------------------------------------------------------------
305
306static unsigned short IntegrateFunc_GlOrder = 0;
307static vector<double> IntegrateFunc_x;
308static vector<double> IntegrateFunc_w;
309
310///////////////////////////////////////////////////////////
311//************** Integration of Functions ***************//
312///////////////////////////////////////////////////////////
313
314
315////////////////////////////////////////////////////////////////////////////////////
316double IntegrateFunc(GenericFunc& func,double xmin,double xmax
317 ,double perc,double dxinc,double dxmax,unsigned short glorder)
318// --- Integration adaptative ---
319// Integrate[func[x], {x,xmin,xmax}]
320// ..xmin,xmax are the integration limits
321// ..dxinc is the searching increment x = xmin+i*dxinc
322// if <0 npt = int(|dxinc|) (si<2 alors npt=100)
323// et dxinc = (xmax-xmin)/npt
324// ..dxmax is the maximum possible increment (if dxmax<=0 no test)
325// ---
326// Partition de [xmin,xmax] en intervalles [x(i),x(i+1)]:
327// on parcourt [xmin,xmax] par pas de "dxinc" : x(i) = xmin + i*dxinc
328// on cree un intervalle [x(i),x(i+1)]
329// - si |f[x(i+1)] - f[x(i)]| > perc*|f[[x(i)]|
330// - si |x(i+1)-x(i)| >= dxmax (si dxmax>0.)
331// Dans un intervalle [x(i),x(i+1)] la fonction est integree
332// avec une methode de gauss-legendre d'ordre "glorder"
333{
334 double signe = 1.;
335 if(xmin>xmax) {double tmp=xmax; xmax=xmin; xmin=tmp; signe=-1.;}
336
337 if(glorder==0) glorder = 4;
338 if(perc<=0.) perc=0.1;
339 if(dxinc<=0.) {int n=int(-dxinc); if(n<2) n=100; dxinc=(xmax-xmin)/n;}
340 if(glorder != IntegrateFunc_GlOrder) {
341 IntegrateFunc_GlOrder = glorder;
342 Compute_GaussLeg(glorder,IntegrateFunc_x,IntegrateFunc_w,0.,1.);
343 }
344
345 // Recherche des intervalles [x(i),x(i+1)]
346 int_4 ninter = 0;
347 double sum = 0., xbas=xmin, fbas=func(xbas), fabsfbas=fabs(fbas);
348 for(double x=xmin+dxinc; x<xmax+dxinc/2.; x += dxinc) {
349 double f = func(x);
350 double dx = x-xbas;
351 bool doit = false;
352 if( x>xmax ) {doit = true; x=xmax;}
353 else if( dxmax>0. && dx>dxmax ) doit = true;
354 else if( fabs(f-fbas)>perc*fabsfbas ) doit = true;
355 if( ! doit ) continue;
356 double s = 0.;
357 for(unsigned short i=0;i<IntegrateFunc_GlOrder;i++)
358 s += IntegrateFunc_w[i]*func(xbas+IntegrateFunc_x[i]*dx);
359 sum += s*dx;
360 xbas = x; fbas = f; fabsfbas=fabs(fbas); ninter++;
361 }
362 //cout<<"Ninter="<<ninter<<endl;
363
364 return sum*signe;
365}
366
367////////////////////////////////////////////////////////////////////////////////////
368double IntegrateFuncLog(GenericFunc& func,double lxmin,double lxmax
369 ,double perc,double dlxinc,double dlxmax,unsigned short glorder)
370// --- Integration adaptative ---
371// Idem IntegrateFunc but integrate on logarithmic (base 10) intervals:
372// Int[ f(x) dx ] = Int[ x*f(x) dlog10(x) ] * log(10)
373// ..lxmin,lxmax are the log10() of the integration limits
374// ..dlxinc is the searching logarithmic (base 10) increment lx = lxmin+i*dlxinc
375// ..dlxmax is the maximum possible logarithmic (base 10) increment (if dlxmax<=0 no test)
376// Remarque: to be use if "x*f(x) versus log10(x)" looks like a polynomial
377// better than "f(x) versus x"
378// ATTENTION: la fonction func que l'on passe en argument
379// est "func(x)" et non pas "func(log10(x))"
380{
381 double signe = 1.;
382 if(lxmin>lxmax) {double tmp=lxmax; lxmax=lxmin; lxmin=tmp; signe=-1.;}
383
384 if(glorder==0) glorder = 4;
385 if(perc<=0.) perc=0.1;
386 if(dlxinc<=0.) {int n=int(-dlxinc); if(n<2) n=100; dlxinc=(lxmax-lxmin)/n;}
387 if(glorder != IntegrateFunc_GlOrder) {
388 IntegrateFunc_GlOrder = glorder;
389 Compute_GaussLeg(glorder,IntegrateFunc_x,IntegrateFunc_w,0.,1.);
390 }
391
392 // Recherche des intervalles [lx(i),lx(i+1)]
393 int_4 ninter = 0;
394 double sum = 0., lxbas=lxmin, fbas=func(pow(10.,lxbas)), fabsfbas=fabs(fbas);
395 for(double lx=lxmin+dlxinc; lx<lxmax+dlxinc/2.; lx += dlxinc) {
396 double f = func(pow(10.,lx));
397 double dlx = lx-lxbas;
398 bool doit = false;
399 if( lx>lxmax ) {doit = true; lx=lxmax;}
400 else if( dlxmax>0. && dlx>dlxmax ) doit = true;
401 else if( fabs(f-fbas)>perc*fabsfbas ) doit = true;
402 if( ! doit ) continue;
403 double s = 0.;
404 for(unsigned short i=0;i<IntegrateFunc_GlOrder;i++) {
405 double y = pow(10.,lxbas+IntegrateFunc_x[i]*dlx);
406 s += IntegrateFunc_w[i]*y*func(y);
407 }
408 sum += s*dlx;
409 lxbas = lx; fbas = f; fabsfbas=fabs(fbas); ninter++;
410 }
411 //cout<<"Ninter="<<ninter<<endl;
412
413 return M_LN10*sum*signe;
414}
415
416////////////////////////////////////////////////////////////////////////////////////
417/*
418Integration des fonctions a une dimension y=f(x) par la Methode de Gauss-Legendre.
419 --> Calcul des coefficients du developpement pour [x1,x2]
420| INPUT:
421| x1,x2 : bornes de l'intervalle (dans nbinteg.h -> x1=-0.5 x2=0.5)
422| glorder = degre n du developpement de Gauss-Legendre
423| OUTPUT:
424| x[] = valeur des abscisses ou l'on calcule (dim=n)
425| w[] = valeur des coefficients associes
426| REMARQUES:
427| - x et w seront dimensionnes a n.
428| - l'integration est rigoureuse si sur l'intervalle d'integration
429| la fonction f(x) peut etre approximee par un polynome
430| de degre 2*m (monome le + haut x**(2*n-1) )
431*/
432void Compute_GaussLeg(unsigned short glorder,vector<double>& x,vector<double>& w,double x1,double x2)
433{
434 if(glorder==0) return;
435 int n = (int)glorder;
436 x.resize(n,0.); w.resize(n,0.);
437
438 int m,j,i;
439 double z1,z,xm,xl,pp,p3,p2,p1;
440
441 m=(n+1)/2;
442 xm=0.5*(x2+x1);
443 xl=0.5*(x2-x1);
444 for (i=1;i<=m;i++) {
445 z=cos(M_PI*(i-0.25)/(n+0.5));
446 do {
447 p1=1.0;
448 p2=0.0;
449 for (j=1;j<=n;j++) {
450 p3=p2;
451 p2=p1;
452 p1=((2.0*j-1.0)*z*p2-(j-1.0)*p3)/j;
453 }
454 pp=n*(z*p1-p2)/(z*z-1.0);
455 z1=z;
456 z=z1-p1/pp;
457 } while (fabs(z-z1) > 3.0e-11); // epsilon
458 x[i-1]=xm-xl*z;
459 x[n-i]=xm+xl*z;
460 w[i-1]=2.0*xl/((1.0-z*z)*pp*pp);
461 w[n-i]=w[i-1];
462 }
463}
464
465} // Fin namespace SOPHYA
Note: See TracBrowser for help on using the repository browser.