| [2615] | 1 | #include "sopnamsp.h"
 | 
|---|
| [220] | 2 | #include "integ.h"
 | 
|---|
 | 3 | #include "generalfit.h"
 | 
|---|
 | 4 | 
 | 
|---|
 | 5 | // A faire :
 | 
|---|
 | 6 | //   1 dans GLInteg, optimiser recalcul en utilisant des poids
 | 
|---|
 | 7 | //     calcules sur [0,1] et chgt de variable
 | 
|---|
 | 8 | //   2 dans GLInteg, distinction nStep = subdivisions de l'intervalle,
 | 
|---|
 | 9 | //     et on applique GL d'ordre "ordre" sur chaque subdivision, en
 | 
|---|
 | 10 | //     utilisant le (1)
 | 
|---|
 | 11 | 
 | 
|---|
| [2808] | 12 | /*!
 | 
|---|
 | 13 |   \ingroup NTools
 | 
|---|
 | 14 |   \class SOPHYA::Integrator
 | 
|---|
| [220] | 15 | 
 | 
|---|
| [2808] | 16 |   \brief Classe abstraite d'intégration numérique 1D.
 | 
|---|
| [220] | 17 | 
 | 
|---|
| [2808] | 18 |   On fournit une fonction double f(double) au constructeur, ou
 | 
|---|
| [3098] | 19 |   une classe-fonction (ClassFuc) ou 
 | 
|---|
| [2808] | 20 |   une GeneralFunction avec des paramètres définis.
 | 
|---|
 | 21 |   L'objet Integrator est convertible en valeur double qui est la valeur
 | 
|---|
 | 22 |   de l'intégrale. Diverses méthodes permettent de choisir des options
 | 
|---|
 | 23 |   de calcul, et ces méthodes retournent une référence sur l'objet, pour
 | 
|---|
 | 24 |   permettre une notation chaînée.
 | 
|---|
 | 25 | 
 | 
|---|
 | 26 |   \sa TrpzInteg GLInteg
 | 
|---|
 | 27 | */
 | 
|---|
 | 28 | 
 | 
|---|
 | 29 | //!  Constructeur par défaut. L'objet n'est pas utilisable en l'état.
 | 
|---|
| [220] | 30 | Integrator::Integrator()
 | 
|---|
| [3098] | 31 | : mFunc(NULL), mClFun(NULL), mGFF(NULL), mGFFParm(NULL),
 | 
|---|
| [220] | 32 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 33 |   mXMin(0), mXMax(1)
 | 
|---|
 | 34 | {}
 | 
|---|
 | 35 | 
 | 
|---|
| [3098] | 36 | //! Constructeur à partir d'une classe-fonction, et des bornes d'intégration.
 | 
|---|
| [3409] | 37 | Integrator::Integrator(ClassFunc & f, double xmin, double xmax)
 | 
|---|
| [3098] | 38 | : mFunc(NULL), mClFun(&f), mGFF(NULL), mGFFParm(NULL),
 | 
|---|
| [220] | 39 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 40 |   mXMin(xmin), mXMax(xmax)
 | 
|---|
 | 41 | {}
 | 
|---|
 | 42 | 
 | 
|---|
| [2808] | 43 | //! Constructeur à partir de la fonction double->double, et des bornes d'intégration.
 | 
|---|
| [3098] | 44 | Integrator::Integrator(FUNC f, double xmin, double xmax)
 | 
|---|
 | 45 | : mFunc(f), mClFun(NULL), mGFF(NULL), mGFFParm(NULL),
 | 
|---|
| [2283] | 46 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 47 |   mXMin(xmin), mXMax(xmax)
 | 
|---|
 | 48 | {}
 | 
|---|
 | 49 | 
 | 
|---|
 | 50 | 
 | 
|---|
| [2808] | 51 | /*!
 | 
|---|
| [3098] | 52 |   Constructeur a partir d'une classe-fonction 
 | 
|---|
 | 53 |   sans spécifier les bornes. Elles sont positionnées
 | 
|---|
| [2808] | 54 |   à [0,1], et on pourra les modifier plus tard.
 | 
|---|
 | 55 | */
 | 
