1 | #include "machdefs.h"
|
---|
2 | #include "rk4cdifeq.h"
|
---|
3 | #include "ctimer.h"
|
---|
4 | #include <iostream>
|
---|
5 |
|
---|
6 | inline double max_dbl(double a, double b) { return ((a < b) ? b : a); }
|
---|
7 |
|
---|
8 | //++
|
---|
9 | // Class RK4CDiffEq
|
---|
10 | // Lib Outils++
|
---|
11 | // include rk4cdifeq.h
|
---|
12 | //
|
---|
13 | // Classe de résolution d'équadif par la méthode de
|
---|
14 | // Runge-Kutta d'ordre 4 adaptatif.
|
---|
15 | // Voir DiffEqSolver et R4KDiffEq pour les méthodes.
|
---|
16 | //
|
---|
17 | // On peut demander une précision relative ou absolue
|
---|
18 | // sur chacune des fonctions du système. Le pas d'intégration
|
---|
19 | // est adapté automatiquement pour fournir au moins
|
---|
20 | // cette précision.
|
---|
21 | //--
|
---|
22 |
|
---|
23 | //++
|
---|
24 | // Links Parents
|
---|
25 | // RK4DiffEq
|
---|
26 | // DiffEqSolver
|
---|
27 | //--
|
---|
28 |
|
---|
29 | //++
|
---|
30 | // Titre Constructeurs
|
---|
31 | // Voir RK4DiffEq et DiffEqSolver
|
---|
32 | //--
|
---|
33 |
|
---|
34 | //++
|
---|
35 | RK4CDiffEq::RK4CDiffEq()
|
---|
36 | //
|
---|
37 | //--
|
---|
38 | : RK4DiffEq(), eps(1.e-4), relAccuracy(true), accScale(10), yTemp(10), ySave(10)
|
---|
39 | {}
|
---|
40 |
|
---|
41 | //++
|
---|
42 | RK4CDiffEq::RK4CDiffEq(DiffEqFunction* f)
|
---|
43 | //
|
---|
44 | //--
|
---|
45 | : RK4DiffEq(f), eps(1.e-4), relAccuracy(true), accScale(f->NFuncReal()), yTemp(f->NFuncReal()), ySave(f->NFuncReal())
|
---|
46 | {}
|
---|
47 |
|
---|
48 | //++
|
---|
49 | RK4CDiffEq::RK4CDiffEq(DIFEQFCN1 f)
|
---|
50 | //
|
---|
51 | //--
|
---|
52 | : RK4DiffEq(f), eps(1.e-4), relAccuracy(true), accScale(1), yTemp(1), ySave(1)
|
---|
53 | {}
|
---|
54 |
|
---|
55 | //++
|
---|
56 | // Titre Méthodes
|
---|
57 | //--
|
---|
58 |
|
---|
59 | //++
|
---|
60 | // RK4CDiffEq& RK4CDiffEq::Accuracy(double eps)
|
---|
61 | // Fixe la précision requise sur les fonctions. Cette précision
|
---|
62 | // est par défaut relative. Elle peut être absolue, auquel cas
|
---|
63 | // il faut fixer, pour chaque fonction, un facteur d'échelle, et
|
---|
64 | // la précision sur chaque fonction est alors scale[i]*eps.
|
---|
65 | //--
|
---|
66 |
|
---|
67 | RK4CDiffEq&
|
---|
68 | RK4CDiffEq::Accuracy(double x)
|
---|
69 | {
|
---|
70 | eps = x;
|
---|
71 | return *this;
|
---|
72 | }
|
---|
73 |
|
---|
74 | //++
|
---|
75 | RK4CDiffEq&
|
---|
76 | RK4CDiffEq::AbsAccuracy(Vector const& vScal)
|
---|
77 | //
|
---|
78 | // On souhaite une précision absolue, et le vecteur vScal contient
|
---|
79 | // le facteur d'échelle pour chaque fonction (voir Accuracy).
|
---|
80 | //--
|
---|
81 | {
|
---|
82 | accScale = vScal;
|
---|
83 | relAccuracy = false;
|
---|
84 | return *this;
|
---|
85 | }
|
---|
86 |
|
---|
87 | //++
|
---|
88 | RK4CDiffEq&
|
---|
89 | RK4CDiffEq::AbsAccuracy(double scal)
|
---|
90 | //
|
---|
91 | // On souhaite une précision absolue, et le vecteur scal contient
|
---|
92 | // le facteur d'échelle à appliquer à toutes les fonctions.
|
---|
93 | // La précision absolue souhaitée est alors eps*scal.
|
---|
94 | //--
|
---|
95 | {
|
---|
96 | for (int i=0; i<accScale.NElts(); i++)
|
---|
97 | accScale(i) = scal;
|
---|
98 | relAccuracy = false;
|
---|
99 | return *this;
|
---|
100 | }
|
---|
101 |
|
---|
102 | //++
|
---|
103 | RK4CDiffEq&
|
---|
104 | RK4CDiffEq::RelAccuracy()
|
---|
105 | //
|
---|
106 | // On souhaite une précision relative. En quelque sorte, le facteur d'échelle
|
---|
107 | // pour chaque fonction est alors la valeur de la fonction.
|
---|
108 | //--
|
---|
109 | {
|
---|
110 | relAccuracy = true;
|
---|
111 | return *this;
|
---|
112 | }
|
---|
113 |
|
---|
114 | static const double pgrow = -0.20;
|
---|
115 | static const double pshrink = -0.25;
|
---|
116 | static const double safety = 0.9;
|
---|
117 | static const double errcon = pow((4/safety),1/pgrow);
|
---|
118 |
|
---|
119 | void
|
---|
120 | RK4CDiffEq::RKCStep(Vector& newY, Vector const& y0, Vector const& yScale,
|
---|
121 | double dttry, double& dtdone, double& dtnext)
|
---|
122 | {
|
---|
123 | double err;
|
---|
124 | ySave = y0;
|
---|
125 | do {
|
---|
126 | // Deux petits pas
|
---|
127 | RKStep(yTemp, ySave, dttry/2.);
|
---|
128 | RKStep(newY, yTemp, dttry/2.);
|
---|
129 |
|
---|
130 | // Un grand pas
|
---|
131 | RKStep(yTemp, ySave, dttry);
|
---|
132 |
|
---|
133 | yTemp -= newY; // l'erreur courante
|
---|
134 |
|
---|
135 | err = 0;
|
---|
136 | for (int i=0; i<yTemp.NElts(); i++)
|
---|
137 | err = max_dbl(fabs(yScale(i) ? (yTemp(i)/yScale(i)) : yTemp(i)), err);
|
---|
138 |
|
---|
139 | err /= eps;
|
---|
140 | if (err > 1)
|
---|
141 | dttry *= safety*pow(err,pshrink);
|
---|
142 | else {
|
---|
143 | dtdone = dttry;
|
---|
144 | if (err > errcon)
|
---|
145 | dtnext = safety*dttry*pow(err,pgrow);
|
---|
146 | else
|
---|
147 | dtnext = dttry*4;
|
---|
148 | }
|
---|
149 | } while (err > 1);
|
---|
150 |
|
---|
151 | // Et on corrige a l'ordre 5
|
---|
152 |
|
---|
153 | newY += yTemp/15.;
|
---|
154 | }
|
---|
155 |
|
---|
156 | void
|
---|
157 | RK4CDiffEq::SolveArr(Matrix& y, double* t, double tf, int n)
|
---|
158 | {
|
---|
159 | //TIMEF;
|
---|
160 | // Les intervalles a stocker dans la matrice des resultats
|
---|
161 |
|
---|
162 | double dxres = (tf - mXStart)/n;
|
---|
163 |
|
---|
164 | Vector yt(mYStart,false);
|
---|
165 |
|
---|
166 | k1.Realloc(mFunc->NFuncReal());
|
---|
167 | k2.Realloc(mFunc->NFuncReal());
|
---|
168 | k3.Realloc(mFunc->NFuncReal());
|
---|
169 | k4.Realloc(mFunc->NFuncReal());
|
---|
170 | yTemp.Realloc(mFunc->NFuncReal());
|
---|
171 | ySave.Realloc(mFunc->NFuncReal());
|
---|
172 |
|
---|
173 | double x = mXStart;
|
---|
174 | double step = mStep;
|
---|
175 | double nstep;
|
---|
176 |
|
---|
177 | for (int i=0; i<n; i++) {
|
---|
178 | double xEndStep = (i+1)*dxres + mXStart;
|
---|
179 | do {
|
---|
180 | if (relAccuracy) accScale = yt;
|
---|
181 | if ((x+step-xEndStep)*(x+step-mXStart)>0) step = xEndStep-x;
|
---|
182 | RKCStep(yt, yt, accScale, step, step, nstep);
|
---|
183 | x += step;
|
---|
184 | //cout << x << " " << step << endl;
|
---|
185 | step = nstep;
|
---|
186 | } while ((x-xEndStep)*(x-mXStart) < 0);
|
---|
187 | for (int j=0; j<mFunc->NFunc(); j++)
|
---|
188 | y(i,j) = yt(j);
|
---|
189 | t[i] = xEndStep;
|
---|
190 | //cout << "fin boucle "<<i<<" "<<x<<" "<<yt(0)<<endl;
|
---|
191 | }
|
---|
192 |
|
---|
193 | }
|
---|