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-2023 PCOpt/NTUA
9  Copyright (C) 2013-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 
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
71  (
72  1,
74  extendedVariableName("Ua")
75  );
76 
77  // Check if cellZones provided include at least one cell
78  checkCellZonesSize(zones_);
79 
80  // Allocate source term to the adjoint momentum equations
81  dJdvPtr_.reset
82  (
83  createZeroFieldPtr<vector>
84  (
85  mesh_,
86  "dJdv" + type(),
88  )
89  );
90  // Allocate term to be added to volume-based sensitivity derivatives
91  divDxDbMultPtr_.reset
92  (
93  new volScalarField
94  (
95  IOobject
96  (
97  "divDxDbMult" + objectiveName_,
98  mesh_.time().timeName(),
99  mesh_,
102  ),
103  mesh_,
106  )
107  );
108 }
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 {
115  J_ = Zero;
116 
117  // References
118  const volVectorField& U = vars_.UInst();
119  const scalarField& V = mesh_.V().field();
120 
121  for (const label zI : zones_)
122  {
123  const cellZone& zoneI = mesh_.cellZones()[zI];
124  scalarField VZone(V, zoneI);
125  vectorField UZone(U.primitiveField(), zoneI);
126  volZone_[zI] = gSum(VZone);
127  UMean_[zI] = gSum(UZone*VZone)/volZone_[zI];
128  UVar_[zI] = gSum(magSqr(UZone - UMean_[zI])*VZone)/volZone_[zI];
129 
130  J_ += 0.5*UVar_[zI];
131  }
132 
133  return J_;
134 }
135 
136 
138 {
139  const volVectorField& U = vars_.U();
140 
141  for (const label zI : zones_)
142  {
143  const cellZone& zoneI = mesh_.cellZones()[zI];
144  for (const label cellI : zoneI)
145  {
146  dJdvPtr_()[cellI] = (U[cellI] - UMean_[zI])/volZone_[zI];
147  }
148  }
149 }
150 
151 
153 {
154  volScalarField& divDxDbMult = divDxDbMultPtr_();
155  const volVectorField& U = vars_.U();
156 
157  for (const label zI : zones_)
158  {
159  const cellZone& zoneI = mesh_.cellZones()[zI];
160  for (const label cellI : zoneI)
161  {
162  divDxDbMult[cellI] =
163  0.5*(magSqr(U[cellI] - UMean_[zI]) - UVar_[zI])/volZone_[zI];
164  }
165  }
166  divDxDbMult.correctBoundaryConditions();
167 }
168 
169 
171 {
172  OFstream& file = objFunctionFilePtr_();
173  for (const label zI : zones_)
174  {
175  const word zoneName(mesh_.cellZones()[zI].name());
176  file<< setw(width_) << word(zoneName + "-" + "UMean") << " ";
177  file<< setw(width_) << word(zoneName + "-" + "UVar") << " ";
178  file<< setw(width_) << word(zoneName + "-" + "UStd") << " ";
179  file<< setw(width_) << word(zoneName + "-" + "Vol") << " ";
180  }
181 }
182 
183 
185 {
186  OFstream& file = objFunctionFilePtr_();
187  forAll(UMean_, zI)
188  {
189  file<< setw(width_) << mag(UMean_[zI]) << " ";
190  file<< setw(width_) << UVar_[zI] << " ";
191  file<< setw(width_) << sqrt(UVar_[zI]) << " ";
192  file<< setw(width_) << volZone_[zI] << " ";
193  }
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace objectives
200 } // End namespace Foam
201 
202 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
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
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
Base solver class.
Definition: solver.H:45
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:129
virtual void addColumnValues() const
Write information to additional columns.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:209
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.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
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
const word objectiveName_
Definition: objective.H:67
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
Definition: objective.H:130
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
virtual void addHeaderColumns() const
Write headers for additional columns.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
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_
defineTypeNameAndDebug(objectivePartialVolume, 1)
virtual void update_dJdv()
Update values to be added to the adjoint outlet velocity.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const fvMesh & mesh_
Definition: objective.H:63
scalar J_
Objective function value and weight.
Definition: objective.H:76
const volVectorField & U() const
Return const reference to velocity.
Istream and Ostream manipulators taking arguments.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
A subset of mesh cells.
Definition: cellZone.H:58
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:194
virtual 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:226
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
const volVectorField & UInst() const
Return const reference to velocity.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
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:127