|---|
| [3409] | 56 | Integrator::Integrator(ClassFunc & f)
 | 
|---|
| [3098] | 57 | : mFunc(NULL), mClFun(&f), mGFF(NULL), mGFFParm(NULL),
 | 
|---|
| [220] | 58 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 59 |   mXMin(0), mXMax(1)
 | 
|---|
 | 60 | {}
 | 
|---|
 | 61 | 
 | 
|---|
| [2808] | 62 | /*!
 | 
|---|
 | 63 |   Constructeur sans spécifier les bornes. Elles sont positionnées
 | 
|---|
 | 64 |   à [0,1], et on pourra les modifier plus tard.
 | 
|---|
 | 65 | */
 | 
|---|
| [3098] | 66 | Integrator::Integrator(FUNC f)
 | 
|---|
 | 67 | : mFunc(f), mClFun(NULL), mGFF(NULL), mGFFParm(NULL),
 | 
|---|
| [2283] | 68 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 69 |   mXMin(0), mXMax(1)
 | 
|---|
 | 70 | {}
 | 
|---|
 | 71 | 
 | 
|---|
| [2808] | 72 | /*!
 | 
|---|
 | 73 |   Constructeur à partir d'une GeneralFunction. La fonction doit être une
 | 
|---|
 | 74 |   fonction de une variable, et on fournit les valeurs des paramètres ainsi
 | 
|---|
 | 75 |   que les bornes d'intégration.
 | 
|---|
 | 76 | */
 | 
|---|
| [220] | 77 | Integrator::Integrator(GeneralFunction* gff, double* par, double xmin, double xmax)
 | 
|---|
| [3098] | 78 | : mFunc(NULL), mClFun(NULL), mGFF(gff), mGFFParm(par),
 | 
|---|
| [220] | 79 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 80 |   mXMin(xmin), mXMax(xmax)
 | 
|---|
 | 81 | {
 | 
|---|
 | 82 |   DBASSERT(gff->NVar() == 1);
 | 
|---|
 | 83 | }
 | 
|---|
 | 84 | 
 | 
|---|
| [2808] | 85 | /*!
 | 
|---|
 | 86 |   Constructeur à partir d'une GeneralFunction. La fonction doit être une
 | 
|---|
 | 87 |   fonction de une variable, et on fournit les valeurs des paramètres. 
 | 
|---|
 | 88 |   On ne spécifie pas les bornes. Elles sont positionnées
 | 
|---|
 | 89 |   à [0,1], et on pourra les modifier plus tard.
 | 
|---|
 | 90 | */
 | 
|---|
| [220] | 91 | Integrator::Integrator(GeneralFunction* gff, double* par)
 | 
|---|
| [3098] | 92 | : mFunc(NULL), mClFun(NULL), mGFF(gff), mGFFParm(par),
 | 
|---|
| [220] | 93 |   mNStep(50), mDX(-1), mReqPrec(-1),
 | 
|---|
 | 94 |   mXMin(0), mXMax(1)
 | 
|---|
 | 95 | {
 | 
|---|
 | 96 |   DBASSERT(gff->NVar() == 1);
 | 
|---|
 | 97 | }
 | 
|---|
 | 98 | 
 | 
|---|
 | 99 | Integrator::~Integrator()
 | 
|---|
| [3098] | 100 | {
 | 
|---|
 | 101 | }
 | 
|---|
| [220] | 102 | 
 | 
|---|
| [3098] | 103 | 
 | 
|---|
| [2808] | 104 | /*! 
 | 
|---|
 | 105 |   \brief Spécifie le nombre de pas pour l'intégration numérique.
 | 
|---|
 | 106 |   
 | 
|---|
 | 107 |   La signification peut dépendre de la méthode d'intégration.
 | 
|---|
 | 108 | */
 | 
|---|
| [220] | 109 | Integrator& 
 | 
|---|
 | 110 | Integrator::NStep(int n)
 | 
