mappedVelocityFluxFixedValueFvPatchField.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) 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
30#include "fvPatchFieldMapper.H"
31#include "mappedPatchBase.H"
32#include "volFields.H"
33#include "surfaceFields.H"
35
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
39Foam::mappedVelocityFluxFixedValueFvPatchField::
40mappedVelocityFluxFixedValueFvPatchField
41(
42 const fvPatch& p,
44)
45:
46 fixedValueFvPatchVectorField(p, iF),
47 phiName_("phi")
48{}
49
50
51Foam::mappedVelocityFluxFixedValueFvPatchField::
52mappedVelocityFluxFixedValueFvPatchField
53(
55 const fvPatch& p,
57 const fvPatchFieldMapper& mapper
58)
59:
60 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
61 phiName_(ptf.phiName_)
62{
63 if (!isA<mappedPatchBase>(this->patch().patch()))
64 {
66 << "Patch type '" << p.type()
67 << "' not type '" << mappedPatchBase::typeName << "'"
68 << " for patch " << p.name()
69 << " of field " << internalField().name()
70 << " in file " << internalField().objectPath()
71 << exit(FatalError);
72 }
73}
74
75
76Foam::mappedVelocityFluxFixedValueFvPatchField::
77mappedVelocityFluxFixedValueFvPatchField
78(
79 const fvPatch& p,
81 const dictionary& dict
82)
83:
84 fixedValueFvPatchVectorField(p, iF, dict),
85 phiName_(dict.getOrDefault<word>("phi", "phi"))
86{
87 if (!isA<mappedPatchBase>(this->patch().patch()))
88 {
90 << "Patch type '" << p.type()
91 << "' not type '" << mappedPatchBase::typeName << "'"
92 << " for patch " << p.name()
93 << " of field " << internalField().name()
94 << " in file " << internalField().objectPath()
95 << exit(FatalError);
96 }
97
98 const mappedPatchBase& mpp = refCast<const mappedPatchBase>
99 (
100 this->patch().patch(),
101 dict
102 );
104 {
106 << "Patch " << p.name()
107 << " of type '" << p.type()
108 << "' can not be used in 'nearestCell' mode"
109 << " of field " << internalField().name()
110 << " in file " << internalField().objectPath()
111 << exit(FatalError);
112 }
113}
114
115
116Foam::mappedVelocityFluxFixedValueFvPatchField::
117mappedVelocityFluxFixedValueFvPatchField
118(
120)
121:
122 fixedValueFvPatchVectorField(ptf),
123 phiName_(ptf.phiName_)
124{}
125
126
127Foam::mappedVelocityFluxFixedValueFvPatchField::
128mappedVelocityFluxFixedValueFvPatchField
129(
132)
133:
134 fixedValueFvPatchVectorField(ptf, iF),
135 phiName_(ptf.phiName_)
136{}
137
138
139// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140
142{
143 if (updated())
144 {
145 return;
146 }
147
148 // Since we're inside initEvaluate/evaluate there might be processor
149 // comms underway. Change the tag we use.
150 int oldTag = UPstream::msgType();
151 UPstream::msgType() = oldTag+1;
152
153 // Get the mappedPatchBase
154 const mappedPatchBase& mpp = refCast<const mappedPatchBase>
155 (
157 );
158 const fvMesh& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
159 const word& fieldName = internalField().name();
160 const volVectorField& UField =
161 nbrMesh.lookupObject<volVectorField>(fieldName);
162
163 const surfaceScalarField& phiField =
164 nbrMesh.lookupObject<surfaceScalarField>(phiName_);
165
166 vectorField newUValues;
167 scalarField newPhiValues;
168
169 switch (mpp.mode())
170 {
172 {
173 vectorField allUValues(nbrMesh.nFaces(), Zero);
174 scalarField allPhiValues(nbrMesh.nFaces(), Zero);
175
176 forAll(UField.boundaryField(), patchi)
177 {
178 const fvPatchVectorField& Upf = UField.boundaryField()[patchi];
179 const scalarField& phipf = phiField.boundaryField()[patchi];
180
181 label faceStart = Upf.patch().start();
182
183 forAll(Upf, facei)
184 {
185 allUValues[faceStart + facei] = Upf[facei];
186 allPhiValues[faceStart + facei] = phipf[facei];
187 }
188 }
189
190 mpp.distribute(allUValues);
191 newUValues.transfer(allUValues);
192
193 mpp.distribute(allPhiValues);
194 newPhiValues.transfer(allPhiValues);
195
196 break;
197 }
200 {
201 const label nbrPatchID =
202 nbrMesh.boundaryMesh().findPatchID(mpp.samplePatch());
203
204 newUValues = UField.boundaryField()[nbrPatchID];
205 mpp.distribute(newUValues);
206
207 newPhiValues = phiField.boundaryField()[nbrPatchID];
208 mpp.distribute(newPhiValues);
209
210 break;
211 }
212 default:
213 {
215 << "patch can only be used in NEARESTPATCHFACE, "
216 << "NEARESTPATCHFACEAMI or NEARESTFACE mode" << nl
217 << abort(FatalError);
218 }
219 }
220
221 operator==(newUValues);
222 const_cast<surfaceScalarField&>
223 (
224 phiField
225 ).boundaryFieldRef()[patch().index()] == newPhiValues;
226
227 // Restore tag
228 UPstream::msgType() = oldTag;
229
230 fixedValueFvPatchVectorField::updateCoeffs();
231}
232
233
235(
236 Ostream& os
237) const
238{
240 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
241 this->writeEntry("value", os);
242}
243
244
245// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246
247namespace Foam
248{
250 (
253 );
254}
255
256
257// ************************************************************************* //
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 Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void transfer(List< T > &list)
Definition: List.C:447
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A FieldMapper for finite-volume patch fields.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:179
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
@ NEARESTCELL
nearest cell containing sample
@ NEARESTPATCHFACE
nearest face on selected patch
@ NEARESTPATCHFACEAMI
nearest patch face + AMI interpolation
@ NEARESTFACE
nearest face
const polyMesh & sampleMesh() const
Get the region mesh.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
sampleMode mode() const noexcept
What to sample.
This boundary condition maps the velocity and flux from a neighbour patch to this patch.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const Type & lookupObject(const word &name, const bool recursive=false) const
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
label nFaces() const noexcept
Number of mesh faces.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Rectangular patch, selectable as patch.
static const char *const typeName
The type name used in ensight case files.
Foam::surfaceFields.