adjointNull.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) 2023 PCOpt/NTUA
9  Copyright (C) 2023 FOSS GP
10 -------------------------------------------------------------------------------
11 License
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 
29 #include "adjointNull.H"
30 #include "findRefCell.H"
31 #include "constrainHbyA.H"
32 #include "adjustPhi.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(adjointNull, 0);
41  (
42  adjointSolver,
43  adjointNull,
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
53  adjointSensitivity_->accumulateIntegrand(scalar(1));
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
59 Foam::adjointNull::adjointNull
60 (
61  fvMesh& mesh,
62  const word& managerType,
63  const dictionary& dict,
64  const word& primalSolverName,
65  const word& solverName
66 )
67 :
69  (
70  mesh,
71  managerType,
72  dict,
73  primalSolverName,
74  solverName
75  )
76 {
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
82 
84 (
85  fvMesh& mesh,
86  const word& managerType,
87  const dictionary& dict,
88  const word& primalSolverName,
89  const word& solverName
90 )
91 {
93  (
94  new adjointNull
95  (
96  mesh,
97  managerType,
98  dict,
99  primalSolverName,
100  solverName
101  )
102  );
103 }
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 {
110  return word("null");
111 }
112 
115 {
116  return pow3(dimLength/dimTime);
117 }
118 
121 {
122  // Does nothing
123 }
124 
127 {
128  // Does nothing
129 }
130 
133 {
134  return false;
135 }
136 
137 
139 {
140  // Update objective function related quantities
141  objectiveManager_.updateAndWrite();
142 }
143 
144 
146 (
147  volTensorField& gradDxDbMult,
148  const scalar dt
149 )
150 {
151  tmp<volTensorField> tsens
152  (
154  (
155  IOobject
156  (
157  "flowTerm",
158  mesh_.time().timeName(),
159  mesh_,
162  ),
163  mesh_,
166  )
167  );
168  volTensorField& sens = tsens.ref();
169 
170  PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
171 
172  for (objective& objI : functions)
173  {
174  if (objI.hasGradDxDbMult())
175  {
176  sens += objI.weight()*objI.gradDxDbMultiplier();
177  }
178  }
180 
181  gradDxDbMult += sens.T()*dt;
182 }
183 
184 
186 (
187  autoPtr<scalarField>& divDxDbMult,
188  const scalar dt
189 )
190 {
191  PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
192  for (objective& func : functions)
193  {
194  if (func.hasDivDxDbMult())
195  {
196  divDxDbMult() +=
197  func.weight()*func.divDxDbMultiplier().primitiveField()*dt;
198  }
199  }
200 }
201 
202 
204 (
205  autoPtr<boundaryVectorField>& dSfdbMult,
206  autoPtr<boundaryVectorField>& dnfdbMult,
207  autoPtr<boundaryVectorField>& dxdbDirectMult,
208  autoPtr<pointBoundaryVectorField>& pointDxDbDirectMult,
209  const labelHashSet& sensitivityPatchIDs,
210  const scalar dt
211 )
212 {
213  PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
214  for (const label patchI : sensitivityPatchIDs)
215  {
216  const scalarField magSfDt(mesh_.boundary()[patchI].magSf()*dt);
217  for (objective& func : functions)
218  {
219  const scalar wei(func.weight());
220  if (func.hasdSdbMult())
221  {
222  dSfdbMult()[patchI] += wei*func.dSdbMultiplier(patchI)*dt;
223  }
224  if (func.hasdndbMult())
225  {
226  dnfdbMult()[patchI] += wei*func.dndbMultiplier(patchI)*magSfDt;
227  }
228  if (func.hasdxdbDirectMult())
229  {
230  dxdbDirectMult()[patchI] +=
231  wei*func.dxdbDirectMultiplier(patchI)*magSfDt;
232  }
233  }
234  }
235 }
236 
237 
239 (
240  scalarField& betaMult,
241  const word& designVariablesName,
242  const scalar dt
243 )
244 {
245  // Terms resulting directly from the objective function
246  PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
247  for (objective& objI : functions)
248  {
249  const scalar weight = objI.weight();
250  if (objI.hasdJdb())
251  {
252  betaMult += weight*objI.dJdb()*dt;
253  }
254 
255  if (objI.hasdJdbField())
256  {
257  SubField<scalar> betaSens(objI.dJdbField(), mesh_.nCells(), 0);
258  betaMult += weight*betaSens*dt;
259  }
260  }
261 }
262 
263 
264 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
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...
Definition: adjointNull.C:197
virtual void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb)
Definition: adjointNull.C:179
dictionary dict
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
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual const word simulationType() const
Return the type of simulation this solver pertains to.
Definition: adjointNull.C:101
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void allocateSensitivities()
Allocate the sensitivity derivatives.
Definition: adjointSolver.C:38
virtual void updatePrimalBasedQuantities()
Update primal based quantities related to the objective functions.
Definition: adjointNull.C:131
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
Definition: adjointNull.C:232
Base class for adjoint solvers.
Definition: adjointSolver.H:53
autoPtr< adjointSensitivity > adjointSensitivity_
Sensitivity Derivatives engine.
Ignore writing from objectRegistry::writeObject()
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic...
Macros for easy insertion into run-time selection tables.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb)
Definition: adjointNull.C:139
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
virtual bool loop()
Looper (advances iters, time step)
Definition: adjointNull.C:125
Dummy adjoint solver. Used to add the derivatives of geometric constraints which do not require the s...
Definition: adjointNull.H:44
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
Definition: adjointNull.C:107
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
virtual void preCalculateSensitivities()
Accumulate the sensitivities integrand before calculating them.
Definition: adjointNull.C:44
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual void solve()
Main control loop.
Definition: adjointNull.C:119
Namespace for OpenFOAM.
virtual void solveIter()
Execute one iteration of the solution algorithm.
Definition: adjointNull.C:113
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127