objectiveUniformityCellZone.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 
30 #include "createZeroField.H"
31 #include "IOmanip.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 namespace objectives
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(objectiveUniformityCellZone, 0);
46 (
47  objectiveIncompressible,
48  objectiveUniformityCellZone,
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  UMean_(zones_.size(), Zero),
66  UVar_(zones_.size(), Zero),
67  volZone_(zones_.size(), Zero)
68 {
69  // Append Ua name to fieldNames
70  /*
71  fieldNames_.setSize
72  (
73  1,
74  mesh_.lookupObject<solver>(adjointSolverName_).
75  extendedVariableName("Ua")
76  );
77  */
78 
79  // Check if cellZones provided include at least one cell
80  checkCellZonesSize(zones_);
81 
82  // Allocate source term to the adjoint momentum equations
83  dJdvPtr_.reset
84  (
85  createZeroFieldPtr<vector>
86  (
87  mesh_,
88  "dJdv" + type(),
90  )
91  );
92  // Allocate term to be added to volume-based sensitivity derivatives
93  divDxDbMultPtr_.reset
94  (
95  createZeroFieldPtr<scalar>
96  (
97  mesh_,
98  ("divDxdbMult" + type()) ,
99  // Dimensions are set in a way that the gradient of this term
100  // matches the source of the adjoint grid displacement PDE
102  )
103  );
104 }
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  J_ = Zero;
112 
113  // References
114  const volVectorField& U = vars_.UInst();
115  const scalarField& V = mesh_.V().field();
116 
117  for (const label zI : zones_)
118  {
119  const cellZone& zoneI = mesh_.cellZones()[zI];
120  scalarField VZone(V, zoneI);
121  vectorField UZone(U.primitiveField(), zoneI);
122  volZone_[zI] = gSum(VZone);
123  UMean_[zI] = gSum(UZone*VZone)/volZone_[zI];
124  UVar_[zI] = gSum(magSqr(UZone - UMean_[zI])*VZone)/volZone_[zI];
125 
126  J_ += 0.5*UVar_[zI];
127  }
128 
129  return J_;
130 }
131 
132 
134 {
135  const volVectorField& U = vars_.U();
136 
137  for (const label zI : zones_)
138  {
139  const cellZone& zoneI = mesh_.cellZones()[zI];
140  for (const label cellI : zoneI)
141  {
142  dJdvPtr_()[cellI] = (U[cellI] - UMean_[zI])/volZone_[zI];
143  }
144  }
145 }
146 
147 
149 {
150  volScalarField& divDxDbMult = divDxDbMultPtr_();
151  const volVectorField& U = vars_.U();
152 
153  for (const label zI : zones_)
154  {
155  const cellZone& zoneI = mesh_.cellZones()[zI];
156  for (const label cellI : zoneI)
157  {
158  divDxDbMult[cellI] =
159  0.5*(magSqr(U[cellI] - UMean_[zI]) - UVar_[zI])/volZone_[zI];
160  }
161  }
162  divDxDbMult.correctBoundaryConditions();
163 }
164 
165 
167 {
168  OFstream& file = objFunctionFilePtr_();
169  for (const label zI : zones_)
170  {
171  const word zoneName(mesh_.cellZones()[zI].name());
172  file<< setw(width_) << word(zoneName + "-" + "UMean") << " ";
173  file<< setw(width_) << word(zoneName + "-" + "UVar") << " ";
174  file<< setw(width_) << word(zoneName + "-" + "UStd") << " ";
175  file<< setw(width_) << word(zoneName + "-" + "Vol") << " ";
176  }
177 }
178 
179 
181 {
182  OFstream& file = objFunctionFilePtr_();
183  forAll(UMean_, zI)
184  {
185  file<< setw(width_) << mag(UMean_[zI]) << " ";
186  file<< setw(width_) << UVar_[zI] << " ";
187  file<< setw(width_) << sqrt(UVar_[zI]) << " ";
188  file<< setw(width_) << volZone_[zI] << " ";
189  }
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace objectives
196 } // End namespace Foam
197 
198 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
virtual void addColumnValues() const
Write information to additional columns.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:173
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Output to file stream, using an OSstream.
Definition: OFstream.H:49
dimensionedScalar sqrt(const dimensionedScalar &ds)
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
defineTypeNameAndDebug(objectiveFlowRate, 0)
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Macros for easy insertion into run-time selection tables.
objectiveUniformityCellZone(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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
virtual void addHeaderColumns() const
Write headers for additional columns.
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_
void update_dJdv()
Update values to be added to the adjoint outlet velocity.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
const fvMesh & mesh_
Definition: objective.H:63
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.
A subset of mesh cells.
Definition: cellZone.H:58
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:158
scalar J()
Return the objective function value.
const incompressibleVars & vars_
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:190
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
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
const volVectorField & UInst() const
Return const reference to velocity.
Abstract base class for objective functions in incompressible flows.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157