Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
12License
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
31#include "createZeroField.H"
32#include "IOmanip.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
40namespace 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_));
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());
106 J_ -= initVol_;
107 J_ /= initVol_;
108 return J_;
109}
110
111
112void objectivePartialVolume::update_dxdbDirectMultiplier()
113{
114 const scalar oneThird(1.0/3.0);
115 for (const label patchi : objectivePatches_)
116 {
117 const fvPatch& patch = mesh_.boundary()[patchi];
118 tmp<vectorField> tnf = patch.nf();
119 const vectorField& nf = tnf();
120 bdxdbDirectMultPtr_()[patchi] = -oneThird*nf/initVol_;
121 }
122}
123
124
125void objectivePartialVolume::update_dSdbMultiplier()
126{
127 const scalar oneThird(1.0/3.0);
128 for (const label patchi : objectivePatches_)
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{
145 objFunctionFilePtr_()
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// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
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...
wordList sortedToc() const
Return the sorted table of contents.
Definition dictionary.C:601
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Abstract base class for objective functions that contain only geometric quantities.
const fvMesh & mesh_
Definition objective.H:63
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition objective.H:161
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 averaged objective for continuation.
Definition objective.C:641
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x,...
Definition objective.H:178
unsigned int width_
Default width of entries when writing in the objective files.
Definition objective.H:226
scalar J_
Objective function value and weight.
Definition objective.H:76
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
virtual bool writeData(Ostream &os) const
Write initial volume for continuation.
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
objectivePartialVolume(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
virtual scalar J()
Return the objective function value.
virtual void addHeaderInfo() const
Write headers for additional columns.
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface and volume-based sensitivity term.
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
static constexpr scalar oneThird
Definition colourTools.C:35
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)