

Protected Attributes | |
These are state variables, state residuals, and partial derivatives | |
| volVectorField & | U_ |
| volVectorField | URes_ |
| volScalarField & | p_ |
| volScalarField | pRes_ |
| volScalarField & | T_ |
| volScalarField | TRes_ |
| surfaceScalarField & | phi_ |
| surfaceScalarField | phiRes_ |
| volVectorField & | fvSource_ |
| fvSource term More... | |
| volScalarField & | fvSourceEnergy_ |
| fvSource term for the energy equation More... | |
| fluidThermo & | thermo_ |
| thermophysical property More... | |
These are intermediate variables | |
| volScalarField & | he_ |
| volScalarField & | rho_ |
| volScalarField & | alphat_ |
| volScalarField & | psi_ |
| volScalarField & | dpdt_ |
| volScalarField & | K_ |
| DATurbulenceModel & | daTurb_ |
| DATurbulenceModel object. More... | |
| pimpleControl | pimple_ |
| pimpleControl object which will be initialized in this class More... | |
Protected Attributes inherited from DAResidual | |
| const fvMesh & | mesh_ |
| fvMesh More... | |
| const DAOption & | daOption_ |
| DAOption object. More... | |
| const DAModel & | daModel_ |
| DAModel object. More... | |
| const DAIndex & | daIndex_ |
| DAIndex. More... | |
| DAField | daField_ |
| DAField object. More... | |
These are constants to update the intermediate variables | |
| scalar | molWeight_ |
| scalar | Cp_ |
| label | hasFvSource_ = 0 |
| whether to have fvSource term More... | |
| TypeName ("DARhoPimpleFoam") | |
| DAResidualRhoPimpleFoam (const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex) | |
| virtual | ~DAResidualRhoPimpleFoam () |
| virtual void | clear () |
| clear the members More... | |
| virtual void | calcResiduals (const dictionary &options) |
| compute residual More... | |
| virtual void | updateIntermediateVariables () |
| update any intermediate variables that are dependent on state variables and are used in calcResiduals More... | |
| virtual void | correctBoundaryConditions () |
| update the boundary condition for all the states in the selected solver More... | |
| virtual void | calcPCMatWithFvMatrix (Mat PCMat) |
| calculating the adjoint preconditioner matrix using fvMatrix More... | |
Additional Inherited Members | |
Public Member Functions inherited from DAResidual | |
| TypeName ("DAResidual") | |
| Runtime type information. More... | |
| declareRunTimeSelectionTable (autoPtr, DAResidual, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex),(modelType, mesh, daOption, daModel, daIndex)) | |
| DAResidual (const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex) | |
| virtual | ~DAResidual () |
| void | masterFunction (const dictionary &options, const Vec xvVec, const Vec wVec, Vec resVec) |
| the master function that compute the residual vector given the state and point vectors More... | |
| bool | writeData (Ostream &os) const |
| virtual function for regIOobject More... | |
Static Public Member Functions inherited from DAResidual | |
| static autoPtr< DAResidual > | New (const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex) |
Definition at line 31 of file DAResidualRhoPimpleFoam.H.
| DAResidualRhoPimpleFoam | ( | const word | modelType, |
| const fvMesh & | mesh, | ||
| const DAOption & | daOption, | ||
| const DAModel & | daModel, | ||
| const DAIndex & | daIndex | ||
| ) |
Definition at line 19 of file DAResidualRhoPimpleFoam.C.
References allOptions, DAResidualRhoPimpleFoam::Cp_, DAResidual::daOption_, DAOption::getAllOptions(), DAOption::getOption(), DAResidualRhoPimpleFoam::hasFvSource_, mesh, and DAResidualRhoPimpleFoam::molWeight_.

|
inlinevirtual |
Definition at line 98 of file DAResidualRhoPimpleFoam.H.
| TypeName | ( | "DARhoPimpleFoam" | ) |
|
virtual |
clear the members
Implements DAResidual.
Definition at line 76 of file DAResidualRhoPimpleFoam.C.
References DAResidualRhoPimpleFoam::phiRes_, DAResidualRhoPimpleFoam::pRes_, DAResidualRhoPimpleFoam::TRes_, and DAResidualRhoPimpleFoam::URes_.
|
virtual |
compute residual
Implements DAResidual.
Definition at line 90 of file DAResidualRhoPimpleFoam.C.
References alphaEff(), DAResidualRhoPimpleFoam::alphat_, DAFvSource::calcFvSource(), DAResidual::daOption_, DAResidualRhoPimpleFoam::daTurb_, DATurbulenceModel::divDevRhoReff(), DAResidualRhoPimpleFoam::dpdt_, EEqn(), DAResidualRhoPimpleFoam::fvSource_, DAResidualRhoPimpleFoam::fvSourceEnergy_, DAOption::getOption(), DAResidualRhoPimpleFoam::hasFvSource_, HbyA, HbyAPtr, DAResidualRhoPimpleFoam::he_, DAResidualRhoPimpleFoam::K_, DAResidual::mesh_, normalizePhiResiduals, normalizeResiduals, DAResidualRhoPimpleFoam::p_, DAResidualRhoPimpleFoam::phi_, phiHbyA, DAResidualRhoPimpleFoam::phiRes_, DAResidualRhoPimpleFoam::pRes_, DAResidualRhoPimpleFoam::psi_, rAU(), DAResidualRhoPimpleFoam::rho_, rhorAUf(), DAResidualRhoPimpleFoam::thermo_, DAResidualRhoPimpleFoam::TRes_, tUEqn(), DAResidualRhoPimpleFoam::U_, UEqn, DAResidualRhoPimpleFoam::URes_, and useConstrainHbyA.

|
virtual |
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Implements DAResidual.
Definition at line 212 of file DAResidualRhoPimpleFoam.C.
References DAResidualRhoPimpleFoam::Cp_, DAResidualRhoPimpleFoam::dpdt_, DAResidualRhoPimpleFoam::he_, DAResidualRhoPimpleFoam::K_, DAResidualRhoPimpleFoam::molWeight_, DAResidualRhoPimpleFoam::p_, DAResidualRhoPimpleFoam::psi_, DAResidualRhoPimpleFoam::rho_, DAResidualRhoPimpleFoam::T_, and DAResidualRhoPimpleFoam::U_.
|
virtual |
update the boundary condition for all the states in the selected solver
Implements DAResidual.
Definition at line 300 of file DAResidualRhoPimpleFoam.C.
References DAResidualRhoPimpleFoam::p_, DAResidualRhoPimpleFoam::T_, and DAResidualRhoPimpleFoam::U_.
|
virtual |
calculating the adjoint preconditioner matrix using fvMatrix
Reimplemented from DAResidual.
Definition at line 312 of file DAResidualRhoPimpleFoam.C.
References assignValueCheckAD, DAFvSource::calcFvSource(), D, DAResidual::daIndex_, DAResidual::daOption_, DAResidualRhoPimpleFoam::daTurb_, DATurbulenceModel::divDevRhoReff(), forAll(), DAResidualRhoPimpleFoam::fvSource_, DAResidualRhoPimpleFoam::fvSourceEnergy_, DAOption::getAllOptions(), DAIndex::getGlobalAdjointStateIndex(), DAOption::getOption(), DAResidualRhoPimpleFoam::hasFvSource_, HbyA, HbyAPtr, DAResidual::mesh_, DAIndex::nLocalInternalFaces, DAResidualRhoPimpleFoam::p_, DAResidualRhoPimpleFoam::phi_, phiHbyA, DAResidualRhoPimpleFoam::psi_, rAU(), DAResidualRhoPimpleFoam::rho_, rhorAUf(), DAResidualRhoPimpleFoam::U_, UEqn, and useConstrainHbyA.

|
protected |
|
protected |
Definition at line 39 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::clear().
|
protected |
|
protected |
Definition at line 42 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::clear().
|
protected |
Definition at line 44 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::correctBoundaryConditions(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
Definition at line 45 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::clear().
|
protected |
Definition at line 47 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), and DAResidualRhoPimpleFoam::calcResiduals().
|
protected |
Definition at line 48 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::clear().
|
protected |
fvSource term
Definition at line 52 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), and DAResidualRhoPimpleFoam::calcResiduals().
|
protected |
fvSource term for the energy equation
Definition at line 55 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), and DAResidualRhoPimpleFoam::calcResiduals().
|
protected |
thermophysical property
Definition at line 58 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals().
|
protected |
Definition at line 62 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
Definition at line 63 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
Definition at line 64 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals().
|
protected |
Definition at line 65 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
Definition at line 66 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
Definition at line 67 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
DATurbulenceModel object.
Definition at line 71 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), and DAResidualRhoPimpleFoam::calcResiduals().
|
protected |
pimpleControl object which will be initialized in this class
Definition at line 74 of file DAResidualRhoPimpleFoam.H.
|
protected |
Definition at line 78 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::DAResidualRhoPimpleFoam(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
Definition at line 79 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::DAResidualRhoPimpleFoam(), and DAResidualRhoPimpleFoam::updateIntermediateVariables().
|
protected |
whether to have fvSource term
Definition at line 83 of file DAResidualRhoPimpleFoam.H.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), DAResidualRhoPimpleFoam::calcResiduals(), and DAResidualRhoPimpleFoam::DAResidualRhoPimpleFoam().