objectiveManagerIncompressible.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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
42(
46);
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const fvMesh& mesh,
53 const dictionary& dict,
54 const word& adjointSolverName,
55 const word& primalSolverName
56)
57:
58 objectiveManager(mesh, dict, adjointSolverName, primalSolverName)
59{}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
65{
66 // Add contributions from objective functions
67 for (objective& obj : objectives_)
68 {
69 auto& icoObj = refCast<objectiveIncompressible>(obj);
70
71 if (icoObj.hasdJdv())
72 {
73 scalar weight = icoObj.weight();
74 UaEqn += weight*icoObj.dJdv();
75 }
76 }
77}
78
79
81{
82 // Add contributions from objective functions
83 for (objective& obj : objectives_)
84 {
85 auto& icoObj = refCast<objectiveIncompressible>(obj);
86
87 if (icoObj.hasdJdp())
88 {
89 scalar weight = icoObj.weight();
90 paEqn += weight*icoObj.dJdp();
91 }
92 }
93}
94
95
97{
98 // Add contributions from objective functions
99 for (objective& obj : objectives_)
100 {
101 auto& icoObj = refCast<objectiveIncompressible>(obj);
102
103 if (icoObj.hasdJdTMVar1())
104 {
105 scalar weight = icoObj.weight();
106 adjTMEqn1 += weight*icoObj.dJdTMvar1();
107 }
108 }
109}
110
111
113{
114 // Add contributions from objective functions
115 for (objective& obj : objectives_)
116 {
117 auto& icoObj = refCast<objectiveIncompressible>(obj);
118
119 if (icoObj.hasdJdTMVar2())
120 {
121 scalar weight = icoObj.weight();
122 adjTMEqn2 += weight*icoObj.dJdTMvar2();
123 }
124 }
125}
126
127
128// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129
130} // End namespace Foam
131
132// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
class for managing incompressible objective functions.
virtual void addTMEqn1Source(fvScalarMatrix &adjTMEqn1)
Add contribution to adjoint turbulence model PDE.
virtual void addTMEqn2Source(fvScalarMatrix &adjTMEqn2)
Add contribution to adjoint turbulence model PDE.
virtual void addUaEqnSource(fvVectorMatrix &UaEqn)
Add contribution to adjoint momentum PDEs.
virtual void addPaEqnSource(fvScalarMatrix &paEqn)
Add contribution to adjoint momentum PDEs.
class for managing incompressible objective functions.
PtrList< objective > objectives_
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
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
Namespace for OpenFOAM.
fvScalarMatrix paEqn(fvm::d2dt2(pa) - sqr(c0) *fvc::laplacian(pa))
dictionary dict