DAHeatTransferFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
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  daFvSourcePtr_(nullptr)
27 {
28 }
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
32 {
33  /*
34  Description:
35  Initialize variables for DASolver
36  */
37 
38  Info << "Initializing fields for DAHeatTransferFoam" << endl;
39  Time& runTime = runTimePtr_();
40  fvMesh& mesh = meshPtr_();
42 #include "createAdjoint.H"
43 
44  // initialize fvSource and compute the source term
45  const dictionary& allOptions = daOptionPtr_->getAllOptions();
46  if (allOptions.subDict("fvSource").toc().size() != 0)
47  {
48  hasFvSource_ = 1;
49  Info << "Initializing DASource" << endl;
50  word sourceName = allOptions.subDict("fvSource").toc()[0];
51  word fvSourceType = allOptions.subDict("fvSource").subDict(sourceName).getWord("type");
53  fvSourceType, mesh, daOptionPtr_(), daModelPtr_(), daIndexPtr_()));
54  }
55 }
56 
58 {
59  /*
60  Description:
61  Call the primal solver to get converged state variables
62 
63  Input:
64  xvVec: a vector that contains all volume mesh coordinates
65 
66  Output:
67  wVec: state variable vector
68  */
69 
70 #include "createRefsHeatTransfer.H"
71 
72  Info << "\nCalculating temperature distribution\n"
73  << endl;
74 
75  // main loop
76  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
77  {
78 
79  if (printToScreen_)
80  {
81  Info << "Time = " << runTime.timeName() << nl << endl;
82  }
83 
84  if (hasFvSource_)
85  {
86  volScalarField& fvSource = fvSourcePtr_();
87  daFvSourcePtr_->calcFvSource(fvSource);
88  }
89 
90  fvScalarMatrix TEqn(
91  fvm::laplacian(k, T)
92  + fvSource);
93 
94  // get the solver performance info such as initial
95  // and final residuals
96  SolverPerformance<scalar> solverT = TEqn.solve();
98 
100 
101  // print run time
103 
104  runTime.write();
105  }
106 
107  // write the mesh to files
108  mesh.write();
109 
110  Info << "End\n"
111  << endl;
112 
113  return 0;
114 }
115 
116 } // End namespace Foam
117 
118 // ************************************************************************* //
Foam::DAHeatTransferFoam::fvSourcePtr_
autoPtr< volScalarField > fvSourcePtr_
heat source
Definition: DAHeatTransferFoam.H:35
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::daModelPtr_
autoPtr< DAModel > daModelPtr_
DAModel pointer.
Definition: DASolver.H:84
Foam::DASolver::daIndexPtr_
autoPtr< DAIndex > daIndexPtr_
DAIndex pointer.
Definition: DASolver.H:87
Foam::DAHeatTransferFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DAHeatTransferFoam.C:31
DAHeatTransferFoam.H
Foam::DASolver
Definition: DASolver.H:54
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAHeatTransferFoam::daFvSourcePtr_
autoPtr< DAFvSource > daFvSourcePtr_
DASource pointer.
Definition: DAHeatTransferFoam.H:41
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:81
fvSource
volScalarField & fvSource
Definition: createRefsHeatTransfer.H:7
createAdjoint.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
solverT
SolverPerformance< scalar > solverT
Definition: TEqnPimpleDyM.H:21
Foam::DAHeatTransferFoam::DAHeatTransferFoam
DAHeatTransferFoam(char *argsAll, PyObject *pyOptions)
Definition: DAHeatTransferFoam.C:19
Foam::DAUtility::primalResidualControl
static void primalResidualControl(const SolverPerformance< scalar > &solverP, const label printToScreen, const word varName, scalar &primalMaxRes)
control when to print the residual and also compute the maxInitRes
Definition: DAUtility.C:734
Foam::DAFvSource::New
static autoPtr< DAFvSource > New(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAFvSource.C:47
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
k
dimensionedScalar & k
Definition: createRefsHeatTransfer.H:6
Foam
Definition: checkGeometry.C:32
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
allOptions
const dictionary & allOptions
Definition: createRefsPimpleDyM.H:14
createRefsHeatTransfer.H
createFieldsHeatTransfer.H
Foam::DAHeatTransferFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DAHeatTransferFoam.H:44
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::DAHeatTransferFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DAHeatTransferFoam.C:57
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DASolver::printElapsedTime
void printElapsedTime(const Time &runTime, const label printToScreen)
Definition: DASolver.H:256
Foam::DASolver::calcAllFunctions
void calcAllFunctions(label print=0)
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:204
Foam::DASolver::daGlobalVarPtr_
autoPtr< DAGlobalVar > daGlobalVarPtr_
DAGlobalVar pointer.
Definition: DASolver.H:114