Loading...
Searching...
No Matches
adjointNull.H
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) 2023 PCOpt/NTUA
9 Copyright (C) 2023 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
27Class
28 Foam::adjointNull
29
30Description
31 Dummy adjoint solver.
32 Used to add the derivatives of geometric constraints which do not require
33 the solution of adjoint PDEs
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef adjointNull_H
38#define adjointNull_H
39
40#include "adjointSolver.H"
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44namespace Foam
45{
46
47/*---------------------------------------------------------------------------*\
48 Class adjointNull Declaration
49\*---------------------------------------------------------------------------*/
50
51class adjointNull
52:
53 public adjointSolver
54{
55
56private:
57
58 // Private Member Functions
59
60 //- No copy construct
61 adjointNull(const adjointNull&) = delete;
62
63 //- No copy assignment
64 void operator=(const adjointNull&) = delete;
65
66
67protected:
68
69 // Protected data
70
71 //- Accumulate the sensitivities integrand before calculating them
72 virtual void preCalculateSensitivities();
73
74
75public:
76
77 // Static Data Members
78
79 //- Run-time type information
80 TypeName("null");
82
83 // Constructors
84
85 //- Construct from mesh and dictionary
86 adjointNull
87 (
88 fvMesh& mesh,
89 const word& managerType,
90 const dictionary& dict,
92 const word& solverName
93 );
94
95
96 // Selectors
97
98 //- Return a reference to the selected turbulence model
100 (
101 fvMesh& mesh,
102 const word& managerType,
103 const dictionary& dict,
104 const word& primalSolverName,
105 const word& solverName
106 );
107
108
109 //- Destructor
110 virtual ~adjointNull() = default;
111
112
113 // Member Functions
114
115 // Access
116
117 //- Return the type of simulation this solver pertains to
118 virtual const word simulationType() const;
119
120 //- Return the dimensions of the adjoint grid displacement variable
121 virtual dimensionSet maDimensions() const;
122
123
124 // Evolution
125
126 //- Execute one iteration of the solution algorithm
127 virtual void solveIter();
128
129 //- Main control loop
130 virtual void solve();
131
132 //- Looper (advances iters, time step)
133 virtual bool loop();
134
135 //- Update primal based quantities related to the objective
136 //- functions.
137 // Also writes the objective function values to files.
138 virtual void updatePrimalBasedQuantities();
139
140 // Functions related to the computation of sensitivity derivatives
141
142 // Shape optimisation
143
144 //- Compute the multiplier for grad(dxdb)
145 // Used in shape sensitivity derivatives, computed with
146 // the FI and E-SI approaches
148 (
149 volTensorField& gradDxDbMult,
150 const scalar dt
151 );
152
153 //- Compute the multiplier for div(dxdb)
154 // Used in shape sensitivity derivatives, computed with
155 // the FI and E-SI approaches
156 virtual void accumulateDivDxDbMultiplier
157 (
158 autoPtr<scalarField>& divDxDbMult,
159 const scalar dt
160 );
161
162 //- Accumulate the multipliers of geometric quantities
163 //- defined at the boundary, usually through an objective
164 //- or constraint function
166 (
169 autoPtr<boundaryVectorField>& dxdbDirectMult,
170 autoPtr<pointBoundaryVectorField>& pointDxDirectDbMult,
171 const labelHashSet& sensitivityPatchIDs,
172 const scalar dt
173 );
174
175
176 // Topology optimisation
177
178 //- Compute the multiplier of beta
179 virtual void topOSensMultiplier
180 (
181 scalarField& betaMult,
182 const word& designVariablesName,
183 const scalar dt
184 );
185};
186
187
188// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189
190} // End namespace Foam
191
192
193// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194
195#endif
196
197// ************************************************************************* //
virtual const word simulationType() const
Return the type of simulation this solver pertains to.
TypeName("null")
Run-time type information.
virtual void accumulateGeometryVariationsMultipliers(autoPtr< boundaryVectorField > &dSfdbMult, autoPtr< boundaryVectorField > &dnfdbMult, autoPtr< boundaryVectorField > &dxdbDirectMult, autoPtr< pointBoundaryVectorField > &pointDxDirectDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Accumulate the multipliers of geometric quantities defined at the boundary, usually through an object...
virtual void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb).
virtual ~adjointNull()=default
Destructor.
virtual void updatePrimalBasedQuantities()
Update primal based quantities related to the objective functions.
virtual void preCalculateSensitivities()
Accumulate the sensitivities integrand before calculating them.
Definition adjointNull.C:44
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb).
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
virtual bool loop()
Looper (advances iters, time step).
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
static autoPtr< adjointNull > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected turbulence model.
Definition adjointNull.C:77
virtual void solveIter()
Execute one iteration of the solution algorithm.
virtual void solve()
Main control loop.
const word & primalSolverName() const
Return the primal solver name.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const word & managerType() const
Return the manager type.
Definition solverI.H:72
const dictionary & dict() const
Return the solver dictionary.
Definition solverI.H:54
const fvMesh & mesh() const
Return the solver mesh.
Definition solverI.H:24
const word & solverName() const
Return the solver name.
Definition solverI.H:30
A class for handling words, derived from Foam::string.
Definition word.H:66
Namespace for OpenFOAM.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68