Loading...
Searching...
No Matches
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-2025 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
30#include "pointMesh.H"
31#include "dictionary.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class Type>
37(
38 const pointPatch& p,
40)
43 internalField_(iF)
44{}
45
46
47template<class Type>
49(
50 const pointPatch& p,
52 const dictionary& dict
53)
56 internalField_(iF)
57{}
58
59
60template<class Type>
62(
63 const pointPatchField<Type>& pfld,
64 const pointPatch& p,
67)
70 internalField_(iF)
71{}
72
73
74template<class Type>
76(
77 const pointPatchField<Type>& pfld,
78 const pointPatch& p,
80 const Type&
81)
84 internalField_(iF)
85{}
86
87
88template<class Type>
90(
91 const pointPatchField<Type>& pfld,
93)
94:
96 internalField_(iF)
97{}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
102template<class Type>
104{
105 os.writeEntry("type", type());
106
107 if (!patchType().empty())
108 {
109 os.writeEntry("patchType", patchType());
110 }
111}
112
113
114template<class Type>
115template<class Type1>
117(
118 const UList<Type1>& internalData,
119 const labelUList& addressing,
120 UList<Type1>& pfld
121) const
122{
123 if (FOAM_UNLIKELY(internalData.size() != primitiveField().size()))
124 {
126 << "Internal field size: " << internalData.size()
127 << " != mesh size: " << primitiveField().size() << nl
128 << abort(FatalError);
129 }
130
131 // For v2412 and earlier this was a field:
132 // const label len = this->size();
133 // pfld.resize_nocopy(len);
134 //
135 // Now uses pre-sized storage
136
137 const label len = pfld.size();
138
139 #ifdef FULLDEBUG
140 if (FOAM_UNLIKELY((addressing.size() < len) || (this->size() < len)))
141 {
143 << "patchField size = " << len
144 << " but patch size = " << this->size()
145 << " and addressing size = " << addressing.size() << nl
146 << abort(FatalError);
147 }
148 #endif
149
150 for (label i = 0; i < len; ++i)
151 {
152 pfld[i] = internalData[addressing[i]];
154}
155
156
157template<class Type>
158template<class Type1>
161(
162 const UList<Type1>& internalData,
163 const labelUList& addressing
164) const
165{
166 auto tpfld = tmp<Field<Type1>>::New(this->size());
167 this->patchInternalField(internalData, addressing, tpfld.ref());
168 return tpfld;
169}
170
171
172template<class Type>
173template<class Type1>
176(
177 const UList<Type1>& internalData
178) const
179{
180 auto tpfld = tmp<Field<Type1>>::New(this->size());
181 this->patchInternalField(internalData, patch().meshPoints(), tpfld.ref());
182 return tpfld;
183}
184
185
186template<class Type>
191}
192
193
194template<class Type>
195template<class Type1>
197(
198 Field<Type1>& iF,
199 const Field<Type1>& pF
200) const
201{
202 if (FOAM_UNLIKELY(iF.size() != primitiveField().size()))
203 {
205 << "Internal field size: " << iF.size()
206 << " != mesh size: " << primitiveField().size() << nl
207 << abort(FatalError);
208 }
209
210 if (FOAM_UNLIKELY(pF.size() != size()))
211 {
213 << "Patch field size: " << pF.size()
214 << " != patch size: " << size() << nl
215 << abort(FatalError);
216 }
217
218 // Get the addressing
219 const labelList& mp = patch().meshPoints();
220
221 forAll(mp, pointi)
222 {
223 iF[mp[pointi]] += pF[pointi];
224 }
225}
226
227
228template<class Type>
229template<class Type1>
231(
232 Field<Type1>& iF,
233 const Field<Type1>& pF,
234 const labelUList& points
235) const
236{
237 if (FOAM_UNLIKELY(iF.size() != primitiveField().size()))
238 {
240 << "Internal field size: " << iF.size()
241 << " != mesh size: " << primitiveField().size() << nl
242 << abort(FatalError);
243 }
244
245 if (FOAM_UNLIKELY(pF.size() != size()))
246 {
248 << "Patch field size: " << pF.size()
249 << " != patch size: " << size() << nl
250 << abort(FatalError);
251 }
252
253 // Get the addressing
254 const labelList& mp = patch().meshPoints();
255
256 forAll(points, i)
257 {
258 label pointi = points[i];
259 iF[mp[pointi]] += pF[pointi];
260 }
261}
262
263
264template<class Type>
265template<class Type1>
267(
268 Field<Type1>& iF,
269 const Field<Type1>& pF,
270 const labelUList& meshPoints
271) const
272{
273 if (FOAM_UNLIKELY(iF.size() != primitiveField().size()))
274 {
276 << "Internal field size: " << iF.size()
277 << " != mesh size: " << primitiveField().size() << nl
278 << abort(FatalError);
279 }
280
281 if (FOAM_UNLIKELY(pF.size() != meshPoints.size()))
282 {
284 << "Patch field size: " << pF.size()
285 << " != meshPoints size: " << meshPoints.size() << nl
286 << abort(FatalError);
287 }
288
289 forAll(meshPoints, pointi)
290 {
291 iF[meshPoints[pointi]] = pF[pointi];
292 }
293}
294
295
296template<class Type>
297template<class Type1>
299(
300 Field<Type1>& iF,
301 const Field<Type1>& pF
302) const
303{
304 setInInternalField(iF, pF, patch().meshPoints());
305}
306
307
308template<class Type>
310{
312}
313
314
315template<class Type>
317{
318 if (!updated())
319 {
320 updateCoeffs();
321 }
322
323 pointPatchFieldBase::setUpdated(false);
324 pointPatchFieldBase::setManipulated(false);
325}
326
327
328// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
329
330template<class Type>
331Foam::Ostream& Foam::operator<<
332(
333 Ostream& os,
334 const pointPatchField<Type>& ptf
335)
336{
337 ptf.write(os);
338
339 os.check(FUNCTION_NAME);
340
341 return os;
342}
343
344
345// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Template invariant parts for pointPatchField.
const pointPatch & patch() const noexcept
Return the patch.
void setManipulated(bool state) noexcept
Set matrix manipulated state. Currently a no-op for pointPatchField.
void setUpdated(bool state) noexcept
Set updated state.
pointPatchFieldBase(const pointPatch &p)
Construct from patch.
const word & patchType() const noexcept
The optional patch type.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
void patchInternalField(const UList< Type1 > &internalData, const labelUList &addressing, UList< Type1 > &pfld) const
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelUList &meshPoints) const
Given the internal field and a patch field, set the patch field in the internal field.
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field, add the patch field to the internal field.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
label size() const noexcept
Return the patch size.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
pointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
A class for managing temporary objects.
Definition tmp.H:75
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
#define FUNCTION_NAME
const dimensionedScalar mp
Proton mass.
const std::string patch
OpenFOAM patch number as a std::string.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define FOAM_UNLIKELY(cond)
Definition stdFoam.H:64