incompressibleAdjointVars.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46(
47 fvMesh& mesh,
48 solverControl& SolverControl,
49 objectiveManager& objManager,
50 incompressibleVars& primalVars
51)
52:
53 incompressibleAdjointMeanFlowVars(mesh, SolverControl, primalVars),
54 objectiveManager_(objManager),
55
56 adjointTurbulence_
57 (
58 incompressibleAdjoint::adjointRASModel::New
59 (
60 primalVars_,
61 *this,
62 objManager
63 )
64 )
65{}
66
67
68// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69
71{
73 {
74 Info<< "Resetting adjoint mean fields to zero" << endl;
75
76 // Reset fields to zero
77 paMeanPtr_() == dimensionedScalar(paPtr_().dimensions(), Zero);
78 UaMeanPtr_() == dimensionedVector(UaPtr_().dimensions(), Zero);
79 phiaMeanPtr_() == dimensionedScalar(phiaPtr_().dimensions(), Zero);
80 adjointTurbulence_().resetMeanFields();
81
82 // Reset averaging iteration index to 0
84 }
85}
86
87
89{
91 {
92 Info<< "Averaging adjoint fields" << endl;
93 label& iAverageIter = solverControl_.averageIter();
94 scalar avIter(iAverageIter);
95 scalar oneOverItP1 = 1./(avIter+1);
96 scalar mult = avIter*oneOverItP1;
97 paMeanPtr_() == paMeanPtr_() *mult + paPtr_() *oneOverItP1;
98 UaMeanPtr_() == UaMeanPtr_() *mult + UaPtr_() *oneOverItP1;
99 phiaMeanPtr_() == phiaMeanPtr_()*mult + phiaPtr_()*oneOverItP1;
100 adjointTurbulence_().computeMeanFields();
101 ++iAverageIter;
102 }
103}
104
105
107{
109 adjointTurbulence_->nullify();
110}
111
112
114{
115 /*
116 // WIP
117 for (fvPatchVectorField& pf : UaInst().boundaryFieldRef())
118 {
119 if (isA<adjointBoundaryCondition<vector>>(pf))
120 {
121 adjointBoundaryCondition<vector>& adjointBC =
122 refCast<adjointBoundaryCondition<vector>>(pf);
123 adjointBC.updatePrimalBasedQuantities();
124 }
125 }
126
127 for (fvPatchScalarField& pf : paInst().boundaryFieldRef())
128 {
129 if (isA<adjointBoundaryCondition<scalar>>(pf))
130 {
131 adjointBoundaryCondition<scalar>& adjointBC =
132 refCast<adjointBoundaryCondition<scalar>>(pf);
133 adjointBC.updatePrimalBasedQuantities();
134 }
135 }
136 */
137}
138
139
140// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141
142} // End namespace Foam
143
144// ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Manages the adjoint mean flow fields and their mean values.
solverControl & solverControl_
Reference to the solverControl of the solver allocating the fields.
virtual void nullify()
Nullify all adjoint fields.
Class including all adjoint fields for incompressible flows.
virtual void nullify()
Nullify all adjoint fields.
virtual void updatePrimalBasedQuantities()
Update primal based quantities of the adjoint boundary.
void computeMeanFields()
Compute mean fields on the fly.
void resetMeanFields()
Reset mean fields to zero.
autoPtr< incompressibleAdjoint::adjointRASModel > adjointTurbulence_
Adjoint to the turbulence model.
Base class for solution control classes.
class for managing incompressible objective functions.
Base class for solver control classes.
Definition: solverControl.H:52
bool doAverageIter() const
label & averageIter()
Return average iteration index reference.
bool average() const
Whether averaging is enabled or not.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.