objectiveTopOVolume.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 
29 #include "objectiveTopOVolume.H"
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(objectiveTopOVolume, 1);
46 (
47  objectiveGeometric,
48  objectiveTopOVolume,
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  objectiveGeometric(mesh, dict, adjointSolverName, primalSolverName),
64  targetPercentage_(Function1<scalar>::New("percentage", dict)),
65  percentInDenom_(dict.getOrDefault<bool>("percentInDenom", true))
66 {
67  // Allocate boundary field pointers
68  dJdbPtr_.reset(createZeroFieldPtr<scalar>(mesh_, "dJdb", dimless));
69 }
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 
76 {
77  J_ = Zero;
78  if (mesh_.foundObject<volScalarField>("beta"))
79  {
82  const scalar time = mesh_.time().timeOutputValue();
83  J_ =
84  scalar(1) - gSum(beta.primitiveField()*V)/gSum(V)
85  - targetPercentage_->value(time);
86  if (percentInDenom_)
87  {
88  J_ /= targetPercentage_->value(time);
89  }
90  }
91  else
92  {
94  << "Beta field not yet registered in database. OK for start-up"
95  << endl;
96  }
97  return J_;
98 }
99 
100 
102 {
103  const scalar time = mesh_.time().timeOutputValue();
104  dJdbPtr_().primitiveFieldRef() = -scalar(1)/gSum(mesh_.V());
105  if (percentInDenom_)
106  {
107  dJdbPtr_().primitiveFieldRef() /= targetPercentage_->value(time);
108  }
109 }
110 
111 
113 {
115  << setw(width_) << "TargetVolume" << " ";
116 }
117 
118 
120 {
121  const scalar time = mesh_.time().timeOutputValue();
123  << setw(width_) << targetPercentage_->value(time) << " ";
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 } // End namespace objectives
130 } // End namespace Foam
131 
132 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
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...
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
Abstract base class for objective functions that contain only geometric quantities.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
virtual void update_dJdb()
Contribution to field sensitivities.
Macros for easy insertion into run-time selection tables.
objectiveTopOVolume(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
From components.
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)
const fvMesh & mesh_
Definition: objective.H:63
scalar J_
Objective function value and weight.
Definition: objective.H:76
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
Istream and Ostream manipulators taking arguments.
virtual void addHeaderColumns() const
Write headers for additional columns.
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:226
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion) ...
Definition: TimeStateI.H:24
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
virtual void addColumnValues() const
Write information to additional columns.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Namespace for OpenFOAM.
virtual scalar J()
Return the objective function value.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127