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>
129 {
130  return patchInternalField(primitiveField());
131 }
132 
133 
134 template<class Type>
135 template<class Type1>
138 (
139  const Field<Type1>& iF,
140  const labelList& meshPoints
141 ) const
142 {
143  // Check size
144  if (iF.size() != primitiveField().size())
145  {
147  << "given internal field does not correspond to the mesh. "
148  << "Field size: " << iF.size()
149  << " mesh size: " << primitiveField().size()
150  << abort(FatalError);
151  }
152 
153  return tmp<Field<Type1>>(new Field<Type1>(iF, meshPoints));
154 }
155 
156 
157 template<class Type>
158 template<class Type1>
161 (
162  const Field<Type1>& iF
163 ) const
164 {
165  return patchInternalField(iF, patch().meshPoints());
166 }
167 
168 
169 template<class Type>
170 template<class Type1>
172 (
173  Field<Type1>& iF,
174  const Field<Type1>& pF
175 ) const
176 {
177  // Check size
178  if (iF.size() != primitiveField().size())
179  {
181  << "given internal field does not correspond to the mesh. "
182  << "Field size: " << iF.size()
183  << " mesh size: " << primitiveField().size()
184  << abort(FatalError);
185  }
186 
187  if (pF.size() != size())
188  {
190  << "given patch field does not correspond to the mesh. "
191  << "Field size: " << pF.size()
192  << " mesh size: " << size()
193  << abort(FatalError);
194  }
195 
196  // Get the addressing
197  const labelList& mp = patch().meshPoints();
198 
199  forAll(mp, pointi)
200  {
201  iF[mp[pointi]] += pF[pointi];
202  }
203 }
204 
205 
206 template<class Type>
207 template<class Type1>
209 (
210  Field<Type1>& iF,
211  const Field<Type1>& pF,
212  const labelList& points
213 ) const
214 {
215  // Check size
216  if (iF.size() != primitiveField().size())
217  {
219  << "given internal field does not correspond to the mesh. "
220  << "Field size: " << iF.size()
221  << " mesh size: " << primitiveField().size()
222  << abort(FatalError);
223  }
224 
225  if (pF.size() != size())
226  {
228  << "given patch field does not correspond to the mesh. "
229  << "Field size: " << pF.size()
230  << " mesh size: " << size()
231  << abort(FatalError);
232  }
233 
234  // Get the addressing
235  const labelList& mp = patch().meshPoints();
236 
237  forAll(points, i)
238  {
239  label pointi = points[i];
240  iF[mp[pointi]] += pF[pointi];
241  }
242 }
243 
244 
245 template<class Type>
246 template<class Type1>
248 (
249  Field<Type1>& iF,
250  const Field<Type1>& pF,
251  const labelList& meshPoints
252 ) const
253 {
254  // Check size
255  if (iF.size() != primitiveField().size())
256  {
258  << "given internal field does not correspond to the mesh. "
259  << "Field size: " << iF.size()
260  << " mesh size: " << primitiveField().size()
261  << abort(FatalError);
262  }
263 
264  if (pF.size() != meshPoints.size())
265  {
267  << "given patch field does not correspond to the meshPoints. "
268  << "Field size: " << pF.size()
269  << " meshPoints size: " << size()
270  << abort(FatalError);
271  }
272 
273  forAll(meshPoints, pointi)
274  {
275  iF[meshPoints[pointi]] = pF[pointi];
276  }
277 }
278 
279 
280 template<class Type>
281 template<class Type1>
283 (
284  Field<Type1>& iF,
285  const Field<Type1>& pF
286 ) const
287 {
288  setInInternalField(iF, pF, patch().meshPoints());
289 }
290 
291 
292 template<class Type>
294 {
295  if (!updated_)
296  {
297  updateCoeffs();
298  }
299 
300  updated_ = false;
301 }
302 
303 
304 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
305 
306 template<class Type>
307 Foam::Ostream& Foam::operator<<
308 (
309  Ostream& os,
310  const pointPatchField<Type>& ptf
311 )
312 {
313  ptf.write(os);
314 
315  os.check("Ostream& operator<<(Ostream&, const pointPatchField<Type>&)");
316 
317  return os;
318 }
319 
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 #include "pointPatchFieldNew.C"
324 
325 // ************************************************************************* //
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field,.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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: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.
Foam::pointPatchFieldMapper.
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelList &meshPoints) const
Given the internal field and a patch field,.
const objectRegistry & db() const
Return local objectRegistry.
Abstract base class for point-mesh patch fields.
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void write(Ostream &) const
Write.
A class for managing temporary objects.
Definition: PtrList.H:54
Registry of regIOobjects.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
const dimensionedScalar mp
Proton mass.