FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

ODEChemistryModel.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2009-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "ODEChemistryModel.H"
00027 #include <chemistryModel/chemistrySolver.H>
00028 #include <reactionThermophysicalModels/reactingMixture.H>
00029 
00030 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00031 
00032 template<class CompType, class ThermoType>
00033 Foam::ODEChemistryModel<CompType, ThermoType>::ODEChemistryModel
00034 (
00035     const fvMesh& mesh,
00036     const word& compTypeName,
00037     const word& thermoTypeName
00038 )
00039 :
00040     CompType(mesh, thermoTypeName),
00041 
00042     ODE(),
00043 
00044     Y_(this->thermo().composition().Y()),
00045 
00046     reactions_
00047     (
00048         dynamic_cast<const reactingMixture<ThermoType>&>(this->thermo())
00049     ),
00050     specieThermo_
00051     (
00052         dynamic_cast<const reactingMixture<ThermoType>&>
00053             (this->thermo()).speciesData()
00054     ),
00055 
00056     nSpecie_(Y_.size()),
00057     nReaction_(reactions_.size()),
00058 
00059     solver_
00060     (
00061         chemistrySolver<CompType, ThermoType>::New
00062         (
00063             *this,
00064             compTypeName,
00065             thermoTypeName
00066         )
00067     ),
00068 
00069     RR_(nSpecie_)
00070 {
00071     // create the fields for the chemistry sources
00072     forAll(RR_, fieldI)
00073     {
00074         RR_.set
00075         (
00076             fieldI,
00077             new scalarField(mesh.nCells(), 0.0)
00078         );
00079     }
00080 
00081     Info<< "ODEChemistryModel: Number of species = " << nSpecie_
00082         << " and reactions = " << nReaction_ << endl;
00083 }
00084 
00085 
00086 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00087 
00088 template<class CompType, class ThermoType>
00089 Foam::ODEChemistryModel<CompType, ThermoType>::~ODEChemistryModel()
00090 {}
00091 
00092 
00093 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00094 
00095 template<class CompType, class ThermoType>
00096 Foam::scalarField Foam::ODEChemistryModel<CompType, ThermoType>::omega
00097 (
00098     const scalarField& c,
00099     const scalar T,
00100     const scalar p
00101 ) const
00102 {
00103     scalar pf, cf, pr, cr;
00104     label lRef, rRef;
00105 
00106     scalarField om(nEqns(), 0.0);
00107 
00108     forAll(reactions_, i)
00109     {
00110         const Reaction<ThermoType>& R = reactions_[i];
00111 
00112         scalar omegai = omega
00113         (
00114             R, c, T, p, pf, cf, lRef, pr, cr, rRef
00115         );
00116 
00117         forAll(R.lhs(), s)
00118         {
00119             label si = R.lhs()[s].index;
00120             scalar sl = R.lhs()[s].stoichCoeff;
00121             om[si] -= sl*omegai;
00122         }
00123 
00124         forAll(R.rhs(), s)
00125         {
00126             label si = R.rhs()[s].index;
00127             scalar sr = R.rhs()[s].stoichCoeff;
00128             om[si] += sr*omegai;
00129         }
00130     }
00131 
00132     return om;
00133 }
00134 
00135 
00136 template<class CompType, class ThermoType>
00137 Foam::scalar Foam::ODEChemistryModel<CompType, ThermoType>::omega
00138 (
00139     const Reaction<ThermoType>& R,
00140     const scalarField& c,
00141     const scalar T,
00142     const scalar p,
00143     scalar& pf,
00144     scalar& cf,
00145     label& lRef,
00146     scalar& pr,
00147     scalar& cr,
00148     label& rRef
00149 ) const
00150 {
00151     scalarField c2(nSpecie_, 0.0);
00152     for (label i=0; i<nSpecie_; i++)
00153     {
00154         c2[i] = max(0.0, c[i]);
00155     }
00156 
00157     scalar kf = R.kf(T, p, c2);
00158     scalar kr = R.kr(kf, T, p, c2);
00159 
00160     pf = 1.0;
00161     pr = 1.0;
00162 
00163     label Nl = R.lhs().size();
00164     label Nr = R.rhs().size();
00165 
00166     label slRef = 0;
00167     lRef = R.lhs()[slRef].index;
00168 
00169     pf = kf;
00170     for (label s=1; s<Nl; s++)
00171     {
00172         label si = R.lhs()[s].index;
00173 
00174         if (c[si] < c[lRef])
00175         {
00176             scalar exp = R.lhs()[slRef].exponent;
00177             pf *= pow(max(0.0, c[lRef]), exp);
00178             lRef = si;
00179             slRef = s;
00180         }
00181         else
00182         {
00183             scalar exp = R.lhs()[s].exponent;
00184             pf *= pow(max(0.0, c[si]), exp);
00185         }
00186     }
00187     cf = max(0.0, c[lRef]);
00188 
00189     {
00190         scalar exp = R.lhs()[slRef].exponent;
00191         if (exp<1.0)
00192         {
00193             if (cf > SMALL)
00194             {
00195                 pf *= pow(cf, exp - 1.0);
00196             }
00197             else
00198             {
00199                 pf = 0.0;
00200             }
00201         }
00202         else
00203         {
00204             pf *= pow(cf, exp - 1.0);
00205         }
00206     }
00207 
00208     label srRef = 0;
00209     rRef = R.rhs()[srRef].index;
00210 
00211     // find the matrix element and element position for the rhs
00212     pr = kr;
00213     for (label s=1; s<Nr; s++)
00214     {
00215         label si = R.rhs()[s].index;
00216         if (c[si] < c[rRef])
00217         {
00218             scalar exp = R.rhs()[srRef].exponent;
00219             pr *= pow(max(0.0, c[rRef]), exp);
00220             rRef = si;
00221             srRef = s;
00222         }
00223         else
00224         {
00225             scalar exp = R.rhs()[s].exponent;
00226             pr *= pow(max(0.0, c[si]), exp);
00227         }
00228     }
00229     cr = max(0.0, c[rRef]);
00230 
00231     {
00232         scalar exp = R.rhs()[srRef].exponent;
00233         if (exp<1.0)
00234         {
00235             if (cr>SMALL)
00236             {
00237                 pr *= pow(cr, exp - 1.0);
00238             }
00239             else
00240             {
00241                 pr = 0.0;
00242             }
00243         }
00244         else
00245         {
00246             pr *= pow(cr, exp - 1.0);
00247         }
00248     }
00249 
00250     return pf*cf - pr*cr;
00251 }
00252 
00253 
00254 template<class CompType, class ThermoType>
00255 void Foam::ODEChemistryModel<CompType, ThermoType>::derivatives
00256 (
00257     const scalar time,
00258     const scalarField &c,
00259     scalarField& dcdt
00260 ) const
00261 {
00262     scalar T = c[nSpecie_];
00263     scalar p = c[nSpecie_ + 1];
00264 
00265     dcdt = omega(c, T, p);
00266 
00267     // constant pressure
00268     // dT/dt = ...
00269     scalar rho = 0.0;
00270     scalar cSum = 0.0;
00271     for (label i=0; i<nSpecie_; i++)
00272     {
00273         scalar W = specieThermo_[i].W();
00274         cSum += c[i];
00275         rho += W*c[i];
00276     }
00277     scalar mw = rho/cSum;
00278     scalar cp = 0.0;
00279     for (label i=0; i<nSpecie_; i++)
00280     {
00281         scalar cpi = specieThermo_[i].cp(T);
00282         scalar Xi = c[i]/rho;
00283         cp += Xi*cpi;
00284     }
00285     cp /= mw;
00286 
00287     scalar dT = 0.0;
00288     for (label i=0; i<nSpecie_; i++)
00289     {
00290         scalar hi = specieThermo_[i].h(T);
00291         dT += hi*dcdt[i];
00292     }
00293     dT /= rho*cp;
00294 
00295     // limit the time-derivative, this is more stable for the ODE
00296     // solver when calculating the allowed time step
00297     scalar dtMag = min(500.0, mag(dT));
00298     dcdt[nSpecie_] = -dT*dtMag/(mag(dT) + 1.0e-10);
00299 
00300     // dp/dt = ...
00301     dcdt[nSpecie_+1] = 0.0;
00302 }
00303 
00304 
00305 template<class CompType, class ThermoType>
00306 void Foam::ODEChemistryModel<CompType, ThermoType>::jacobian
00307 (
00308     const scalar t,
00309     const scalarField& c,
00310     scalarField& dcdt,
00311     scalarSquareMatrix& dfdc
00312 ) const
00313 {
00314     scalar T = c[nSpecie_];
00315     scalar p = c[nSpecie_ + 1];
00316 
00317     scalarField c2(nSpecie_, 0.0);
00318     for (label i=0; i<nSpecie_; i++)
00319     {
00320         c2[i] = max(c[i], 0.0);
00321     }
00322 
00323     for (label i=0; i<nEqns(); i++)
00324     {
00325         for (label j=0; j<nEqns(); j++)
00326         {
00327             dfdc[i][j] = 0.0;
00328         }
00329     }
00330 
00331     // length of the first argument must be nSpecie()
00332     dcdt = omega(c2, T, p);
00333 
00334     for (label ri=0; ri<reactions_.size(); ri++)
00335     {
00336         const Reaction<ThermoType>& R = reactions_[ri];
00337 
00338         scalar kf0 = R.kf(T, p, c2);
00339         scalar kr0 = R.kr(T, p, c2);
00340 
00341         forAll(R.lhs(), j)
00342         {
00343             label sj = R.lhs()[j].index;
00344             scalar kf = kf0;
00345             forAll(R.lhs(), i)
00346             {
00347                 label si = R.lhs()[i].index;
00348                 scalar el = R.lhs()[i].exponent;
00349                 if (i == j)
00350                 {
00351                     if (el < 1.0)
00352                     {
00353                         if (c2[si]>SMALL)
00354                         {
00355                             kf *= el*pow(c2[si] + VSMALL, el - 1.0);
00356                         }
00357                         else
00358                         {
00359                             kf = 0.0;
00360                         }
00361                     }
00362                     else
00363                     {
00364                         kf *= el*pow(c2[si], el - 1.0);
00365                     }
00366                 }
00367                 else
00368                 {
00369                     kf *= pow(c2[si], el);
00370                 }
00371             }
00372 
00373             forAll(R.lhs(), i)
00374             {
00375                 label si = R.lhs()[i].index;
00376                 scalar sl = R.lhs()[i].stoichCoeff;
00377                 dfdc[si][sj] -= sl*kf;
00378             }
00379             forAll(R.rhs(), i)
00380             {
00381                 label si = R.rhs()[i].index;
00382                 scalar sr = R.rhs()[i].stoichCoeff;
00383                 dfdc[si][sj] += sr*kf;
00384             }
00385         }
00386 
00387         forAll(R.rhs(), j)
00388         {
00389             label sj = R.rhs()[j].index;
00390             scalar kr = kr0;
00391             forAll(R.rhs(), i)
00392             {
00393                 label si = R.rhs()[i].index;
00394                 scalar er = R.rhs()[i].exponent;
00395                 if (i==j)
00396                 {
00397                     if (er<1.0)
00398                     {
00399                         if (c2[si]>SMALL)
00400                         {
00401                             kr *= er*pow(c2[si] + VSMALL, er - 1.0);
00402                         }
00403                         else
00404                         {
00405                             kr = 0.0;
00406                         }
00407                     }
00408                     else
00409                     {
00410                         kr *= er*pow(c2[si], er - 1.0);
00411                     }
00412                 }
00413                 else
00414                 {
00415                     kr *= pow(c2[si], er);
00416                 }
00417             }
00418 
00419             forAll(R.lhs(), i)
00420             {
00421                 label si = R.lhs()[i].index;
00422                 scalar sl = R.lhs()[i].stoichCoeff;
00423                 dfdc[si][sj] += sl*kr;
00424             }
00425             forAll(R.rhs(), i)
00426             {
00427                 label si = R.rhs()[i].index;
00428                 scalar sr = R.rhs()[i].stoichCoeff;
00429                 dfdc[si][sj] -= sr*kr;
00430             }
00431         }
00432     }
00433 
00434     // calculate the dcdT elements numerically
00435     scalar delta = 1.0e-8;
00436     scalarField dcdT0 = omega(c2, T - delta, p);
00437     scalarField dcdT1 = omega(c2, T + delta, p);
00438 
00439     for (label i=0; i<nEqns(); i++)
00440     {
00441         dfdc[i][nSpecie()] = 0.5*(dcdT1[i] - dcdT0[i])/delta;
00442     }
00443 
00444 }
00445 
00446 
00447 template<class CompType, class ThermoType>
00448 Foam::tmp<Foam::volScalarField>
00449 Foam::ODEChemistryModel<CompType, ThermoType>::tc() const
00450 {
00451     scalar pf, cf, pr, cr;
00452     label lRef, rRef;
00453 
00454     const volScalarField rho
00455     (
00456         IOobject
00457         (
00458             "rho",
00459             this->time().timeName(),
00460             this->mesh(),
00461             IOobject::NO_READ,
00462             IOobject::NO_WRITE,
00463             false
00464         ),
00465         this->thermo().rho()
00466     );
00467 
00468     tmp<volScalarField> ttc
00469     (
00470         new volScalarField
00471         (
00472             IOobject
00473             (
00474                 "tc",
00475                 this->time().timeName(),
00476                 this->mesh(),
00477                 IOobject::NO_READ,
00478                 IOobject::NO_WRITE
00479             ),
00480             this->mesh(),
00481             dimensionedScalar("zero", dimTime, SMALL),
00482             zeroGradientFvPatchScalarField::typeName
00483         )
00484     );
00485 
00486     scalarField& tc = ttc();
00487 
00488     label nReaction = reactions_.size();
00489 
00490     if (this->chemistry_)
00491     {
00492         forAll(rho, celli)
00493         {
00494             scalar rhoi = rho[celli];
00495             scalar Ti = this->thermo().T()[celli];
00496             scalar pi = this->thermo().p()[celli];
00497             scalarField c(nSpecie_);
00498             scalar cSum = 0.0;
00499 
00500             for (label i=0; i<nSpecie_; i++)
00501             {
00502                 scalar Yi = Y_[i][celli];
00503                 c[i] = rhoi*Yi/specieThermo_[i].W();
00504                 cSum += c[i];
00505             }
00506 
00507             forAll(reactions_, i)
00508             {
00509                 const Reaction<ThermoType>& R = reactions_[i];
00510 
00511                 omega
00512                 (
00513                     R, c, Ti, pi, pf, cf, lRef, pr, cr, rRef
00514                 );
00515 
00516                 forAll(R.rhs(), s)
00517                 {
00518                     scalar sr = R.rhs()[s].stoichCoeff;
00519                     tc[celli] += sr*pf*cf;
00520                 }
00521             }
00522             tc[celli] = nReaction*cSum/tc[celli];
00523         }
00524     }
00525 
00526 
00527     ttc().correctBoundaryConditions();
00528 
00529     return ttc;
00530 }
00531 
00532 
00533 template<class CompType, class ThermoType>
00534 Foam::tmp<Foam::volScalarField>
00535 Foam::ODEChemistryModel<CompType, ThermoType>::Sh() const
00536 {
00537     tmp<volScalarField> tSh
00538     (
00539         new volScalarField
00540         (
00541             IOobject
00542             (
00543                 "Sh",
00544                 this->mesh_.time().timeName(),
00545                 this->mesh_,
00546                 IOobject::NO_READ,
00547                 IOobject::NO_WRITE,
00548                 false
00549             ),
00550             this->mesh_,
00551             dimensionedScalar("zero", dimEnergy/dimTime/dimVolume, 0.0),
00552             zeroGradientFvPatchScalarField::typeName
00553         )
00554     );
00555 
00556     if (this->chemistry_)
00557     {
00558         scalarField& Sh = tSh();
00559 
00560         forAll(Y_, i)
00561         {
00562             forAll(Sh, cellI)
00563             {
00564                 scalar hi = specieThermo_[i].Hc();
00565                 Sh[cellI] -= hi*RR_[i][cellI];
00566             }
00567         }
00568     }
00569 
00570     return tSh;
00571 }
00572 
00573 
00574 template<class CompType, class ThermoType>
00575 Foam::tmp<Foam::volScalarField>
00576 Foam::ODEChemistryModel<CompType, ThermoType>::dQ() const
00577 {
00578     tmp<volScalarField> tdQ
00579     (
00580         new volScalarField
00581         (
00582             IOobject
00583             (
00584                 "dQ",
00585                 this->mesh_.time().timeName(),
00586                 this->mesh_,
00587                 IOobject::NO_READ,
00588                 IOobject::NO_WRITE,
00589                 false
00590             ),
00591             this->mesh_,
00592             dimensionedScalar("dQ", dimEnergy/dimTime, 0.0),
00593             zeroGradientFvPatchScalarField::typeName
00594         )
00595     );
00596 
00597     if (this->chemistry_)
00598     {
00599         volScalarField& dQ = tdQ();
00600         dQ.dimensionedInternalField() = this->mesh_.V()*Sh()();
00601     }
00602 
00603     return tdQ;
00604 }
00605 
00606 
00607 template<class CompType, class ThermoType>
00608 Foam::label Foam::ODEChemistryModel<CompType, ThermoType>::nEqns() const
00609 {
00610     // nEqns = number of species + temperature + pressure
00611     return nSpecie_ + 2;
00612 }
00613 
00614 
00615 template<class CompType, class ThermoType>
00616 void Foam::ODEChemistryModel<CompType, ThermoType>::calculate()
00617 {
00618     const volScalarField rho
00619     (
00620         IOobject
00621         (
00622             "rho",
00623             this->time().timeName(),
00624             this->mesh(),
00625             IOobject::NO_READ,
00626             IOobject::NO_WRITE,
00627             false
00628         ),
00629         this->thermo().rho()
00630     );
00631 
00632     for (label i=0; i<nSpecie_; i++)
00633     {
00634         RR_[i].setSize(rho.size());
00635     }
00636 
00637     if (this->chemistry_)
00638     {
00639         forAll(rho, celli)
00640         {
00641             for (label i=0; i<nSpecie_; i++)
00642             {
00643                 RR_[i][celli] = 0.0;
00644             }
00645 
00646             scalar rhoi = rho[celli];
00647             scalar Ti = this->thermo().T()[celli];
00648             scalar pi = this->thermo().p()[celli];
00649 
00650             scalarField c(nSpecie_);
00651             scalarField dcdt(nEqns(), 0.0);
00652 
00653             for (label i=0; i<nSpecie_; i++)
00654             {
00655                 scalar Yi = Y_[i][celli];
00656                 c[i] = rhoi*Yi/specieThermo_[i].W();
00657             }
00658 
00659             dcdt = omega(c, Ti, pi);
00660 
00661             for (label i=0; i<nSpecie_; i++)
00662             {
00663                 RR_[i][celli] = dcdt[i]*specieThermo_[i].W();
00664             }
00665         }
00666     }
00667 }
00668 
00669 
00670 template<class CompType, class ThermoType>
00671 Foam::scalar Foam::ODEChemistryModel<CompType, ThermoType>::solve
00672 (
00673     const scalar t0,
00674     const scalar deltaT
00675 )
00676 {
00677     const volScalarField rho
00678     (
00679         IOobject
00680         (
00681             "rho",
00682             this->time().timeName(),
00683             this->mesh(),
00684             IOobject::NO_READ,
00685             IOobject::NO_WRITE,
00686             false
00687         ),
00688         this->thermo().rho()
00689     );
00690 
00691     for (label i=0; i<nSpecie_; i++)
00692     {
00693         RR_[i].setSize(rho.size());
00694     }
00695 
00696     if (!this->chemistry_)
00697     {
00698         return GREAT;
00699     }
00700 
00701     scalar deltaTMin = GREAT;
00702 
00703     tmp<volScalarField> thc = this->thermo().hc();
00704     const scalarField& hc = thc();
00705 
00706     forAll(rho, celli)
00707     {
00708         for (label i=0; i<nSpecie_; i++)
00709         {
00710             RR_[i][celli] = 0.0;
00711         }
00712 
00713         scalar rhoi = rho[celli];
00714         scalar Ti = this->thermo().T()[celli];
00715         scalar hi = this->thermo().hs()[celli] + hc[celli];
00716         scalar pi = this->thermo().p()[celli];
00717 
00718         scalarField c(nSpecie_);
00719         scalarField c0(nSpecie_);
00720         scalarField dc(nSpecie_, 0.0);
00721 
00722         for (label i=0; i<nSpecie_; i++)
00723         {
00724             c[i] = rhoi*Y_[i][celli]/specieThermo_[i].W();
00725         }
00726         c0 = c;
00727 
00728         scalar t = t0;
00729         scalar tauC = this->deltaTChem_[celli];
00730         scalar dt = min(deltaT, tauC);
00731         scalar timeLeft = deltaT;
00732 
00733         // calculate the chemical source terms
00734         scalar cTot = 0.0;
00735 
00736         while (timeLeft > SMALL)
00737         {
00738             tauC = solver().solve(c, Ti, pi, t, dt);
00739             t += dt;
00740 
00741             // update the temperature
00742             cTot = sum(c);
00743             ThermoType mixture(0.0*specieThermo_[0]);
00744             for (label i=0; i<nSpecie_; i++)
00745             {
00746                 mixture += (c[i]/cTot)*specieThermo_[i];
00747             }
00748             Ti = mixture.TH(hi, Ti);
00749 
00750             timeLeft -= dt;
00751             this->deltaTChem_[celli] = tauC;
00752             dt = min(timeLeft, tauC);
00753             dt = max(dt, SMALL);
00754         }
00755         deltaTMin = min(tauC, deltaTMin);
00756 
00757         dc = c - c0;
00758         scalar WTot = 0.0;
00759         for (label i=0; i<nSpecie_; i++)
00760         {
00761             WTot += c[i]*specieThermo_[i].W();
00762         }
00763         WTot /= cTot;
00764 
00765         for (label i=0; i<nSpecie_; i++)
00766         {
00767             RR_[i][celli] = dc[i]*specieThermo_[i].W()/deltaT;
00768         }
00769     }
00770 
00771     // Don't allow the time-step to change more than a factor of 2
00772     deltaTMin = min(deltaTMin, 2*deltaT);
00773 
00774     return deltaTMin;
00775 }
00776 
00777 
00778 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines