objectivePowerDissipation.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 -------------------------------------------------------------------------------
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 
31 #include "createZeroField.H"
32 #include "IOmanip.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 namespace objectives
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(objectivePowerDissipation, 0);
47 (
48  objectiveIncompressible,
49  objectivePowerDissipation,
51 );
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const fvMesh& mesh,
59  const dictionary& dict,
60  const word& adjointSolverName,
61  const word& primalSolverName
62 )
63 :
64  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
65  zones_(mesh_.cellZones().indices(dict.get<wordRes>("zones")))
66 {
67  // Append Ua name to fieldNames
68  /*
69  fieldNames_.setSize
70  (
71  1,
72  mesh_.lookupObject<solver>(adjointSolverName_).
73  extendedVariableName("Ua")
74  );
75  */
76 
77  // Check if cellZones provided include at least one cell
78  checkCellZonesSize(zones_);
79 
80  // Allocate dJdTMvar1Ptr_ and dJdTMvar2Ptr_ if needed
82 
83  // Allocate source term to the adjoint momentum equations
84  dJdvPtr_.reset
85  (
86  createZeroFieldPtr<vector>
87  (
88  mesh_,
89  "dJdv" + type(),
91  )
92  );
93  // Allocate terms to be added to volume-based sensitivity derivatives
94  divDxDbMultPtr_.reset
95  (
96  createZeroFieldPtr<scalar>
97  (
98  mesh_,
99  ("divDxdbMult" + type()),
100  // Dimensions are set in a way that the gradient of this term
101  // matches the source of the adjoint grid displacement PDE
103  )
104  );
105  gradDxDbMultPtr_.reset
106  (
107  createZeroFieldPtr<tensor>
108  (
109  mesh_,
110  ("gradDxdbMult" + type()),
112  )
113  );
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  J_ = Zero;
122 
123  // References
124  const volVectorField& U = vars_.UInst();
126  const scalarField& V = mesh_.V().field();
127 
128  volScalarField integrand(turb->nuEff()*magSqr(twoSymm(fvc::grad(U))));
129 
130  for (const label zI : zones_)
131  {
132  const cellZone& zoneI = mesh_.cellZones()[zI];
133  scalarField VZone(V, zoneI);
134  scalarField integrandZone(integrand.primitiveField(), zoneI);
135 
136  J_ += 0.5*gSum(integrandZone*VZone);
137  }
138 
139  return J_;
140 }
141 
142 
144 {
145  dJdvPtr_().primitiveFieldRef() = Zero;
146 
147  const volVectorField& U = vars_.U();
150 
151  // Add source from possible dependencies of nut on U
153  {
154  const incompressibleAdjointSolver& adjSolver =
158  adjSolver.getAdjointVars().adjointTurbulence();
159  tmp<volScalarField> dnutdUMult = 0.5*magSqr(S());
160  tmp<volVectorField> dnutdU = adjointRAS->nutJacobianU(dnutdUMult);
161  if (dnutdU)
162  {
163  for (const label zI : zones_)
164  {
165  const cellZone& zoneI = mesh_.cellZones()[zI];
166  for (const label cellI : zoneI)
167  {
168  dJdvPtr_()[cellI] = dnutdU()[cellI];
169  }
170  }
171  }
172  }
173 
174  // Add source from the strain rate magnitude
175  volVectorField integrand(-2.0*fvc::div(turb->nuEff()*S));
176  for (const label zI : zones_)
177  {
178  const cellZone& zoneI = mesh_.cellZones()[zI];
179  for (const label cellI : zoneI)
180  {
181  dJdvPtr_()[cellI] += integrand[cellI];
182  }
183  }
184 }
185 
186 
188 {
189  const volVectorField& U = vars_.U();
190  volScalarField JacobianMultiplier(0.5*magSqr(twoSymm(fvc::grad(U))));
192  (
195  JacobianMultiplier,
196  zones_
197  );
198 }
199 
200 
202 {
203  const volVectorField& U = vars_.U();
204  volScalarField JacobianMultiplier(0.5*magSqr(twoSymm(fvc::grad(U))));
206  (
209  JacobianMultiplier,
210  zones_
211  );
212 }
213 
214 
216 {
217  // References
218  volScalarField& divDxDbMult = divDxDbMultPtr_();
219  const volVectorField& U = vars_.U();
221  volScalarField integrand(0.5*turb->nuEff()*magSqr(twoSymm(fvc::grad(U))));
222 
223  for (const label zI : zones_)
224  {
225  const cellZone& zoneI = mesh_.cellZones()[zI];
226  for (const label cellI : zoneI)
227  {
228  divDxDbMult[cellI] = integrand[cellI];
229  }
230  }
231  divDxDbMult.correctBoundaryConditions();
232 }
233 
234 
236 {
237  // References
238  volTensorField& gradDxDbMult = gradDxDbMultPtr_();
239  const volVectorField& U = vars_.U();
240  const autoPtr<incompressible::turbulenceModel>& turb = vars_.turbulence();
241  tmp<volTensorField> gradU = fvc::grad(U);
242  volTensorField integrand(-2.0*turb->nuEff()*(gradU() & twoSymm(gradU())));
243 
244  for (const label zI : zones_)
245  {
246  const cellZone& zoneI = mesh_.cellZones()[zI];
247  for (const label cellI : zoneI)
248  {
249  gradDxDbMult[cellI] = integrand[cellI];
250  }
251  }
252  gradDxDbMult.correctBoundaryConditions();
253  // Missing contribution from gradU in nut
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 } // End namespace objectives
260 } // End namespace Foam
261 
262 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
objectivePowerDissipation(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const word adjointSolverName_
Definition: objective.H:65
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
compressible::turbulenceModel & turb
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
defineTypeNameAndDebug(objectiveFlowRate, 0)
void update_dJdv()
Update values to be added to the adjoint outlet velocity.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Base class for incompressibleAdjoint solvers.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar J()
Return the objective function value.
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
autoPtr< volVectorField > dJdvPtr_
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const fvMesh & mesh_
Definition: objective.H:63
void update_dJdTMvar1()
Update field to be added to the first adjoint turbulence model PDE.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
scalar J_
Objective function value and weight.
Definition: objective.H:75
const volVectorField & U() const
Return const reference to velocity.
Istream and Ostream manipulators taking arguments.
virtual void update_divDxDbMultiplier()
Update div(dx/db multiplier). Volume-based sensitivity term.
void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
A subset of mesh cells.
Definition: cellZone.H:58
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:158
virtual void update_gradDxDbMultiplier()
Update grad(dx/db multiplier). Volume-based sensitivity term.
const incompressibleVars & vars_
dimensionedScalar pow3(const dimensionedScalar &ds)
addToRunTimeSelectionTable(objectiveIncompressible, objectiveFlowRate, dictionary)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:654
void correctBoundaryConditions()
Correct boundary field.
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:163
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const volVectorField & UInst() const
Return const reference to velocity.
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
Abstract base class for objective functions in incompressible flows.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133