DimensionedField.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) 2015-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 
30 #include "dimensionedType.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 // Check that both fields use the same mesh
35 #undef checkField
36 #define checkField(fld1, fld2, op) \
37 if (&(fld1).mesh() != &(fld2).mesh()) \
38 { \
39  FatalErrorInFunction \
40  << "Different mesh for fields " \
41  << (fld1).name() << " and " << (fld2).name() \
42  << " during operation " << op \
43  << abort(FatalError); \
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 template<class Type, class GeoMesh>
51 {
52  const label fieldSize = this->size();
53  if (fieldSize)
54  {
55  const label meshSize = GeoMesh::size(this->mesh_);
56  if (fieldSize != meshSize)
57  {
59  << "size of field = " << fieldSize
60  << " is not the same as the size of mesh = "
61  << meshSize
62  << abort(FatalError);
63  }
64  }
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
70 template<class Type, class GeoMesh>
72 (
73  const IOobject& io,
74  const Mesh& mesh,
75  const dimensionSet& dims,
76  const Field<Type>& field
77 )
78 :
79  regIOobject(io),
80  Field<Type>(field),
81  mesh_(mesh),
82  dimensions_(dims),
83  oriented_()
84 {
85  checkFieldSize();
86 }
87 
88 
89 template<class Type, class GeoMesh>
91 (
92  const IOobject& io,
93  const Mesh& mesh,
94  const dimensionSet& dims,
96 )
97 :
98  regIOobject(io),
99  Field<Type>(std::move(field)),
100  mesh_(mesh),
101  dimensions_(dims)
102 {
103  checkFieldSize();
104 }
105 
106 
107 template<class Type, class GeoMesh>
109 (
110  const IOobject& io,
111  const Mesh& mesh,
112  const dimensionSet& dims,
113  List<Type>&& field
114 )
115 :
116  regIOobject(io),
117  Field<Type>(std::move(field)),
118  mesh_(mesh),
119  dimensions_(dims)
120 {
121  checkFieldSize();
122 }
123 
124 
125 template<class Type, class GeoMesh>
127 (
128  const IOobject& io,
129  const Mesh& mesh,
130  const dimensionSet& dims,
131  const tmp<Field<Type>>& tfield
132 )
133 :
134  regIOobject(io),
135  Field<Type>(tfield.constCast(), tfield.movable()),
136  mesh_(mesh),
137  dimensions_(dims),
138  oriented_()
139 {
140  tfield.clear();
141  checkFieldSize();
142 }
143 
144 
145 template<class Type, class GeoMesh>
147 (
148  const IOobject& io,
149  const Mesh& mesh,
150  const dimensionSet& dims,
151  const bool checkIOFlags
152 )
153 :
154  regIOobject(io),
155  Field<Type>(GeoMesh::size(mesh)),
156  mesh_(mesh),
157  dimensions_(dims),
158  oriented_()
159 {
160  if (checkIOFlags)
161  {
163  }
164 }
165 
166 
167 template<class Type, class GeoMesh>
169 (
170  const IOobject& io,
171  const Mesh& mesh,
172  const Type& value,
173  const dimensionSet& dims,
174  const bool checkIOFlags
175 )
176 :
177  regIOobject(io),
178  Field<Type>(GeoMesh::size(mesh), value),
179  mesh_(mesh),
180  dimensions_(dims),
181  oriented_()
182 {
183  if (checkIOFlags)
184  {
186  }
187 }
188 
189 
190 template<class Type, class GeoMesh>
192 (
193  const IOobject& io,
194  const Mesh& mesh,
195  const dimensioned<Type>& dt,
196  const bool checkIOFlags
197 )
198 :
199  DimensionedField<Type, GeoMesh>
200  (
201  io,
202  mesh,
203  dt.value(),
204  dt.dimensions(),
205  checkIOFlags
206  )
207 {}
208 
209 
210 template<class Type, class GeoMesh>
212 (
214 )
215 :
216  regIOobject(df),
217  Field<Type>(df),
218  mesh_(df.mesh_),
219  dimensions_(df.dimensions_),
220  oriented_(df.oriented_)
221 {}
222 
223 
224 template<class Type, class GeoMesh>
226 (
228 )
229 :
230  DimensionedField<Type, GeoMesh>(df, true)
231 {}
232 
233 
234 template<class Type, class GeoMesh>
236 (
238  bool reuse
239 )
240 :
241  regIOobject(df, reuse),
242  Field<Type>(df, reuse),
243  mesh_(df.mesh_),
244  dimensions_(df.dimensions_),
245  oriented_(df.oriented_)
246 {}
247 
248 
249 template<class Type, class GeoMesh>
251 (
253 )
254 :
255  DimensionedField<Type, GeoMesh>(tdf.constCast(), tdf.movable())
256 {
257  tdf.clear();
258 }
259 
260 
261 template<class Type, class GeoMesh>
263 (
264  const IOobject& io,
266 )
267 :
268  regIOobject(io),
269  Field<Type>(df),
270  mesh_(df.mesh_),
271  dimensions_(df.dimensions_),
272  oriented_(df.oriented_)
273 {}
274 
275 
276 template<class Type, class GeoMesh>
278 (
279  const IOobject& io,
281 )
282 :
283  DimensionedField<Type, GeoMesh>(io, df, true)
284 {}
285 
286 
287 template<class Type, class GeoMesh>
289 (
290  const IOobject& io,
292  bool reuse
293 )
294 :
295  regIOobject(io, df),
296  Field<Type>(df, reuse),
297  mesh_(df.mesh_),
298  dimensions_(df.dimensions_),
299  oriented_(df.oriented_)
300 {}
301 
302 
303 template<class Type, class GeoMesh>
305 (
306  const IOobject& io,
308 )
309 :
310  DimensionedField<Type, GeoMesh>(io, tdf.constCast(), tdf.movable())
311 {
312  tdf.clear();
313 }
314 
315 
316 template<class Type, class GeoMesh>
318 (
319  const word& newName,
321 )
322 :
323  regIOobject(newName, df, newName != df.name()),
324  Field<Type>(df),
325  mesh_(df.mesh_),
326  dimensions_(df.dimensions_),
327  oriented_(df.oriented_)
328 {}
329 
330 
331 template<class Type, class GeoMesh>
333 (
334  const word& newName,
336 )
337 :
338  DimensionedField<Type, GeoMesh>(newName, df, true)
339 {}
340 
341 
342 template<class Type, class GeoMesh>
344 (
345  const word& newName,
347  bool reuse
348 )
349 :
350  regIOobject(newName, df, true),
351  Field<Type>(df, reuse),
352  mesh_(df.mesh_),
353  dimensions_(df.dimensions_),
354  oriented_(df.oriented_)
355 {}
356 
357 
358 template<class Type, class GeoMesh>
360 (
361  const word& newName,
363 )
364 :
365  DimensionedField<Type, GeoMesh>(newName, tdf.constCast(), tdf.movable())
366 {
367  tdf.clear();
368 }
369 
370 
371 template<class Type, class GeoMesh>
374 {
376 }
377 
378 
379 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
380 
381 template<class Type, class GeoMesh>
383 {
384  // FUTURE: register cache field info
385  // // this->db().cacheTemporaryObject(*this);
386 }
387 
388 
389 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
390 
391 template<class Type, class GeoMesh>
392 Foam::tmp
393 <
395  <
397  >
398 >
400 (
401  const direction d
402 ) const
403 {
405  (
406  name() + ".component(" + ::Foam::name(d) + ')',
407  mesh_,
408  dimensions_
409  );
410 
411  Foam::component(tresult.ref(), *this, d);
412 
413  return tresult;
414 }
415 
416 
417 template<class Type, class GeoMesh>
419 (
420  const direction d,
421  const DimensionedField
422  <
423  typename DimensionedField<Type, GeoMesh>::cmptType, GeoMesh
424  >& df
425 )
426 {
427  Field<Type>::replace(d, df);
428 }
429 
430 
431 template<class Type, class GeoMesh>
433 (
434  const direction d,
435  const tmp
436  <
437  DimensionedField
438  <
439  typename DimensionedField<Type, GeoMesh>::cmptType, GeoMesh
440  >
441  >& tdf
442 )
443 {
444  replace(d, tdf());
445  tdf.clear();
446 }
447 
448 
449 template<class Type, class GeoMesh>
452 {
454  (
455  name() + ".T()",
456  mesh_,
457  dimensions_
458  );
459 
460  Foam::T(tresult.ref(), *this);
461 
462  return tresult;
463 }
464 
465 
466 template<class Type, class GeoMesh>
468 {
469  return
471  (
472  this->name() + ".average()",
473  this->dimensions(),
475  );
476 }
477 
478 
479 template<class Type, class GeoMesh>
481 (
482  const DimensionedField<scalar, GeoMesh>& weightField
483 ) const
484 {
485  return
487  (
488  this->name() + ".weightedAverage(weights)",
489  this->dimensions(),
490  gSum(weightField*field())/gSum(weightField)
491  );
492 }
493 
494 
495 template<class Type, class GeoMesh>
497 (
498  const tmp<DimensionedField<scalar, GeoMesh>>& tweightField
499 ) const
500 {
501  dimensioned<Type> result = weightedAverage(tweightField());
502  tweightField.clear();
503  return result;
504 }
505 
506 
507 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
508 
509 template<class Type, class GeoMesh>
511 (
513 )
514 {
515  if (this == &df)
516  {
517  return; // Self-assignment is a no-op
518  }
519 
520  checkField(*this, df, "=");
521 
522  dimensions_ = df.dimensions();
523  oriented_ = df.oriented();
525 }
526 
527 
528 template<class Type, class GeoMesh>
530 (
532 )
533 {
534  auto& df = tdf.constCast();
535 
536  if (this == &df)
537  {
538  return; // Self-assignment is a no-op
539  }
540 
541  checkField(*this, df, "=");
542 
543  dimensions_ = df.dimensions();
544  oriented_ = df.oriented();
545  this->transfer(df);
546  tdf.clear();
547 }
548 
549 
550 template<class Type, class GeoMesh>
552 (
553  const dimensioned<Type>& dt
554 )
555 {
556  dimensions_ = dt.dimensions();
557  Field<Type>::operator=(dt.value());
558 }
559 
560 
561 #define COMPUTED_ASSIGNMENT(TYPE, op) \
562  \
563 template<class Type, class GeoMesh> \
564 void Foam::DimensionedField<Type, GeoMesh>::operator op \
565 ( \
566  const DimensionedField<TYPE, GeoMesh>& df \
567 ) \
568 { \
569  checkField(*this, df, #op); \
570  \
571  dimensions_ op df.dimensions(); \
572  oriented_ op df.oriented(); \
573  Field<Type>::operator op(df); \
574 } \
575  \
576 template<class Type, class GeoMesh> \
577 void Foam::DimensionedField<Type, GeoMesh>::operator op \
578 ( \
579  const tmp<DimensionedField<TYPE, GeoMesh>>& tdf \
580 ) \
581 { \
582  operator op(tdf()); \
583  tdf.clear(); \
584 } \
585  \
586 template<class Type, class GeoMesh> \
587 void Foam::DimensionedField<Type, GeoMesh>::operator op \
588 ( \
589  const dimensioned<TYPE>& dt \
590 ) \
591 { \
592  dimensions_ op dt.dimensions(); \
593  Field<Type>::operator op(dt.value()); \
594 }
595 
596 COMPUTED_ASSIGNMENT(Type, +=)
597 COMPUTED_ASSIGNMENT(Type, -=)
598 COMPUTED_ASSIGNMENT(scalar, *=)
599 COMPUTED_ASSIGNMENT(scalar, /=)
600 
601 #undef COMPUTED_ASSIGNMENT
602 
603 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
604 
605 #undef checkField
606 
607 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
608 
609 #include "DimensionedFieldIO.C"
610 #include "DimensionedFieldNew.C"
612 
613 // ************************************************************************* //
tmp< DimensionedField< Type, GeoMesh > > clone() const
Clone.
const Type & value() const noexcept
Return const reference to value.
dimensioned< Type > average() const
Calculate and return arithmetic average.
rDeltaTY field()
uint8_t direction
Definition: direction.H:46
void replace(const direction d, const DimensionedField< cmptType, GeoMesh > &df)
Replace a component field of the field.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
DimensionedField(const IOobject &io, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &field)
Construct from components, copy initial field content.
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.
Generic dimensioned Type class.
Field< Type >::cmptType cmptType
Component type of the field elements.
tmp< DimensionedField< cmptType, GeoMesh > > component(const direction d) const
Return a component field of the field.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Generic templated field type.
Definition: Field.H:62
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< DimensionedField< Type, GeoMesh > > T() const
Return the field transpose (only defined for second rank tensors)
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
virtual ~DimensionedField()
Destructor.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define COMPUTED_ASSIGNMENT(TYPE, op)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define checkField(fld1, fld2, op)
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
meshDefDict readIfPresent("polyMeshPatches", polyPatchNames)
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)