seulex.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "seulex.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37
39
40 const scalar
41 seulex::stepFactor1_ = 0.6,
42 seulex::stepFactor2_ = 0.93,
43 seulex::stepFactor3_ = 0.1,
44 seulex::stepFactor4_ = 4,
45 seulex::stepFactor5_ = 0.5,
46 seulex::kFactor1_ = 0.7,
47 seulex::kFactor2_ = 0.9;
48}
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
54:
56 jacRedo_(min(1e-4, min(relTol_))),
57 nSeq_(iMaxx_),
58 cpu_(iMaxx_),
59 coeff_(iMaxx_, iMaxx_),
60 theta_(2*jacRedo_),
61 table_(kMaxx_, n_),
62 dfdx_(n_),
63 dfdy_(n_),
64 a_(n_),
65 pivotIndices_(n_),
66 dxOpt_(iMaxx_),
67 temp_(iMaxx_),
68 y0_(n_),
69 ySequence_(n_),
70 scale_(n_),
71 dy_(n_),
72 yTemp_(n_),
73 dydx_(n_)
74{
75 // The CPU time factors for the major parts of the algorithm
76 const scalar cpuFunc = 1, cpuJac = 5, cpuLU = 1, cpuSolve = 1;
77
78 nSeq_[0] = 2;
79 nSeq_[1] = 3;
80
81 for (int i=2; i<iMaxx_; i++)
82 {
83 nSeq_[i] = 2*nSeq_[i-2];
84 }
85 cpu_[0] = cpuJac + cpuLU + nSeq_[0]*(cpuFunc + cpuSolve);
86
87 for (int k=0; k<kMaxx_; k++)
88 {
89 cpu_[k+1] = cpu_[k] + (nSeq_[k+1]-1)*(cpuFunc + cpuSolve) + cpuLU;
90 }
91
92 // Set the extrapolation coefficients array
93 for (int k=0; k<iMaxx_; k++)
94 {
95 for (int l=0; l<k; l++)
96 {
97 scalar ratio = scalar(nSeq_[k])/nSeq_[l];
98 coeff_(k, l) = 1/(ratio - 1);
99 }
100 }
101}
102
103
104// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105
106bool Foam::seulex::seul
107(
108 const scalar x0,
109 const scalarField& y0,
110 const scalar dxTot,
111 const label k,
112 scalarField& y,
113 const scalarField& scale
114) const
115{
116 label nSteps = nSeq_[k];
117 scalar dx = dxTot/nSteps;
118
119 for (label i=0; i<n_; i++)
120 {
121 for (label j=0; j<n_; j++)
122 {
123 a_(i, j) = -dfdy_(i, j);
124 }
125
126 a_(i, i) += 1/dx;
127 }
128
129 LUDecompose(a_, pivotIndices_);
130
131 scalar xnew = x0 + dx;
132 odes_.derivatives(xnew, y0, dy_);
133 LUBacksubstitute(a_, pivotIndices_, dy_);
134
135 yTemp_ = y0;
136
137 for (label nn=1; nn<nSteps; nn++)
138 {
139 yTemp_ += dy_;
140 xnew += dx;
141
142 if (nn == 1 && k<=1)
143 {
144 scalar dy1 = 0;
145 for (label i=0; i<n_; i++)
146 {
147 dy1 += sqr(dy_[i]/scale[i]);
148 }
149 dy1 = sqrt(dy1);
150
151 odes_.derivatives(x0 + dx, yTemp_, dydx_);
152 for (label i=0; i<n_; i++)
153 {
154 dy_[i] = dydx_[i] - dy_[i]/dx;
155 }
156
157 LUBacksubstitute(a_, pivotIndices_, dy_);
158
159 const scalar denom = min(1, dy1 + SMALL);
160 scalar dy2 = 0;
161 for (label i=0; i<n_; i++)
162 {
163 // Test of dy_[i] to avoid overflow
164 if (mag(dy_[i]) > scale[i]*denom)
165 {
166 theta_ = 1;
167 return false;
168 }
169
170 dy2 += sqr(dy_[i]/scale[i]);
171 }
172 dy2 = sqrt(dy2);
173 theta_ = dy2/denom;
174
175 if (theta_ > 1)
176 {
177 return false;
178 }
179 }
180
181 odes_.derivatives(xnew, yTemp_, dy_);
182 LUBacksubstitute(a_, pivotIndices_, dy_);
183 }
184
185 for (label i=0; i<n_; i++)
186 {
187 y[i] = yTemp_[i] + dy_[i];
188 }
189
190 return true;
191}
192
193
194void Foam::seulex::extrapolate
195(
196 const label k,
199) const
200{
201 for (int j=k-1; j>0; j--)
202 {
203 for (label i=0; i<n_; i++)
204 {
205 table[j-1][i] =
206 table(j, i) + coeff_(k, j)*(table(j, i) - table[j-1][i]);
207 }
208 }
209
210 for (int i=0; i<n_; i++)
211 {
212 y[i] = table(0, i) + coeff_(k, 0)*(table(0, i) - y[i]);
213 }
214}
215
216
218{
219 if (ODESolver::resize())
220 {
221 table_.shallowResize(kMaxx_, n_);
222 resizeField(dfdx_);
223 resizeMatrix(dfdy_);
224 resizeMatrix(a_);
225 resizeField(pivotIndices_);
226 resizeField(y0_);
227 resizeField(ySequence_);
228 resizeField(scale_);
229 resizeField(dy_);
230 resizeField(yTemp_);
231 resizeField(dydx_);
232
233 return true;
234 }
235
236 return false;
237}
238
239
241(
242 scalar& x,
243 scalarField& y,
244 stepState& step
245) const
246{
247 temp_[0] = GREAT;
248 scalar dx = step.dxTry;
249 y0_ = y;
250 dxOpt_[0] = mag(0.1*dx);
251
252 if (step.first || step.prevReject)
253 {
254 theta_ = 2*jacRedo_;
255 }
256
257 if (step.first)
258 {
259 // NOTE: the first element of relTol_ and absTol_ are used here.
260 scalar logTol = -log10(relTol_[0] + absTol_[0])*0.6 + 0.5;
261 kTarg_ = max(1, min(kMaxx_ - 1, int(logTol)));
262 }
263
264 forAll(scale_, i)
265 {
266 scale_[i] = absTol_[i] + relTol_[i]*mag(y[i]);
267 }
268
269 bool jacUpdated = false;
270
271 if (theta_ > jacRedo_)
272 {
273 odes_.jacobian(x, y, dfdx_, dfdy_);
274 jacUpdated = true;
275 }
276
277 int k;
278 scalar dxNew = mag(dx);
279 bool firstk = true;
280
281 while (firstk || step.reject)
282 {
283 dx = step.forward ? dxNew : -dxNew;
284 firstk = false;
285 step.reject = false;
286
287 if (mag(dx) <= mag(x)*sqr(SMALL))
288 {
290 << "step size underflow :" << dx << endl;
291 }
292
293 scalar errOld = 0;
294
295 for (k=0; k<=kTarg_+1; k++)
296 {
297 bool success = seul(x, y0_, dx, k, ySequence_, scale_);
298
299 if (!success)
300 {
301 step.reject = true;
302 dxNew = mag(dx)*stepFactor5_;
303 break;
304 }
305
306 if (k == 0)
307 {
308 y = ySequence_;
309 }
310 else
311 {
312 forAll(ySequence_, i)
313 {
314 table_[k-1][i] = ySequence_[i];
315 }
316 }
317
318 if (k != 0)
319 {
320 extrapolate(k, table_, y);
321 scalar err = 0;
322 forAll(scale_, i)
323 {
324 scale_[i] = absTol_[i] + relTol_[i]*mag(y0_[i]);
325 err += sqr((y[i] - table_(0, i))/scale_[i]);
326 }
327 err = sqrt(err/n_);
328 if (err > 1/SMALL || (k > 1 && err >= errOld))
329 {
330 step.reject = true;
331 dxNew = mag(dx)*stepFactor5_;
332 break;
333 }
334 errOld = min(4*err, 1);
335 scalar expo = 1.0/(k + 1);
336 scalar facmin = pow(stepFactor3_, expo);
337 scalar fac;
338 if (err == 0)
339 {
340 fac = 1/facmin;
341 }
342 else
343 {
344 fac = stepFactor2_/pow(err/stepFactor1_, expo);
345 fac = max(facmin/stepFactor4_, min(1/facmin, fac));
346 }
347 dxOpt_[k] = mag(dx*fac);
348 temp_[k] = cpu_[k]/dxOpt_[k];
349
350 if ((step.first || step.last) && err <= 1)
351 {
352 break;
353 }
354
355 if
356 (
357 k == kTarg_ - 1
358 && !step.prevReject
359 && !step.first && !step.last
360 )
361 {
362 if (err <= 1)
363 {
364 break;
365 }
366 else if (err > nSeq_[kTarg_]*nSeq_[kTarg_ + 1]*4)
367 {
368 step.reject = true;
369 kTarg_ = k;
370 if (kTarg_>1 && temp_[k-1] < kFactor1_*temp_[k])
371 {
372 kTarg_--;
373 }
374 dxNew = dxOpt_[kTarg_];
375 break;
376 }
377 }
378
379 if (k == kTarg_)
380 {
381 if (err <= 1)
382 {
383 break;
384 }
385 else if (err > nSeq_[k + 1]*2)
386 {
387 step.reject = true;
388 if (kTarg_>1 && temp_[k-1] < kFactor1_*temp_[k])
389 {
390 kTarg_--;
391 }
392 dxNew = dxOpt_[kTarg_];
393 break;
394 }
395 }
396
397 if (k == kTarg_+1)
398 {
399 if (err > 1)
400 {
401 step.reject = true;
402 if
403 (
404 kTarg_ > 1
405 && temp_[kTarg_-1] < kFactor1_*temp_[kTarg_]
406 )
407 {
408 kTarg_--;
409 }
410 dxNew = dxOpt_[kTarg_];
411 }
412 break;
413 }
414 }
415 }
416 if (step.reject)
417 {
418 step.prevReject = true;
419 if (!jacUpdated)
420 {
421 theta_ = 2*jacRedo_;
422
423 if (theta_ > jacRedo_ && !jacUpdated)
424 {
425 odes_.jacobian(x, y, dfdx_, dfdy_);
426 jacUpdated = true;
427 }
428 }
429 }
430 }
431
432 jacUpdated = false;
433
434 step.dxDid = dx;
435 x += dx;
436
437 label kopt;
438 if (k == 1)
439 {
440 kopt = 2;
441 }
442 else if (k <= kTarg_)
443 {
444 kopt=k;
445 if (temp_[k-1] < kFactor1_*temp_[k])
446 {
447 kopt = k - 1;
448 }
449 else if (temp_[k] < kFactor2_*temp_[k - 1])
450 {
451 kopt = min(k + 1, kMaxx_ - 1);
452 }
453 }
454 else
455 {
456 kopt = k - 1;
457 if (k > 2 && temp_[k-2] < kFactor1_*temp_[k - 1])
458 {
459 kopt = k - 2;
460 }
461 if (temp_[k] < kFactor2_*temp_[kopt])
462 {
463 kopt = min(k, kMaxx_ - 1);
464 }
465 }
466
467 if (step.prevReject)
468 {
469 kTarg_ = min(kopt, k);
470 dxNew = min(mag(dx), dxOpt_[kTarg_]);
471 step.prevReject = false;
472 }
473 else
474 {
475 if (kopt <= k)
476 {
477 dxNew = dxOpt_[kopt];
478 }
479 else
480 {
481 if (k < kTarg_ && temp_[k] < kFactor2_*temp_[k - 1])
482 {
483 dxNew = dxOpt_[k]*cpu_[kopt + 1]/cpu_[k];
484 }
485 else
486 {
487 dxNew = dxOpt_[k]*cpu_[kopt]/cpu_[k];
488 }
489 }
490 kTarg_ = kopt;
491 }
492
493 step.dxTry = step.forward ? dxNew : -dxNew;
494}
495
496
497// ************************************************************************* //
bool success
scalar y
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Abstract base-class for ODE system solvers.
Definition: ODESolver.H:57
virtual bool resize()=0
Resize the ODE solver.
Definition: ODESolver.C:92
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
An ODE solver for chemistry.
Definition: ode.H:55
An extrapolation-algorithm, based on the linearly implicit Euler method with step size control and or...
Definition: seulex.H:68
virtual bool resize()
Resize the ODE solver.
Definition: seulex.C:217
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
RectangularMatrix< scalar > scalarRectangularMatrix
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar y0(const dimensionedScalar &ds)
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
dimensionedScalar log10(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Calculate the second temporal derivative.
dictionary dict
CEqn solve()
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333