adjointSensitivityIncompressible.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-2022 PCOpt/NTUA
9 Copyright (C) 2013-2022 FOSS GP
10 Copyright (C) 2019-2021 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
34#include "wallFvPatch.H"
35#include "fvOptions.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41
42namespace incompressible
43{
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
53(
54 const fvMesh& mesh,
55 const dictionary& dict,
57)
58:
60 derivatives_(0),
61 adjointSolver_(adjointSolver),
62 primalVars_(adjointSolver.getPrimalVars()),
63 adjointVars_(adjointSolver.getAdjointVars()),
64 objectiveManager_(adjointSolver.getObjectiveManager())
65{}
66
67
68// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
69
71(
72 const fvMesh& mesh,
73 const dictionary& dict,
75)
76{
77 const word modelType(dict.get<word>("type"));
78
79 Info<< "adjointSensitivity type : " << modelType << endl;
80
81 auto* ctorPtr = dictionaryConstructorTable(modelType);
82
83 if (!ctorPtr)
84 {
86 (
87 dict,
88 "adjointSensitivity",
89 modelType,
90 *dictionaryConstructorTablePtr_
91 ) << exit(FatalIOError);
92 }
93
95 (
96 ctorPtr(mesh, dict, adjointSolver)
97 );
98}
99
100
101// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
102
104{
106 write(type());
107 return derivatives_;
108}
109
110
112{
113 return derivatives_;
114}
115
116
118{
119 derivatives_ = scalar(0);
120 if (fieldSensPtr_)
121 {
122 fieldSensPtr_().primitiveFieldRef() = scalar(0);
123 }
124}
125
126
127void adjointSensitivity::write(const word& baseName)
128{
129 sensitivity::write(baseName);
130}
131
132
134{
136}
137
138
140{
142 volTensorField& gradDxDbMult = tgradDxDbMult.ref();
143
144 tmp<volVectorField> tadjointMeshMovementSource
145 (
147 (
149 (
150 "adjointMeshMovementSource",
151 mesh_.time().timeName(),
152 mesh_,
155 ),
156 mesh_,
158 )
159 );
160
161 volVectorField& source = tadjointMeshMovementSource.ref();
162
163 source -= fvc::div(gradDxDbMult.T());
164
165 // Terms from fvOptions
166 fv::options::New(this->mesh_).postProcessSens
167 (
169 );
170
171 return (tadjointMeshMovementSource);
172}
173
174
175// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176
177} // End namespace incompressible
178} // End namespace Foam
179
180// ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Base class for adjoint solvers.
Definition: adjointSolver.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Base class for incompressibleAdjoint solvers.
virtual tmp< volTensorField > computeGradDxDbMultiplier()
Compute the multiplier for grad(dxdb)
Abstract base class for adjoint-based sensitivities in incompressible flows.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
virtual const scalarField & calculateSensitivities()
Calculates and returns sensitivity fields.
virtual void assembleSensitivities()=0
Assemble sensitivities.
tmp< volVectorField > adjointMeshMovementSource()
Compute source term for adjoint mesh movement equation.
const scalarField & getSensitivities() const
Returns the sensitivity fields.
Abstract base class for adjoint sensitivities.
Definition: sensitivity.H:64
const fvMesh & mesh_
Definition: sensitivity.H:69
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:57
autoPtr< volScalarField > fieldSensPtr_
Definition: sensitivity.H:74
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Namespace for OpenFOAM.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
runTime write()
Macros to ease declaration of run-time selection tables.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict