DATurboFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DATurboFoam.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DATurboFoam, 0);
16 addToRunTimeSelectionTable(DASolver, DATurboFoam, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  char* argsAll,
21  PyObject* pyOptions)
22  : DASolver(argsAll, pyOptions),
23  simplePtr_(nullptr),
24  pThermoPtr_(nullptr),
25  pPtr_(nullptr),
26  rhoPtr_(nullptr),
27  UPtr_(nullptr),
28  URelPtr_(nullptr),
29  phiPtr_(nullptr),
30  pressureControlPtr_(nullptr),
31  turbulencePtr_(nullptr),
32  daTurbulenceModelPtr_(nullptr),
33  MRFPtr_(nullptr),
34  initialMass_(dimensionedScalar("initialMass", dimensionSet(1, 0, 0, 0, 0, 0, 0), 0.0))
35 {
36 }
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
40 {
41  /*
42  Description:
43  Initialize variables for DASolver
44  */
45 
46  Info << "Initializing fields for DATurboFoam" << endl;
47  Time& runTime = runTimePtr_();
48  fvMesh& mesh = meshPtr_();
49  argList& args = argsPtr_();
51 #include "createFieldsTurbo.H"
52 
53  // read the RAS model from constant/turbulenceProperties
54  const word turbModelName(
55  IOdictionary(
56  IOobject(
57  "turbulenceProperties",
58  mesh.time().constant(),
59  mesh,
60  IOobject::MUST_READ,
61  IOobject::NO_WRITE,
62  false))
63  .subDict("RAS")
64  .lookup("RASModel"));
66 
67 #include "createAdjoint.H"
68 }
69 
71 {
72  /*
73  Description:
74  Call the primal solver to get converged state variables
75  */
76 
77 #include "createRefsTurbo.H"
78 
79  // call correctNut, this is equivalent to turbulence->validate();
80  daTurbulenceModelPtr_->updateIntermediateVariables();
81 
82  Info << "\nStarting time loop\n"
83  << endl;
84 
85  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
86  {
87  if (printToScreen_)
88  {
89  Info << "Time = " << runTime.timeName() << nl << endl;
90  }
91 
92  p.storePrevIter();
93  rho.storePrevIter();
94 
95  // Pressure-velocity SIMPLE corrector
96 #include "UEqnTurbo.H"
97 #include "pEqnTurbo.H"
98 #include "EEqnTurbo.H"
99 
101 
102  // calculate all functions
104  // calculate yPlus
106  // compute the regression model and print the feature
107  regModelFail_ = daRegressionPtr_->compute();
108  daRegressionPtr_->printInputInfo(printToScreen_);
109  // print run time
111 
112  runTime.write();
113  }
114 
115  // write the mesh to files
116  mesh.write();
117 
118  Info << "End\n"
119  << endl;
120 
121  return this->checkPrimalFailure();
122 }
123 
124 } // End namespace Foam
125 
126 // ************************************************************************* //
Foam::DASolver::checkPrimalFailure
label checkPrimalFailure()
check whether the primal fails based on residual and regression fail flag
Definition: DASolver.C:2472
Foam::DASolver::loop
label loop(Time &runTime)
return whether to loop the primal solution, similar to runTime::loop() except we don't do file IO
Definition: DASolver.C:147
Foam::DASolver::printToScreen_
label printToScreen_
whether to print primal information to the screen
Definition: DASolver.H:150
Foam::DASolver::daRegressionPtr_
autoPtr< DARegression > daRegressionPtr_
DARegression pointer.
Definition: DASolver.H:111
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
Foam::DASolver
Definition: DASolver.H:54
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:81
UEqnTurbo.H
Foam::DATurbulenceModel::New
static autoPtr< DATurbulenceModel > New(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DATurbulenceModel.C:161
createAdjoint.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DATurboFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DATurboFoam.H:65
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
createSimpleControlPython.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
pEqnTurbo.H
rho
volScalarField & rho
Definition: createRefsRhoPimple.H:8
Foam::DASolver::regModelFail_
label regModelFail_
whether the regModel compute fails
Definition: DASolver.H:153
Foam::DATurboFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DATurboFoam.C:70
Foam::DASolver::runTimePtr_
autoPtr< Time > runTimePtr_
runTime pointer
Definition: DASolver.H:75
Foam::DASolver::meshPtr_
autoPtr< fvMesh > meshPtr_
fvMesh pointer
Definition: DASolver.H:78
Foam::DASolver::argsPtr_
autoPtr< argList > argsPtr_
args pointer
Definition: DASolver.H:72
DATurboFoam.H
Foam::DATurboFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DATurboFoam.C:39
args
Foam::argList & args
Definition: setRootCasePython.H:42
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DATurboFoam::DATurboFoam
DATurboFoam(char *argsAll, PyObject *pyOptions)
Definition: DATurboFoam.C:19
createFieldsTurbo.H
createRefsTurbo.H
Foam::DASolver::printElapsedTime
void printElapsedTime(const Time &runTime, const label printToScreen)
Definition: DASolver.H:256
EEqnTurbo.H
Foam::DASolver::calcAllFunctions
void calcAllFunctions(label print=0)
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:204