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 -------------------------------------------------------------------------------
10 License
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 \*---------------------------------------------------------------------------*/
27 
29 #include "IndirectList.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const fvPatch& p,
39 )
40 :
42 {}
43 
44 
45 template<class Type>
48 (
50  const fvPatch& p,
52  const fvPatchFieldMapper& mapper
53 )
54 :
55  mappedFixedValueFvPatchField<Type>(ptf, p, iF, mapper)
56 {}
57 
58 
59 template<class Type>
62 (
63  const fvPatch& p,
65  const dictionary& dict
66 )
67 :
69 {}
70 
71 
72 template<class Type>
75 (
77 )
78 :
80 {}
81 
82 
83 template<class Type>
86 (
89 )
90 :
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
97 template<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
112  mappedFixedValueFvPatchField<Type>::updateCoeffs();
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  {
123  case mappedPatchBase::NEARESTCELL:
124  {
126  << "Cannot apply "
127  << mappedPatchBase::sampleModeNames_
128  [
129  mappedPatchBase::NEARESTCELL
130  ]
131  << " mapping mode for patch " << this->patch().name()
132  << exit(FatalError);
133 
134  break;
135  }
136  case mappedPatchBase::NEARESTPATCHFACE:
137  case mappedPatchBase::NEARESTPATCHFACEAMI:
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  }
147  case mappedPatchBase::NEARESTFACE:
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 
187 template<class Type>
189 (
190  Ostream& os
191 ) const
192 {
194 }
195 
196 
197 // ************************************************************************* //
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition maps the value at a set of cells or patch faces back to *this.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Generic templated field type.
Definition: Field.H:62
A FieldMapper for finite-volume patch fields.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
mappedFixedInternalValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
This boundary condition maps the boundary and internal values of a neighbour patch field to the bound...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127