objectiveNutSqr.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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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 "objectiveNutSqr.H"
31 #include "createZeroField.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 namespace objectives
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(objectiveNutSqr, 0);
46 (
47  objectiveIncompressible,
48  objectiveNutSqr,
50 );
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const fvMesh& mesh,
58  const dictionary& dict,
59  const word& adjointSolverName,
60  const word& primalSolverName
61 )
62 :
63  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
64  zones_(mesh_.cellZones().indices(dict.get<wordRes>("zones")))
65 {
66  // Check if cellZones provided include at least one cell
67  checkCellZonesSize(zones_);
68  // Allocate dJdTMvar1Ptr_ and dJdTMvar2Ptr_ if needed
70  // Allocate term to be added to volume-based sensitivity derivatives
71  divDxDbMultPtr_.reset
72  (
73  createZeroFieldPtr<scalar>
74  (
75  mesh_,
76  ("divDxdbMult"+type()) ,
77  // Dimensions are set in a way that the gradient of this term
78  // matches the source of the adjoint grid displacement PDE
80  )
81  );
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 scalar objectiveNutSqr::J()
88 {
89  J_ = Zero;
90 
92  turbVars = vars_.RASModelVariables();
93  const volScalarField& nut = turbVars->nutRefInst();
94 
95  //scalar zoneVol(Zero);
96  for (const label zI : zones_)
97  {
98  const cellZone& zoneI = mesh_.cellZones()[zI];
99  for (const label cellI : zoneI)
100  {
101  J_ += sqr(nut[cellI])*(mesh_.V()[cellI]);
102  //zoneVol += mesh_.V()[cellI];
103  }
104  }
106  //reduce(zoneVol, sumOp<scalar>());
107 
108  return J_;
109 }
110 
111 
113 {
114  // Add source from possible dependencies of nut on U
116  {
117  const incompressibleAdjointSolver& adjSolver =
121  adjSolver.getAdjointVars().adjointTurbulence();
124  tmp<volScalarField> dnutdUMult = 2*turbVars->nutRef();
125  tmp<volVectorField> dnutdU = adjointRAS->nutJacobianU(dnutdUMult);
126  if (dnutdU)
127  {
128  if (!dJdvPtr_)
129  {
130  dJdvPtr_.reset
131  (
132  createZeroFieldPtr<vector>
133  (
134  mesh_,
135  "dJdv" + type(),
137  )
138  );
139  }
140  for (const label zI : zones_)
141  {
142  const cellZone& zoneI = mesh_.cellZones()[zI];
143  for (const label cellI : zoneI)
144  {
145  dJdvPtr_()[cellI] = dnutdU()[cellI];
146  }
147  }
148  }
149  }
150 }
151 
152 
154 {
155  const autoPtr<incompressible::RASModelVariables>&
156  turbVars = vars_.RASModelVariables();
157  const volScalarField& nut = turbVars->nutRef();
158  volScalarField JacobianMultiplier(2*nut);
159 
161  (
164  JacobianMultiplier,
165  zones_
166  );
167 }
168 
169 
171 {
173  turbVars = vars_.RASModelVariables();
174  const volScalarField& nut = turbVars->nutRef();
175  volScalarField JacobianMultiplier(2*nut);
176 
178  (
181  JacobianMultiplier,
182  zones_
183  );
184 }
185 
186 
188 {
190  turbVars = vars_.RASModelVariables();
191  const volScalarField& nut = turbVars->nutRef();
192 
193  volScalarField& divDxDbMult = divDxDbMultPtr_();
194 
195  for (const label zI : zones_)
196  {
197  const cellZone& zoneI = mesh_.cellZones()[zI];
198  for (const label cellI : zoneI)
199  {
200  divDxDbMult[cellI] = sqr(nut[cellI]);
201  }
202  }
203  divDxDbMult.correctBoundaryConditions();
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 } // End namespace objectives
210 } // End namespace Foam
211 
212 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
objectiveNutSqr(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
dictionary dict
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
void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
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.
void update_dJdv()
Update values to be added to the adjoint outlet velocity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
defineTypeNameAndDebug(objectiveFlowRate, 0)
scalar J()
Return the objective function value.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
void update_divDxDbMultiplier()
Update field to be added to be added to volume-based sensitivity derivatives, emerging from delta ( d...
Macros for easy insertion into run-time selection tables.
Base class for incompressibleAdjoint solvers.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
dynamicFvMesh & mesh
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:47
const fvMesh & mesh_
Definition: objective.H:63
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
scalar J_
Objective function value and weight.
Definition: objective.H:75
void update_dJdTMvar1()
Update field to be added to the first 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.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:158
const incompressibleVars & vars_
dimensionedScalar pow3(const dimensionedScalar &ds)
addToRunTimeSelectionTable(objectiveIncompressible, objectiveFlowRate, dictionary)
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 reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
void correctBoundaryConditions()
Correct boundary field.
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
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.
scalar nut
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:157