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 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 :
42  patch_(p),
43  internalField_(iF),
44  updated_(false),
45  patchType_(word::null)
46 {}
47 
48 
49 template<class Type>
51 (
52  const pointPatch& p,
54  const dictionary& dict
55 )
56 :
57  patch_(p),
58  internalField_(iF),
59  updated_(false),
60  patchType_(dict.getOrDefault<word>("patchType", word::null))
61 {}
62 
63 
64 template<class Type>
66 (
67  const pointPatchField<Type>& ptf,
68  const pointPatch& p,
71 )
72 :
73  patch_(p),
74  internalField_(iF),
75  updated_(false),
76  patchType_(ptf.patchType_)
77 {}
78 
79 
80 template<class Type>
82 (
83  const pointPatchField<Type>& ptf
84 )
85 :
86  patch_(ptf.patch_),
87  internalField_(ptf.internalField_),
88  updated_(false),
89  patchType_(ptf.patchType_)
90 {}
91 
92 
93 template<class Type>
95 (
96  const pointPatchField<Type>& ptf,
98 )
99 :
100  patch_(ptf.patch_),
101  internalField_(iF),
102  updated_(false),
103  patchType_(ptf.patchType_)
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
109 template<class Type>
111 {
112  return patch_.boundaryMesh().mesh()();
113 }
114 
115 
116 template<class Type>
118 {
119  os.writeEntry("type", type());
120 
121  if (patchType_.size())
122  {
123  os.writeEntry("patchType", patchType_);
124  }
125 }
126 
127 
128 template<class Type>
131 {
132  return patchInternalField(primitiveField());
133 }
134 
135 
136 template<class Type>
137 template<class Type1>
140 (
141  const Field<Type1>& iF,
142  const labelList& meshPoints
143 ) const
144 {
145  // Check size
146  if (iF.size() != primitiveField().size())
147  {
149  << "given internal field does not correspond to the mesh. "
150  << "Field size: " << iF.size()
151  << " mesh size: " << primitiveField().size()
152  << abort(FatalError);
153  }
154 
155  return tmp<Field<Type1>>::New(iF, meshPoints);
156 }
157 
158 
159 template<class Type>
160 template<class Type1>
163 (
164  const Field<Type1>& iF
165 ) const
166 {
167  return patchInternalField(iF, patch().meshPoints());
168 }
169 
170 
171 template<class Type>
172 template<class Type1>
174 (
175  Field<Type1>& iF,
176  const Field<Type1>& pF
177 ) const
178 {
179  // Check size
180  if (iF.size() != primitiveField().size())
181  {
183  << "given internal field does not correspond to the mesh. "
184  << "Field size: " << iF.size()
185  << " mesh size: " << primitiveField().size()
186  << abort(FatalError);
187  }
188 
189  if (pF.size() != size())
190  {
192  << "given patch field does not correspond to the mesh. "
193  << "Field size: " << pF.size()
194  << " mesh size: " << size()
195  << abort(FatalError);
196  }
197 
198  // Get the addressing
199  const labelList& mp = patch().meshPoints();
200 
201  forAll(mp, pointi)
202  {
203  iF[mp[pointi]] += pF[pointi];
204  }
205 }
206 
207 
208 template<class Type>
209 template<class Type1>
211 (
212  Field<Type1>& iF,
213  const Field<Type1>& pF,
214  const labelList& points
215 ) const
216 {
217  // Check size
218  if (iF.size() != primitiveField().size())
219  {
221  << "given internal field does not correspond to the mesh. "
222  << "Field size: " << iF.size()
223  << " mesh size: " << primitiveField().size()
224  << abort(FatalError);
225  }
226 
227  if (pF.size() != size())
228  {
230  << "given patch field does not correspond to the mesh. "
231  << "Field size: " << pF.size()
232  << " mesh size: " << size()
233  << abort(FatalError);
234  }
235 
236  // Get the addressing
237  const labelList& mp = patch().meshPoints();
238 
239  forAll(points, i)
240  {
241  label pointi = points[i];
242  iF[mp[pointi]] += pF[pointi];
243  }
244 }
245 
246 
247 template<class Type>
248 template<class Type1>
250 (
251  Field<Type1>& iF,
252  const Field<Type1>& pF,
253  const labelList& meshPoints
254 ) const
255 {
256  // Check size
257  if (iF.size() != primitiveField().size())
258  {
260  << "given internal field does not correspond to the mesh. "
261  << "Field size: " << iF.size()
262  << " mesh size: " << primitiveField().size()
263  << abort(FatalError);
264  }
265 
266  if (pF.size() != meshPoints.size())
267  {
269  << "given patch field does not correspond to the meshPoints. "
270  << "Field size: " << pF.size()
271  << " meshPoints size: " << size()
272  << abort(FatalError);
273  }
274 
275  forAll(meshPoints, pointi)
276  {
277  iF[meshPoints[pointi]] = pF[pointi];
278  }
279 }
280 
281 
282 template<class Type>
283 template<class Type1>
285 (
286  Field<Type1>& iF,
287  const Field<Type1>& pF
288 ) const
289 {
290  setInInternalField(iF, pF, patch().meshPoints());
291 }
292 
293 
294 template<class Type>
296 {
297  if (!updated_)
298  {
299  updateCoeffs();
300  }
301 
302  updated_ = false;
303 }
304 
305 
306 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
307 
308 template<class Type>
309 Foam::Ostream& Foam::operator<<
310 (
311  Ostream& os,
312  const pointPatchField<Type>& ptf
313 )
314 {
315  ptf.write(os);
316 
317  os.check(FUNCTION_NAME);
318 
319  return os;
320 }
321 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 #include "pointPatchFieldNew.C"
326 
327 // ************************************************************************* //
Foam::pointPatchField::addToInternalField
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field,.
Definition: pointPatchField.C:174
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
pointPatchField.H
Foam::pointPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
Definition: pointPatchField.C:295
Foam::pointPatchField::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: pointPatchField.C:110
Foam::pointPatchField::pointPatchField
pointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
Definition: pointPatchField.C:37
pointPatchFieldNew.C
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:58
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:60
Foam::pointPatchField::setInInternalField
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelList &meshPoints) const
Given the internal field and a patch field,.
Definition: pointPatchField.C:250
Foam::pointPatchFieldMapper
Foam::pointPatchFieldMapper.
Definition: pointPatchFieldMapper.H:48
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::pointPatchField::write
virtual void write(Ostream &) const
Write.
Definition: pointPatchField.C:117
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
dictionary.H
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::pointPatchField::patchInternalField
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
Definition: pointPatchField.C:130
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
pointMesh.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54