35 Foam::scalar Foam::chemPointISAT::tolerance_;
39 void Foam::chemPointISAT::qrDecompose
52 for (
label i=
k; i<nCols; i++)
62 for (
label i=
k; i<nCols; i++)
67 for (
label i=
k; i<nCols; i++)
75 for (
label j=
k+1; j<nCols; j++)
78 for (
label i=
k; i<nCols; i++)
82 scalar tau =
sum/
c[
k];
83 for (
label i=
k; i<nCols; i++)
85 R(i, j) -= tau*
R(i,
k);
91 d[nCols-1] =
R(nCols-1, nCols-1);
94 for (
label i=0; i<nCols; i++)
97 for (
label j=0; j<i; j++)
105 void Foam::chemPointISAT::qrUpdate
116 for (
k=
n-1;
k>=0;
k--)
129 for (
label i=
k-1; i>=0; i--)
131 rotate(
R, i, w[i], -w[i+1],
n);
136 else if (
mag(w[i]) >
mag(w[i+1]))
138 w[i] =
mag(w[i])*
sqrt(1.0 +
sqr(w[i+1]/w[i]));
142 w[i] =
mag(w[i+1])*
sqrt(1.0 +
sqr(w[i]/w[i+1]));
148 R(0, i) += w[0]*v[i];
153 rotate(
R, i,
R(i, i), -
R(i+1, i),
n);
158 void Foam::chemPointISAT::rotate
167 scalar
c, fact,
s, w,
y;
172 s = (
b >= 0 ? 1 : -1);
206 const scalar tolerance,
207 const label completeSpaceSize,
217 scaleFactor_(scaleFactor),
219 completeSpaceSize_(completeSpaceSize),
222 timeTag_(table.timeSteps()),
223 lastTimeUsed_(table.timeSteps()),
225 maxNumNewDim_(coeffsDict.lookupOrDefault(
"maxNumNewDim",0)),
226 printProportion_(coeffsDict.lookupOrDefault(
"printProportion",false)),
233 scaleFactor_[iddeltaT_] *= phi_[iddeltaT_]/tolerance_;
240 simplifiedToCompleteIndex_.
setSize(nActive_);
243 forAll(completeToSimplifiedIndex_, i)
248 forAll(simplifiedToCompleteIndex_, i)
254 const label reduOrCompDim =
264 for (
label i=0; i<reduOrCompDim; i++)
266 D[i] =
max(
S[i], 0.5);
275 for (
label i=0; i<reduOrCompDim; i++)
277 for (
label j=0; j<reduOrCompDim; j++)
299 qrDecompose(reduOrCompDim, LT_);
313 scaleFactor_(
p.scaleFactor()),
315 completeSpaceSize_(
p.completeSpaceSize()),
316 nGrowth_(
p.nGrowth()),
317 nActive_(
p.nActive()),
318 simplifiedToCompleteIndex_(
p.simplifiedToCompleteIndex()),
319 timeTag_(
p.timeTag()),
320 lastTimeUsed_(
p.lastTimeUsed()),
321 toRemove_(
p.toRemove()),
322 maxNumNewDim_(
p.maxNumNewDim()),
325 completeToSimplifiedIndex_(
p.completeToSimplifiedIndex())
327 tolerance_ =
p.tolerance();
342 table_.reduction() ? nActive_ : completeSpaceSize() - 3;
347 for (
label i=0; i<completeSpaceSize()-3; i++)
356 !(table_.reduction())
357 ||(table_.reduction() && completeToSimplifiedIndex_[i] != -1)
362 ? completeToSimplifiedIndex_[i]
365 for (
label j=si; j<dim; j++)
369 ? simplifiedToCompleteIndex_[j]
372 temp += LT_(si, j)*dphi[sj];
375 temp += LT_(si, dim)*dphi[idT_];
376 temp += LT_(si, dim+1)*dphi[idp_];
377 temp += LT_(si, dim+2)*dphi[iddeltaT_];
381 temp = dphi[i]/(tolerance_*scaleFactor_[i]);
384 epsTemp +=
sqr(temp);
386 if (printProportion_)
396 LT_(dim, dim)*dphi[idT_]
397 +LT_(dim, dim+1)*dphi[idp_]
398 +LT_(dim, dim+2)*dphi[iddeltaT_]
405 LT_(dim+1, dim+1)*dphi[idp_]
406 +LT_(dim+1, dim+2)*dphi[iddeltaT_]
409 epsTemp +=
sqr(LT_[dim+2][dim+2]*dphi[iddeltaT_]);
411 if (printProportion_)
415 LT_(dim, dim)*dphi[idT_]
416 + LT_(dim, dim+1)*dphi[idp_]
420 sqr(LT_(dim+1, dim+1)*dphi[idp_]);
423 sqr(LT_[dim+2][dim+2]*dphi[iddeltaT_]);
426 if (
sqrt(epsTemp) > 1 + tolerance_)
428 if (printProportion_)
432 for (
label i=0; i<completeSpaceSize(); i++)
441 if (maxIndex >= completeSpaceSize() - 3)
443 if (maxIndex == idT_)
447 else if (maxIndex == idp_)
451 else if (maxIndex == iddeltaT_)
458 propName = table_.chemistry().Y()[maxIndex].member();
461 Info<<
"Direction maximum impact to error in ellipsoid: "
464 Info<<
"Proportion to the total error on the retrieve: "
465 <<
max/(epsTemp+small) <<
endl;
490 table_.reduction() ? nActive_ : completeSpaceSize() - 2;
494 for (
label i=0; i<completeSpaceSize()-3; i++)
497 if (table_.reduction())
499 const label si = completeToSimplifiedIndex_[i];
504 for (
label j=0; j<dim; j++)
506 const label sj = simplifiedToCompleteIndex_[j];
507 dRl += Avar(si, j)*dphi[sj];
509 dRl += Avar(si, nActive_)*dphi[idT_];
510 dRl += Avar(si, nActive_+1)*dphi[idp_];
511 dRl += Avar(si, nActive_+2)*dphi[iddeltaT_];
520 for (
label j=0; j<completeSpaceSize(); j++)
522 dRl += Avar(i, j)*dphi[j];
525 eps2 +=
sqr((dR[i]-dRl)/scaleFactorV[i]);
529 if (eps2 > tolerance())
544 const label initNActiveSpecies(nActive_);
546 if (table_.reduction())
548 label activeAdded(0);
553 for (
label i=0; i<completeSpaceSize()-3; i++)
559 completeToSimplifiedIndex_[i] == -1
560 && table_.chemistry().cTos(i) != -1
570 completeToSimplifiedIndex_[i] != -1
571 && table_.chemistry().cTos(i) == -1
582 completeToSimplifiedIndex_[i] == -1
583 && table_.chemistry().cTos(i) == -1
593 if (activeAdded > maxNumNewDim_)
599 nActive_ += dimToAdd.
size();
600 simplifiedToCompleteIndex_.setSize(nActive_);
603 label si = nActive_ - dimToAdd.
size() + i;
605 simplifiedToCompleteIndex_[si] = dimToAdd[i];
606 completeToSimplifiedIndex_[dimToAdd[i]] = si;
615 if (nActive_ > initNActiveSpecies)
623 for (
label i=0; i<initNActiveSpecies; i++)
625 for (
label j=0; j<initNActiveSpecies; j++)
627 LT_(i, j) = LTvar(i, j);
628 A_(i, j) = Avar(i, j);
633 for (
label i=0; i<initNActiveSpecies; i++)
635 for (
label j=1; j>=0; j--)
637 LT_(i, nActive_+j)=LTvar(i, initNActiveSpecies+j);
638 A_(i, nActive_+j)=Avar(i, initNActiveSpecies+j);
639 LT_(nActive_+j, i)=LTvar(initNActiveSpecies+j, i);
640 A_(nActive_+j, i)=Avar(initNActiveSpecies+j, i);
644 LT_(nActive_, nActive_)=
645 LTvar(initNActiveSpecies, initNActiveSpecies);
646 A_(nActive_, nActive_)=
647 Avar(initNActiveSpecies, initNActiveSpecies);
648 LT_(nActive_+1, nActive_+1)=
649 LTvar(initNActiveSpecies+1, initNActiveSpecies+1);
650 A_(nActive_+1, nActive_+1)=
651 Avar(initNActiveSpecies+1, initNActiveSpecies+1);
652 LT_(nActive_+2, nActive_+2)=
653 LTvar(initNActiveSpecies+2, initNActiveSpecies+2);
654 A_(nActive_+2, nActive_+2)=
655 Avar(initNActiveSpecies+2, initNActiveSpecies+2);
657 for (
label i=initNActiveSpecies; i<nActive_;i++)
661 /(tolerance_*scaleFactor_[simplifiedToCompleteIndex_[i]]);
668 table_.reduction() ? nActive_ + 3 : completeSpaceSize();
672 scalar normPhiTilde = 0;
675 for (
label i=0; i<dim; i++)
677 for (
label j=i; j<dim-3; j++)
681 ? simplifiedToCompleteIndex_[j]
684 phiTilde[i] += LT_(i, j)*dphi[sj];
687 phiTilde[i] += LT_(i, dim-3)*dphi[idT_];
688 phiTilde[i] += LT_(i, dim-3+1)*dphi[idp_];
689 phiTilde[i] += LT_(i, dim-3+2)*dphi[iddeltaT_];
691 normPhiTilde +=
sqr(phiTilde[i]);
694 const scalar invSqrNormPhiTilde = 1.0/normPhiTilde;
695 normPhiTilde =
sqrt(normPhiTilde);
698 const scalar gamma = (1/normPhiTilde - 1)*invSqrNormPhiTilde;
702 for (
label i=0; i<dim; i++)
704 for (
label j=0; j<=i;j++)
706 v[i] += phiTilde[j]*LT_(j, i);
710 qrUpdate(LT_,dim, u, v);
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
#define R(A, B, C, D, E, F, K, M)
#define forAll(list, i)
Loop across all elements in list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
Form T() const
Return the transpose of the matrix.
Singular value decomposition of a rectangular matrix.
const scalarRectangularMatrix & V() const
Return the square matrix V.
const scalarRectangularMatrix & U() const
Return U.
const scalarDiagonalMatrix & S() const
Return the singular values.
Leaf of the binary tree. The chemPoint stores the composition 'phi', the mapping of this composition ...
bool inEOA(const scalarField &phiq)
To RETRIEVE the mapping from the stored chemPoint phi, the query.
chemPointISAT(chemistryTabulationMethods::ISAT &table, const scalarField &phi, const scalarField &Rphi, const scalarSquareMatrix &A, const scalarField &scaleFactor, const scalar tolerance, const label completeSpaceSize, const label nActive, const dictionary &coeffsDict, binaryNode *node=nullptr)
Construct from components.
const scalarField & scaleFactor() const
bool grow(const scalarField &phiq)
More details about the minimum-volume ellipsoid covering an.
bool checkSolution(const scalarField &phiq, const scalarField &Rphiq)
If phiq is not in the EOA, then the mapping is computed.
const scalarSquareMatrix & A() const
const scalar & tolerance() const
label & completeSpaceSize()
const List< label > & simplifiedToCompleteIndex() const
Implementation of the ISAT (In-situ adaptive tabulation), for chemistry calculation.
const odeChemistryModel & chemistry()
bool reduction() const
Return true if reduction is applied to the state variables.
A list of keyword definitions, which are a keyword followed by any number of values (e....
label sToc(const label si) const
Return the index in the complete set of species.
label cTos(const label ci) const
Return the index in the simplified set of species.
A class for handling words, derived from string.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar c
Speed of light in a vacuum.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
SquareMatrix< scalar > scalarSquareMatrix
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)