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

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

1./ grosses modifs de structure pour cmvdefsurv.cc
2./ PoissRandLimit enleve de geneutils.{h,cc} et mis dans srandgen.{h,cc}

cmv 29/10/2007

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