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 | Copyright (C) 2011-2015 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 "fvMatrices.H"
28 #include "DimensionedField.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  option,
43  meanVelocityForce,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 (
54  const scalar gradP
55 ) const
56 {
57  // Only write on output time
58  if (mesh_.time().outputTime())
59  {
61  (
62  IOobject
63  (
64  name_ + "Properties",
65  mesh_.time().timeName(),
66  "uniform",
67  mesh_,
68  IOobject::NO_READ,
69  IOobject::NO_WRITE
70  )
71  );
72  propsDict.add("gradient", gradP);
73  propsDict.regIOobject::write();
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const word& sourceName,
83  const word& modelType,
84  const dictionary& dict,
85  const fvMesh& mesh
86 )
87 :
88  cellSetOption(sourceName, modelType, dict, mesh),
89  Ubar_(coeffs_.lookup("Ubar")),
90  gradP0_(0.0),
91  dGradP_(0.0),
92  flowDir_(Ubar_/mag(Ubar_)),
93  relaxation_(coeffs_.lookupOrDefault<scalar>("relaxation", 1.0)),
94  rAPtr_(NULL)
95 {
96  coeffs_.lookup("fieldNames") >> fieldNames_;
97 
98  if (fieldNames_.size() != 1)
99  {
101  << "settings are:" << fieldNames_ << exit(FatalError);
102  }
103 
104  applied_.setSize(fieldNames_.size(), false);
105 
106  // Read the initial pressure gradient from file if it exists
107  IFstream propsFile
108  (
109  mesh_.time().timePath()/"uniform"/(name_ + "Properties")
110  );
111 
112  if (propsFile.good())
113  {
114  Info<< " Reading pressure gradient from file" << endl;
115  dictionary propsDict(dictionary::null, propsFile);
116  propsDict.lookup("gradient") >> gradP0_;
117  }
118 
119  Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 (
127  const volVectorField& U
128 ) const
129 {
130  scalar magUbarAve = 0.0;
131 
132  const scalarField& cv = mesh_.V();
133  forAll(cells_, i)
134  {
135  label cellI = cells_[i];
136  scalar volCell = cv[cellI];
137  magUbarAve += (flowDir_ & U[cellI])*volCell;
138  }
139 
140  reduce(magUbarAve, sumOp<scalar>());
141 
142  magUbarAve /= V_;
143 
144  return magUbarAve;
145 }
146 
147 
149 {
150  const scalarField& rAU = rAPtr_().internalField();
151 
152  // Integrate flow variables over cell set
153  scalar rAUave = 0.0;
154  const scalarField& cv = mesh_.V();
155  forAll(cells_, i)
156  {
157  label cellI = cells_[i];
158  scalar volCell = cv[cellI];
159  rAUave += rAU[cellI]*volCell;
160  }
161 
162  // Collect across all processors
163  reduce(rAUave, sumOp<scalar>());
164 
165  // Volume averages
166  rAUave /= V_;
167 
168  scalar magUbarAve = this->magUbarAve(U);
169 
170  // Calculate the pressure gradient increment needed to adjust the average
171  // flow-rate to the desired value
172  dGradP_ = relaxation_*(mag(Ubar_) - magUbarAve)/rAUave;
173 
174  // Apply correction to velocity field
175  forAll(cells_, i)
176  {
177  label cellI = cells_[i];
178  U[cellI] += flowDir_*rAU[cellI]*dGradP_;
179  }
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  IOobject
199  (
200  name_ + fieldNames_[fieldI] + "Sup",
201  mesh_.time().timeName(),
202  mesh_,
205  ),
206  mesh_,
208  );
209 
210  scalar gradP = gradP0_ + dGradP_;
211 
212  UIndirectList<vector>(Su, cells_) = flowDir_*gradP;
213 
214  eqn += Su;
215 }
216 
217 
219 (
220  const volScalarField& rho,
221  fvMatrix<vector>& eqn,
222  const label fieldI
223 )
224 {
225  this->addSup(eqn, fieldI);
226 }
227 
228 
230 (
231  fvMatrix<vector>& eqn,
232  const label
233 )
234 {
235  if (rAPtr_.empty())
236  {
237  rAPtr_.reset
238  (
239  new volScalarField
240  (
241  IOobject
242  (
243  name_ + ":rA",
244  mesh_.time().timeName(),
245  mesh_,
248  ),
249  1.0/eqn.A()
250  )
251  );
252  }
253  else
254  {
255  rAPtr_() = 1.0/eqn.A();
256  }
257 
258  gradP0_ += dGradP_;
259  dGradP_ = 0.0;
260 }
261 
262 
263 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fv::meanVelocityForce::addSup
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add explicit contribution to momentum equation.
Definition: meanVelocityForce.C:191
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fv::cellSetOption
Cell-set options abtract base class. Provides a base set of controls, e.g.
Definition: cellSetOption.H:71
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::IFstream
Input from file stream.
Definition: IFstream.H:81
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:286
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:251
Foam::fv::meanVelocityForce::writeProps
void writeProps(const scalar gradP) const
Write the pressure gradient to file (for restarts etc)
Definition: meanVelocityForce.C:53
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:78
U
U
Definition: pEqn.H:46
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
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:43
Foam::fv::meanVelocityForce::Ubar_
vector Ubar_
Average velocity.
Definition: meanVelocityForce.H:79
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
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))
IFstream.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::fv::meanVelocityForce::constrain
virtual void constrain(fvMatrix< vector > &eqn, const label fieldI)
Set 1/A coefficient.
Definition: meanVelocityForce.C:230
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:734
gradP
dimensionedVector gradP("gradP", dimensionSet(0, 1, -2, 0, 0), vector::zero)
Foam::fv::meanVelocityForce::meanVelocityForce
meanVelocityForce(const meanVelocityForce &)
Disallow default bitwise copy construct.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
rho
rho
Definition: pEqn.H:3
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::sumOp
Definition: ops.H:162
Foam::fv::meanVelocityForce::magUbarAve
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity.
Definition: meanVelocityForce.C:126
Foam::fv::meanVelocityForce::gradP0_
scalar gradP0_
Pressure gradient before correction.
Definition: meanVelocityForce.H:82
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
meanVelocityForce.H
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::fv::meanVelocityForce::correct
virtual void correct(volVectorField &U)
Correct the pressure gradient.
Definition: meanVelocityForce.C:148