|---|
 | 111 | {
 | 
|---|
 | 112 |   mNStep = n;
 | 
|---|
 | 113 |   mDX = mReqPrec = -1;
 | 
|---|
 | 114 |   StepsChanged();
 | 
|---|
 | 115 |   return *this;
 | 
|---|
 | 116 | }
 | 
|---|
 | 117 | 
 | 
|---|
| [2808] | 118 | /*! 
 | 
|---|
 | 119 |   \brief Spécifie le nombre de pas pour l'intégration numérique.
 | 
|---|
 | 120 |   
 | 
|---|
 | 121 |   La signification peut dépendre de la méthode d'intégration.
 | 
|---|
 | 122 | */
 | 
|---|
| [220] | 123 | Integrator& 
 | 
|---|
 | 124 | Integrator::DX(double d)
 | 
|---|
 | 125 | {
 | 
|---|
 | 126 |   mDX = d;
 | 
|---|
| [3098] | 127 |   mNStep = -1;
 | 
|---|
 | 128 |   mReqPrec = -1.;
 | 
|---|
| [220] | 129 |   StepsChanged();
 | 
|---|
 | 130 |   return *this;
 | 
|---|
 | 131 | }
 | 
|---|
 | 132 | 
 | 
|---|
| [2808] | 133 | /*!
 | 
|---|
 | 134 |   \brief Spécifie la précision souhaitée. 
 | 
|---|
 | 135 | 
 | 
|---|
 | 136 |   La signification peut dépendre de la méthode d'intégration.
 | 
|---|
 | 137 |   Non disponible dans toutes les méthodes d'intégration.
 | 
|---|
 | 138 | */
 | 
|---|
| [220] | 139 | Integrator& 
 | 
|---|
 | 140 | Integrator::ReqPrec(double p)
 | 
|---|
 | 141 | {
 | 
|---|
 | 142 |   DBASSERT( !"Pas encore implemente !");
 | 
|---|
 | 143 |   mReqPrec = p;
 | 
|---|
 | 144 |   mDX = mNStep = -1;
 | 
|---|
 | 145 |   StepsChanged();
 | 
|---|
 | 146 |   return *this;
 | 
|---|
 | 147 | }
 | 
|---|
 | 148 | 
 | 
|---|
| [2808] | 149 | //! Spécifie les bornes de l'intégrale.
 | 
|---|
| [220] | 150 | Integrator& 
 | 
|---|
 | 151 | Integrator::Limits(double xmin, double xmax)
 | 
|---|
 | 152 | {
 | 
|---|
 | 153 |   mXMin = xmin;
 | 
|---|
 | 154 |   mXMax = xmax;
 | 
|---|
 | 155 |   LimitsChanged();
 | 
|---|
 | 156 |   return *this;
 | 
|---|
 | 157 | }
 | 
|---|
 | 158 | 
 | 
|---|
| [2808] | 159 | //! Spécifie la fonction à intégrer, sous forme double f(double).
 | 
|---|
| [220] | 160 | Integrator& 
 | 
|---|
| [3409] | 161 | Integrator::SetFunc(ClassFunc & f)
 | 
|---|
| [220] | 162 | {
 | 
|---|
| [3098] | 163 |   mFunc = NULL;
 | 
|---|
 | 164 |   mClFun = &f;
 | 
|---|
 | 165 |   mGFFParm = NULL;
 | 
|---|
 | 166 |   mFunc = NULL;
 | 
|---|
 | 167 |   FuncChanged();
 | 
|---|
 | 168 |   return *this;
 | 
|---|
 | 169 | }
 | 
|---|
 | 170 | 
 | 
|---|
 | 171 | //! Spécifie la fonction à intégrer, sous forme double f(double).
 | 
|---|
 | 172 | Integrator& 
 | 
|---|
 | 173 | Integrator::SetFunc(FUNC f)
 | 
|---|
 | 174 | {
 | 
|---|
| [220] | 175 |   mFunc = f;
 | 
|---|
| [3098] | 176 |   mClFun = NULL;
 | 
|---|
| [220] | 177 |   mGFFParm = NULL;
 | 
|---|
 | 178 |   mFunc = NULL;
 | 
|---|
 | 179 |   FuncChanged();
 | 
|---|
 | 180 |   return *this;
 | 
|---|
 | 181 | }
 | 
|---|
 | 182 | 
 | 
|---|
| [2808] | 183 | /*!
 | 
|---|
 | 184 |   \brief Spécifie la fonction à intégrer, sous forme de GeneralFunction 
 | 
|---|
 | 185 |   à une variable, et les paramètres sont fournis.
 | 
|---|
 | 186 | */
 | 
|---|
| [220] | 187 | Integrator& 
 | 
|---|
| [3098] | 188 | Integrator::SetFunc(GeneralFunction* gff, double* par)
 | 
|---|
| [220] | 189 | {
 | 
|---|
 | 190 |   mGFF = gff;
 | 
|---|
 | 191 |   mGFFParm = par;
 | 
|---|
 | 192 |   mFunc = NULL;
 | 
|---|
 | 193 |   DBASSERT(gff->NVar() == 1);
 | 
|---|
 | 194 |   FuncChanged();
 | 
|---|
 | 195 |   return *this;
 | 
|---|
 | 196 | }
 | 
|---|
 | 197 | 
 | 
|---|
 | 198 | double
 | 
|---|
 | 199 | Integrator::FVal(double x) const
 | 
|---|
 | 200 | {
 | 
|---|
| [3098] | 201 |   DBASSERT( mFunc || mClFun || (mGFF && mGFFParm) );
 | 
|---|
 | 202 |   if (mFunc) return mFunc(x);
 | 
|---|
 | 203 |   else return mClFun ? (*mClFun)(x) : mGFF->Value(&x, mGFFParm); 
 | 
|---|
 | 204 | }
 | 
|---|
| [220] | 205 | 
 | 
|---|
 | 206 | double
 | 
|---|
 | 207 | Integrator::ValueBetween(double xmin, double xmax)
 | 
|---|
 | 208 | {
 | 
|---|
 | 209 |   Limits(xmin,xmax);
 | 
|---|
 | 210 |   return Value();
 | 
|---|
 | 211 | }
 | 
|---|
 | 212 | 
 | 
|---|
| [2423] | 213 | void
 | 
|---|
 | 214 | Integrator::Print(int lp)
 | 
|---|
 | 215 | {
 | 
|---|
 | 216 |  cout<<"Integrator between "<<mXMin<<" and "<<mXMax
 | 
|---|
 | 217 |      <<" in "<<mNStep<<" steps"<<endl;
 | 
|---|
| [3098] | 218 |  if(lp>1) cout<<"mFunc="<<mFunc<<"mClFun="<<mClFun<<"mGFF="<<mGFF
 | 
|---|
| [2423] | 219 |               <<"mReqPrec="<<mReqPrec<<"mDX="<<mDX<<endl;
 | 
|---|
 | 220 | }
 | 
|---|
 | 221 | 
 | 
|---|
| [2808] | 222 | /*!
 | 
|---|
 | 223 |   \ingroup NTools
 | 
|---|
 | 224 |   \class SOPHYA::TrpzInteg
 | 
|---|
| [220] | 225 | 
 | 
|---|
| [2808] | 226 |   \brief Implementation de Integrator par la methode des trapezes.
 | 
|---|
 | 227 |   
 | 
|---|
 | 228 |   Classe d'intégration par la méthode des trapèzes.
 | 
|---|
 | 229 |   Voir Integrator pour les méthodes. Le nombre de pas
 | 
|---|
 | 230 |   est le nombre de trapèze, le pas d'intégration est
 | 
|---|
 | 231 |   la largeur des trapèzez. Impossible de demander une
 | 
|---|
 | 232 |   précision.
 | 
|---|
 | 233 |   
 | 
|---|
 | 234 |   \sa SOPHYA::Integrator
 | 
|---|
 | 235 | */
 | 
|---|
| [220] | 236 | 
 | 
|---|
| [2808] | 237 | 
 | 
|---|
| [220] | 238 | TrpzInteg::TrpzInteg()
 | 
|---|
 | 239 | {}
 | 
|---|
 | 240 | 
 | 
|---|
| [3409] | 241 | TrpzInteg::TrpzInteg(ClassFunc & f, double xmin, double xmax)
 | 
|---|
| [220] | 242 | : Integrator(f, xmin, xmax)
 | 
|---|
 | 243 | {}
 | 
|---|
| [2808] | 244 | 
 | 
|---|
| [3098] | 245 | TrpzInteg::TrpzInteg(FUNC f, double xmin, double xmax)
 | 
|---|
| [2283] | 246 | : Integrator(f, xmin, xmax)
 | 
|---|
 | 247 | {}
 | 
|---|
| [220] | 248 | 
 | 
|---|
| [3409] | 249 | TrpzInteg::TrpzInteg(ClassFunc & f)
 | 
|---|
| [220] | 250 | : Integrator(f)
 | 
|---|
 | 251 | {}
 | 
|---|
 | 252 | 
 | 
|---|
| [3098] | 253 | TrpzInteg::TrpzInteg(FUNC f)
 | 
|---|
| [2283] | 254 | : Integrator(f)
 | 
|---|
 | 255 | {}
 | 
|---|
 | 256 | 
 | 
|---|
| [220] | 257 | TrpzInteg::TrpzInteg(GeneralFunction* gff, double* par, double xmin, double xmax)
 | 
|---|
 | 258 | : Integrator(gff, par, xmin, xmax)
 | 
|---|
 | 259 | {}
 | 
|---|
 | 260 | 
 | 
|---|
 | 261 | TrpzInteg::TrpzInteg(GeneralFunction* gff, double* par)
 | 
|---|
 | 262 | : Integrator(gff, par)
 | 
|---|
 | 263 | {}
 | 
|---|
 | 264 | 
 | 
|---|
 | 265 | TrpzInteg::~TrpzInteg()
 | 
|---|
 | 266 | {}
 | 
|---|
 | 267 | 
 | 
|---|
| [2808] | 268 | //! Return the value of the integral.
 | 
|---|
| [220] | 269 | double
 | 
|---|
 | 270 | TrpzInteg::Value()
 | 
|---|
 | 271 | {
 | 
|---|
 | 272 |   double dx = mDX;
 | 
|---|
 | 273 |   double nstep = mNStep;
 | 
|---|
 | 274 | 
 | 
|---|
 | 275 |   if (dx > 0)
 | 
|---|
 | 276 |     nstep = (mXMax - mXMin)/dx;
 | 
|---|
 | 277 | 
 | 
|---|
 | 278 |   if (nstep <= 0)
 | 
|---|
 | 279 |     nstep = 10;
 | 
|---|
 | 280 | 
 | 
|---|
 | 281 |   dx = (mXMax - mXMin) / nstep;
 | 
|---|
 | 282 | 
 | 
|---|
 | 283 |   double s = (FVal(mXMin) + FVal(mXMax))/2;
 | 
|---|
 | 284 |   double x = dx + mXMin;
 | 
|---|
 | 285 |   for (int i=1; i<nstep; i++, x += dx)
 | 
|---|
 | 286 |     s += FVal(x);
 | 
|---|
 | 287 |     
 | 
|---|
 | 288 |   return s * dx;
 | 
|---|
 | 289 | }
 | 
|---|
 | 290 | 
 | 
|---|
 | 291 | 
 | 
|---|
 | 292 | 
 | 
