forcing.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "forcing.H"
27 #include "fvMatrix.H"
28 #include "Function1Evaluate.H"
29 #include "fvcGrad.H"
30 #include "fvcVolumeIntegrate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
44 
46 {
47  lambda_ =
49  (
50  lambda_.name(),
52  coeffs().lookup(lambda_.name())
53  );
54 
57  (
60  coeffs().lookupOrDefault(lambdaBoundary_.name(), 0.0)
61  );
62 }
63 
64 
66 {
67  writeForceFields_ = coeffs().lookupOrDefault("writeForceFields", false);
68 
69  const bool foundScale = coeffs().found("scale");
70  const bool foundOgn = coeffs().found("origin");
71  const bool foundDir = coeffs().found("direction");
72  const bool foundOgns = coeffs().found("origins");
73  const bool foundDirs = coeffs().found("directions");
74  const bool foundAll =
75  foundScale
76  && (
77  (foundOgn && foundDir && !foundOgns && !foundDirs)
78  || (!foundOgn && !foundDir && foundOgns && foundDirs)
79  );
80  const bool foundAny =
81  foundScale || foundOgn || foundDir || foundOgns || foundDirs;
82 
83  if (!foundAll)
84  {
85  scale_ = autoPtr<Function1<scalar>>();
86  origins_.clear();
87  directions_.clear();
88  }
89 
90  if (foundAll)
91  {
92  scale_ = Function1<scalar>::New("scale", dimLength, dimless, coeffs());
93  if (foundOgn)
94  {
95  origins_.setSize(1);
96  directions_.setSize(1);
97  coeffs().lookup("origin") >> origins_.last();
98  coeffs().lookup("direction") >> directions_.last();
99  }
100  else
101  {
102  coeffs().lookup("origins") >> origins_;
103  coeffs().lookup("directions") >> directions_;
104 
105  if
106  (
107  origins_.size() == 0
108  || directions_.size() == 0
109  || origins_.size() != directions_.size()
110  )
111  {
113  << "The same, non-zero number of origins and "
114  << "directions must be provided" << exit(FatalError);
115  }
116  }
117  forAll(directions_, i)
118  {
119  directions_[i] /= mag(directions_[i]);
120  }
121  }
122 
123  if (!foundAll && foundAny)
124  {
126  << "The scaling specification is incomplete. \"scale\", "
127  << "\"origin\" and \"direction\" (or \"origins\" and "
128  << "\"directions\"), must all be specified in order to scale "
129  << "the forcing. The forcing will be applied uniformly across "
130  << "the cell set." << endl << endl;
131  }
132 }
133 
134 
136 {
137  dimensionedScalar vs("vs", dimVolume, 0);
138  dimensionedScalar vgrads("vs", dimArea, 0);
139 
140  forAll(origins_, i)
141  {
142  const volScalarField x
143  (
144  (mesh().C() - dimensionedVector(dimLength, origins_[i]))
145  & directions_[i]
146  );
147 
148  const volScalarField scale(evaluate(*scale_, dimless, x));
149 
150  vs += fvc::domainIntegrate(scale);
151  vgrads += fvc::domainIntegrate(directions_[i] & fvc::grad(scale));
152  }
153 
154  return vs/vgrads;
155 }
156 
157 
158 
160 {
162  (
164  (
165  typedName("scale"),
166  mesh(),
168  )
169  );
170 
171  scalarField& scale = tscale.ref();
172 
173  forAll(origins_, i)
174  {
175  const vectorField& c = mesh().cellCentres();
176  const scalarField x((c - origins_[i]) & directions_[i]);
177  scale = max(scale, scale_->value(x));
178  }
179 
180  return tscale;
181 }
182 
183 
185 {
186  tmp<volScalarField::Internal> tscale(this->scale());
187  const volScalarField::Internal& scale = tscale();
188 
190  (
191  volScalarField::Internal::New(typedName("forceCoeff"), lambda_*scale)
192  );
193 
194  // Damp the cells adjacent to the boundary with lambdaBoundary if specified
195  if (lambdaBoundary_.value() > 0)
196  {
197  const fvBoundaryMesh& bm = mesh().boundary();
198 
199  forAll(bm, patchi)
200  {
201  if (!bm[patchi].coupled())
202  {
203  UIndirectList<scalar>(tforceCoeff.ref(), bm[patchi].faceCells())
204  = lambdaBoundary_.value()
205  *Field<scalar>(scale, bm[patchi].faceCells());
206  }
207  }
208  }
209 
210  return tforceCoeff;
211 }
212 
213 
215 {
216  if (writeForceFields_)
217  {
218  Info<< " Writing forcing fields: forcing:scale, forcing:forceCoeff"
219  << endl;
220 
221  scale()().write();
222  forceCoeff()().write();
223  }
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
228 
230 (
231  const word& name,
232  const word& modelType,
233  const fvMesh& mesh,
234  const dictionary& dict
235 )
236 :
237  fvModel(name, modelType, mesh, dict),
238  writeForceFields_(false),
239  lambda_("lambda", dimless/dimTime, NaN),
240  lambdaBoundary_("lambdaBoundary", dimless/dimTime, 0.0),
241  scale_(nullptr),
242  origins_(),
243  directions_()
244 {
245  readCoeffs();
246 }
247 
248 
249 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
250 
252 {
253  if (fvModel::read(dict))
254  {
255  readCoeffs();
256  return true;
257  }
258  else
259  {
260  return false;
261  }
262 }
263 
264 
265 // ************************************************************************* //
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void clear()
Delete object (if the pointer is valid) and set pointer to.
Definition: autoPtrI.H:126
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dimensionSet & dimensions() const
Return const reference to dimensions.
const word & name() const
Return const reference to name.
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
Base fvModel for forcing functions.
Definition: forcing.H:59
void readCoeffs()
Non-virtual read.
Definition: forcing.C:65
void readLambda()
Read the forcing coefficients.
Definition: forcing.C:45
void writeForceFields() const
Optionally write the forcing fields:
Definition: forcing.C:214
dimensionedScalar lambda_
Forcing coefficient [1/s].
Definition: forcing.H:68
dimensionedScalar lambdaBoundary_
Optional boundary forcing coefficient [1/s].
Definition: forcing.H:71
forcing(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: forcing.C:230
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: forcing.C:251
virtual tmp< volScalarField::Internal > forceCoeff() const
Return the force coefficient.
Definition: forcing.C:184
dimensionedScalar regionLength() const
Calculate and return the volume average forcing region length.
Definition: forcing.C:135
virtual tmp< volScalarField::Internal > scale() const
Return the scale distribution.
Definition: forcing.C:159
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the gradient of the given field.
Volume integrate volField creating a volField.
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
defineTypeNameAndDebug(bound, 0)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
const dimensionSet dimTime
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict