Loading...
Searching...
No Matches
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
41(
42 const fvPatch& p,
45:
46 fixedValueFvPatchVectorField(p, iF),
47 phiName_("phi")
48{}
49
50
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
78(
79 const fvPatch& p,
80 const DimensionedField<vector, volMesh>& iF,
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 );
103 if (mpp.mode() == mappedPolyPatch::NEARESTCELL)
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
118(
119 const mappedVelocityFluxFixedValueFvPatchField& ptf
121:
122 fixedValueFvPatchVectorField(ptf),
123 phiName_(ptf.phiName_)
124{}
125
126
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 const int oldTag = UPstream::incrMsgType();
151
152 // Get the mappedPatchBase
153 const mappedPatchBase& mpp = refCast<const mappedPatchBase>
154 (
155 mappedVelocityFluxFixedValueFvPatchField::patch().patch()
156 );
157 const fvMesh& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
158 const word& fieldName = internalField().name();
159 const volVectorField& UField =
160 nbrMesh.lookupObject<volVectorField>(fieldName);
161
162 const surfaceScalarField& phiField =
163 nbrMesh.lookupObject<surfaceScalarField>(phiName_);
164
165 vectorField newUValues;
166 scalarField newPhiValues;
167
168 switch (mpp.mode())
169 {
171 {
172 vectorField allUValues(nbrMesh.nFaces(), Zero);
173 scalarField allPhiValues(nbrMesh.nFaces(), Zero);
174
175 forAll(UField.boundaryField(), patchi)
176 {
177 const fvPatchVectorField& Upf = UField.boundaryField()[patchi];
178 const scalarField& phipf = phiField.boundaryField()[patchi];
179
180 label faceStart = Upf.patch().start();
181
182 forAll(Upf, facei)
183 {
184 allUValues[faceStart + facei] = Upf[facei];
185 allPhiValues[faceStart + facei] = phipf[facei];
186 }
187 }
188
189 mpp.distribute(allUValues);
190 newUValues.transfer(allUValues);
191
192 mpp.distribute(allPhiValues);
193 newPhiValues.transfer(allPhiValues);
194
195 break;
196 }
199 {
200 const label nbrPatchID =
201 nbrMesh.boundaryMesh().findPatchID(mpp.samplePatch());
202
203 newUValues = UField.boundaryField()[nbrPatchID];
204 mpp.distribute(newUValues);
205
206 newPhiValues = phiField.boundaryField()[nbrPatchID];
207 mpp.distribute(newPhiValues);
208
209 break;
210 }
211 default:
212 {
214 << "patch can only be used in NEARESTPATCHFACE, "
215 << "NEARESTPATCHFACEAMI or NEARESTFACE mode" << nl
216 << abort(FatalError);
217 }
218 }
219
220 operator==(newUValues);
221 phiField.constCast().boundaryFieldRef()[patch().index()] == newPhiValues;
223 UPstream::msgType(oldTag); // Restore tag
224
225 fixedValueFvPatchVectorField::updateCoeffs();
226}
227
228
230(
231 Ostream& os
232) const
233{
235 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
238
239
240// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241
242namespace Foam
243{
245 (
247 mappedVelocityFluxFixedValueFvPatchField
248 );
249}
250
251
252// ************************************************************************* //
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 noexcept
Return const-reference to the boundary field.
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvPatch & patch() const noexcept
Return the patch.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
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
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
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
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.
mappedVelocityFluxFixedValueFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
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 noexcept
Return boundary mesh.
Definition polyMesh.H:609
label nFaces() const noexcept
Number of mesh faces.
A class for handling words, derived from Foam::string.
Definition word.H:66
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,...
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
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...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< vector > fvPatchVectorField
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
Foam::surfaceFields.