|---|
| [2808] | 293 | /*!
 | 
|---|
 | 294 |   \ingroup NTools
 | 
|---|
 | 295 |   \class SOPHYA::GLInteg
 | 
|---|
 | 296 |   
 | 
|---|
 | 297 |   \brief Implementation de Integrator par la methode de Gauss-Legendre.
 | 
|---|
| [220] | 298 | 
 | 
|---|
| [2808] | 299 |   Classe d'intégration par la méthode de Gauss-Legendre.
 | 
|---|
 | 300 |   Voir Integrator pour les méthodes. 
 | 
|---|
 | 301 |   Pour le moment, nstep est l'ordre de la méthode.
 | 
|---|
 | 302 |   Il est prévu un jour de spécifier l'ordre, et que NStep
 | 
|---|
 | 303 |   découpe en intervalles sur chacun desquels on applique GL.
 | 
|---|
 | 304 |   Le principe de la méthode est de calculer les valeurs de la 
 | 
|---|
 | 305 |   fonction aux zéros des polynomes de Legendre. Avec les poids
 | 
|---|
 | 306 |   qui vont bien, GL d'ordre n est exacte pour des polynomes de 
 | 
|---|
 | 307 |   degré <= 2n+1 (monome le + haut x^(2*n-1).
 | 
|---|
 | 308 |   Impossible de demander une précision donnée.
 | 
|---|
| [220] | 309 | 
 | 
|---|
| [2808] | 310 |   \sa SOPHYA::Integrator
 | 
|---|
| [220] | 311 | 
 | 
|---|
| [2808] | 312 |   \warning statut EXPERIMENTAL , NON TESTE
 | 
|---|
 | 313 | */
 | 
|---|
| [220] | 314 | 
 | 
|---|
 | 315 | 
 | 
|---|
| [2808] | 316 | 
 | 
|---|
| [220] | 317 | GLInteg::GLInteg()
 | 
|---|
 | 318 | : mXPos(NULL), mWeights(NULL)
 | 
|---|
 | 319 | {}
 | 
|---|
 | 320 | 
 | 
|---|
 | 321 | 
 | 
|---|
| [3409] | 322 | GLInteg::GLInteg(ClassFunc & f, double xmin, double xmax)
 | 
|---|
| [220] | 323 | : Integrator(f, xmin, xmax), mXPos(NULL), mWeights(NULL), mOrder(8)
 | 
|---|
 | 324 | {
 | 
|---|
 | 325 |   NStep(1);
 | 
|---|
 | 326 | }
 | 
|---|
| [2808] | 327 | 
 | 
|---|
| [3098] | 328 | GLInteg::GLInteg(FUNC f, double xmin, double xmax)
 | 
|---|
| [2283] | 329 | : Integrator(f, xmin, xmax), mXPos(NULL), mWeights(NULL), mOrder(8)
 | 
|---|
 | 330 | {
 | 
|---|
 | 331 |   NStep(1);
 | 
|---|
 | 332 | }
 | 
|---|
| [2808] | 333 | 
 | 
|---|
| [3409] | 334 | GLInteg::GLInteg(ClassFunc & f)
 | 
|---|
| [220] | 335 | : Integrator(f), mXPos(NULL), mWeights(NULL), mOrder(8)
 | 
|---|
 | 336 | {
 | 
|---|
 | 337 |   NStep(1);
 | 
|---|
 | 338 | }
 | 
|---|
| [2808] | 339 | 
 | 
|---|
| [3098] | 340 | GLInteg::GLInteg(FUNC f)
 | 
|---|
| [2283] | 341 | : Integrator(f), mXPos(NULL), mWeights(NULL), mOrder(8)
 | 
|---|
 | 342 | {
 | 
|---|
 | 343 |   NStep(1);
 | 
|---|
 | 344 | }
 | 
|---|
| [220] | 345 | 
 | 
|---|
 | 346 | GLInteg::GLInteg(GeneralFunction* gff, double* par, double xmin, double xmax)
 | 
|---|
 | 347 | : Integrator(gff, par, xmin, xmax), mXPos(NULL), mWeights(NULL), mOrder(8)
 | 
|---|
 | 348 | {
 | 
|---|
 | 349 |   NStep(1);
 | 
|---|
 | 350 | }
 | 
|---|
 | 351 | 
 | 
|---|
 | 352 | GLInteg::GLInteg(GeneralFunction* gff, double* par)
 | 
|---|
 | 353 | : Integrator(gff, par), mXPos(NULL), mWeights(NULL), mOrder(8)
 | 
|---|
 | 354 | {
 | 
|---|
 | 355 |   NStep(1);
 | 
|---|
 | 356 | }
 | 
|---|
 | 357 | 
 | 
|---|
 | 358 | GLInteg::~GLInteg()
 | 
