TDACChemistryModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "TDACChemistryModel.H"
27 #include "UniformField.H"
28 #include "localEulerDdtScheme.H"
29 #include "clockTime.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class ThermoType>
35 (
36  const fluidReactionThermo& thermo
37 )
38 :
40  timeSteps_(0),
41  NsDAC_(this->nSpecie_),
42  completeC_(this->nSpecie_, 0),
43  reactionsDisabled_(this->reactions_.size(), false),
44  completeToSimplifiedIndex_(this->nSpecie_, -1),
45  simplifiedToCompleteIndex_(this->nSpecie_),
46  tabulationResults_
47  (
48  IOobject
49  (
50  thermo.phasePropertyName("TabulationResults"),
51  this->time().timeName(),
52  this->mesh(),
53  IOobject::NO_READ,
54  IOobject::AUTO_WRITE
55  ),
56  this->mesh(),
57  scalar(0)
58  )
59 {
60  const basicSpecieMixture& composition = this->thermo().composition();
61 
63  (
64  *this,
65  *this
66  );
67 
68  // When the mechanism reduction method is used, the 'active' flag for every
69  // species should be initialised (by default 'active' is true)
70  if (mechRed_->active())
71  {
72  forAll(this->Y(), i)
73  {
74  IOobject header
75  (
76  this->Y()[i].name(),
77  this->mesh().time().timeName(),
78  this->mesh(),
79  IOobject::NO_READ
80  );
81 
82  // Check if the species file is provided, if not set inactive
83  // and NO_WRITE
84  if (!header.typeHeaderOk<volScalarField>(true))
85  {
86  composition.setInactive(i);
87  }
88  }
89  }
90 
92  (
93  *this,
94  *this
95  );
96 
97  if (mechRed_->log())
98  {
99  cpuReduceFile_ = logFile("cpu_reduce.out");
100  nActiveSpeciesFile_ = logFile("nActiveSpecies.out");
101  }
102 
103  if (tabulation_->log())
104  {
105  cpuAddFile_ = logFile("cpu_add.out");
106  cpuGrowFile_ = logFile("cpu_grow.out");
107  cpuRetrieveFile_ = logFile("cpu_retrieve.out");
108  }
109 
110  if (mechRed_->log() || tabulation_->log())
111  {
112  cpuSolveFile_ = logFile("cpu_solve.out");
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
118 
119 template<class ThermoType>
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
126 template<class ThermoType>
128 (
129  const scalar p,
130  const scalar T,
131  const scalarField& c, // Contains all species even when mechRed is active
132  const label li,
133  scalarField& dcdt
134 ) const
135 {
136  const bool reduced = mechRed_->active();
137 
138  scalar omegaf, omegar;
139 
140  dcdt = Zero;
141 
142  forAll(this->reactions_, i)
143  {
144  if (!reactionsDisabled_[i])
145  {
146  const Reaction<ThermoType>& R = this->reactions_[i];
147  const scalar omegaI = R.omega(p, T, c, li, omegaf, omegar);
148 
149  forAll(R.lhs(), s)
150  {
151  const label si =
152  reduced
153  ? completeToSimplifiedIndex_[R.lhs()[s].index]
154  : R.lhs()[s].index;
155  const scalar sl = R.lhs()[s].stoichCoeff;
156  dcdt[si] -= sl*omegaI;
157  }
158 
159  forAll(R.rhs(), s)
160  {
161  const label si =
162  reduced
163  ? completeToSimplifiedIndex_[R.rhs()[s].index]
164  : R.rhs()[s].index;
165  const scalar sr = R.rhs()[s].stoichCoeff;
166  dcdt[si] += sr*omegaI;
167  }
168  }
169  }
170 }
171 
172 
173 template<class ThermoType>
175 (
176  const scalar time,
177  const scalarField& c,
178  const label li,
179  scalarField& dcdt
180 ) const
181 {
182  const bool reduced = mechRed_->active();
183 
184  if (reduced)
185  {
186  // When using DAC, the ODE solver submit a reduced set of species the
187  // complete set is used and only the species in the simplified mechanism
188  // are updated
189  this->c_ = completeC_;
190 
191  // Update the concentration of the species in the simplified mechanism
192  // the other species remain the same and are used only for third-body
193  // efficiencies
194  for (label i=0; i<NsDAC_; i++)
195  {
196  this->c_[simplifiedToCompleteIndex_[i]] = max(c[i], 0);
197  }
198  }
199  else
200  {
201  for (label i=0; i<this->nSpecie(); i++)
202  {
203  this->c_[i] = max(c[i], 0);
204  }
205  }
206 
207  const scalar T = c[this->nSpecie_];
208  const scalar p = c[this->nSpecie_ + 1];
209 
210  dcdt = Zero;
211 
212  // Evaluate contributions from reactions
213  omega(p, T, this->c_, li, dcdt);
214 
215  // Evaluate the effect on the thermodynamic system ...
216 
217  // c*Cp
218  scalar ccp = 0;
219  for (label i=0; i<this->c_.size(); i++)
220  {
221  ccp += this->c_[i]*this->specieThermos_[i].cp(p, T);
222  }
223 
224  // dT/dt
225  scalar& dTdt = dcdt[this->nSpecie_];
226  for (label i=0; i<this->nSpecie_; i++)
227  {
228  const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
229  dTdt -= dcdt[i]*this->specieThermos_[si].ha(p, T);
230  }
231  dTdt /= ccp;
232 
233  // dp/dt = 0 (pressure is assumed constant)
234 }
235 
236 
237 template<class ThermoType>
239 (
240  const scalar t,
241  const scalarField& c,
242  const label li,
243  scalarField& dcdt,
245 ) const
246 {
247  const bool reduced = mechRed_->active();
248 
249  // If the mechanism reduction is active, the computed Jacobian
250  // is compact (size of the reduced set of species)
251  // but according to the information of the complete set
252  // (i.e. for the third-body efficiencies)
253 
254  if (reduced)
255  {
256  this->c_ = completeC_;
257 
258  for (label i=0; i<NsDAC_; i++)
259  {
260  this->c_[simplifiedToCompleteIndex_[i]] = max(c[i], 0);
261  }
262  }
263  else
264  {
265  forAll(this->c_, i)
266  {
267  this->c_[i] = max(c[i], 0);
268  }
269  }
270 
271  const scalar T = c[this->nSpecie_];
272  const scalar p = c[this->nSpecie_ + 1];
273 
274  dcdt = Zero;
275  J = Zero;
276 
277  // Evaluate contributions from reactions
278  forAll(this->reactions_, ri)
279  {
280  if (!reactionsDisabled_[ri])
281  {
282  const Reaction<ThermoType>& R = this->reactions_[ri];
283  scalar omegaI, kfwd, kbwd;
284  R.dwdc
285  (
286  p,
287  T,
288  this->c_,
289  li,
290  J,
291  dcdt,
292  omegaI,
293  kfwd,
294  kbwd,
295  reduced,
296  completeToSimplifiedIndex_
297  );
298  R.dwdT
299  (
300  p,
301  T,
302  this->c_,
303  li,
304  omegaI,
305  kfwd,
306  kbwd,
307  J,
308  reduced,
309  completeToSimplifiedIndex_,
310  this->nSpecie_
311  );
312  }
313  }
314 
315  // Evaluate the effect on the thermodynamic system ...
316 
317  // c*Cp
318  scalar ccp = 0, dccpdT = 0;
319  forAll(this->c_, i)
320  {
321  ccp += this->c_[i]*this->specieThermos_[i].cp(p, T);
322  dccpdT += this->c_[i]*this->specieThermos_[i].dcpdT(p, T);
323  }
324 
325  // dT/dt
326  scalar& dTdt = dcdt[this->nSpecie_];
327  for (label i=0; i<this->nSpecie_; i++)
328  {
329  const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
330  dTdt -= dcdt[i]*this->specieThermos_[si].ha(p, T);
331  }
332  dTdt /= ccp;
333 
334  // dp/dt = 0 (pressure is assumed constant)
335 
336  // d(dTdt)/dc
337  for (label i = 0; i < this->nSpecie_; i++)
338  {
339  scalar& d2Tdtdci = J(this->nSpecie_, i);
340  for (label j = 0; j < this->nSpecie_; j++)
341  {
342  const scalar d2cjdtdci = J(j, i);
343  const label sj = reduced ? simplifiedToCompleteIndex_[j] : j;
344  d2Tdtdci -= d2cjdtdci*this->specieThermos_[sj].ha(p, T);
345  }
346  const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
347  d2Tdtdci -= this->specieThermos_[si].cp(p, T)*dTdt;
348  d2Tdtdci /= ccp;
349  }
350 
351  // d(dTdt)/dT
352  scalar& d2TdtdT = J(this->nSpecie_, this->nSpecie_);
353  for (label i = 0; i < this->nSpecie_; i++)
354  {
355  const scalar d2cidtdT = J(i, this->nSpecie_);
356  const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
357  d2TdtdT -=
358  dcdt[i]*this->specieThermos_[si].cp(p, T)
359  + d2cidtdT*this->specieThermos_[si].ha(p, T);
360  }
361  d2TdtdT -= dTdt*dccpdT;
362  d2TdtdT /= ccp;
363 
364  // d(dpdt)/dc = 0 (pressure is assumed constant)
365 
366  // d(dpdt)/dT = 0 (pressure is assumed constant)
367 }
368 
369 
370 template<class ThermoType>
371 template<class DeltaTType>
373 (
374  const DeltaTType& deltaT
375 )
376 {
377  // Increment counter of time-step
378  timeSteps_++;
379 
380  const bool reduced = mechRed_->active();
381 
382  const basicSpecieMixture& composition = this->thermo().composition();
383 
384  // CPU time analysis
385  const clockTime clockTime_= clockTime();
386  clockTime_.timeIncrement();
387  scalar reduceMechCpuTime_ = 0;
388  scalar addNewLeafCpuTime_ = 0;
389  scalar growCpuTime_ = 0;
390  scalar solveChemistryCpuTime_ = 0;
391  scalar searchISATCpuTime_ = 0;
392 
393  this->resetTabulationResults();
394 
395  // Average number of active species
396  scalar nActiveSpecies = 0;
397  scalar nAvg = 0;
398 
400 
401  scalar deltaTMin = great;
402 
403  if (!this->chemistry_)
404  {
405  return deltaTMin;
406  }
407 
408  tmp<volScalarField> trho(this->thermo().rho0());
409  const scalarField& rho = trho();
410 
411  const scalarField& T = this->thermo().T().oldTime();
412  const scalarField& p = this->thermo().p().oldTime();
413 
414  scalarField c(this->nSpecie_);
415  scalarField c0(this->nSpecie_);
416 
417  // Composition vector (Yi, T, p, deltaT)
418  scalarField phiq(this->nEqns() + 1);
419  scalarField Rphiq(this->nEqns() + 1);
420 
421  forAll(rho, celli)
422  {
423  const scalar rhoi = rho[celli];
424  scalar pi = p[celli];
425  scalar Ti = T[celli];
426 
427  for (label i=0; i<this->nSpecie_; i++)
428  {
429  c[i] = c0[i] =
430  rhoi*this->Y_[i].oldTime()[celli]/this->specieThermos_[i].W();
431  phiq[i] = this->Y()[i].oldTime()[celli];
432  }
433  phiq[this->nSpecie()] = Ti;
434  phiq[this->nSpecie() + 1] = pi;
435  phiq[this->nSpecie() + 2] = deltaT[celli];
436 
437  // Initialise time progress
438  scalar timeLeft = deltaT[celli];
439 
440  // Not sure if this is necessary
441  Rphiq = Zero;
442 
443  clockTime_.timeIncrement();
444 
445  // When tabulation is active (short-circuit evaluation for retrieve)
446  // It first tries to retrieve the solution of the system with the
447  // information stored through the tabulation method
448  if (tabulation_->active() && tabulation_->retrieve(phiq, Rphiq))
449  {
450  // Retrieved solution stored in Rphiq
451  for (label i=0; i<this->nSpecie(); i++)
452  {
453  c[i] = rhoi*Rphiq[i]/this->specieThermos_[i].W();
454  }
455 
456  searchISATCpuTime_ += clockTime_.timeIncrement();
457  }
458  // This position is reached when tabulation is not used OR
459  // if the solution is not retrieved.
460  // In the latter case, it adds the information to the tabulation
461  // (it will either expand the current data or add a new stored point).
462  else
463  {
464  // Reset the time
465  clockTime_.timeIncrement();
466 
467  if (reduced)
468  {
469  // Reduce mechanism change the number of species (only active)
470  mechRed_->reduceMechanism(pi, Ti, c, celli);
471  nActiveSpecies += mechRed_->NsSimp();
472  nAvg++;
473  reduceMechCpuTime_ += clockTime_.timeIncrement();
474  }
475 
476  // Calculate the chemical source terms
477  while (timeLeft > small)
478  {
479  scalar dt = timeLeft;
480  if (reduced)
481  {
482  // completeC_ used in the overridden ODE methods
483  // to update only the active species
484  completeC_ = c;
485 
486  // Solve the reduced set of ODE
487  this->solve
488  (
489  pi,
490  Ti,
491  simplifiedC_,
492  celli,
493  dt,
494  this->deltaTChem_[celli]
495  );
496 
497  for (label i=0; i<NsDAC_; i++)
498  {
499  c[simplifiedToCompleteIndex_[i]] = simplifiedC_[i];
500  }
501  }
502  else
503  {
504  this->solve(pi, Ti, c, celli, dt, this->deltaTChem_[celli]);
505  }
506  timeLeft -= dt;
507  }
508 
509  {
510  solveChemistryCpuTime_ += clockTime_.timeIncrement();
511  }
512 
513  // If tabulation is used, we add the information computed here to
514  // the stored points (either expand or add)
515  if (tabulation_->active())
516  {
517  forAll(c, i)
518  {
519  Rphiq[i] = c[i]/rhoi*this->specieThermos_[i].W();
520  }
521  Rphiq[Rphiq.size()-3] = Ti;
522  Rphiq[Rphiq.size()-2] = pi;
523  Rphiq[Rphiq.size()-1] = deltaT[celli];
524 
525  label growOrAdd =
526  tabulation_->add(phiq, Rphiq, celli, rhoi, deltaT[celli]);
527  if (growOrAdd)
528  {
529  this->setTabulationResultsAdd(celli);
530  addNewLeafCpuTime_ += clockTime_.timeIncrement();
531  }
532  else
533  {
534  this->setTabulationResultsGrow(celli);
535  growCpuTime_ += clockTime_.timeIncrement();
536  }
537  }
538 
539  // When operations are done and if mechanism reduction is active,
540  // the number of species (which also affects nEqns) is set back
541  // to the total number of species (stored in the mechRed object)
542  if (reduced)
543  {
544  this->nSpecie_ = mechRed_->nSpecie();
545  }
546  deltaTMin = min(this->deltaTChem_[celli], deltaTMin);
547 
548  this->deltaTChem_[celli] =
549  min(this->deltaTChem_[celli], this->deltaTChemMax_);
550  }
551 
552  // Set the RR vector (used in the solver)
553  for (label i=0; i<this->nSpecie_; i++)
554  {
555  this->RR_[i][celli] =
556  (c[i] - c0[i])*this->specieThermos_[i].W()/deltaT[celli];
557  }
558  }
559 
560  if (mechRed_->log() || tabulation_->log())
561  {
562  cpuSolveFile_()
563  << this->time().timeOutputValue()
564  << " " << solveChemistryCpuTime_ << endl;
565  }
566 
567  if (mechRed_->log())
568  {
569  cpuReduceFile_()
570  << this->time().timeOutputValue()
571  << " " << reduceMechCpuTime_ << endl;
572  }
573 
574  if (tabulation_->active())
575  {
576  // Every time-step, look if the tabulation should be updated
577  tabulation_->update();
578 
579  // Write the performance of the tabulation
580  tabulation_->writePerformance();
581 
582  if (tabulation_->log())
583  {
584  cpuRetrieveFile_()
585  << this->time().timeOutputValue()
586  << " " << searchISATCpuTime_ << endl;
587 
588  cpuGrowFile_()
589  << this->time().timeOutputValue()
590  << " " << growCpuTime_ << endl;
591 
592  cpuAddFile_()
593  << this->time().timeOutputValue()
594  << " " << addNewLeafCpuTime_ << endl;
595  }
596  }
597 
598  if (reduced && nAvg && mechRed_->log())
599  {
600  // Write average number of species
601  nActiveSpeciesFile_()
602  << this->time().timeOutputValue()
603  << " " << nActiveSpecies/nAvg << endl;
604  }
605 
606  if (Pstream::parRun())
607  {
608  List<bool> active(composition.active());
609  Pstream::listCombineGather(active, orEqOp<bool>());
610  Pstream::listCombineScatter(active);
611 
612  forAll(active, i)
613  {
614  if (active[i])
615  {
616  composition.setActive(i);
617  }
618  }
619  }
620 
621  return deltaTMin;
622 }
623 
624 
625 template<class ThermoType>
627 (
628  const scalar deltaT
629 )
630 {
631  // Don't allow the time-step to change more than a factor of 2
632  return min
633  (
635  2*deltaT
636  );
637 }
638 
639 
640 template<class ThermoType>
642 (
643  const scalarField& deltaT
644 )
645 {
646  return this->solve<scalarField>(deltaT);
647 }
648 
649 
650 template<class ThermoType>
653 (
654  const label celli
655 )
656 {
657  tabulationResults_[celli] = 0;
658 }
659 
660 
661 template<class ThermoType>
664 {
665  tabulationResults_[celli] = 1;
666 }
667 
668 
669 template<class ThermoType>
672 (
673  const label celli
674 )
675 {
676  tabulationResults_[celli] = 2;
677 }
678 
679 
680 // ************************************************************************* //
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:47
void setInactive(label speciei) const
Set speciei inactive.
Extends standardChemistryModel by adding the TDAC method.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
fluidReactionThermo & thermo
Definition: createFields.H:28
void setTabulationResultsGrow(const label celli)
static word phasePropertyName(const word &name, const word &phaseName)
Return the name of a property for a given phase.
Definition: basicThermo.H:156
void setTabulationResultsRetrieve(const label celli)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool active(label speciei) const
Return true for active species.
basicSpecieMixture & composition
void dwdc(const scalar p, const scalar T, const scalarField &c, const label li, scalarSquareMatrix &J, scalarField &dcdt, scalar &omegaI, scalar &kfwd, scalar &kbwd, const bool reduced, const List< label > &c2s) const
Derivative of the net reaction rate for each species involved.
Definition: Reaction.C:352
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< volScalarField > trho
virtual ~TDACChemistryModel()
Destructor.
Base-class for multi-component fluid thermodynamic properties.
scalar rho0
const label nSpecie
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: reactionI.H:36
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
const dimensionedScalar c
Speed of light in a vacuum.
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
Net reaction rate for individual species.
Definition: Reaction.C:290
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
virtual void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
dc/dt = omega, rate of change in concentration, for each species
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
TDACChemistryModel(const fluidReactionThermo &thermo)
Construct from thermo.
void setTabulationResultsAdd(const label celli)
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
bool active(const label i) const
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: reactionI.H:42
rhoEqn solve()
#define R(A, B, C, D, E, F, K, M)
const scalarList W(::W(thermo))
PtrList< volScalarField > & Y
virtual void derivatives(const scalar t, const scalarField &c, const label li, scalarField &dcdt) const
Calculate the derivatives in dydx.
double timeIncrement() const
Return time (in seconds) since last call to timeIncrement()
Definition: clockTime.C:62
void setActive(label speciei) const
Set speciei active.
void dwdT(const scalar p, const scalar T, const scalarField &c, const label li, const scalar omegaI, const scalar kfwd, const scalar kbwd, scalarSquareMatrix &J, const bool reduced, const List< label > &c2s, const label indexT) const
Derivative of the net reaction rate for each species involved.
Definition: Reaction.C:499
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void jacobian(const scalar t, const scalarField &c, const label li, scalarField &dcdt, scalarSquareMatrix &J) const
Calculate the Jacobian of the system.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Starts timing (using rtc) and returns elapsed time from start. Better resolution (2uSec instead of ~2...
Definition: clockTime.H:49