meanVelocityForce.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2023 OpenCFD Ltd.
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 "meanVelocityForce.H"
30 #include "fvMatrices.H"
31 #include "DimensionedField.H"
32 #include "IFstream.H"
34 
35 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fv
40 {
41  defineTypeNameAndDebug(meanVelocityForce, 0);
42  addToRunTimeSelectionTable(option, meanVelocityForce, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 (
51  const scalar gradP
52 ) const
53 {
54  // Only write on output time
55  if (mesh_.time().writeTime())
56  {
58  (
59  IOobject
60  (
61  name_ + "Properties",
62  mesh_.time().timeName(),
63  "uniform",
64  mesh_,
68  )
69  );
70  propsDict.add("gradient", gradP);
71  propsDict.regIOobject::write();
72  }
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
79 (
80  const word& sourceName,
81  const word& modelType,
82  const dictionary& dict,
83  const fvMesh& mesh
84 )
85 :
86  fv::cellSetOption(sourceName, modelType, dict, mesh),
87  Ubar_(coeffs_.get<vector>("Ubar")),
88  gradP0_(0.0),
89  dGradP_(0.0),
90  flowDir_(Ubar_/mag(Ubar_)),
91  relaxation_(coeffs_.getOrDefault<scalar>("relaxation", 1)),
92  rAPtr_(nullptr)
93 {
94  coeffs_.readEntry("fields", fieldNames_);
95 
96  if (fieldNames_.size() != 1)
97  {
99  << "settings are:" << fieldNames_ << exit(FatalError);
100  }
101 
103 
104  // Read the initial pressure gradient from file if it exists
105  IFstream propsFile
106  (
107  mesh_.time().timePath()/"uniform"/(name_ + "Properties")
108  );
109 
110  if (propsFile.good())
111  {
112  Info<< " Reading pressure gradient from file" << endl;
113  dictionary propsDict(propsFile);
114  propsDict.readEntry("gradient", gradP0_);
115  }
116 
117  Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
118 }
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 (
125  const volVectorField& U
126 ) const
127 {
128  scalar magUbarAve = 0.0;
129 
130  const scalarField& cv = mesh_.V();
131  forAll(cells_, i)
132  {
133  label celli = cells_[i];
134  scalar volCell = cv[celli];
135  magUbarAve += (flowDir_ & U[celli])*volCell;
136  }
137 
138  reduce(magUbarAve, sumOp<scalar>());
140  magUbarAve /= V_;
141 
142  return magUbarAve;
143 }
144 
145 
147 {
148  const scalarField& rAU = rAPtr_();
149 
150  // Integrate flow variables over cell set
151  scalar rAUave = 0.0;
152  const scalarField& cv = mesh_.V();
153  forAll(cells_, i)
154  {
155  label celli = cells_[i];
156  scalar volCell = cv[celli];
157  rAUave += rAU[celli]*volCell;
158  }
159 
160  // Collect across all processors
161  reduce(rAUave, sumOp<scalar>());
162 
163  // Volume averages
164  rAUave /= V_;
165 
166  scalar magUbarAve = this->magUbarAve(U);
167 
168  // Calculate the pressure gradient increment needed to adjust the average
169  // flow-rate to the desired value
170  dGradP_ = relaxation_*(mag(Ubar_) - magUbarAve)/rAUave;
171 
172  // Apply correction to velocity field
173  forAll(cells_, i)
174  {
175  label celli = cells_[i];
176  U[celli] += flowDir_*rAU[celli]*dGradP_;
177  }
178 
179  U.correctBoundaryConditions();
180 
181  scalar gradP = gradP0_ + dGradP_;
182 
183  Info<< "Pressure gradient source: uncorrected Ubar = " << magUbarAve
184  << ", pressure gradient = " << gradP << endl;
185 
186  writeProps(gradP);
187 }
188 
189 
191 (
192  fvMatrix<vector>& eqn,
193  const label fieldi
194 )
195 {
197  (
198  name_ + fieldNames_[fieldi] + "Sup",
200  mesh_,
202  );
203  auto& Su = tSu.ref();
204 
205  scalar gradP = gradP0_ + dGradP_;
206 
207  UIndirectList<vector>(Su, cells_) = flowDir_*gradP;
208 
209  eqn += tSu;
210 }
211 
212 
214 (
215  const volScalarField& rho,
216  fvMatrix<vector>& eqn,
217  const label fieldi
218 )
219 {
220  this->addSup(eqn, fieldi);
221 }
222 
223 
225 (
226  fvMatrix<vector>& eqn,
227  const label
228 )
229 {
230  if (!rAPtr_)
231  {
232  rAPtr_.reset
233  (
234  new volScalarField
235  (
236  mesh_.newIOobject(IOobject::scopedName(name_, "rA")),
237  1.0/eqn.A()
238  )
239  );
240  }
241  else
242  {
243  rAPtr_() = 1.0/eqn.A();
244  }
245 
246  gradP0_ += dGradP_;
247  dGradP_ = 0.0;
248 }
249 
250 
252 {
254 
255  return false;
256 }
257 
258 
259 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity averaged over the selected cellSet...
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual void correct(volVectorField &U)
Correct the pressure gradient.
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
fileName timePath() const
Return current time path = path/timeName.
Definition: Time.H:520
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual bool read(const dictionary &dict)
Read source dictionary.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
bool writeTime() const noexcept
True if this is a write interval.
Definition: TimeStateI.H:73
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Macros for easy insertion into run-time selection tables.
IOdictionary propsDict(dictIO)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
volVectorField gradP(fvc::grad(p))
void writeProps(const scalar gradP) const
Write the pressure gradient to file (for restarts etc)
meanVelocityForce(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
dynamicFvMesh & mesh
const word name_
Source name.
Definition: fvOption.H:132
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
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
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1307
U
Definition: pEqn.H:72
const dimensionSet & dimensions() const noexcept
Definition: fvMatrix.H:528
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field...
Nothing to be read.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
Information stream (stdout output on master, null elsewhere)
scalar gradP0_
Pressure gradient before correction.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:696
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
Do not request registration (bool: false)
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Namespace for OpenFOAM.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Set 1/A coefficient.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127