Loading...
Searching...
No Matches
fixedFluxPressureFvPatchScalarField.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) 2015-2020 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
31#include "volFields.H"
32#include "surfaceFields.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
38(
39 const fvPatch& p,
42:
43 fixedGradientFvPatchScalarField(p, iF),
44 curTimeIndex_(-1)
45{}
46
47
49(
50 const fvPatch& p,
52 const dictionary& dict
53)
54:
55 fixedGradientFvPatchScalarField(p, iF), // Bypass dictionary constructor
56 curTimeIndex_(-1)
57{
59
60 if (!this->readGradientEntry(dict) || !this->readValueEntry(dict))
61 {
62 extrapolateInternal();
63 gradient() = Zero;
64 }
65}
66
67
69(
70 const fixedFluxPressureFvPatchScalarField& ptf,
71 const fvPatch& p,
72 const DimensionedField<scalar, volMesh>& iF,
73 const fvPatchFieldMapper& mapper
74)
75:
76 fixedGradientFvPatchScalarField(p, iF),
77 curTimeIndex_(-1)
78{
79 patchType() = ptf.patchType();
80
81 // Map gradient. Set unmapped values and overwrite with mapped ptf
82 gradient() = Zero;
83 gradient().map(ptf.gradient(), mapper);
84
85 // Evaluate the value field from the gradient if the internal field is valid
86 if (notNull(iF))
87 {
88 if (iF.size())
89 {
90 // Note: cannot ask for nf() if zero faces
91
92 scalarField::operator=
93 (
94 //patchInternalField() + gradient()/patch().deltaCoeffs()
95 // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
96 // which fails for AMI patches for some mapping operations
97 patchInternalField()
98 + gradient()*(patch().nf() & patch().delta())
99 );
100 }
101 }
102 else
104 // Enforce mapping of values so we have a valid starting value
105 this->map(ptf, mapper);
106 }
107}
108
109
111(
112 const fixedFluxPressureFvPatchScalarField& wbppsf
114:
115 fixedGradientFvPatchScalarField(wbppsf),
116 curTimeIndex_(-1)
117{}
118
119
121(
124)
125:
126 fixedGradientFvPatchScalarField(wbppsf, iF),
127 curTimeIndex_(-1)
128{}
129
130
131// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132
134(
135 const scalarField& snGradp
136)
137{
138 if (updated())
139 {
140 return;
141 }
142
143 curTimeIndex_ = this->db().time().timeIndex();
144
145 gradient() = snGradp;
146 fixedGradientFvPatchScalarField::updateCoeffs();
147}
148
149
151{
152 if (updated())
153 {
154 return;
155 }
156
157 if (curTimeIndex_ != this->db().time().timeIndex())
158 {
160 << "updateCoeffs(const scalarField& snGradp) MUST be called before"
161 " updateCoeffs() or evaluate() to set the boundary gradient."
162 << exit(FatalError);
163 }
164}
165
166
168{
172
173
174// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175
176namespace Foam
177{
179 (
181 fixedFluxPressureFvPatchScalarField
182 );
183}
184
185
186// ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
fixedFluxPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the patch pressure gradient field.
virtual void updateSnGrad(const scalarField &snGradp)
Update the patch pressure gradient field from the given snGradp.
virtual void write(Ostream &) const
Write.
virtual void readDict(const dictionary &dict)
Read dictionary entries.
A FieldMapper for finite-volume patch fields.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition nullObject.H:267
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
label timeIndex
dictionary dict
Foam::surfaceFields.