pointPatchField.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-2019 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 
26 #include "pointPatchField.H"
27 #include "pointMesh.H"
28 #include "dictionary.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
34 (
35  const pointPatch& p,
37 )
38 :
39  patch_(p),
40  internalField_(iF),
41  updated_(false),
42  patchType_(word::null)
43 {}
44 
45 
46 template<class Type>
48 (
49  const pointPatch& p,
51  const dictionary& dict
52 )
53 :
54  patch_(p),
55  internalField_(iF),
56  updated_(false),
57  patchType_(dict.lookupOrDefault<word>("patchType", word::null))
58 {}
59 
60 
61 template<class Type>
63 (
64  const pointPatchField<Type>& ptf,
65  const pointPatch& p,
68 )
69 :
70  patch_(p),
71  internalField_(iF),
72  updated_(false),
73  patchType_(ptf.patchType_)
74 {}
75 
76 
77 template<class Type>
79 (
80  const pointPatchField<Type>& ptf
81 )
82 :
83  patch_(ptf.patch_),
84  internalField_(ptf.internalField_),
85  updated_(false),
86  patchType_(ptf.patchType_)
87 {}
88 
89 
90 template<class Type>
92 (
93  const pointPatchField<Type>& ptf,
95 )
96 :
97  patch_(ptf.patch_),
98  internalField_(iF),
99  updated_(false),
100  patchType_(ptf.patchType_)
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class Type>
108 {
109  return patch_.boundaryMesh().mesh()();
110 }
111 
112 
113 template<class Type>
115 {
116  writeEntry(os, "type", type());
117 
118  if (patchType_.size())
119  {
120  writeEntry(os, "patchType", patchType_);
121  }
122 }
123 
124 
125 template<class Type>
128 {
129  return patchInternalField(primitiveField());
130 }
131 
132 
133 template<class Type>
134 template<class Type1>
137 (
138  const Field<Type1>& iF,
139  const labelList& meshPoints
140 ) const
141 {
142  // Check size
143  if (iF.size() != primitiveField().size())
144  {
146  << "given internal field does not correspond to the mesh. "
147  << "Field size: " << iF.size()
148  << " mesh size: " << primitiveField().size()
149  << abort(FatalError);
150  }
151 
152  return tmp<Field<Type1>>(new Field<Type1>(iF, meshPoints));
153 }
154 
155 
156 template<class Type>
157 template<class Type1>
160 (
161  const Field<Type1>& iF
162 ) const
163 {
164  return patchInternalField(iF, patch().meshPoints());
165 }
166 
167 
168 template<class Type>
169 template<class Type1>
171 (
172  Field<Type1>& iF,
173  const Field<Type1>& pF
174 ) const
175 {
176  // Check size
177  if (iF.size() != primitiveField().size())
178  {
180  << "given internal field does not correspond to the mesh. "
181  << "Field size: " << iF.size()
182  << " mesh size: " << primitiveField().size()
183  << abort(FatalError);
184  }
185 
186  if (pF.size() != size())
187  {
189  << "given patch field does not correspond to the mesh. "
190  << "Field size: " << pF.size()
191  << " mesh size: " << size()
192  << abort(FatalError);
193  }
194 
195  // Get the addressing
196  const labelList& mp = patch().meshPoints();
197 
198  forAll(mp, pointi)
199  {
200  iF[mp[pointi]] += pF[pointi];
201  }
202 }
203 
204 
205 template<class Type>
206 template<class Type1>
208 (
209  Field<Type1>& iF,
210  const Field<Type1>& pF,
211  const labelList& points
212 ) const
213 {
214  // Check size
215  if (iF.size() != primitiveField().size())
216  {
218  << "given internal field does not correspond to the mesh. "
219  << "Field size: " << iF.size()
220  << " mesh size: " << primitiveField().size()
221  << abort(FatalError);
222  }
223 
224  if (pF.size() != size())
225  {
227  << "given patch field does not correspond to the mesh. "
228  << "Field size: " << pF.size()
229  << " mesh size: " << size()
230  << abort(FatalError);
231  }
232 
233  // Get the addressing
234  const labelList& mp = patch().meshPoints();
235 
236  forAll(points, i)
237  {
238  label pointi = points[i];
239  iF[mp[pointi]] += pF[pointi];
240  }
241 }
242 
243 
244 template<class Type>
245 template<class Type1>
247 (
248  Field<Type1>& iF,
249  const Field<Type1>& pF,
250  const labelList& meshPoints
251 ) const
252 {
253  // Check size
254  if (iF.size() != primitiveField().size())
255  {
257  << "given internal field does not correspond to the mesh. "
258  << "Field size: " << iF.size()
259  << " mesh size: " << primitiveField().size()
260  << abort(FatalError);
261  }
262 
263  if (pF.size() != meshPoints.size())
264  {
266  << "given patch field does not correspond to the meshPoints. "
267  << "Field size: " << pF.size()
268  << " meshPoints size: " << size()
269  << abort(FatalError);
270  }
271 
272  forAll(meshPoints, pointi)
273  {
274  iF[meshPoints[pointi]] = pF[pointi];
275  }
276 }
277 
278 
279 template<class Type>
280 template<class Type1>
282 (
283  Field<Type1>& iF,
284  const Field<Type1>& pF
285 ) const
286 {
287  setInInternalField(iF, pF, patch().meshPoints());
288 }
289 
290 
291 template<class Type>
293 {
294  if (!updated_)
295  {
296  updateCoeffs();
297  }
298 
299  updated_ = false;
300 }
301 
302 
303 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
304 
305 template<class Type>
306 Foam::Ostream& Foam::operator<<
307 (
308  Ostream& os,
309  const pointPatchField<Type>& ptf
310 )
311 {
312  ptf.write(os);
313 
314  os.check("Ostream& operator<<(Ostream&, const pointPatchField<Type>&)");
315 
316  return os;
317 }
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 #include "pointPatchFieldNew.C"
323 
324 // ************************************************************************* //
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
pointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
commsTypes
Types of communications.
Definition: UPstream.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
const objectRegistry & db() const
Return local objectRegistry.
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelList &meshPoints) const
Given the internal field and a patch field,.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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:53
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void write(Ostream &) const
Write.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field,.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
const dimensionedScalar mp
Proton mass.