Loading...
Searching...
No Matches
mappedFixedInternalValueFvPatchField.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
29#include "IndirectList.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
36(
37 const fvPatch& p,
50 const fvPatch& p,
52 const fvPatchFieldMapper& mapper
63 const fvPatch& p,
65 const dictionary& dict
67:
69{}
70
71
72template<class Type>
75(
86(
89)
90:
92{}
93
94
95// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96
97template<class Type>
99{
101
102 if (this->updated())
103 {
104 return;
105 }
106
107 // Since we're inside initEvaluate/evaluate there might be processor
108 // comms underway. Change the tag we use.
109 const int oldTag = UPstream::incrMsgType();
110
111 // Retrieve the neighbour values and assign to this patch boundary field
113
114 // Get the coupling information from the mappedPatchBase
115 const mappedPatchBase& mpp =
116 refCast<const mappedPatchBase>(this->patch().patch());
117 const fvMesh& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
118
119 Field<Type> nbrIntFld;
120
121 switch (mpp.mode())
122 {
124 {
126 << "Cannot apply "
128 [
130 ]
131 << " mapping mode for patch " << this->patch().name()
132 << exit(FatalError);
133
134 break;
135 }
138 {
139 const label samplePatchi = mpp.samplePolyPatch().index();
140 const fvPatchField<Type>& nbrPatchField =
141 this->sampleField().boundaryField()[samplePatchi];
142 nbrIntFld = nbrPatchField.patchInternalField();
143 mpp.distribute(nbrIntFld);
144
145 break;
146 }
148 {
149 Field<Type> allValues(nbrMesh.nFaces(), Zero);
150
151 const FieldType& nbrField = this->sampleField();
152
153 forAll(nbrField.boundaryField(), patchi)
154 {
155 const fvPatchField<Type>& pf = nbrField.boundaryField()[patchi];
156 const Field<Type> pif(pf.patchInternalField());
157
158 label faceStart = pf.patch().start();
159
160 forAll(pf, facei)
161 {
162 allValues[faceStart++] = pif[facei];
163 }
164 }
165
166 mpp.distribute(allValues);
167 nbrIntFld.transfer(allValues);
168
169 break;
170 }
171 default:
172 {
174 << "Unknown sampling mode: " << mpp.mode()
175 << abort(FatalError);
176 }
177 }
178
179 UPstream::msgType(oldTag); // Restore tag
180
181 // Assign to (this) patch internal field its neighbour values
182 Field<Type>& intFld = const_cast<Field<Type>&>(this->primitiveField());
183 UIndirectList<Type>(intFld, this->patch().faceCells()) = nbrIntFld;
184}
185
186
187template<class Type>
189(
190 Ostream& os
191) const
192{
194}
195
196
197// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
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
A List with indirect addressing. Like IndirectList but does not store addressing.
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
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
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
This boundary condition maps the boundary and internal values of a neighbour patch field to the bound...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
mappedFixedInternalValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
mappedFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
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
static const Enum< sampleMode > sampleModeNames_
const polyMesh & sampleMesh() const
Get the region mesh.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
sampleMode mode() const noexcept
What to sample.
const GeometricField< T, fvPatchField, volMesh > & sampleField(const word &fieldName) const
Field to sample. Either on my or nbr mesh.
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
label index() const noexcept
The index of this patch in the boundaryMesh.
label nFaces() const noexcept
Number of mesh faces.
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
errorManip< error > abort(error &err)
Definition errorManip.H:139
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
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299