shapeSensitivitiesIncompressible.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) 2020 PCOpt/NTUA
9 Copyright (C) 2020 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37namespace incompressible
38{
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
43
44
45// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46
48{
49 // Accumulate direct sensitivities
51 for (const label patchI : sensitivityPatchIDs_)
52 {
53 const scalarField magSfDt(mesh_.boundary()[patchI].magSf()*dt);
54 for (objective& func : functions)
55 {
56 const scalar wei(func.weight());
57 dSfdbMult_()[patchI] += wei*func.dSdbMultiplier(patchI)*dt;
58 dnfdbMult_()[patchI] += wei*func.dndbMultiplier(patchI)*magSfDt;
59 dxdbDirectMult_()[patchI] +=
60 wei*func.dxdbDirectMultiplier(patchI)*magSfDt;
61 }
62 }
63}
64
65
67{
68 // Avoid updating the event number to keep consistency with cases caching
69 // gradUa
70 auto& UaBoundary = adjointVars_.Ua().boundaryFieldRef(false);
72
73 // Accumulate sensitivities due to boundary conditions
74 for (const label patchI : sensitivityPatchIDs_)
75 {
76 const scalarField magSfDt(mesh_.boundary()[patchI].magSf()*dt);
77 fvPatchVectorField& Uab = UaBoundary[patchI];
78 if (isA<adjointVectorBoundaryCondition>(Uab))
79 {
80 bcDxDbMult_()[patchI] +=
81 (
82 DvDbMult()[patchI]
83 & refCast<adjointVectorBoundaryCondition>(Uab).dxdbMult()
84 )*magSfDt;
85 }
86 }
87}
88
89
91{
93 tres(createZeroBoundaryPtr<vector>(meshShape_).ptr());
94 boundaryVectorField& res = tres.ref();
95
96 // Grab references
97 const volScalarField& pa = adjointVars_.pa();
98 const volVectorField& Ua = adjointVars_.Ua();
101
102 // Fields needed to calculate adjoint sensitivities
104 turbVars = primalVars_.RASModelVariables();
106 volScalarField nuEff(lamTransp.nu() + turbVars->nutRef());
107 tmp<volTensorField> tgradUa = fvc::grad(Ua);
108 const volTensorField::Boundary& gradUabf = tgradUa.cref().boundaryField();
109
110 for (const label patchI : sensitivityPatchIDs_)
111 {
112 const fvPatch& patch = meshShape_.boundary()[patchI];
113 tmp<vectorField> tnf = patch.nf();
114 const vectorField& nf = tnf();
115
116 res[patchI] =
117 (
118 nuEff.boundaryField()[patchI]
119 * (
120 Ua.boundaryField()[patchI].snGrad()
121 + (gradUabf[patchI] & nf)
122 )
123 )
124 - (nf*pa.boundaryField()[patchI])
125 + adjointTurbulence().adjointMomentumBCSource()[patchI];
126 }
127
128 return tres;
129}
130
131
132// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133
135(
136 const fvMesh& mesh,
137 const dictionary& dict,
139)
140:
143 dSfdbMult_(createZeroBoundaryPtr<vector>(mesh_)),
144 dnfdbMult_(createZeroBoundaryPtr<vector>(mesh_)),
145 dxdbDirectMult_(createZeroBoundaryPtr<vector>(mesh_)),
146 bcDxDbMult_(createZeroBoundaryPtr<vector>(mesh_))
147{}
148
149
150// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151
153{
158
161}
162
163
164void shapeSensitivities::write(const word& baseName)
165{
168}
169
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173} // End namespace incompressible
174} // End namespace Foam
175
176// ************************************************************************* //
Generic GeometricBoundaryField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
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
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const volScalarField & pa() const
Return const reference to pressure.
const volVectorField & Ua() const
Return const reference to velocity.
Base class for incompressibleAdjoint solvers.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Abstract base class for adjoint-based sensitivities in incompressible flows.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
virtual void accumulateBCSensitivityIntegrand(const scalar dt)
Accumulate sensitivities enamating from the boundary conditions.
autoPtr< boundaryVectorField > dSfdbMult_
Fields related to direct sensitivities.
virtual void accumulateDirectSensitivityIntegrand(const scalar dt)
Accumulate direct sensitivities.
tmp< boundaryVectorField > dvdbMult() const
Compute multiplier of dv_i/db.
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
const fvMesh & mesh_
Definition: sensitivity.H:69
void clearSensitivities()
Zero sensitivity fields and their constituents.
void write()
Write sensitivity fields.
A simple single-phase transport model based on viscosityModel.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
T & ref() const
Definition: tmpI.H:227
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
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
void func(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
dictionary dict
A non-counting (dummy) refCount.
Definition: refCount.H:59