38 seulex::stepFactor1_ = 0.6,
39 seulex::stepFactor2_ = 0.93,
40 seulex::stepFactor3_ = 0.1,
41 seulex::stepFactor4_ = 4,
42 seulex::stepFactor5_ = 0.5,
43 seulex::kFactor1_ = 0.7,
44 seulex::kFactor2_ = 0.9;
53 jacRedo_(
min(1
e-4,
min(relTol_))),
56 coeff_(iMaxx_, iMaxx_),
73 const scalar cpuFunc = 1, cpuJac = 5, cpuLU = 1, cpuSolve = 1;
78 for (
int i=2; i<iMaxx_; i++)
80 nSeq_[i] = 2*nSeq_[i-2];
82 cpu_[0] = cpuJac + cpuLU + nSeq_[0]*(cpuFunc + cpuSolve);
84 for (
int k=0; k<kMaxx_; k++)
86 cpu_[k+1] = cpu_[
k] + (nSeq_[k+1]-1)*(cpuFunc + cpuSolve) + cpuLU;
90 for (
int k=0; k<iMaxx_; k++)
92 for (
int l=0; l<
k; l++)
94 scalar ratio = scalar(nSeq_[k])/nSeq_[l];
95 coeff_(k, l) = 1/(ratio - 1);
103 bool Foam::seulex::seul
115 scalar dx = dxTot/nSteps;
121 a_(i, j) = -dfdy_(i, j);
129 scalar xnew = x0 + dx;
135 for (
label nn=1; nn<nSteps; nn++)
145 dy1 +=
sqr(dy_[i]/scale[i]);
152 dy_[i] = dydx_[i] - dy_[i]/dx;
162 const scalar denom =
min(1, dy1 + small);
168 if (
mag(dy_[i]) > scale[i]*denom)
174 dy2 +=
sqr(dy_[i]/scale[i]);
191 y[i] = yTemp_[i] + dy_[i];
198 void Foam::seulex::extrapolate
205 for (
int j=k-1; j>0; j--)
210 table(j, i) + coeff_(k, j)*(table(j, i) - table[j-1][i]);
214 for (
int i=0; i<
n_; i++)
216 y[i] = table(0, i) + coeff_(k, 0)*(table(0, i) - y[i]);
255 scalar dx = step.
dxTry;
257 dxOpt_[0] =
mag(0.1*dx);
268 kTarg_ =
max(1,
min(kMaxx_ - 1,
int(logTol)));
276 bool jacUpdated =
false;
278 if (theta_ > jacRedo_)
285 scalar dxNew =
mag(dx);
288 while (firstk || step.
reject)
290 dx = step.
forward ? dxNew : -dxNew;
297 <<
"step size underflow :" << dx <<
endl;
302 for (k=0; k<=kTarg_+1; k++)
304 bool success = seul(x, y0_, li, dx, k, ySequence_, scale_);
309 dxNew =
mag(dx)*stepFactor5_;
321 table_[k-1][i] = ySequence_[i];
327 extrapolate(k, table_, y);
332 err +=
sqr((y[i] - table_(0, i))/scale_[i]);
335 if (err > 1/small || (k > 1 && err >= errOld))
338 dxNew =
mag(dx)*stepFactor5_;
341 errOld =
min(4*err, 1);
342 scalar expo = 1.0/(k + 1);
343 scalar facmin =
pow(stepFactor3_, expo);
351 fac = stepFactor2_/
pow(err/stepFactor1_, expo);
352 fac =
max(facmin/stepFactor4_,
min(1/facmin, fac));
354 dxOpt_[
k] =
mag(dx*fac);
355 temp_[
k] = cpu_[
k]/dxOpt_[
k];
357 if ((step.
first || step.
last) && err <= 1)
373 else if (err > nSeq_[kTarg_]*nSeq_[kTarg_ + 1]*4)
377 if (kTarg_>1 && temp_[k-1] < kFactor1_*temp_[k])
381 dxNew = dxOpt_[kTarg_];
392 else if (err > nSeq_[k + 1]*2)
395 if (kTarg_>1 && temp_[k-1] < kFactor1_*temp_[k])
399 dxNew = dxOpt_[kTarg_];
412 && temp_[kTarg_-1] < kFactor1_*temp_[kTarg_]
417 dxNew = dxOpt_[kTarg_];
430 if (theta_ > jacRedo_ && !jacUpdated)
449 else if (k <= kTarg_)
452 if (temp_[k-1] < kFactor1_*temp_[k])
456 else if (temp_[k] < kFactor2_*temp_[k - 1])
458 kopt =
min(k + 1, kMaxx_ - 1);
464 if (k > 2 && temp_[k-2] < kFactor1_*temp_[k - 1])
468 if (temp_[k] < kFactor2_*temp_[kopt])
470 kopt =
min(k, kMaxx_ - 1);
476 kTarg_ =
min(kopt, k);
477 dxNew =
min(
mag(dx), dxOpt_[kTarg_]);
484 dxNew = dxOpt_[kopt];
488 if (k < kTarg_ && temp_[k] < kFactor2_*temp_[k - 1])
490 dxNew = dxOpt_[
k]*cpu_[kopt + 1]/cpu_[
k];
494 dxNew = dxOpt_[
k]*cpu_[kopt]/cpu_[
k];
void shallowResize(const label m, const label n)
Resize the matrix without reallocating storage (unsafe)
#define forAll(list, i)
Loop across all elements in list.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Abstract base class for the systems of ordinary differential equations.
virtual bool resize()=0
Resize the ODE solver.
A list of keyword definitions, which are a keyword followed by any number of values (e...
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual void derivatives(const scalar x, const scalarField &y, const label li, scalarField &dydx) const =0
Calculate the derivatives in dydx.
dimensionedScalar y0(const dimensionedScalar &ds)
virtual void solve(scalar &x, scalarField &y, const label li, stepState &step) const
Solve the ODE system and the update the state.
label k
Boltzmann constant.
An ODE solver for chemistry.
Macros for easy insertion into run-time selection tables.
virtual bool resize()
Resize the ODE solver.
void resizeMatrix(scalarSquareMatrix &m) const
const ODESystem & odes_
Reference to ODESystem.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define WarningInFunction
Report a warning using Foam::Warning.
scalarField relTol_
Relative convergence tolerance per step.
Abstract base-class for ODE system solvers.
virtual void jacobian(const scalar x, const scalarField &y, const label li, scalarField &dfdx, scalarSquareMatrix &dfdy) const =0
Calculate the Jacobian of the system.
scalarField absTol_
Absolute convergence tolerance per step.
seulex(const ODESystem &ode, const dictionary &dict)
Construct from ODESystem.
dimensioned< scalar > mag(const dimensioned< Type > &)
static void resizeField(UList< Type > &f, const label n)
const doubleScalar e
Elementary charge.
label n_
Size of the ODESystem (adjustable)
dimensionedScalar log10(const dimensionedScalar &ds)
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.