objectivePartialVolume.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  Copyright (C) 2019-2020 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "objectivePartialVolume.H"
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 
47 (
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  objectiveGeometric(mesh, dict, adjointSolverName, primalSolverName),
65  initVol_(Zero),
66  objectivePatches_
67  (
68  mesh_.boundaryMesh().patchSet
69  (
70  dict.get<wordRes>("patches")
71  ).sortedToc()
72  )
73 {
74  // Read target volume if present. Else use the current one as a target
75  if
76  (
77  !objective::readIfPresent("initialVolume", initVol_)
78  && !dict.readIfPresent("initialVolume", initVol_)
79  )
80  {
81  const scalar oneThird(1.0/3.0);
82  for (const label patchi : objectivePatches_)
83  {
84  const fvPatch& patch = mesh_.boundary()[patchi];
85  initVol_ -= oneThird*gSum(patch.Sf() & patch.Cf());
86  }
87  }
88  // Allocate boundary field pointers
89  bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
90  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 
98 {
99  J_ = Zero;
100  const scalar oneThird(1.0/3.0);
101  for (const label patchi : objectivePatches_)
102  {
103  const fvPatch& patch = mesh_.boundary()[patchi];
104  J_ -= oneThird*gSum(patch.Sf() & patch.Cf());
105  }
106  J_ -= initVol_;
107  J_ /= initVol_;
108  return J_;
109 }
110 
111 
113 {
114  const scalar oneThird(1.0/3.0);
115  for (const label patchi : objectivePatches_)
116  {
117  const fvPatch& patch = mesh_.boundary()[patchi];
119  const vectorField& nf = tnf();
120  bdxdbDirectMultPtr_()[patchi] = -oneThird*nf/initVol_;
121  }
122 }
123 
124 
126 {
127  const scalar oneThird(1.0/3.0);
128  for (const label patchi : objectivePatches_)
129  {
130  const fvPatch& patch = mesh_.boundary()[patchi];
131  bdSdbMultPtr_()[patchi] = -oneThird*patch.Cf()/initVol_;
132  }
133 }
134 
135 
137 {
138  os.writeEntry("initialVolume", initVol_);
139  return objective::writeData(os);
140 }
141 
142 
144 {
146  << setw(4) << "#" << " "
147  << setw(width_) << "VInit" << " "
148  << setw(width_) << initVol_ << endl;
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 } // End namespace objectives
155 } // End namespace Foam
156 
157 // ************************************************************************* //
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
objectivePartialVolume(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:91
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:209
virtual bool writeData(Ostream &os) const
Write initial volume for continuation.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Abstract base class for objective functions that contain only geometric quantities.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
static constexpr scalar oneThird
Definition: colourTools.C:35
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
defineTypeNameAndDebug(objectivePartialVolume, 1)
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
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
Istream and Ostream manipulators taking arguments.
virtual void addHeaderInfo() const
Write headers for additional columns.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
virtual bool writeData(Ostream &os) const
Write averaged objective for continuation.
Definition: objective.C:638
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface and volume-based sensitivity term.
OBJstream os(runTime.globalPath()/outputName)
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:226
virtual scalar J()
Return the objective function value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const std::string patch
OpenFOAM patch number as a std::string.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x...
Definition: objective.H:178
A class for managing temporary objects.
Definition: HashPtrTable.H:50
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:161
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127