Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
35// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace fv
40{
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 (
60 (
61 name_ + "Properties",
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 {
98 FatalErrorInFunction
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 }
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
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_;
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 (
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();
245
246 gradP0_ += dGradP_;
247 dGradP_ = 0.0;
248}
249
250
252{
254
255 return false;
256}
257
258
259// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
IOdictionary propsDict(dictIO)
if(patchID !=-1)
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
Definition IFstream.H:55
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
bool writeTime() const noexcept
True if this is a write interval.
Definition TimeStateI.H:73
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
fileName timePath() const
Return current time path = path/timeName.
Definition Time.H:512
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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,...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
tmp< volScalarField > A() const
Return the central coefficient.
Definition fvMatrix.C:1306
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
labelList cells_
Set of cells to apply source to.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
scalar V_
Sum of cell volumes.
Applies the force within a specified region to maintain the specified mean velocity for incompressibl...
vector Ubar_
Desired mean velocity.
meanVelocityForce(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
scalar dGradP_
Change in pressure gradient.
autoPtr< volScalarField > rAPtr_
Matrix 1/A coefficients field pointer.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Set 1/A coefficient.
scalar gradP0_
Pressure gradient before correction.
scalar relaxation_
Relaxation factor.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
void writeProps(const scalar gradP) const
Write the pressure gradient to file (for restarts etc).
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity averaged over the selected cellSet.
virtual void correct(volVectorField &U)
Correct the pressure gradient.
vector flowDir_
Flow direction.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volVectorField gradP(fvc::grad(p))
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
A special matrix type and solver, designed for finite volume solutions of scalar equations.
tmp< volScalarField > rAU
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
Vector< scalar > vector
Definition vector.H:57
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299