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-------------------------------------------------------------------------------
11License
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
35template<class Type>
37(
38 const pointPatch& p,
40)
41:
42 patch_(p),
43 internalField_(iF),
44 updated_(false),
45 patchType_()
46{}
47
48
49template<class Type>
51(
52 const pointPatch& p,
54 const dictionary& dict
55)
56:
57 patch_(p),
58 internalField_(iF),
59 updated_(false),
60 patchType_()
61{
62 dict.readIfPresent("patchType", patchType_, keyType::LITERAL);
63}
64
65
66template<class Type>
68(
69 const pointPatchField<Type>& ptf,
70 const pointPatch& p,
73)
74:
75 patch_(p),
76 internalField_(iF),
77 updated_(false),
78 patchType_(ptf.patchType_)
79{}
80
81
82template<class Type>
84(
85 const pointPatchField<Type>& ptf
86)
87:
88 patch_(ptf.patch_),
89 internalField_(ptf.internalField_),
90 updated_(false),
91 patchType_(ptf.patchType_)
92{}
93
94
95template<class Type>
97(
98 const pointPatchField<Type>& ptf,
100)
101:
102 patch_(ptf.patch_),
103 internalField_(iF),
104 updated_(false),
105 patchType_(ptf.patchType_)
106{}
107
108
109// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110
111template<class Type>
113{
114 return patch_.boundaryMesh().mesh()();
115}
116
117
118template<class Type>
120{
121 os.writeEntry("type", type());
122
123 if (!patchType_.empty())
124 {
125 os.writeEntry("patchType", patchType_);
126 }
127}
128
129
130template<class Type>
133{
134 return patchInternalField(primitiveField());
135}
136
137
138template<class Type>
139template<class Type1>
142(
143 const Field<Type1>& iF,
144 const labelList& meshPoints
145) const
146{
147 // Check size
148 if (iF.size() != primitiveField().size())
149 {
151 << "given internal field does not correspond to the mesh. "
152 << "Field size: " << iF.size()
153 << " mesh size: " << primitiveField().size()
154 << abort(FatalError);
155 }
156
157 return tmp<Field<Type1>>::New(iF, meshPoints);
158}
159
160
161template<class Type>
162template<class Type1>
165(
166 const Field<Type1>& iF
167) const
168{
169 return patchInternalField(iF, patch().meshPoints());
170}
171
172
173template<class Type>
174template<class Type1>
176(
178 const Field<Type1>& pF
179) const
180{
181 // Check size
182 if (iF.size() != primitiveField().size())
183 {
185 << "given internal field does not correspond to the mesh. "
186 << "Field size: " << iF.size()
187 << " mesh size: " << primitiveField().size()
188 << abort(FatalError);
189 }
190
191 if (pF.size() != size())
194 << "given patch field does not correspond to the mesh. "
195 << "Field size: " << pF.size()
196 << " mesh size: " << size()
197 << abort(FatalError);
198 }
199
200 // Get the addressing
201 const labelList& mp = patch().meshPoints();
202
203 forAll(mp, pointi)
204 {
205 iF[mp[pointi]] += pF[pointi];
206 }
207}
208
209
210template<class Type>
211template<class Type1>
213(
214 Field<Type1>& iF,
215 const Field<Type1>& pF,
216 const labelList& points
217) const
218{
219 // Check size
220 if (iF.size() != primitiveField().size())
221 {
223 << "given internal field does not correspond to the mesh. "
224 << "Field size: " << iF.size()
225 << " mesh size: " << primitiveField().size()
226 << abort(FatalError);
227 }
228
229 if (pF.size() != size())
230 {
232 << "given patch field does not correspond to the mesh. "
233 << "Field size: " << pF.size()
234 << " mesh size: " << size()
235 << abort(FatalError);
236 }
237
238 // Get the addressing
239 const labelList& mp = patch().meshPoints();
240
241 forAll(points, i)
242 {
243 label pointi = points[i];
244 iF[mp[pointi]] += pF[pointi];
245 }
246}
247
248
249template<class Type>
250template<class Type1>
252(
253 Field<Type1>& iF,
254 const Field<Type1>& pF,
255 const labelList& meshPoints
256) const
257{
258 // Check size
259 if (iF.size() != primitiveField().size())
260 {
262 << "given internal field does not correspond to the mesh. "
263 << "Field size: " << iF.size()
264 << " mesh size: " << primitiveField().size()
265 << abort(FatalError);
266 }
268 if (pF.size() != meshPoints.size())
269 {
271 << "given patch field does not correspond to the meshPoints. "
272 << "Field size: " << pF.size()
273 << " meshPoints size: " << size()
274 << abort(FatalError);
275 }
276
277 forAll(meshPoints, pointi)
278 {
279 iF[meshPoints[pointi]] = pF[pointi];
280 }
281}
282
283
284template<class Type>
285template<class Type1>
287(
288 Field<Type1>& iF,
289 const Field<Type1>& pF
290) const
291{
292 setInInternalField(iF, pF, patch().meshPoints());
293}
294
295
296template<class Type>
298{
299 if (!updated_)
300 {
301 updateCoeffs();
302 }
303
304 updated_ = false;
305}
306
307
308// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
309
310template<class Type>
311Foam::Ostream& Foam::operator<<
312(
313 Ostream& os,
314 const pointPatchField<Type>& ptf
315)
316{
317 ptf.write(os);
318
320
321 return os;
322}
323
324
325// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326
327#include "pointPatchFieldNew.C"
328
329// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:82
void evaluate()
Evaluate boundary conditions.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
virtual bool write()
Write the output fields.
@ LITERAL
String literal.
Definition: keyType.H:81
Registry of regIOobjects.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
const objectRegistry & db() const
Return local objectRegistry.
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field,.
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelList &meshPoints) const
Given the internal field and a patch field,.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
#define FUNCTION_NAME
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
errorManip< error > abort(error &err)
Definition: errorManip.H:144
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333