mappedVelocityFluxFvPatchField.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) 2011-2023 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 
27 #include "fvPatchFieldMapper.H"
28 #include "mappedPatchBase.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF, dict),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
44 {
46  (
47  *this,
48  iF,
49  dict,
51  );
52 }
53 
54 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
63  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
64  phiName_(ptf.phiName_)
65 {}
66 
67 
69 (
72 )
73 :
74  fixedValueFvPatchVectorField(ptf, iF),
75  phiName_(ptf.phiName_)
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  if (updated())
84  {
85  return;
86  }
87 
88  // Since we're inside initEvaluate/evaluate there might be processor
89  // comms underway. Change the tag we use.
90  int oldTag = UPstream::msgType();
91  UPstream::msgType() = oldTag+1;
92 
93  const mappedPatchBase& mapper = mappedPatchBase::getMap(patch().patch());
94  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper.nbrMesh());
95  const label nbrPatchi = mapper.nbrPolyPatch().index();
96 
97  const volVectorField& UField =
98  refCast<const volVectorField>(internalField());
99  surfaceScalarField& phiField =
100  nbrMesh.lookupObjectRef<surfaceScalarField>(phiName_);
101 
102  operator==(mapper.fromNeighbour(UField.boundaryField()[nbrPatchi]));
103  phiField.boundaryFieldRef()[patch().index()] ==
104  mapper.fromNeighbour(phiField.boundaryField()[nbrPatchi]);
105 
106  // Restore tag
107  UPstream::msgType() = oldTag;
108 
109  fixedValueFvPatchVectorField::updateCoeffs();
110 }
111 
112 
114 (
115  Ostream& os
116 ) const
117 {
119  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
120  writeEntry(os, "value", *this);
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 namespace Foam
127 {
129  (
132  );
133 }
134 
135 
136 // ************************************************************************* //
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...
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Engine which provides mapping between two patches.
const polyPatch & nbrPolyPatch() const
Get the patch to map from.
static const mappedPatchBase & getMap(const polyPatch &patch)
Cast the given polyPatch to a mappedPatchBase. Handle errors.
const polyMesh & nbrMesh() const
Get the mesh for the region to map from.
static void validateMapForField(const PatchFieldType &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
tmp< Field< Type > > fromNeighbour(const Field< Type > &nbrFld) const
Map/interpolate the neighbour patch field to this patch.
This boundary condition maps the velocity and flux from a neighbouring patch to this patch.
virtual void write(Ostream &) const
Write.
mappedVelocityFluxFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
label index() const
Return the index of this patch in the boundaryMesh.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
static const label differentPatch
Foam::surfaceFields.