Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
36namespace Foam
39namespace objectives
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46(
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
75{
76 J_ = Zero;
77
78 if (const auto* betap = mesh_.cfindObject<volScalarField>("beta"))
79 {
80 const auto& beta = *betap;
81 const scalar time = mesh_.time().timeOutputValue();
82
83 J_ =
84 scalar(1) - gWeightedAverage(mesh_.V(), beta.primitiveField())
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
105 dJdbPtr_().primitiveFieldRef() = -scalar(1)/gSum(mesh_.V());
106 if (percentInDenom_)
107 {
108 dJdbPtr_().primitiveFieldRef() /= targetPercentage_->value(time);
109 }
110}
111
112
114{
116 << setw(width_) << "TargetVolume" << " ";
117}
118
119
121{
122 const scalar time = mesh_.time().timeOutputValue();
124 << setw(width_) << targetPercentage_->value(time) << " ";
125}
126
127
128// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129
130} // End namespace objectives
131} // End namespace Foam
132
133// ************************************************************************* //
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.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion).
Definition TimeStateI.H:24
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Abstract base class for objective functions that contain only geometric quantities.
static autoPtr< objectiveGeometric > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
const fvMesh & mesh_
Definition objective.H:63
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition objective.H:209
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
Objective quantifying the difference between the volume occupied by fluid in topology optimisation an...
virtual void update_dJdb()
Contribution to field sensitivities.
virtual void addHeaderColumns() const
Write headers for additional columns.
virtual scalar J()
Return the objective function value.
objectiveTopOVolume(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
From components.
virtual void addColumnValues() const
Write information to additional columns.
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
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
autoPtr< GeometricField< Type, fvPatchField, volMesh > > createZeroFieldPtr(const fvMesh &mesh, const word &name, const dimensionSet dims, bool printAllocation=false)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)