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-2023 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 :
43  internalField_(iF)
44 {}
45 
46 
47 template<class Type>
49 (
50  const pointPatch& p,
52  const dictionary& dict
53 )
54 :
56  internalField_(iF)
57 {}
58 
59 
60 template<class Type>
62 (
63  const pointPatchField<Type>& ptf,
64  const pointPatch& p,
67 )
68 :
70  internalField_(iF)
71 {}
72 
73 
74 template<class Type>
76 (
77  const pointPatchField<Type>& ptf
78 )
79 :
81  internalField_(ptf.internalField_)
82 {}
83 
84 
85 template<class Type>
87 (
88  const pointPatchField<Type>& ptf,
90 )
91 :
93  internalField_(iF)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class Type>
101 {
102  os.writeEntry("type", type());
103 
104  if (!patchType().empty())
105  {
106  os.writeEntry("patchType", patchType());
107  }
108 }
109 
110 
111 template<class Type>
112 template<class Type1>
114 (
115  const UList<Type1>& internalData,
116  const labelUList& addressing,
117  Field<Type1>& pfld
118 ) const
119 {
120  // Check size
121  if (internalData.size() != primitiveField().size())
122  {
124  << "Internal field size: " << internalData.size()
125  << " != mesh size: " << primitiveField().size() << nl
126  << abort(FatalError);
127  }
128 
129  const label len = this->size();
130 
131  pfld.resize_nocopy(len);
132 
133  for (label i = 0; i < len; ++i)
134  {
135  pfld[i] = internalData[addressing[i]];
136  }
137 }
138 
139 
140 template<class Type>
141 template<class Type1>
144 (
145  const UList<Type1>& internalData,
146  const labelUList& addressing
147 ) const
148 {
149  auto tpfld = tmp<Field<Type1>>::New();
150  this->patchInternalField(internalData, addressing, tpfld.ref());
151  return tpfld;
152 }
153 
154 
155 template<class Type>
156 template<class Type1>
159 (
160  const UList<Type1>& internalData
161 ) const
162 {
163  auto tpfld = tmp<Field<Type1>>::New();
164  this->patchInternalField(internalData, patch().meshPoints(), tpfld.ref());
165  return tpfld;
166 }
167 
168 
169 template<class Type>
172 {
173  return patchInternalField(primitiveField());
174 }
175 
176 
177 template<class Type>
178 template<class Type1>
180 (
181  Field<Type1>& iF,
182  const Field<Type1>& pF
183 ) const
184 {
185  // Check size
186  if (iF.size() != primitiveField().size())
187  {
189  << "Internal field size: " << iF.size()
190  << " != mesh size: " << primitiveField().size() << nl
191  << abort(FatalError);
192  }
193 
194  if (pF.size() != size())
195  {
197  << "Patch field size: " << pF.size()
198  << " != patch size: " << size() << nl
199  << abort(FatalError);
200  }
201 
202  // Get the addressing
203  const labelList& mp = patch().meshPoints();
204 
205  forAll(mp, pointi)
206  {
207  iF[mp[pointi]] += pF[pointi];
208  }
209 }
210 
211 
212 template<class Type>
213 template<class Type1>
215 (
216  Field<Type1>& iF,
217  const Field<Type1>& pF,
218  const labelUList& points
219 ) const
220 {
221  // Check size
222  if (iF.size() != primitiveField().size())
223  {
225  << "Internal field size: " << iF.size()
226  << " != mesh size: " << primitiveField().size() << nl
227  << abort(FatalError);
228  }
229 
230  if (pF.size() != size())
231  {
233  << "Patch field size: " << pF.size()
234  << " != patch size: " << size() << nl
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 
249 template<class Type>
250 template<class Type1>
252 (
253  Field<Type1>& iF,
254  const Field<Type1>& pF,
255  const labelUList& meshPoints
256 ) const
257 {
258  // Check size
259  if (iF.size() != primitiveField().size())
260  {
262  << "Internal field size: " << iF.size()
263  << " != mesh size: " << primitiveField().size() << nl
264  << abort(FatalError);
265  }
266 
267  if (pF.size() != meshPoints.size())
268  {
270  << "Patch field size: " << pF.size()
271  << " != meshPoints size: " << meshPoints.size() << nl
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  pointPatchFieldBase::setUpdated(true);
298 }
299 
300 
301 template<class Type>
303 {
304  if (!updated())
305  {
306  updateCoeffs();
307  }
308 
309  pointPatchFieldBase::setUpdated(false);
310 }
311 
312 
313 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
314 
315 template<class Type>
316 Foam::Ostream& Foam::operator<<
317 (
318  Ostream& os,
319  const pointPatchField<Type>& ptf
320 )
321 {
322  ptf.write(os);
323 
324  os.check(FUNCTION_NAME);
325 
326  return os;
327 }
328 
329 
330 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Template invariant parts for pointPatchField.
type
Types of root.
Definition: Roots.H:52
pointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
commsTypes
Communications types.
Definition: UPstream.H:77
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Foam::pointPatchFieldMapper.
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.
virtual void write(Ostream &os) const
Write.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Abstract base class for point-mesh patch fields.
const pointField & points
Generic templated field type.
Definition: Field.H:62
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
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:56
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
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.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const std::string patch
OpenFOAM patch number as a std::string.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar mp
Proton mass.