topOSource.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) 2020-2023 PCOpt/NTUA
9  Copyright (C) 2020-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 "topOSource.H"
30 #include "fvMatrices.H"
31 #include "fvmSup.H"
32 #include "topOVariablesBase.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace fv
40  {
41  defineTypeNameAndDebug(topOSource, 1);
43  (
44  option,
45  topOSource,
46  dictionary
47  );
48  }
49 }
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
56 {
57  auto tinterpolant
58  (
60  (
61  IOobject
62  (
63  "source",
64  mesh_.time().timeName(),
65  mesh_,
68  ),
69  mesh_,
72  )
73  );
74  DimensionedField<scalar, volMesh>& interpolant = tinterpolant.ref();
75 
76  if (mesh_.foundObject<topOVariablesBase>("topOVars"))
77  {
78  const topOVariablesBase& vars =
80  vars.sourceTerm
82 
83  if (darcyFlow_)
84  {
85  interpolant.field() += betaMax_*Da_();
86  }
87  }
88 
89  return tinterpolant;
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
95 Foam::fv::topOSource::topOSource
96 (
97  const word& name,
98  const word& modelType,
99  const dictionary& dict,
100  const fvMesh& mesh
101 )
102 :
103  option(name, modelType, dict, mesh),
104  interpolation_(topOInterpolationFunction::New(mesh, dict)),
105  interpolationFieldName_(word::null),
106  betaMax_(0),
107  darcyFlow_(false),
108  Da_(nullptr)
109 {
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 (
118  fvMatrix<vector>& eqn,
119  const label fieldi
120 )
121 {
122  DebugInfo
123  << "Adding Brinkman source to " << eqn.psi().name() << endl;
124 
125  eqn -= fvm::Sp(getSource(), eqn.psi());
126 }
127 
128 
130 (
131  fvMatrix<scalar>& eqn,
132  const label fieldi
133 )
134 {
135  DebugInfo
136  << "Adding Brinkman source to " << eqn.psi().name() << endl;
137 
138  eqn -= fvm::Sp(getSource(), eqn.psi());
139 }
140 
141 
143 (
144  const volScalarField& rho,
145  fvMatrix<vector>& eqn,
146  const label fieldi
147 )
148 {
149  DebugInfo
150  << "Adding Brinkman source to " << eqn.psi().name() << endl;
151 
152  eqn -= fvm::Sp(rho*getSource(), eqn.psi());
153 }
154 
155 
157 (
158  const volScalarField& rho,
159  fvMatrix<scalar>& eqn,
160  const label fieldi
161 )
162 {
163  DebugInfo
164  << "Adding Brinkman source to " << eqn.psi().name() << endl;
165 
166  eqn -= fvm::Sp(rho*getSource(), eqn.psi());
167 }
168 
169 
171 (
172  scalarField& sens,
173  const word& fieldName,
174  const word& designVariablesName
175 )
176 {
177  const label fieldi = applyToField(fieldName);
178  if
179  (
180  fieldi != -1
181  && mesh_.foundObject<topOVariablesBase>("topOVars")
182  )
183  {
184  DebugInfo
185  << "Postprocessing Brinkman sensitivities for field "
186  << fieldName << endl;
187  const topOVariablesBase& vars =
188  mesh_.lookupObject<topOVariablesBase>("topOVars");
190  (
191  sens,
192  interpolation_(),
193  betaMax_,
194  designVariablesName,
195  interpolationFieldName_
196  );
197  }
198 }
199 
200 
202 {
203  if (option::read(dict))
204  {
205  fieldNames_ = coeffs_.get<wordList>("names");
206  interpolationFieldName_ = coeffs_.get<word>("interpolationField");
207  applied_.setSize(fieldNames_.size(), false);
208  if (mesh_.foundObject<topOVariablesBase>("topOVars"))
209  {
210  const topOVariablesBase& vars =
211  mesh_.lookupObject<topOVariablesBase>("topOVars");
212  betaMax_ =
213  coeffs_.getOrDefault<scalar>("betaMax", vars.getBetaMax());
214  }
215 
216  darcyFlow_ = coeffs_.getOrDefault<bool>("darcyFlow", false);
217  if (darcyFlow_)
218  {
219  Da_.reset(new scalar(coeffs_.getOrDefault<scalar>("Da", 1.e-5)));
220  }
221 
222  return true;
223  }
224 
225  return false;
226 }
227 
228 
229 // ************************************************************************* //
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...
bool darcyFlow_
Does this option apply to a Darcy flow model.
Definition: topOSource.H:86
autoPtr< topOInterpolationFunction > interpolation_
Interpolation function.
Definition: topOSource.H:69
scalar getBetaMax() const
Get betaMax value.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
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.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition: UList.H:791
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: topOSource.C:194
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add implicit contribution to momentum equation.
Definition: topOSource.C:110
Macros for easy insertion into run-time selection tables.
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
Definition: fvOption.C:75
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual void sourceTerm(DimensionedField< scalar, volMesh > &field, const topOInterpolationFunction &interpolationFunc, const scalar betaMax, const word &interpolationFieldName="beta") const
Populate source terms for the flow equations.
word interpolationFieldName_
Interpolation field name.
Definition: topOSource.H:74
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
autoPtr< scalar > Da_
Dimensionless Darcy number.
Definition: topOSource.H:91
virtual void sourceTermSensitivities(scalarField &sens, const topOInterpolationFunction &interpolationFunc, const scalar betaMax, const word &designVariablesName, const word &interpolationFieldName="beta") const
Post-processing sensitivities due to interpolations based on the indicator fields.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:48
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
#define DebugInfo
Report an information message using Foam::Info.
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
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
label nCells() const noexcept
Number of mesh cells.
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...
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
scalar betaMax_
Optional betaMax.
Definition: topOSource.H:81
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
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...
virtual void postProcessSens(scalarField &sensField, const word &fieldName=word::null, const word &designValue=word::null)
Multiply sensitivities with the derivative of the interpolation function.
Definition: topOSource.C:164
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Base class for all design variables related to topology optimisation (topO). Provides the lookup func...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
virtual tmp< DimensionedField< scalar, volMesh > > getSource()
Compute the source term based on the indicator field.
Definition: topOSource.C:48
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127