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 | Copyright (C) 2011-2016 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  os.writeKeyword("type") << type() << token::END_STATEMENT << nl;
117 
118  if (patchType_.size())
119  {
120  os.writeKeyword("patchType") << patchType_
121  << token::END_STATEMENT << nl;
122  }
123 }
124 
125 
126 template<class Type>
127 template<class EntryType>
129 (
130  Ostream& os,
131  const word& entryName,
132  const EntryType& value1,
133  const EntryType& value2
134 ) const
135 {
136  if (value1 != value2)
137  {
138  os.writeKeyword(entryName) << value2 << token::END_STATEMENT << nl;
139  }
140 }
141 
142 
143 template<class Type>
146 {
147  return patchInternalField(primitiveField());
148 }
149 
150 
151 template<class Type>
152 template<class Type1>
155 (
156  const Field<Type1>& iF,
157  const labelList& meshPoints
158 ) const
159 {
160  // Check size
161  if (iF.size() != primitiveField().size())
162  {
164  << "given internal field does not correspond to the mesh. "
165  << "Field size: " << iF.size()
166  << " mesh size: " << primitiveField().size()
167  << abort(FatalError);
168  }
169 
170  return tmp<Field<Type1>>(new Field<Type1>(iF, meshPoints));
171 }
172 
173 
174 template<class Type>
175 template<class Type1>
178 (
179  const Field<Type1>& iF
180 ) const
181 {
182  return patchInternalField(iF, patch().meshPoints());
183 }
184 
185 
186 template<class Type>
187 template<class Type1>
189 (
190  Field<Type1>& iF,
191  const Field<Type1>& pF
192 ) const
193 {
194  // Check size
195  if (iF.size() != primitiveField().size())
196  {
198  << "given internal field does not correspond to the mesh. "
199  << "Field size: " << iF.size()
200  << " mesh size: " << primitiveField().size()
201  << abort(FatalError);
202  }
203 
204  if (pF.size() != size())
205  {
207  << "given patch field does not correspond to the mesh. "
208  << "Field size: " << pF.size()
209  << " mesh size: " << size()
210  << abort(FatalError);
211  }
212 
213  // Get the addressing
214  const labelList& mp = patch().meshPoints();
215 
216  forAll(mp, pointi)
217  {
218  iF[mp[pointi]] += pF[pointi];
219  }
220 }
221 
222 
223 template<class Type>
224 template<class Type1>
226 (
227  Field<Type1>& iF,
228  const Field<Type1>& pF,
229  const labelList& points
230 ) const
231 {
232  // Check size
233  if (iF.size() != primitiveField().size())
234  {
236  << "given internal field does not correspond to the mesh. "
237  << "Field size: " << iF.size()
238  << " mesh size: " << primitiveField().size()
239  << abort(FatalError);
240  }
241 
242  if (pF.size() != size())
243  {
245  << "given patch field does not correspond to the mesh. "
246  << "Field size: " << pF.size()
247  << " mesh size: " << size()
248  << abort(FatalError);
249  }
250 
251  // Get the addressing
252  const labelList& mp = patch().meshPoints();
253 
254  forAll(points, i)
255  {
256  label pointi = points[i];
257  iF[mp[pointi]] += pF[pointi];
258  }
259 }
260 
261 
262 template<class Type>
263 template<class Type1>
265 (
266  Field<Type1>& iF,
267  const Field<Type1>& pF,
268  const labelList& meshPoints
269 ) const
270 {
271  // Check size
272  if (iF.size() != primitiveField().size())
273  {
275  << "given internal field does not correspond to the mesh. "
276  << "Field size: " << iF.size()
277  << " mesh size: " << primitiveField().size()
278  << abort(FatalError);
279  }
280 
281  if (pF.size() != meshPoints.size())
282  {
284  << "given patch field does not correspond to the meshPoints. "
285  << "Field size: " << pF.size()
286  << " meshPoints size: " << size()
287  << abort(FatalError);
288  }
289 
290  forAll(meshPoints, pointi)
291  {
292  iF[meshPoints[pointi]] = pF[pointi];
293  }
294 }
295 
296 
297 template<class Type>
298 template<class Type1>
300 (
301  Field<Type1>& iF,
302  const Field<Type1>& pF
303 ) const
304 {
305  setInInternalField(iF, pF, patch().meshPoints());
306 }
307 
308 
309 template<class Type>
311 {
312  if (!updated_)
313  {
314  updateCoeffs();
315  }
316 
317  updated_ = false;
318 }
319 
320 
321 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
322 
323 template<class Type>
324 Foam::Ostream& Foam::operator<<
325 (
326  Ostream& os,
327  const pointPatchField<Type>& ptf
328 )
329 {
330  ptf.write(os);
331 
332  os.check("Ostream& operator<<(Ostream&, const pointPatchField<Type>&)");
333 
334  return os;
335 }
336 
337 
338 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
339 
340 #include "pointPatchFieldNew.C"
341 
342 // ************************************************************************* //
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:428
type
Types of root.
Definition: Roots.H:52
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
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2) const
Helper function to write the keyword and entry only if the.
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:137
#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:57
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
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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,.
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.