DimensionedField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "DimensionedField.H"
27 #include "dimensionedType.H"
28 #include "Time.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 #define checkFieldAssignment(df1, df2) \
38  \
39  if \
40  ( \
41  static_cast<const regIOobject*>(&df1) \
42  == static_cast<const regIOobject*>(&df2) \
43  ) \
44  { \
45  FatalErrorInFunction \
46  << "attempted assignment to self for field " \
47  << (df1).name() << abort(FatalError); \
48  }
49 
50 
51 #define checkFieldOperation(df1, df2, op) \
52  \
53  if (&(df1).mesh() != &(df2).mesh()) \
54  { \
55  FatalErrorInFunction \
56  << "different mesh for fields " \
57  << (df1).name() << " and " << (df2).name() \
58  << " during operation " << op \
59  << abort(FatalError); \
60  }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 template<class Type, class GeoMesh, template<class> class PrimitiveField>
67 (
68  const IOobject& io,
69  const Mesh& mesh,
70  const dimensionSet& dims,
71  const PrimitiveField<Type>& field
72 )
73 :
74  regIOobject(io),
75  PrimitiveField<Type>(field),
76  OldTimeField<DimensionedField>(this->time().timeIndex()),
77  mesh_(mesh),
78  dimensions_(dims)
79 {
80  if (field.size() && field.size() != GeoMesh::size(mesh))
81  {
83  << "size of field = " << field.size()
84  << " is not the same as the size of mesh = "
85  << GeoMesh::size(mesh)
86  << abort(FatalError);
87  }
88 }
89 
90 
91 template<class Type, class GeoMesh, template<class> class PrimitiveField>
93 (
94  const IOobject& io,
95  const Mesh& mesh,
96  const dimensionSet& dims,
97  const bool checkIOFlags
98 )
99 :
100  regIOobject(io),
101  PrimitiveField<Type>(GeoMesh::size(mesh)),
102  OldTimeField<DimensionedField>(this->time().timeIndex()),
103  mesh_(mesh),
104  dimensions_(dims)
105 {
106  // Expand dynamic primitive fields to their full size
107  PrimitiveField<Type>::setSize(GeoMesh::size(mesh));
108 
109  if (checkIOFlags)
110  {
111  readIfPresent();
112  }
113 }
114 
115 
116 template<class Type, class GeoMesh, template<class> class PrimitiveField>
118 (
119  const IOobject& io,
120  const Mesh& mesh,
121  const dimensioned<Type>& dt,
122  const bool checkIOFlags
123 )
124 :
125  regIOobject(io),
126  PrimitiveField<Type>(GeoMesh::size(mesh), dt.value()),
127  OldTimeField<DimensionedField>(this->time().timeIndex()),
128  mesh_(mesh),
129  dimensions_(dt.dimensions())
130 {
131  if (checkIOFlags)
132  {
133  readIfPresent();
134  }
135 }
136 
137 
138 template<class Type, class GeoMesh, template<class> class PrimitiveField>
140 (
142 )
143 :
144  regIOobject(df),
145  PrimitiveField<Type>(df),
147  mesh_(df.mesh_),
148  dimensions_(df.dimensions_)
149 {}
150 
151 
152 template<class Type, class GeoMesh, template<class> class PrimitiveField>
154 (
156 )
157 :
158  regIOobject(move(df)),
159  PrimitiveField<Type>(move(df)),
160  OldTimeField<DimensionedField>(move(df)),
161  mesh_(df.mesh_),
162  dimensions_(move(df.dimensions_))
163 {}
164 
165 
166 template<class Type, class GeoMesh, template<class> class PrimitiveField>
167 template<template<class> class PrimitiveField2>
169 (
171 )
172 :
173  regIOobject(df),
174  PrimitiveField<Type>(df),
175  OldTimeField<DimensionedField>(this->time().timeIndex()),
176  mesh_(df.mesh()),
177  dimensions_(df.dimensions())
178 {}
179 
180 
181 template<class Type, class GeoMesh, template<class> class PrimitiveField>
182 template<template<class> class PrimitiveField2>
184 (
186  const bool reuse
187 )
188 :
189  regIOobject(df, reuse),
190  PrimitiveField<Type>(df, reuse),
191  OldTimeField<DimensionedField>(this->time().timeIndex()),
192  mesh_(df.mesh()),
193  dimensions_(df.dimensions())
194 {}
195 
196 
197 template<class Type, class GeoMesh, template<class> class PrimitiveField>
199 (
201 )
202 :
203  regIOobject(tdf(), tdf.isTmp()),
204  PrimitiveField<Type>
205  (
206  const_cast<DimensionedField<Type, GeoMesh, PrimitiveField>&>(tdf()),
207  tdf.isTmp()
208  ),
209  OldTimeField<DimensionedField>(this->time().timeIndex()),
210  mesh_(tdf().mesh_),
211  dimensions_(tdf().dimensions_)
212 {
213  tdf.clear();
214 }
215 
216 
217 template<class Type, class GeoMesh, template<class> class PrimitiveField>
218 template<template<class> class PrimitiveField2>
220 (
221  const IOobject& io,
223  const bool checkIOFlags
224 )
225 :
227  PrimitiveField<Type>(df),
228  OldTimeField<DimensionedField>(this->time().timeIndex()),
229  mesh_(df.mesh_),
230  dimensions_(df.dimensions_)
231 {
232  if (!checkIOFlags || !readIfPresent())
233  {
234  copyOldTimes(io, df);
235  }
236 }
237 
238 
239 template<class Type, class GeoMesh, template<class> class PrimitiveField>
240 template<template<class> class PrimitiveField2>
242 (
243  const IOobject& io,
245  const bool reuse,
246  const bool checkIOFlags
247 )
248 :
249  regIOobject(io, df),
250  PrimitiveField<Type>(df, reuse),
251  OldTimeField<DimensionedField>(this->time().timeIndex()),
252  mesh_(df.mesh_),
253  dimensions_(df.dimensions_)
254 {
255  if (checkIOFlags)
256  {
257  readIfPresent();
258  }
259 }
260 
261 
262 template<class Type, class GeoMesh, template<class> class PrimitiveField>
264 (
265  const IOobject& io,
267  const bool checkIOFlags
268 )
269 :
270  regIOobject(io),
271  PrimitiveField<Type>
272  (
273  const_cast<DimensionedField<Type, GeoMesh, PrimitiveField>&>(tdf()),
274  tdf.isTmp()
275  ),
276  OldTimeField<DimensionedField>(this->time().timeIndex()),
277  mesh_(tdf().mesh_),
278  dimensions_(tdf().dimensions_)
279 {
280  tdf.clear();
281 
282  if (checkIOFlags)
283  {
284  readIfPresent();
285  }
286 }
287 
288 
289 template<class Type, class GeoMesh, template<class> class PrimitiveField>
290 template<template<class> class PrimitiveField2>
292 (
293  const word& newName,
295 )
296 :
297  regIOobject(newName, df, newName != df.name()),
298  PrimitiveField<Type>(df),
299  OldTimeField<DimensionedField>(this->time().timeIndex()),
300  mesh_(df.mesh_),
301  dimensions_(df.dimensions_)
302 {
303  copyOldTimes(newName, df);
304 }
305 
306 
307 template<class Type, class GeoMesh, template<class> class PrimitiveField>
308 template<template<class> class PrimitiveField2>
310 (
311  const word& newName,
313  const bool reuse
314 )
315 :
316  regIOobject(newName, df, true),
317  PrimitiveField<Type>(df, reuse),
318  OldTimeField<DimensionedField>(this->time().timeIndex()),
319  mesh_(df.mesh_),
320  dimensions_(df.dimensions_)
321 {}
322 
323 
324 template<class Type, class GeoMesh, template<class> class PrimitiveField>
326 (
327  const word& newName,
329 )
330 :
331  regIOobject(newName, tdf(), true),
332  PrimitiveField<Type>
333  (
334  const_cast<DimensionedField<Type, GeoMesh, PrimitiveField>&>(tdf()),
335  tdf.isTmp()
336  ),
337  OldTimeField<DimensionedField>(this->time().timeIndex()),
338  mesh_(tdf().mesh_),
339  dimensions_(tdf().dimensions_)
340 {
341  tdf.clear();
342 }
343 
344 
345 template<class Type, class GeoMesh, template<class> class PrimitiveField>
348 {
350  (
352  );
353 }
354 
355 
356 template<class Type, class GeoMesh, template<class> class PrimitiveField>
359 (
360  const word& name,
361  const Mesh& mesh,
362  const dimensionSet& ds,
363  const PrimitiveField<Type>& field
364 )
365 {
366  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
367 
369  (
371  (
372  IOobject
373  (
374  name,
375  mesh.thisDb().time().name(),
376  mesh.thisDb(),
377  IOobject::NO_READ,
378  IOobject::NO_WRITE,
379  cacheTmp
380  ),
381  mesh,
382  ds,
383  field
384  ),
385  cacheTmp
386  );
387 }
388 
389 
390 template<class Type, class GeoMesh, template<class> class PrimitiveField>
393 (
394  const word& name,
395  const Mesh& mesh,
396  const dimensionSet& ds
397 )
398 {
399  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
400 
402  (
404  (
405  IOobject
406  (
407  name,
408  mesh.thisDb().time().name(),
409  mesh.thisDb(),
410  IOobject::NO_READ,
411  IOobject::NO_WRITE,
412  cacheTmp
413  ),
415  ds,
416  false
417  ),
418  cacheTmp
419  );
420 }
421 
422 
423 template<class Type, class GeoMesh, template<class> class PrimitiveField>
426 (
427  const word& name,
428  const Mesh& mesh,
430 )
431 {
432  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
433 
435  (
437  (
438  IOobject
439  (
440  name,
443  IOobject::NO_READ,
444  IOobject::NO_WRITE,
445  cacheTmp
446  ),
447  mesh,
448  dt,
449  false
450  ),
451  cacheTmp
452  );
453 }
454 
455 
456 template<class Type, class GeoMesh, template<class> class PrimitiveField>
457 template<template<class> class PrimitiveField2>
460 (
461  const word& newName,
463 )
464 {
465  const bool cacheTmp = df.db().cacheTemporaryObject(newName);
466 
468  (
470  (
471  IOobject
472  (
473  newName,
474  df.instance(),
475  df.local(),
476  df.db(),
477  IOobject::NO_READ,
478  IOobject::NO_WRITE,
479  cacheTmp
480  ),
481  df
482  ),
483  cacheTmp
484  );
485 }
486 
487 
488 template<class Type, class GeoMesh, template<class> class PrimitiveField>
491 (
492  const word& newName,
494 )
495 {
496  const bool cacheTmp = tdf().db().cacheTemporaryObject(newName);
497 
499  (
501  (
502  IOobject
503  (
504  newName,
505  tdf().instance(),
506  tdf().local(),
507  tdf().db(),
508  IOobject::NO_READ,
509  IOobject::NO_WRITE,
510  cacheTmp
511  ),
512  tdf
513  ),
514  cacheTmp
515  );
516 }
517 
518 
519 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
520 
521 template<class Type, class GeoMesh, template<class> class PrimitiveField>
523 {
524  db().cacheTemporaryObject(*this);
525 }
526 
527 
528 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
529 
530 template<class Type, class GeoMesh, template<class> class PrimitiveField>
531 PrimitiveField<Type>&
533 {
534  this->setUpToDate();
535  storeOldTimes();
536  return *this;
537 }
538 
539 
540 template<class Type, class GeoMesh, template<class> class PrimitiveField>
541 tmp
542 <
544  <
546  GeoMesh,
547  Field
548  >
549 >
551 (
552  const direction d
553 ) const
554 {
556  (
558  (
559  name() + ".component(" + ::Foam::name(d) + ')',
560  mesh_,
561  dimensions_
562  )
563  );
564 
565  Foam::component(result.ref(), *this, d);
566 
567  return result;
568 }
569 
570 
571 template<class Type, class GeoMesh, template<class> class PrimitiveField>
572 template<template<class> class PrimitiveField2>
574 (
575  const direction d,
576  const DimensionedField
577  <
579  GeoMesh,
580  PrimitiveField2
581  >& df
582 )
583 {
584  PrimitiveField<Type>::replace(d, df);
585 }
586 
587 
588 template<class Type, class GeoMesh, template<class> class PrimitiveField>
589 template<template<class> class PrimitiveField2>
591 (
592  const direction d,
593  const tmp
594  <
596  <
598  GeoMesh,
599  PrimitiveField2
600  >
601  >& tdf
602 )
603 {
604  replace(d, tdf());
605  tdf.clear();
606 }
607 
608 
609 template<class Type, class GeoMesh, template<class> class PrimitiveField>
612 {
614  (
616  (
617  name() + ".T()",
618  mesh_,
619  dimensions_
620  )
621  );
622 
623  Foam::T(result.ref(), *this);
624 
625  return result;
626 }
627 
628 
629 template<class Type, class GeoMesh, template<class> class PrimitiveField>
632 {
633  dimensioned<Type> Average
634  (
635  this->name() + ".average()",
636  this->dimensions(),
637  gAverage(primitiveField())
638  );
639 
640  return Average;
641 }
642 
643 
644 template<class Type, class GeoMesh, template<class> class PrimitiveField>
645 template<template<class> class PrimitiveField2>
648 (
650 ) const
651 {
652  return
653  (
655  (
656  this->name() + ".weightedAverage(weights)",
657  this->dimensions(),
658  gSum(weightField*primitiveField())/gSum(weightField)
659  )
660  );
661 }
662 
663 
664 template<class Type, class GeoMesh, template<class> class PrimitiveField>
665 template<template<class> class PrimitiveField2>
668 (
670 ) const
671 {
672  dimensioned<Type> wa = weightedAverage(tweightField());
673  tweightField.clear();
674  return wa;
675 }
676 
677 
678 template<class Type, class GeoMesh, template<class> class PrimitiveField>
679 template<template<class> class PrimitiveField2>
681 (
683 )
684 {
685  checkFieldAssignment(*this, df);
686 
687  dimensions_ = df.dimensions();
688  PrimitiveField<Type>::operator=(df.primitiveField());
689 }
690 
691 
692 template<class Type, class GeoMesh, template<class> class PrimitiveField>
694 (
696 )
697 {
699 
700  checkFieldAssignment(*this, df);
701 
702  dimensions_ = df.dimensions();
703 
704  if (tdf.isTmp())
705  {
706  PrimitiveField<Type>::transfer(tdf.ref());
707  }
708  else
709  {
710  PrimitiveField<Type>::operator=(df.primitiveField());
711  }
712 
713  tdf.clear();
714 }
715 
716 
717 template<class Type, class GeoMesh, template<class> class PrimitiveField>
718 template<template<class> class PrimitiveField2>
720 (
722 )
723 {
724  reset(tdf());
725 
726  tdf.clear();
727 }
728 
729 
730 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
731 
732 template<class Type, class GeoMesh, template<class> class PrimitiveField>
734 (
736 )
737 {
738  checkFieldAssignment(*this, df);
739  checkFieldOperation(*this, df, "=");
740 
741  dimensions_ = df.dimensions();
742  PrimitiveField<Type>::operator=(df);
743 }
744 
745 
746 template<class Type, class GeoMesh, template<class> class PrimitiveField>
748 (
750 )
751 {
752  checkFieldAssignment(*this, df);
753  checkFieldOperation(*this, df, "=");
754 
755  dimensions_ = move(df.dimensions());
756  PrimitiveField<Type>::operator=(move(df));
757 }
758 
759 
760 template<class Type, class GeoMesh, template<class> class PrimitiveField>
762 (
764 )
765 {
767 
768  checkFieldAssignment(*this, df);
769  checkFieldOperation(*this, df, "=");
770 
771  dimensions_ = df.dimensions();
772 
773  if (tdf.isTmp())
774  {
775  primitiveFieldRef().transfer(tdf.ref());
776  }
777  else
778  {
780  }
781 
782  tdf.clear();
783 }
784 
785 
786 template<class Type, class GeoMesh, template<class> class PrimitiveField>
787 template<template<class> class PrimitiveField2>
789 (
791 )
792 {
793  checkFieldOperation(*this, df, "=");
794 
795  dimensions_ = df.dimensions();
796  PrimitiveField<Type>::operator=(df);
797 }
798 
799 
800 template<class Type, class GeoMesh, template<class> class PrimitiveField>
801 template<template<class> class PrimitiveField2>
803 (
805 )
806 {
808 
809  checkFieldOperation(*this, df, "=");
810 
811  dimensions_ = df.dimensions();
812  PrimitiveField<Type>::operator=(df);
813  tdf.clear();
814 }
815 
816 
817 template<class Type, class GeoMesh, template<class> class PrimitiveField>
819 (
820  const dimensioned<Type>& dt
821 )
822 {
823  dimensions_ = dt.dimensions();
824  PrimitiveField<Type>::operator=(dt.value());
825 }
826 
827 
828 template<class Type, class GeoMesh, template<class> class PrimitiveField>
830 {
831  PrimitiveField<Type>::operator=(Zero);
832 }
833 
834 
835 template<class Type, class GeoMesh, template<class> class PrimitiveField>
836 template<template<class> class PrimitiveField2>
838 (
840 )
841 {
842  this->operator=(df);
843 }
844 
845 
846 template<class Type, class GeoMesh, template<class> class PrimitiveField>
847 template<template<class> class PrimitiveField2>
849 (
851 )
852 {
853  this->operator=(tdf);
854 }
855 
856 
857 template<class Type, class GeoMesh, template<class> class PrimitiveField>
859 (
860  const dimensioned<Type>& dt
861 )
862 {
863  this->operator=(dt);
864 }
865 
866 
867 template<class Type, class GeoMesh, template<class> class PrimitiveField>
869 {
870  this->operator=(Zero);
871 }
872 
873 
874 #define COMPUTED_ASSIGNMENT(TYPE, op) \
875  \
876 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
877 template<template<class> class PrimitiveField2> \
878 void DimensionedField<Type, GeoMesh, PrimitiveField>::operator op \
879 ( \
880  const DimensionedField<TYPE, GeoMesh, PrimitiveField2>& df \
881 ) \
882 { \
883  checkFieldOperation(*this, df, #op); \
884  \
885  dimensions_ op df.dimensions(); \
886  PrimitiveField<Type>::operator op(df); \
887 } \
888  \
889 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
890 template<template<class> class PrimitiveField2> \
891 void DimensionedField<Type, GeoMesh, PrimitiveField>::operator op \
892 ( \
893  const tmp<DimensionedField<TYPE, GeoMesh, PrimitiveField2>>& tdf \
894 ) \
895 { \
896  operator op(tdf()); \
897  tdf.clear(); \
898 } \
899  \
900 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
901 void DimensionedField<Type, GeoMesh, PrimitiveField>::operator op \
902 ( \
903  const dimensioned<TYPE>& dt \
904 ) \
905 { \
906  dimensions_ op dt.dimensions(); \
907  PrimitiveField<Type>::operator op(dt.value()); \
908 }
909 
914 
915 #undef COMPUTED_ASSIGNMENT
916 
917 
918 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
919 
920 #undef checkFieldAssignment
921 #undef checkFieldOperation
922 
923 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
924 
925 } // End namespace Foam
926 
927 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
928 
929 #include "DimensionedFieldIO.C"
931 
932 // ************************************************************************* //
#define checkFieldOperation(df1, df2, op)
#define COMPUTED_ASSIGNMENT(TYPE, op)
#define checkFieldAssignment(df1, df2)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
PrimitiveField< Type >::cmptType cmptType
Component type of the elements of the field.
friend class DimensionedField
Declare friendship with other dimensioned fields.
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const dimensionSet & dimensions() const
Return dimensions.
PrimitiveField< Type > & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const fileName & local() const
Definition: IOobject.H:400
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:352
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:309
Class to add into field types to provide old-time storage and retrieval.
Definition: OldTimeField.H:113
void copyOldTimes(const IOobject &io, const OtherOldTime< OtherPrimitiveField > &)
Copy the old-times from the given field.
Definition: OldTimeField.C:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:426
const Time & time() const
Return time.
bool cacheTemporaryObject(const word &name) const
Return true if given name is in the cacheTemporaryObjects set.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
T clone(const T &t)
Definition: List.H:55
Type gAverage(const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
uint8_t direction
Definition: direction.H:45
label timeIndex
Definition: getTimeIndex.H:4
points setSize(newPointi)
conserve primitiveFieldRef()+