|---|
 | 359 | {
 | 
|---|
| [2423] | 360 |  if(mXPos)    delete[] mXPos;     mXPos    = NULL;
 | 
|---|
 | 361 |  if(mWeights) delete[] mWeights;  mWeights = NULL;
 | 
|---|
| [220] | 362 | }
 | 
|---|
 | 363 | 
 | 
|---|
 | 364 | void
 | 
|---|
 | 365 | GLInteg::InvalWeights()
 | 
|---|
 | 366 | {
 | 
|---|
| [2423] | 367 |  if(mXPos)    delete[] mXPos;     mXPos    = NULL;
 | 
|---|
 | 368 |  if(mWeights) delete[] mWeights;  mWeights = NULL;
 | 
|---|
| [220] | 369 | }
 | 
|---|
 | 370 | 
 | 
|---|
 | 371 | 
 | 
|---|
 | 372 | void
 | 
|---|
 | 373 | GLInteg::LimitsChanged()
 | 
|---|
 | 374 | {
 | 
|---|
| [2423] | 375 |  ComputeBounds();
 | 
|---|
| [220] | 376 | }
 | 
|---|
 | 377 | 
 | 
|---|
 | 378 | void
 | 
|---|
 | 379 | GLInteg::StepsChanged()
 | 
|---|
 | 380 | {
 | 
|---|
| [2423] | 381 |  ComputeBounds();
 | 
|---|
| [220] | 382 | }
 | 
|---|
 | 383 | 
 | 
|---|
 | 384 | GLInteg&
 | 
|---|
 | 385 | GLInteg::AddBound(double x)
 | 
|---|
 | 386 | {
 | 
|---|
| [2506] | 387 |  if (x<=mXMin || x>=mXMax) throw RangeCheckError("GLInteg::AddBound()") ; 
 | 
|---|
| [2423] | 388 |  // On introduira les classes d'exections apres reflexion et de maniere systematique (Rz+cmv)
 | 
|---|
 | 389 |  // if (x<=mXMin || x>=mXMax) throw range_error("GLInteg::AddBound  bound outside interval");
 | 
|---|
 | 390 |  mBounds.insert(x);
 | 
|---|
 | 391 |  return *this;
 | 
|---|
| [220] | 392 | }
 | 
|---|
 | 393 | 
 | 
|---|
 | 394 | 
 | 
|---|
| [2808] | 395 | //! Definit l'ordre de la methode d'integration Gauus-Legendre
 | 
|---|
| [220] | 396 | GLInteg&
 | 
|---|
 | 397 | GLInteg::SetOrder(int order)
 | 
|---|
 | 398 | {
 | 
|---|
| [2423] | 399 |  mOrder = order;
 | 
|---|
 | 400 |  InvalWeights();
 | 
|---|
 | 401 |  return *this;
 | 
|---|
| [220] | 402 | }
 | 
|---|
 | 403 | 
 | 
|---|
 | 404 | 
 | 
|---|
| [2808] | 405 | //! Retourne la valeur de l'integrale
 | 
|---|
| [220] | 406 | double
 | 
|---|
 | 407 | GLInteg::Value()
 | 
|---|
 | 408 | {
 | 
|---|
| [2423] | 409 |  if (!mXPos) ComputeWeights();
 | 
|---|
 | 410 |  double s = 0;
 | 
|---|
 | 411 |  set<double>::iterator i=mBounds.begin();
 | 
|---|
 | 412 |  set<double>::iterator j=i; j++;
 | 
|---|
 | 413 |  while (j != mBounds.end()) {
 | 
|---|
 | 414 |    double s1 = 0;
 | 
|---|
 | 415 |    double x1 = *i;
 | 
|---|
 | 416 |    double x2 = *j;
 | 
|---|
 | 417 |    for(int k=0; k<mOrder; k++)
 | 
|---|
 | 418 |       s1 += mWeights[k] * FVal(x1 + (x2-x1)*mXPos[k]);
 | 
|---|
 | 419 |    s += s1*(x2-x1);
 | 
|---|
 | 420 |    i++; j++;
 | 
|---|
 | 421 |  }
 | 
|---|
 | 422 |  return s;
 | 
|---|
| [220] | 423 | }
 | 
