pointPatchField.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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
29 #include "pointPatchField.H"
30 #include "pointMesh.H"
31 #include "dictionary.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const pointPatch& p,
40 )
41 :
43  internalField_(iF)
44 {}
45 
46 
47 template<class Type>
49 (
50  const pointPatch& p,
52  const dictionary& dict
53 )
54 :
56  internalField_(iF)
57 {}
58 
59 
60 template<class Type>
62 (
63  const pointPatchField<Type>& ptf,
64  const pointPatch& p,
67 )
68 :
70  internalField_(iF)
71 {}
72 
73 
74 template<class Type>
76 (
77  const pointPatchField<Type>& ptf
78 )
79 :
81  internalField_(ptf.internalField_)
82 {}
83 
84 
85 template<class Type>
87 (
88  const pointPatchField<Type>& ptf,
90 )
91 :
93  internalField_(iF)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class Type>
101 {
102  os.writeEntry("type", type());
103 
104  if (!patchType().empty())
105  {
106  os.writeEntry("patchType", patchType());
107  }
108 }
109 
110 
111 template<class Type>
114 {
115  return patchInternalField(primitiveField());
116 }
117 
118 
119 template<class Type>
120 template<class Type1>
123 (
124  const Field<Type1>& iF,
125  const labelUList& meshPoints
126 ) const
127 {
128  // Check size
129  if (iF.size() != primitiveField().size())
130  {
132  << "given internal field does not correspond to the mesh. "
133  << "Field size: " << iF.size()
134  << " mesh size: " << primitiveField().size()
135  << abort(FatalError);
136  }
137 
138  return tmp<Field<Type1>>::New(iF, meshPoints);
139 }
140 
141 
142 template<class Type>
143 template<class Type1>
146 (
147  const Field<Type1>& iF
148 ) const
149 {
150  return patchInternalField(iF, patch().meshPoints());
151 }
152 
153 
154 template<class Type>
155 template<class Type1>
157 (
158  Field<Type1>& iF,
159  const Field<Type1>& pF
160 ) const
161 {
162  // Check size
163  if (iF.size() != primitiveField().size())
164  {
166  << "given internal field does not correspond to the mesh. "
167  << "Field size: " << iF.size()
168  << " mesh size: " << primitiveField().size()
169  << abort(FatalError);
170  }
171 
172  if (pF.size() != size())
173  {
175  << "given patch field does not correspond to the mesh. "
176  << "Field size: " << pF.size()
177  << " mesh size: " << size()
178  << abort(FatalError);
179  }
180 
181  // Get the addressing
182  const labelList& mp = patch().meshPoints();
183 
184  forAll(mp, pointi)
185  {
186  iF[mp[pointi]] += pF[pointi];
187  }
188 }
189 
190 
191 template<class Type>
192 template<class Type1>
194 (
195  Field<Type1>& iF,
196  const Field<Type1>& pF,
197  const labelUList& points
198 ) const
199 {
200  // Check size
201  if (iF.size() != primitiveField().size())
202  {
204  << "given internal field does not correspond to the mesh. "
205  << "Field size: " << iF.size()
206  << " mesh size: " << primitiveField().size()
207  << abort(FatalError);
208  }
209 
210  if (pF.size() != size())
211  {
213  << "given patch field does not correspond to the mesh. "
214  << "Field size: " << pF.size()
215  << " mesh size: " << size()
216  << abort(FatalError);
217  }
218 
219  // Get the addressing
220  const labelList& mp = patch().meshPoints();
221 
222  forAll(points, i)
223  {
224  label pointi = points[i];
225  iF[mp[pointi]] += pF[pointi];
226  }
227 }
228 
229 
230 template<class Type>
231 template<class Type1>
233 (
234  Field<Type1>& iF,
235  const Field<Type1>& pF,
236  const labelUList& meshPoints
237 ) const
238 {
239  // Check size
240  if (iF.size() != primitiveField().size())
241  {
243  << "given internal field does not correspond to the mesh. "
244  << "Field size: " << iF.size()
245  << " mesh size: " << primitiveField().size()
246  << abort(FatalError);
247  }
248 
249  if (pF.size() != meshPoints.size())
250  {
252  << "given patch field does not correspond to the meshPoints. "
253  << "Field size: " << pF.size()
254  << " meshPoints size: " << size()
255  << abort(FatalError);
256  }
257 
258  forAll(meshPoints, pointi)
259  {
260  iF[meshPoints[pointi]] = pF[pointi];
261  }
262 }
263 
264 
265 template<class Type>
266 template<class Type1>
268 (
269  Field<Type1>& iF,
270  const Field<Type1>& pF
271 ) const
272 {
273  setInInternalField(iF, pF, patch().meshPoints());
274 }
275 
276 
277 template<class Type>
279 {
280  pointPatchFieldBase::setUpdated(true);
281 }
282 
283 
284 template<class Type>
286 {
287  if (!updated())
288  {
289  updateCoeffs();
290  }
291 
292  pointPatchFieldBase::setUpdated(false);
293 }
294 
295 
296 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
297 
298 template<class Type>
299 Foam::Ostream& Foam::operator<<
300 (
301  Ostream& os,
302  const pointPatchField<Type>& ptf
303 )
304 {
305  ptf.write(os);
306 
307  os.check(FUNCTION_NAME);
308 
309  return os;
310 }
311 
312 
313 // ************************************************************************* //
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Template invariant parts for pointPatchField.
type
Types of root.
Definition: Roots.H:52
pointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
commsTypes
Communications types.
Definition: UPstream.H:72
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
Foam::pointPatchFieldMapper.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void write(Ostream &os) const
Write.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Abstract base class for point-mesh patch fields.
const pointField & points
Generic templated field type.
Definition: Field.H:62
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelUList &meshPoints) const
Given the internal field and a patch field, set the patch field in the internal field.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field, add the patch field to the internal 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.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar mp
Proton mass.