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;
56 jacRedo_(
min(1
e-4,
min(relTol_))),
59 coeff_(iMaxx_, iMaxx_),
76 const scalar cpuFunc = 1, cpuJac = 5, cpuLU = 1, cpuSolve = 1;
81 for (
int i=2; i<iMaxx_; i++)
83 nSeq_[i] = 2*nSeq_[i-2];
85 cpu_[0] = cpuJac + cpuLU + nSeq_[0]*(cpuFunc + cpuSolve);
87 for (
int k=0;
k<kMaxx_;
k++)
89 cpu_[
k+1] = cpu_[
k] + (nSeq_[
k+1]-1)*(cpuFunc + cpuSolve) + cpuLU;
93 for (
int k=0;
k<iMaxx_;
k++)
95 for (
int l=0; l<
k; l++)
97 scalar ratio = scalar(nSeq_[
k])/nSeq_[l];
98 coeff_(
k, l) = 1/(ratio - 1);
106bool Foam::seulex::seul
116 label nSteps = nSeq_[
k];
117 scalar dx = dxTot/nSteps;
119 for (label i=0; i<n_; i++)
121 for (label j=0; j<n_; j++)
123 a_(i, j) = -dfdy_(i, j);
131 scalar xnew = x0 + dx;
132 odes_.derivatives(xnew,
y0, dy_);
137 for (label nn=1; nn<nSteps; nn++)
145 for (label i=0; i<n_; i++)
147 dy1 +=
sqr(dy_[i]/scale[i]);
151 odes_.derivatives(x0 + dx, yTemp_, dydx_);
152 for (label i=0; i<n_; i++)
154 dy_[i] = dydx_[i] - dy_[i]/dx;
159 const scalar denom =
min(1, dy1 + SMALL);
161 for (label i=0; i<n_; i++)
164 if (
mag(dy_[i]) > scale[i]*denom)
170 dy2 +=
sqr(dy_[i]/scale[i]);
181 odes_.derivatives(xnew, yTemp_, dy_);
185 for (label i=0; i<n_; i++)
187 y[i] = yTemp_[i] + dy_[i];
194void Foam::seulex::extrapolate
201 for (
int j=
k-1; j>0; j--)
203 for (label i=0; i<n_; i++)
206 table(j, i) + coeff_(
k, j)*(table(j, i) - table[j-1][i]);
210 for (
int i=0; i<n_; i++)
212 y[i] = table(0, i) + coeff_(
k, 0)*(table(0, i) -
y[i]);
221 table_.shallowResize(kMaxx_, n_);
225 resizeField(pivotIndices_);
227 resizeField(ySequence_);
248 scalar dx = step.
dxTry;
250 dxOpt_[0] =
mag(0.1*dx);
260 scalar logTol = -
log10(relTol_[0] + absTol_[0])*0.6 + 0.5;
261 kTarg_ =
max(1,
min(kMaxx_ - 1,
int(logTol)));
266 scale_[i] = absTol_[i] + relTol_[i]*
mag(
y[i]);
269 bool jacUpdated =
false;
271 if (theta_ > jacRedo_)
273 odes_.jacobian(
x,
y, dfdx_, dfdy_);
278 scalar dxNew =
mag(dx);
281 while (firstk || step.
reject)
283 dx = step.
forward ? dxNew : -dxNew;
290 <<
"step size underflow :" << dx <<
endl;
295 for (
k=0;
k<=kTarg_+1;
k++)
297 bool success = seul(
x, y0_, dx,
k, ySequence_, scale_);
302 dxNew =
mag(dx)*stepFactor5_;
314 table_[
k-1][i] = ySequence_[i];
320 extrapolate(
k, table_,
y);
324 scale_[i] = absTol_[i] + relTol_[i]*
mag(y0_[i]);
325 err +=
sqr((
y[i] - table_(0, i))/scale_[i]);
328 if (err > 1/SMALL || (
k > 1 && err >= errOld))
331 dxNew =
mag(dx)*stepFactor5_;
334 errOld =
min(4*err, 1);
335 scalar expo = 1.0/(
k + 1);
336 scalar facmin =
pow(stepFactor3_, expo);
344 fac = stepFactor2_/
pow(err/stepFactor1_, expo);
348 temp_[
k] = cpu_[
k]/dxOpt_[
k];
350 if ((step.
first || step.
last) && err <= 1)
366 else if (err > nSeq_[kTarg_]*nSeq_[kTarg_ + 1]*4)
370 if (kTarg_>1 && temp_[
k-1] < kFactor1_*temp_[
k])
374 dxNew = dxOpt_[kTarg_];
385 else if (err > nSeq_[
k + 1]*2)
388 if (kTarg_>1 && temp_[
k-1] < kFactor1_*temp_[
k])
392 dxNew = dxOpt_[kTarg_];
405 && temp_[kTarg_-1] < kFactor1_*temp_[kTarg_]
410 dxNew = dxOpt_[kTarg_];
423 if (theta_ > jacRedo_ && !jacUpdated)
425 odes_.jacobian(
x,
y, dfdx_, dfdy_);
442 else if (
k <= kTarg_)
445 if (temp_[
k-1] < kFactor1_*temp_[
k])
449 else if (temp_[
k] < kFactor2_*temp_[
k - 1])
451 kopt =
min(
k + 1, kMaxx_ - 1);
457 if (
k > 2 && temp_[
k-2] < kFactor1_*temp_[
k - 1])
461 if (temp_[
k] < kFactor2_*temp_[kopt])
463 kopt =
min(
k, kMaxx_ - 1);
469 kTarg_ =
min(kopt,
k);
470 dxNew =
min(
mag(dx), dxOpt_[kTarg_]);
477 dxNew = dxOpt_[kopt];
481 if (
k < kTarg_ && temp_[
k] < kFactor2_*temp_[
k - 1])
483 dxNew = dxOpt_[
k]*cpu_[kopt + 1]/cpu_[
k];
487 dxNew = dxOpt_[
k]*cpu_[kopt]/cpu_[
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.
virtual bool resize()=0
Resize the ODE solver.
Abstract base class for the systems of ordinary differential equations.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
An ODE solver for chemistry.
An extrapolation-algorithm, based on the linearly implicit Euler method with step size control and or...
virtual bool resize()
Resize the ODE solver.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define WarningInFunction
Report a warning using Foam::Warning.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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.
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.
Calculate the second temporal derivative.
#define forAll(list, i)
Loop across all elements in list.