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 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 
44  (
45  option,
46  meanVelocityForce,
47  dictionary
48  );
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
56 (
57  const scalar gradP
58 ) const
59 {
60  // Only write on output time
61  if (mesh_.time().writeTime())
62  {
64  (
65  IOobject
66  (
67  name_ + "Properties",
68  mesh_.time().timeName(),
69  "uniform",
70  mesh_,
71  IOobject::NO_READ,
72  IOobject::NO_WRITE
73  )
74  );
75  propsDict.add("gradient", gradP);
76  propsDict.regIOobject::write();
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
83 Foam::fv::meanVelocityForce::meanVelocityForce
84 (
85  const word& sourceName,
86  const word& modelType,
87  const dictionary& dict,
88  const fvMesh& mesh
89 )
90 :
91  cellSetOption(sourceName, modelType, dict, mesh),
92  Ubar_(coeffs_.get<vector>("Ubar")),
93  gradP0_(0.0),
94  dGradP_(0.0),
95  flowDir_(Ubar_/mag(Ubar_)),
96  relaxation_(coeffs_.lookupOrDefault<scalar>("relaxation", 1.0)),
97  rAPtr_(nullptr)
98 {
99  coeffs_.readEntry("fields", fieldNames_);
100 
101  if (fieldNames_.size() != 1)
102  {
104  << "settings are:" << fieldNames_ << exit(FatalError);
105  }
106 
107  applied_.setSize(fieldNames_.size(), false);
108 
109  // Read the initial pressure gradient from file if it exists
110  IFstream propsFile
111  (
112  mesh_.time().timePath()/"uniform"/(name_ + "Properties")
113  );
114 
115  if (propsFile.good())
116  {
117  Info<< " Reading pressure gradient from file" << endl;
118  dictionary propsDict(dictionary::null, propsFile);
119  propsDict.readEntry("gradient", gradP0_);
120  }
121 
122  Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 (
130  const volVectorField& U
131 ) const
132 {
133  scalar magUbarAve = 0.0;
134 
135  const scalarField& cv = mesh_.V();
136  forAll(cells_, i)
137  {
138  label celli = cells_[i];
139  scalar volCell = cv[celli];
140  magUbarAve += (flowDir_ & U[celli])*volCell;
141  }
142 
143  reduce(magUbarAve, sumOp<scalar>());
144 
145  magUbarAve /= V_;
146 
147  return magUbarAve;
148 }
149 
150 
152 {
153  const scalarField& rAU = rAPtr_();
154 
155  // Integrate flow variables over cell set
156  scalar rAUave = 0.0;
157  const scalarField& cv = mesh_.V();
158  forAll(cells_, i)
159  {
160  label celli = cells_[i];
161  scalar volCell = cv[celli];
162  rAUave += rAU[celli]*volCell;
163  }
164 
165  // Collect across all processors
166  reduce(rAUave, sumOp<scalar>());
167 
168  // Volume averages
169  rAUave /= V_;
170 
171  scalar magUbarAve = this->magUbarAve(U);
172 
173  // Calculate the pressure gradient increment needed to adjust the average
174  // flow-rate to the desired value
175  dGradP_ = relaxation_*(mag(Ubar_) - magUbarAve)/rAUave;
176 
177  // Apply correction to velocity field
178  forAll(cells_, i)
179  {
180  label celli = cells_[i];
181  U[celli] += flowDir_*rAU[celli]*dGradP_;
182  }
183 
184  U.correctBoundaryConditions();
185 
186  scalar gradP = gradP0_ + dGradP_;
187 
188  Info<< "Pressure gradient source: uncorrected Ubar = " << magUbarAve
189  << ", pressure gradient = " << gradP << endl;
190 
191  writeProps(gradP);
192 }
193 
194 
196 (
197  fvMatrix<vector>& eqn,
198  const label fieldi
199 )
200 {
202  (
203  IOobject
204  (
205  name_ + fieldNames_[fieldi] + "Sup",
206  mesh_.time().timeName(),
207  mesh_,
210  ),
211  mesh_,
213  );
214 
215  scalar gradP = gradP0_ + dGradP_;
216 
217  UIndirectList<vector>(Su, cells_) = flowDir_*gradP;
218 
219  eqn += Su;
220 }
221 
222 
224 (
225  const volScalarField& rho,
226  fvMatrix<vector>& eqn,
227  const label fieldi
228 )
229 {
230  this->addSup(eqn, fieldi);
231 }
232 
233 
235 (
236  fvMatrix<vector>& eqn,
237  const label
238 )
239 {
240  if (rAPtr_.empty())
241  {
242  rAPtr_.reset
243  (
244  new volScalarField
245  (
246  IOobject
247  (
248  name_ + ":rA",
249  mesh_.time().timeName(),
250  mesh_,
253  ),
254  1.0/eqn.A()
255  )
256  );
257  }
258  else
259  {
260  rAPtr_() = 1.0/eqn.A();
261  }
262 
263  gradP0_ += dGradP_;
264  dGradP_ = 0.0;
265 }
266 
267 
269 {
271 
272  return false;
273 }
274 
275 
276 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
gradP
volVectorField gradP(fvc::grad(p))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Cell-set options abtract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:72
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::IFstream
Input from file stream, using an ISstream.
Definition: IFstream.H:97
DimensionedField.H
Foam::fv::cellSetOption::V_
scalar V_
Sum of cell volumes.
Definition: cellSetOption.H:116
Foam::fv::meanVelocityForce::rAPtr_
autoPtr< volScalarField > rAPtr_
Matrix 1/A coefficients field pointer.
Definition: meanVelocityForce.H:94
Foam::fvMatrix::dimensions
const dimensionSet & dimensions() const
Definition: fvMatrix.H:290
Foam::fv::meanVelocityForce::relaxation_
scalar relaxation_
Relaxation factor.
Definition: meanVelocityForce.H:91
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::fv::meanVelocityForce::writeProps
void writeProps(const scalar gradP) const
Write the pressure gradient to file (for restarts etc)
Definition: meanVelocityForce.C:56
rho
rho
Definition: readInitialConditions.H:96
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:82
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Su
zeroField Su
Definition: alphaSuSp.H:1
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:419
Foam::fv::meanVelocityForce::flowDir_
vector flowDir_
Flow direction.
Definition: meanVelocityForce.H:88
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::fv::meanVelocityForce::Ubar_
vector Ubar_
Average velocity.
Definition: meanVelocityForce.H:79
Foam::label
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:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fv::meanVelocityForce::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: meanVelocityForce.C:268
Foam::fv::meanVelocityForce::dGradP_
scalar dGradP_
Change in pressure gradient.
Definition: meanVelocityForce.H:85
Foam::fv::cellSetOption::cells_
labelList cells_
Set of cells to apply source to.
Definition: cellSetOption.H:113
propsDict
IOdictionary propsDict(IOobject("particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED))
Foam::fv::meanVelocityForce::addSup
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
Definition: meanVelocityForce.C:196
IFstream.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:788
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::fv::meanVelocityForce::constrain
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Set 1/A coefficient.
Definition: meanVelocityForce.C:235
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::Vector< scalar >
Foam::fv::meanVelocityForce::magUbarAve
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity.
Definition: meanVelocityForce.C:129
Foam::fv::meanVelocityForce::gradP0_
scalar gradP0_
Pressure gradient before correction.
Definition: meanVelocityForce.H:82
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
meanVelocityForce.H
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:216
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fv::meanVelocityForce::correct
virtual void correct(volVectorField &U)
Correct the pressure gradient.
Definition: meanVelocityForce.C:151
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:190