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-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 
35 namespace 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<< "Restoring adjoint field to initial ones" << endl;
75  paInst() == dimensionedScalar(paInst().dimensions(), Zero);
76  UaInst() == dimensionedVector(UaInst().dimensions(), Zero);
77  phiaInst() == dimensionedScalar(phiaInst().dimensions(), Zero);
78  adjointTurbulence_().restoreInitValues();
79  }
80 }
81 
82 
84 {
85  if (solverControl_.average())
86  {
87  Info<< "Resetting adjoint mean fields to zero" << endl;
88 
89  // Reset fields to zero
90  paMeanPtr_() == dimensionedScalar(paPtr_().dimensions(), Zero);
91  UaMeanPtr_() == dimensionedVector(UaPtr_().dimensions(), Zero);
92  phiaMeanPtr_() == dimensionedScalar(phiaPtr_().dimensions(), Zero);
93  adjointTurbulence_().resetMeanFields();
94 
95  // Reset averaging iteration index to 0
97  }
98 }
99 
100 
102 {
104  {
105  Info<< "Averaging adjoint fields" << endl;
106  label& iAverageIter = solverControl_.averageIter();
107  scalar avIter(iAverageIter);
108  scalar oneOverItP1 = 1./(avIter+1);
109  scalar mult = avIter*oneOverItP1;
110  paMeanPtr_() == paMeanPtr_() *mult + paPtr_() *oneOverItP1;
111  UaMeanPtr_() == UaMeanPtr_() *mult + UaPtr_() *oneOverItP1;
112  phiaMeanPtr_() == phiaMeanPtr_()*mult + phiaPtr_()*oneOverItP1;
113  adjointTurbulence_().computeMeanFields();
114  ++iAverageIter;
115  }
116 }
117 
118 
120 {
122  adjointTurbulence_->nullify();
123 }
124 
125 
127 {
128  /*
129  // WIP
130  for (fvPatchVectorField& pf : UaInst().boundaryFieldRef())
131  {
132  if (isA<adjointBoundaryCondition<vector>>(pf))
133  {
134  adjointBoundaryCondition<vector>& adjointBC =
135  refCast<adjointBoundaryCondition<vector>>(pf);
136  adjointBC.updatePrimalBasedQuantities();
137  }
138  }
139 
140  for (fvPatchScalarField& pf : paInst().boundaryFieldRef())
141  {
142  if (isA<adjointBoundaryCondition<scalar>>(pf))
143  {
144  adjointBoundaryCondition<scalar>& adjointBC =
145  refCast<adjointBoundaryCondition<scalar>>(pf);
146  adjointBC.updatePrimalBasedQuantities();
147  }
148  }
149  */
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 } // End namespace Foam
156 
157 // ************************************************************************* //
virtual void nullify()
Nullify all adjoint fields.
const volScalarField & paInst() const
Return const reference to pressure.
Class for managing objective functions.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
autoPtr< volScalarField > paPtr_
Fields involved in the solution of the incompressible adjoint NS equations.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
incompressibleAdjointVars(const incompressibleAdjointVars &)=delete
No copy construct.
Class including all adjoint fields for incompressible flows.
virtual void nullify()
Nullify all adjoint fields.
Base class for solution control classes.
dynamicFvMesh & mesh
Base class for solver control classes.
Definition: solverControl.H:45
autoPtr< volScalarField > paMeanPtr_
Mean Adjoint Fields. Actual averaging is done in the incompressibleAdjointVars class to take care of ...
defineTypeNameAndDebug(combustionModel, 0)
bool doAverageIter() const
Whether or not to add fields of the current iteration to the average fields.
Manages the adjoint mean flow fields and their mean values.
const volVectorField & UaInst() const
Return const reference to velocity.
solverControl & solverControl_
Reference to the solverControl of the solver allocating the fields.
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
bool average() const
Whether averaging is enabled or not.
void computeMeanFields()
Compute mean fields on the fly.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void updatePrimalBasedQuantities()
Update primal based quantities of the adjoint boundary.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label & averageIter()
Return average iteration index reference.
bool storeInitValues() const
Re-initialize.
void restoreInitValues()
Restore field values to the initial ones.
void resetMeanFields()
Reset mean fields to zero.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
autoPtr< incompressibleAdjoint::adjointRASModel > adjointTurbulence_
Adjoint to the turbulence model.