|---|
 | 424 | 
 | 
|---|
 | 425 | void
 | 
|---|
 | 426 | GLInteg::ComputeBounds()
 | 
|---|
 | 427 | {
 | 
|---|
 | 428 |   mBounds.erase(mBounds.begin(), mBounds.end());
 | 
|---|
 | 429 |   for (int i=0; i<=mNStep; i++) 
 | 
|---|
 | 430 |     mBounds.insert(mXMin + (mXMax-mXMin)*i/mNStep);
 | 
|---|
 | 431 | }
 | 
|---|
 | 432 | 
 | 
|---|
 | 433 | void
 | 
|---|
 | 434 | GLInteg::ComputeWeights()
 | 
|---|
 | 435 | {
 | 
|---|
| [2423] | 436 |    const double EPS_gauleg = 5.0e-12;
 | 
|---|
| [220] | 437 | 
 | 
|---|
 | 438 |    int m=(mOrder+1)/2;
 | 
|---|
 | 439 |    const double xxMin = 0;
 | 
|---|
 | 440 |    const double xxMax = 1;
 | 
|---|
 | 441 |    
 | 
|---|
 | 442 |    mXPos    = new double[mOrder];
 | 
|---|
 | 443 |    mWeights = new double[mOrder];
 | 
|---|
 | 444 |    
 | 
|---|
 | 445 |    double xm=0.5*(xxMax+xxMin);
 | 
|---|
 | 446 |    double xl=0.5*(xxMax-xxMin);
 | 
|---|
 | 447 |    for (int i=1;i<=m;i++)  {
 | 
|---|
 | 448 |       double z=cos(3.141592654*(i-0.25)/(mOrder+0.5));
 | 
|---|
 | 449 |       double z1, pp;
 | 
|---|
 | 450 |       do {
 | 
|---|
 | 451 |          double p1=1.0;
 | 
|---|
 | 452 |          double p2=0.0;
 | 
|---|
 | 453 |          for (int j=1;j<=mOrder;j++) {
 | 
|---|
 | 454 |             double p3=p2;
 | 
|---|
 | 455 |             p2=p1;
 | 
|---|
 | 456 |             p1=((2.0*j-1.0)*z*p2-(j-1.0)*p3)/j;
 | 
|---|
 | 457 |          }
 | 
|---|
 | 458 |          pp=mOrder*(z*p1-p2)/(z*z-1.0);
 | 
|---|
 | 459 |          z1=z;
 | 
|---|
 | 460 |          z=z1-p1/pp;
 | 
|---|
 | 461 |       } while (fabs(z-z1) > EPS_gauleg);
 | 
|---|
 | 462 |       mXPos[i-1]         = xm-xl*z;
 | 
|---|
 | 463 |       mXPos[mOrder-i]    = xm+xl*z;
 | 
|---|
 | 464 |       mWeights[i-1]      = 2.0*xl/((1.0-z*z)*pp*pp);
 | 
|---|
 | 465 |       mWeights[mOrder-i] = mWeights[i-1];
 | 
|---|
 | 466 |    }
 | 
|---|
 | 467 | 
 | 
|---|
 | 468 | }
 | 
|---|
 | 469 | 
 | 
|---|
| [2808] | 470 | //! Imprime l'ordre et la valeur des poids sur cout
 | 
|---|
| [2423] | 471 | void
 | 
|---|
 | 472 | GLInteg::Print(int lp)
 | 
|---|
 | 473 | {
 | 
|---|
 | 474 |  Integrator::Print(lp);
 | 
|---|
 | 475 |  cout<<"GLInteg order="<<mOrder<<endl;
 | 
|---|
 | 476 |  if(lp>0 && mOrder>0) {
 | 
|---|
 | 477 |    for(int i=0;i<mOrder;i++)
 | 
|---|
 | 478 |      cout<<" ("<<mXPos[i]<<","<<mWeights[i]<<")";
 | 
|---|
 | 479 |    cout<<endl;
 | 
|---|
 | 480 |  }
 | 
|---|
 | 481 | }
 | 
|---|