prghCyclicPressureFvPatchScalarField.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) 2024 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 
28 #include "fieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "fvcSnGrad.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  jumpCyclicFvPatchScalarField(p, iF, dict),
43  rhoName_
44  (
45  cyclicPatch().owner()
46  ? dict.lookupOrDefault<word>("rho", "rho")
47  : word::null
48  ),
49  rhoInf_
50  (
51  cyclicPatch().owner()
52  ? dict.lookup<scalar>("rhoInf", dimDensity)
53  : NaN
54  ),
55  jump_(p.size(), Zero)
56 {
57  if (dict.found("jump"))
58  {
59  jump_ = scalarField("jump", iF.dimensions(), dict, p.size());
60  }
61 
62  evaluateNoUpdateCoeffs();
63 }
64 
65 
67 (
69  const fvPatch& p,
71  const fieldMapper& mapper
72 )
73 :
74  jumpCyclicFvPatchScalarField(ptf, p, iF, mapper),
75  rhoName_(ptf.rhoName_),
76  rhoInf_(ptf.rhoInf_),
77  jump_(mapper(ptf.jump_))
78 {}
79 
80 
82 (
85 )
86 :
87  jumpCyclicFvPatchScalarField(ptf, iF),
88  rhoName_(ptf.rhoName_),
89  rhoInf_(ptf.rhoInf_),
90  jump_(ptf.jump_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
98 {
99  return jump_;
100 }
101 
102 
104 (
105  const fvPatchScalarField& ptf,
106  const fieldMapper& mapper
107 )
108 {
109  jumpCyclicFvPatchScalarField::map(ptf, mapper);
110 
112  refCast<const prghCyclicPressureFvPatchScalarField>(ptf);
113 
114  mapper(jump_, tiptf.jump_);
115 }
116 
117 
119 (
120  const fvPatchScalarField& ptf
121 )
122 {
123  jumpCyclicFvPatchScalarField::reset(ptf);
124 
126  refCast<const prghCyclicPressureFvPatchScalarField>(ptf);
127 
128  jump_.reset(tiptf.jump_);
129 }
130 
131 
133 {
134  if (updated())
135  {
136  return;
137  }
138 
139  const volScalarField& vf =
140  static_cast<const volScalarField&>(internalField());
141 
143  refCast<const prghCyclicPressureFvPatchScalarField>(nbrPatchField());
144 
145  const label patchi = patch().index(), nbrPatchi = nbrPf.patch().index();
146 
147  // Buoyancy fields
148  const volScalarField& rhoVf =
149  db().lookupObject<volScalarField>
150  (cyclicPatch().owner() ? rhoName_ : nbrPf.rhoName_);
151  const volScalarField::Boundary& rhoBf = rhoVf.boundaryField();
152  const surfaceScalarField::Boundary& ghfBf =
153  db().lookupObject<surfaceScalarField>("ghf").boundaryField();
154 
155  // Pressure solution fields
156  const surfaceScalarField::Boundary& rAUfBf =
157  db().lookupObject<surfaceScalarField>("rAUf").boundaryField();
158  const surfaceScalarField::Boundary& phiHbyABf =
159  db().lookupObject<surfaceScalarField>("phiHbyA").boundaryField();
160 
161  // Delta coefficients for this field
162  const tmp<surfaceScalarField> deltaCoeffsSf =
164  (
165  vf.mesh(),
166  vf.mesh().schemes().snGrad(internalField().name())
167  )->deltaCoeffs(vf);
168  const scalarField& deltaCoeffsPf =
169  deltaCoeffsSf->boundaryField()[patch().index()];
170 
171  // Calculate the jump
172  jump_ =
173  (rhoBf[patchi] - (cyclicPatch().owner() ? rhoInf_ : nbrPf.rhoInf_))
174  *(ghfBf[nbrPatchi] - ghfBf[patchi])
175  + (
176  phiHbyABf[patchi]/rAUfBf[patchi]/patch().magSf()
177  + phiHbyABf[nbrPatchi]/rAUfBf[nbrPatchi]/nbrPf.patch().magSf()
178  )*patch().weights()/deltaCoeffsPf;
179 
180  jumpCyclicFvPatchScalarField::updateCoeffs();
181 }
182 
183 
185 (
186  Ostream& os
187 ) const
188 {
190 
191  if (cyclicPatch().owner())
192  {
193  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
194  writeEntry(os, "rhoInf", rhoInf_);
195  }
196 
197  writeEntry(os, "jump", jump_);
198 }
199 
200 
201 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
202 
203 namespace Foam
204 {
206  (
209  );
210 }
211 
212 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
static tmp< snGradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Definition: snGradScheme.C:46
This boundary condition provides a cyclic condition for p_rgh. It applies corrections to the value an...
virtual tmp< scalarField > jump() const
Return the "jump".
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
prghCyclicPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the patch pressure gradient field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the snGrad of the given volField.
label patchi
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.