meanVelocityForce.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) 2011-2021 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 "meanVelocityForce.H"
27 #include "volFields.H"
28 #include "fvMatrices.H"
29 #include "IFstream.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
38  defineTypeNameAndDebug(meanVelocityForce, 0);
39 
41  (
42  fvConstraint,
43  meanVelocityForce,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::fv::meanVelocityForce::readCoeffs()
53 {
54  UName_ = coeffs().lookupOrDefault<word>("U", "U");
55  Ubar_ = coeffs().lookup<vector>("Ubar");
56  relaxation_ = coeffs().lookupOrDefault<scalar>("relaxation", 1);
57 }
58 
59 
60 void Foam::fv::meanVelocityForce::writeProps
61 (
62  const scalar gradP
63 ) const
64 {
65  // Only write on output time
66  if (mesh().time().writeTime())
67  {
68  IOdictionary propsDict
69  (
70  IOobject
71  (
72  name() + "Properties",
73  mesh().time().timeName(),
74  "uniform",
75  mesh(),
78  )
79  );
80  propsDict.add("gradient", gradP);
81  propsDict.regIOobject::write();
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const word& name,
91  const word& modelType,
92  const dictionary& dict,
93  const fvMesh& mesh
94 )
95 :
96  fvConstraint(name, modelType, dict, mesh),
97  set_(coeffs(), mesh),
98  UName_(word::null),
99  Ubar_(vector::uniform(NaN)),
100  relaxation_(NaN),
101  gradP0_(0),
102  dGradP_(0),
103  rAPtr_(nullptr)
104 {
105  readCoeffs();
106 
107  // Read the initial pressure gradient from file if it exists
108  IFstream propsFile
109  (
110  mesh.time().timePath()/"uniform"/(this->name() + "Properties")
111  );
112  if (propsFile.good())
113  {
114  Info<< " Reading pressure gradient from file" << endl;
116  propsDict.lookup("gradient") >> gradP0_;
117  }
118 
119  Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 {
127  return wordList(1, UName_);
128 }
129 
130 
131 Foam::scalar Foam::fv::meanVelocityForce::magUbarAve
132 (
133  const volVectorField& U
134 ) const
135 {
136  const labelList& cells = set_.cells();
137  const scalarField& cv = mesh().V();
138 
139  scalar magUbarAve = 0;
140  forAll(cells, i)
141  {
142  const label celli = cells[i];
143  magUbarAve += (normalised(Ubar_) & U[celli])*cv[celli];
144  }
145  reduce(magUbarAve, sumOp<scalar>());
146  magUbarAve /= set_.V();
147 
148  return magUbarAve;
149 }
150 
151 
153 (
154  fvMatrix<vector>& eqn,
155  const word& fieldName
156 ) const
157 {
159  (
160  IOobject
161  (
162  name() + fieldName + "Sup",
163  mesh().time().timeName(),
164  mesh(),
167  ),
168  mesh(),
170  );
171 
172  const scalar gradP = gradP0_ + dGradP_;
173 
174  UIndirectList<vector>(Su, set_.cells()) = normalised(Ubar_)*gradP;
175 
176  eqn -= Su;
177 
178  if (rAPtr_.empty())
179  {
180  rAPtr_.reset
181  (
182  new volScalarField
183  (
184  IOobject
185  (
186  name() + ":rA",
187  mesh().time().timeName(),
188  mesh(),
191  false
192  ),
193  1/eqn.A()
194  )
195  );
196  }
197  else
198  {
199  rAPtr_() = 1/eqn.A();
200  }
201 
202  gradP0_ += dGradP_;
203  dGradP_ = 0;
204 
205  return true;
206 }
207 
208 
210 {
211  const scalarField& rAU = rAPtr_();
212 
213  const labelList& cells = set_.cells();
214  const scalarField& cv = mesh().V();
215 
216  // Average rAU over the cell set
217  scalar rAUave = 0;
218  forAll(cells, i)
219  {
220  const label celli = cells[i];
221  rAUave += rAU[celli]*cv[celli];
222  }
223  reduce(rAUave, sumOp<scalar>());
224  rAUave /= set_.V();
225 
226  const scalar magUbarAve = this->magUbarAve(U);
227 
228  // Calculate the pressure gradient increment needed to adjust the average
229  // flow-rate to the desired value
230  dGradP_ = relaxation_*(mag(Ubar_) - magUbarAve)/rAUave;
231 
232  // Apply correction to velocity field
233  forAll(cells, i)
234  {
235  label celli = cells[i];
236  U[celli] += normalised(Ubar_)*rAU[celli]*dGradP_;
237  }
238 
239  const scalar gradP = gradP0_ + dGradP_;
240 
241  Info<< "Pressure gradient source: uncorrected Ubar = " << magUbarAve
242  << ", pressure gradient = " << gradP << endl;
243 
244  writeProps(gradP);
245 
246  return true;
247 }
248 
249 
251 {
252  set_.updateMesh(mpm);
253 }
254 
255 
257 {
258  if (fvConstraint::read(dict))
259  {
260  set_.read(coeffs());
261  return true;
262  }
263  else
264  {
265  return false;
266  }
267 }
268 
269 
270 // ************************************************************************* //
IOdictionary propsDict(systemDict("particleTracksDict", args, runTime))
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
fileName timePath() const
Return current time path.
Definition: Time.H:279
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual bool constrain(fvMatrix< vector > &eqn, const word &fieldName) const
Add the momentum source and set the 1/A coefficient.
dimensionedVector gradP("gradP", dimensionSet(0, 1, -2, 0, 0), Zero)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
meanVelocityForce(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
dynamicFvMesh & mesh
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
const cellShapeList & cells
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static const dictionary null
Null dictionary.
Definition: dictionary.H:242
static const word null
An empty word.
Definition: word.H:77
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:72
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
virtual void updateMesh(const mapPolyMesh &)
Update for mesh changes.
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
static const char nl
Definition: Ostream.H:260
Input from file stream.
Definition: IFstream.H:81
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:738
tmp< volScalarField > rAU
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet dimVolume
messageStream Info
static Vector< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
dimensioned< scalar > mag(const dimensioned< Type > &)
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:143
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Finite volume options abstract base class.
Definition: fvConstraint.H:54
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
const dimensionSet & dimensions() const
Definition: fvMatrix.H:287