DAHeatTransferFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAHeatTransferFoam.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAHeatTransferFoam, 0);
16 addToRunTimeSelectionTable(DASolver, DAHeatTransferFoam, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  char* argsAll,
21  PyObject* pyOptions)
22  : DASolver(argsAll, pyOptions),
23  TPtr_(nullptr),
24  fvSourcePtr_(nullptr),
25  kPtr_(nullptr)
26 {
27 }
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
31 {
32  /*
33  Description:
34  Initialize variables for DASolver
35  */
36 
37  Info << "Initializing fields for DAHeatTransferFoam" << endl;
38  Time& runTime = runTimePtr_();
39  fvMesh& mesh = meshPtr_();
41 #include "createAdjointSolid.H"
42  // initialize checkMesh
44 
46 
47  this->setDAObjFuncList();
48 }
49 
51  const Vec xvVec,
52  Vec wVec)
53 {
54  /*
55  Description:
56  Call the primal solver to get converged state variables
57 
58  Input:
59  xvVec: a vector that contains all volume mesh coordinates
60 
61  Output:
62  wVec: state variable vector
63  */
64 
65 #include "createRefsHeatTransfer.H"
66 
67  // change the run status
68  daOptionPtr_->setOption<word>("runStatus", "solvePrimal");
69 
70  Info << "\nCalculating temperature distribution\n"
71  << endl;
72 
73  // deform the mesh based on the xvVec
74  this->pointVec2OFMesh(xvVec);
75 
76  // check mesh quality
77  label meshOK = this->checkMesh();
78 
79  if (!meshOK)
80  {
81  this->writeFailedMesh();
82  return 1;
83  }
84 
85  primalMinRes_ = 1e10;
86  label printInterval = daOptionPtr_->getOption<label>("printIntervalUnsteady");
87  label printToScreen = 0;
88  // main loop
89  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
90  {
91  printToScreen = this->isPrintTime(runTime, printInterval);
92 
93  if (printToScreen)
94  {
95  Info << "Time = " << runTime.timeName() << nl << endl;
96  }
97 
98  fvScalarMatrix TEqn(
99  fvm::laplacian(k, T)
100  + fvSource);
101 
102  // get the solver performance info such as initial
103  // and final residuals
104  SolverPerformance<scalar> solverT = TEqn.solve();
105  this->primalResidualControl<scalar>(solverT, printToScreen, printInterval, "T");
106 
107  if (printToScreen)
108  {
109 
110  this->printAllObjFuncs();
111 
112  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
113  << " ClockTime = " << runTime.elapsedClockTime() << " s"
114  << nl << endl;
115  }
116 
117  runTime.write();
118 
119  }
120 
121  this->calcPrimalResidualStatistics("print");
122 
123  // primal converged, assign the OpenFoam fields to the state vec wVec
124  this->ofField2StateVec(wVec);
125 
126  // write the mesh to files
127  mesh.write();
128 
129  Info << "End\n"
130  << endl;
131 
132  return 0;
133 }
134 
135 } // End namespace Foam
136 
137 // ************************************************************************* //
Foam::DASolver::writeFailedMesh
void writeFailedMesh()
write the failed mesh to disk
Definition: DASolver.C:7955
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:112
Foam::DASolver::daCheckMeshPtr_
autoPtr< DACheckMesh > daCheckMeshPtr_
DACheckMesh object pointer.
Definition: DASolver.H:91
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DACheckMesh
Definition: DACheckMesh.H:29
Foam::DAHeatTransferFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DAHeatTransferFoam.C:50
Foam::DAHeatTransferFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DAHeatTransferFoam.C:30
DAHeatTransferFoam.H
Foam::DASolver
Definition: DASolver.H:49
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
fvSource
volScalarField & fvSource
Definition: createRefsHeatTransfer.H:7
Foam::DASolver::calcPrimalResidualStatistics
void calcPrimalResidualStatistics(const word mode, const label writeRes=0)
calculate the norms of all residuals and print to screen
Definition: DASolver.C:2559
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASolver::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DASolver.C:7460
primalResidualControl< scalar >
this primalResidualControl< scalar >(solverE, printToScreen, printInterval, "he")
Foam::DASolver::daLinearEqnPtr_
autoPtr< DALinearEqn > daLinearEqnPtr_
DALinearEqn pointer.
Definition: DASolver.H:94
Foam::DAHeatTransferFoam::DAHeatTransferFoam
DAHeatTransferFoam(char *argsAll, PyObject *pyOptions)
Definition: DAHeatTransferFoam.C:19
createAdjointSolid.H
k
dimensionedScalar & k
Definition: createRefsHeatTransfer.H:6
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASolver::setDAObjFuncList
void setDAObjFuncList()
initialize DASolver::daObjFuncPtrList_ one needs to call this before calling printAllObjFuncs
Definition: DASolver.C:235
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
Foam::DASolver::printAllObjFuncs
void printAllObjFuncs()
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:160
createRefsHeatTransfer.H
createFieldsHeatTransfer.H
Foam::DASolver::primalMinRes_
scalar primalMinRes_
the maximal residual for primal solution
Definition: DASolver.H:139
Foam::DASolver::runTimePtr_
autoPtr< Time > runTimePtr_
runTime pointer
Definition: DASolver.H:70
Foam::DASolver::meshPtr_
autoPtr< fvMesh > meshPtr_
fvMesh pointer
Definition: DASolver.H:73
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DASolver::ofField2StateVec
void ofField2StateVec(Vec stateVec) const
set the state vector based on the latest fields in OpenFOAM
Definition: DASolver.H:596
Foam::DASolver::checkMesh
label checkMesh() const
check the mesh quality and return meshOK
Definition: DASolver.H:710
Foam::DALinearEqn
Definition: DALinearEqn.H:31
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DASolver::pointVec2OFMesh
void pointVec2OFMesh(const Vec xvVec) const
assign the points in fvMesh of OpenFOAM based on the point vector
Definition: DASolver.H:608
solverT
SolverPerformance< scalar > solverT
Definition: TEqnSimpleT.H:16