GeometricField.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 "GeometricField.H"
27 #include "Time.H"
28 #include "demandDrivenData.H"
29 #include "dictionary.H"
30 #include "localIOdictionary.H"
31 #include "solutionControl.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 #define checkFieldAssignment(gf1, gf2) \
36  \
37  if \
38  ( \
39  static_cast<const regIOobject*>(&gf1) \
40  == static_cast<const regIOobject*>(&gf2) \
41  ) \
42  { \
43  FatalErrorInFunction \
44  << "attempted assignment to self for field " \
45  << (gf1).name() << abort(FatalError); \
46  }
47 
48 
49 #define checkFieldOperation(gf1, gf2, op) \
50  \
51  if ((gf1).mesh() != (gf2).mesh()) \
52  { \
53  FatalErrorInFunction \
54  << "different mesh for fields " \
55  << (gf1).name() << " and " << (gf2).name() \
56  << " during operation " << op \
57  << abort(FatalError); \
58  }
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
62 
63 template<class Type, class GeoMesh, template<class> class PrimitiveField>
65 (
66  const dictionary& dict
67 )
68 {
69  Internal::readField(dict, "internalField");
70 
71  boundaryField_.readField(*this, dict.subDict("boundaryField"));
72 
73  // Don't use subOrEmptyDict here, or line numbers will be lost from any IO
74  // error messages. Use an actual sub-dict reference here if possible.
75  if (dict.found("sources"))
76  {
77  sources_.readField(*this, dict.subDict("sources"));
78  }
79  else
80  {
81  sources_.readField(*this, dictionary(dict.name()/"sources", dict));
82  }
83 
84  if (dict.found("referenceLevel"))
85  {
86  Type fieldAverage(pTraits<Type>(dict.lookup("referenceLevel")));
87 
88  Field<Type>::operator+=(fieldAverage);
89 
90  forAll(boundaryField_, patchi)
91  {
92  boundaryField_[patchi] == boundaryField_[patchi] + fieldAverage;
93  }
94  }
95 }
96 
97 
98 template<class Type, class GeoMesh, template<class> class PrimitiveField>
100 {
101  const localIOdictionary dict
102  (
103  IOobject
104  (
105  this->name(),
106  this->instance(),
107  this->local(),
108  this->db(),
109  IOobject::MUST_READ,
110  IOobject::NO_WRITE,
111  false
112  ),
113  typeName
114  );
115 
116  this->close();
117 
118  readFields(dict);
119 }
120 
121 
122 template<class Type, class GeoMesh, template<class> class PrimitiveField>
124 {
125  if
126  (
127  this->readOpt() == IOobject::MUST_READ
128  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
129  )
130  {
132  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
133  << " suggests that a read constructor for field " << this->name()
134  << " would be more appropriate." << endl;
135  }
136  else if
137  (
138  this->readOpt() == IOobject::READ_IF_PRESENT
139  && this->headerOk()
140  )
141  {
142  readFields();
143 
144  // Check compatibility between field and mesh
145  if (this->size() != GeoMesh::size(this->mesh()))
146  {
147  FatalIOErrorInFunction(this->readStream(typeName))
148  << " number of field elements = " << this->size()
149  << " number of mesh elements = "
150  << GeoMesh::size(this->mesh())
151  << exit(FatalIOError);
152  }
153 
154  readOldTimeIfPresent();
155 
156  return true;
157  }
158 
159  return false;
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
164 
165 template<class Type, class GeoMesh, template<class> class PrimitiveField>
167 (
168  const IOobject& io,
169  const Mesh& mesh,
170  const dimensionSet& ds,
171  const word& patchFieldType
172 )
173 :
174  Internal(io, mesh, ds, false),
175  OldTimeField<GeometricField>(this->time().timeIndex()),
176  fieldPrevIterPtr_(nullptr),
177  boundaryField_(mesh.boundary(), *this, patchFieldType),
178  sources_()
179 {
180  if (debug)
181  {
182  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
183  }
184 
185  readIfPresent();
186 }
187 
188 
189 template<class Type, class GeoMesh, template<class> class PrimitiveField>
191 (
192  const IOobject& io,
193  const Mesh& mesh,
194  const dimensionSet& ds,
195  const wordList& patchFieldTypes,
196  const wordList& actualPatchTypes,
197  const HashTable<word>& fieldSourceTypes
198 )
199 :
200  Internal(io, mesh, ds, false),
201  OldTimeField<GeometricField>(this->time().timeIndex()),
202  fieldPrevIterPtr_(nullptr),
203  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes),
204  sources_(*this, fieldSourceTypes)
205 {
206  if (debug)
207  {
208  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
209  }
210 
211  readIfPresent();
212 }
213 
214 
215 template<class Type, class GeoMesh, template<class> class PrimitiveField>
217 (
218  const IOobject& io,
219  const Mesh& mesh,
220  const dimensioned<Type>& dt,
221  const word& patchFieldType
222 )
223 :
224  Internal(io, mesh, dt, false),
225  OldTimeField<GeometricField>(this->time().timeIndex()),
226  fieldPrevIterPtr_(nullptr),
227  boundaryField_(mesh.boundary(), *this, patchFieldType),
228  sources_()
229 {
230  if (debug)
231  {
232  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
233  }
234 
235  boundaryField_ == dt.value();
236 
237  readIfPresent();
238 }
239 
240 
241 template<class Type, class GeoMesh, template<class> class PrimitiveField>
243 (
244  const IOobject& io,
245  const Mesh& mesh,
246  const dimensioned<Type>& dt,
247  const wordList& patchFieldTypes,
248  const wordList& actualPatchTypes,
249  const HashTable<word>& fieldSourceTypes
250 )
251 :
252  Internal(io, mesh, dt, false),
253  OldTimeField<GeometricField>(this->time().timeIndex()),
254  fieldPrevIterPtr_(nullptr),
255  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes),
256  sources_(*this, fieldSourceTypes)
257 {
258  if (debug)
259  {
260  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
261  }
262 
263  boundaryField_ == dt.value();
264 
265  readIfPresent();
266 }
267 
268 
269 template<class Type, class GeoMesh, template<class> class PrimitiveField>
271 (
272  const IOobject& io,
273  const Internal& diField,
274  const PtrList<Patch>& ptfl,
275  const HashPtrTable<Source>& stft
276 )
277 :
278  Internal(io, diField, false),
280  fieldPrevIterPtr_(nullptr),
281  boundaryField_(this->mesh().boundary(), *this, ptfl),
282  sources_(*this, stft)
283 {
284  if (debug)
285  {
287  << "Constructing from components" << endl << this->info() << endl;
288  }
289 
290  readIfPresent();
291 }
292 
293 
294 template<class Type, class GeoMesh, template<class> class PrimitiveField>
296 (
297  const IOobject& io,
298  const Mesh& mesh,
299  const dimensionSet& ds,
300  const PrimitiveField<Type>& iField,
301  const PtrList<Patch>& ptfl,
302  const HashPtrTable<Source>& stft
303 )
304 :
305  Internal(io, mesh, ds, iField),
306  OldTimeField<GeometricField>(this->time().timeIndex()),
307  fieldPrevIterPtr_(nullptr),
308  boundaryField_(mesh.boundary(), *this, ptfl),
309  sources_(*this, stft)
310 {
311  if (debug)
312  {
314  << "Constructing from components" << endl << this->info() << endl;
315  }
316 }
317 
318 
319 template<class Type, class GeoMesh, template<class> class PrimitiveField>
321 (
322  const IOobject& io,
323  const Mesh& mesh
324 )
325 :
326  Internal(io, mesh, dimless, false),
327  OldTimeField<GeometricField>(this->time().timeIndex()),
328  fieldPrevIterPtr_(nullptr),
329  boundaryField_(mesh.boundary()),
330  sources_()
331 {
332  readFields();
333 
334  // Check compatibility between field and mesh
335 
336  if (this->size() != GeoMesh::size(this->mesh()))
337  {
338  FatalIOErrorInFunction(this->readStream(typeName))
339  << " number of field elements = " << this->size()
340  << " number of mesh elements = " << GeoMesh::size(this->mesh())
341  << exit(FatalIOError);
342  }
343 
344  readOldTimeIfPresent();
345 
346  if (debug)
347  {
349  << "Finishing read-construction of" << endl << this->info() << endl;
350  }
351 }
352 
353 
354 template<class Type, class GeoMesh, template<class> class PrimitiveField>
356 (
357  const IOobject& io,
358  const Mesh& mesh,
359  const dictionary& dict
360 )
361 :
362  Internal(io, mesh, dimless, false),
363  OldTimeField<GeometricField>(this->time().timeIndex()),
364  fieldPrevIterPtr_(nullptr),
365  boundaryField_(mesh.boundary()),
366  sources_()
367 {
368  readFields(dict);
369 
370  // Check compatibility between field and mesh
371 
372  if (this->size() != GeoMesh::size(this->mesh()))
373  {
375  << " number of field elements = " << this->size()
376  << " number of mesh elements = " << GeoMesh::size(this->mesh())
377  << exit(FatalIOError);
378  }
379 
380  if (debug)
381  {
383  << "Finishing dictionary-construct of "
384  << endl << this->info() << endl;
385  }
386 }
387 
388 
389 template<class Type, class GeoMesh, template<class> class PrimitiveField>
391 (
393 )
394 :
395  Internal(gf),
397  fieldPrevIterPtr_(nullptr),
398  boundaryField_(*this, gf.boundaryField_),
399  sources_(*this, gf.sources_)
400 {
401  if (debug)
402  {
404  << "Constructing as copy" << endl << this->info() << endl;
405  }
406 
407  this->writeOpt() = IOobject::NO_WRITE;
408 }
409 
410 
411 template<class Type, class GeoMesh, template<class> class PrimitiveField>
412 template<template<class> class PrimitiveField2>
414 (
416 )
417 :
418  Internal(gf),
420  fieldPrevIterPtr_(nullptr),
421  boundaryField_(*this, gf.boundaryField_),
422  sources_(*this, gf.sources_)
423 {
424  if (debug)
425  {
427  << "Constructing as copy" << endl << this->info() << endl;
428  }
429 
430  this->writeOpt() = IOobject::NO_WRITE;
431 }
432 
433 
434 template<class Type, class GeoMesh, template<class> class PrimitiveField>
436 (
438 )
439 :
440  Internal(move(gf)),
441  OldTimeField<GeometricField>(move(gf)),
442  fieldPrevIterPtr_(nullptr),
443  boundaryField_(*this, gf.boundaryField_),
444  sources_(*this, gf.sources_)
445 {
446  if (debug)
447  {
449  << "Constructing by moving" << endl << this->info() << endl;
450  }
451 
452  this->writeOpt() = IOobject::NO_WRITE;
453 }
454 
455 
456 template<class Type, class GeoMesh, template<class> class PrimitiveField>
458 (
460 )
461 :
462  Internal
463  (
464  const_cast<GeometricField<Type, GeoMesh, PrimitiveField>&>(tgf()),
465  tgf.isTmp()
466  ),
468  fieldPrevIterPtr_(nullptr),
469  boundaryField_(*this, tgf().boundaryField_),
470  sources_(*this, tgf().sources_)
471 {
472  if (debug)
473  {
475  << "Constructing from tmp" << endl << this->info() << endl;
476  }
477 
478  this->writeOpt() = IOobject::NO_WRITE;
479 
480  tgf.clear();
481 }
482 
483 
484 template<class Type, class GeoMesh, template<class> class PrimitiveField>
485 template<template<class> class PrimitiveField2>
487 (
488  const IOobject& io,
490 )
491 :
492  Internal(io, gf, false),
494  fieldPrevIterPtr_(nullptr),
495  boundaryField_(*this, gf.boundaryField_),
496  sources_(*this, gf.sources_)
497 {
498  if (debug)
499  {
501  << "Constructing as copy resetting IO params"
502  << endl << this->info() << endl;
503  }
504 
505  if (!readIfPresent())
506  {
507  copyOldTimes(io, gf);
508  }
509 }
510 
511 
512 template<class Type, class GeoMesh, template<class> class PrimitiveField>
514 (
515  const IOobject& io,
517 )
518 :
519  Internal
520  (
521  io,
523  tgf.isTmp(),
524  false
525  ),
527  fieldPrevIterPtr_(nullptr),
528  boundaryField_(*this, tgf().boundaryField_),
529  sources_(*this, tgf().sources_)
530 {
531  if (debug)
532  {
534  << "Constructing from tmp resetting IO params"
535  << endl << this->info() << endl;
536  }
537 
538  tgf.clear();
539 
540  readIfPresent();
541 }
542 
543 
544 template<class Type, class GeoMesh, template<class> class PrimitiveField>
545 template<template<class> class PrimitiveField2>
547 (
548  const word& newName,
550 )
551 :
552  Internal(newName, gf),
554  fieldPrevIterPtr_(nullptr),
555  boundaryField_(*this, gf.boundaryField_),
556  sources_(*this, gf.sources_)
557 {
558  if (debug)
559  {
561  << "Constructing as copy resetting name"
562  << endl << this->info() << endl;
563  }
564 
565  copyOldTimes(newName, gf);
566 }
567 
568 
569 template<class Type, class GeoMesh, template<class> class PrimitiveField>
571 (
572  const word& newName,
574 )
575 :
576  Internal
577  (
578  newName,
579  const_cast<GeometricField<Type, GeoMesh, PrimitiveField>&>(tgf()),
580  tgf.isTmp()
581  ),
583  fieldPrevIterPtr_(nullptr),
584  boundaryField_(*this, tgf().boundaryField_),
585  sources_(*this, tgf().sources_)
586 {
587  if (debug)
588  {
590  << "Constructing from tmp resetting name"
591  << endl << this->info() << endl;
592  }
593 
594  tgf.clear();
595 }
596 
597 
598 template<class Type, class GeoMesh, template<class> class PrimitiveField>
599 template<template<class> class PrimitiveField2>
601 (
602  const IOobject& io,
604  const word& patchFieldType
605 )
606 :
607  Internal(io, gf, false),
609  fieldPrevIterPtr_(nullptr),
610  boundaryField_(this->mesh().boundary(), *this, patchFieldType),
611  sources_(*this, gf.sources_)
612 {
613  if (debug)
614  {
616  << "Constructing as copy resetting IO params"
617  << endl << this->info() << endl;
618  }
619 
620  boundaryField_ == gf.boundaryField_;
621 
622  if (!readIfPresent())
623  {
624  copyOldTimes(io, gf);
625  }
626 }
627 
628 
629 template<class Type, class GeoMesh, template<class> class PrimitiveField>
631 (
632  const IOobject& io,
634  const word& patchFieldType
635 )
636 :
637  Internal
638  (
639  io,
640  const_cast<GeometricField<Type, GeoMesh, PrimitiveField>&>(tgf()),
641  tgf.isTmp(),
642  false
643  ),
645  fieldPrevIterPtr_(nullptr),
646  boundaryField_(this->mesh().boundary(), *this, patchFieldType),
647  sources_(*this, tgf().sources_)
648 {
649  if (debug)
650  {
652  << "Constructing as copy resetting IO params"
653  << endl << this->info() << endl;
654  }
655 
656  boundaryField_ == tgf().boundaryField_;
657 
658  tgf.clear();
659 
660  readIfPresent();
661 }
662 
663 
664 template<class Type, class GeoMesh, template<class> class PrimitiveField>
665 template<template<class> class PrimitiveField2>
667 (
668  const IOobject& io,
670  const word& patchFieldType
671 )
672 :
673  Internal(io, df, false),
674  OldTimeField<GeometricField>(this->time().timeIndex()),
675  fieldPrevIterPtr_(nullptr),
676  boundaryField_(this->mesh().boundary(), *this, patchFieldType),
677  sources_()
678 {
679  if (debug)
680  {
682  << "Constructing from components" << endl << this->info() << endl;
683  }
684 
685  if (!readIfPresent())
686  {
687  boundaryField_.evaluate();
688  }
689 }
690 
691 
692 template<class Type, class GeoMesh, template<class> class PrimitiveField>
694 (
695  const IOobject& io,
696  const tmp<Internal>& tdf,
697  const word& patchFieldType
698 )
699 :
700  Internal(io, const_cast<Internal&>(tdf()), tdf.isTmp(), false),
702  fieldPrevIterPtr_(nullptr),
703  boundaryField_(this->mesh().boundary(), *this, patchFieldType),
704  sources_()
705 {
706  if (debug)
707  {
709  << "Constructing from components" << endl << this->info() << endl;
710  }
711 
712  if (!readIfPresent())
713  {
714  boundaryField_.evaluate();
715  }
716 }
717 
718 
719 template<class Type, class GeoMesh, template<class> class PrimitiveField>
720 template<template<class> class PrimitiveField2>
722 (
723  const IOobject& io,
725  const wordList& patchFieldTypes,
726  const wordList& actualPatchTypes,
727  const HashTable<word>& fieldSourceTypes
728 )
729 :
730  Internal(io, gf, false),
732  fieldPrevIterPtr_(nullptr),
733  boundaryField_
734  (
735  this->mesh().boundary(),
736  *this,
737  patchFieldTypes,
738  actualPatchTypes
739  ),
740  sources_(*this, fieldSourceTypes)
741 {
742  if (debug)
743  {
745  << "Constructing as copy resetting IO params and patch types"
746  << endl << this->info() << endl;
747  }
748 
749  boundaryField_ == gf.boundaryField_;
750 
751  if (!readIfPresent())
752  {
753  copyOldTimes(io, gf);
754  }
755 }
756 
757 
758 template<class Type, class GeoMesh, template<class> class PrimitiveField>
760 (
761  const IOobject& io,
763  const wordList& patchFieldTypes,
764  const wordList& actualPatchTypes,
765  const HashTable<word>& fieldSourceTypes
766 )
767 :
768  Internal
769  (
770  io,
771  const_cast<GeometricField<Type, GeoMesh, PrimitiveField>&>(tgf()),
772  tgf.isTmp(),
773  false
774  ),
776  fieldPrevIterPtr_(nullptr),
777  boundaryField_
778  (
779  this->mesh().boundary(),
780  *this,
781  patchFieldTypes,
782  actualPatchTypes
783  ),
784  sources_(*this, fieldSourceTypes)
785 {
786  if (debug)
787  {
789  << "Constructing from tmp resetting IO params and patch types"
790  << endl << this->info() << endl;
791  }
792 
793  boundaryField_ == tgf().boundaryField_;
794 
795  tgf.clear();
796 
797  readIfPresent();
798 }
799 
800 
801 template<class Type, class GeoMesh, template<class> class PrimitiveField>
802 template<template<class> class PrimitiveField2>
804 (
805  const IOobject& io,
807  const wordList& patchFieldTypes,
808  const wordList& actualPatchTypes,
809  const HashTable<word>& fieldSourceTypes
810 )
811 :
812  Internal(io, df, false),
814  fieldPrevIterPtr_(nullptr),
815  boundaryField_
816  (
817  this->mesh().boundary(),
818  *this,
819  patchFieldTypes,
820  actualPatchTypes
821  ),
822  sources_(*this, fieldSourceTypes)
823 {
824  if (debug)
825  {
827  << "Constructing from internal field and patch types"
828  << endl << this->info() << endl;
829  }
830 
831  readIfPresent();
832 }
833 
834 
835 template<class Type, class GeoMesh, template<class> class PrimitiveField>
837 (
838  const IOobject& io,
839  const tmp<Internal>& tdf,
840  const wordList& patchFieldTypes,
841  const wordList& actualPatchTypes,
842  const HashTable<word>& fieldSourceTypes
843 )
844 :
845  Internal(io, const_cast<Internal&>(tdf()), tdf.isTmp(), false),
847  fieldPrevIterPtr_(nullptr),
848  boundaryField_
849  (
850  this->mesh().boundary(),
851  *this,
852  patchFieldTypes,
853  actualPatchTypes
854  ),
855  sources_(*this, fieldSourceTypes)
856 {
857  if (debug)
858  {
860  << "Constructing from tmp internal field and patch types"
861  << endl << this->info() << endl;
862  }
863 
864  tdf.clear();
865 
866  readIfPresent();
867 }
868 
869 
870 template<class Type, class GeoMesh, template<class> class PrimitiveField>
873 {
875  (
877  );
878 }
879 
880 
881 template<class Type, class GeoMesh, template<class> class PrimitiveField>
884 {
886  (
888  (
889  IOobject
890  (
891  this->name(),
892  this->mesh().thisDb().time().name(),
893  this->mesh().thisDb(),
896  false
897  ),
898  *this,
899  Patch::calculatedType()
900  )
901  );
902 }
903 
904 
905 template<class Type, class GeoMesh, template<class> class PrimitiveField>
908 (
909  const word& name,
910  const Internal& diField,
911  const PtrList<Patch>& ptfl,
912  const HashPtrTable<Source>& stft
913 )
914 {
915  const bool cacheTmp = diField.mesh().thisDb().cacheTemporaryObject(name);
916 
918  (
920  (
921  IOobject
922  (
923  name,
924  diField.mesh().thisDb().time().name(),
925  diField.mesh().thisDb(),
928  cacheTmp
929  ),
930  diField,
931  ptfl,
932  stft
933  ),
934  cacheTmp
935  );
936 }
937 
938 
939 template<class Type, class GeoMesh, template<class> class PrimitiveField>
942 (
943  const word& name,
944  const Mesh& mesh,
945  const dimensionSet& ds,
946  const word& patchFieldType
947 )
948 {
949  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
950 
952  (
954  (
955  IOobject
956  (
957  name,
958  mesh.thisDb().time().name(),
959  mesh.thisDb(),
962  cacheTmp
963  ),
964  mesh,
965  ds,
967  ),
968  cacheTmp
969  );
970 }
971 
972 
973 template<class Type, class GeoMesh, template<class> class PrimitiveField>
976 (
977  const word& name,
978  const Mesh& mesh,
979  const dimensioned<Type>& dt,
980  const word& patchFieldType
981 )
982 {
983  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
984 
986  (
988  (
989  IOobject
990  (
991  name,
992  mesh.thisDb().time().name(),
993  mesh.thisDb(),
996  cacheTmp
997  ),
998  mesh,
999  dt,
1001  ),
1002  cacheTmp
1003  );
1004 }
1005 
1006 
1007 
1008 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1011 (
1012  const word& name,
1013  const Mesh& mesh,
1014  const dimensioned<Type>& dt,
1015  const wordList& patchFieldTypes,
1016  const wordList& actualPatchTypes,
1017  const HashTable<word>& fieldSourceTypes
1018 )
1019 {
1020  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
1021 
1023  (
1025  (
1026  IOobject
1027  (
1028  name,
1029  mesh.thisDb().time().name(),
1030  mesh.thisDb(),
1033  cacheTmp
1034  ),
1035  mesh,
1036  dt,
1037  patchFieldTypes,
1038  actualPatchTypes,
1039  fieldSourceTypes
1040  ),
1041  cacheTmp
1042  );
1043 }
1044 
1045 
1046 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1049 (
1050  const word& newName,
1052 )
1053 {
1054  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
1055 
1057  (
1059  (
1060  IOobject
1061  (
1062  newName,
1063  tgf().instance(),
1064  tgf().local(),
1065  tgf().db(),
1068  cacheTmp
1069  ),
1070  tgf
1071  ),
1072  cacheTmp
1073  );
1074 }
1075 
1076 
1077 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1078 template<template<class> class PrimitiveField2>
1081 (
1082  const word& newName,
1084  const word& patchFieldType
1085 )
1086 {
1087  const bool cacheTmp = gf.db().cacheTemporaryObject(newName);
1088 
1090  (
1092  (
1093  IOobject
1094  (
1095  newName,
1096  gf.instance(),
1097  gf.local(),
1098  gf.db(),
1101  cacheTmp
1102  ),
1103  gf,
1105  ),
1106  cacheTmp
1107  );
1108 }
1109 
1110 
1111 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1114 (
1115  const word& newName,
1117  const word& patchFieldType
1118 )
1119 {
1120  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
1121 
1123  (
1125  (
1126  IOobject
1127  (
1128  newName,
1129  tgf().instance(),
1130  tgf().local(),
1131  tgf().db(),
1134  cacheTmp
1135  ),
1136  tgf,
1138  ),
1139  cacheTmp
1140  );
1141 }
1142 
1143 
1144 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1145 template<template<class> class PrimitiveField2>
1148 (
1149  const word& newName,
1151  const word& patchFieldType
1152 )
1153 {
1154  const bool cacheTmp = df.db().cacheTemporaryObject(newName);
1155 
1157  (
1159  (
1160  IOobject
1161  (
1162  newName,
1163  df.instance(),
1164  df.local(),
1165  df.db(),
1168  cacheTmp
1169  ),
1170  df,
1172  ),
1173  cacheTmp
1174  );
1175 }
1176 
1177 
1178 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1181 (
1182  const word& newName,
1183  const tmp<Internal>& tdf,
1184  const word& patchFieldType
1185 )
1186 {
1187  const bool cacheTmp = tdf().db().cacheTemporaryObject(newName);
1188 
1190  (
1192  (
1193  IOobject
1194  (
1195  newName,
1196  tdf().instance(),
1197  tdf().local(),
1198  tdf().db(),
1201  cacheTmp
1202  ),
1203  tdf,
1205  ),
1206  cacheTmp
1207  );
1208 }
1209 
1210 
1211 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1212 template<template<class> class PrimitiveField2>
1215 (
1216  const word& newName,
1218  const wordList& patchFieldTypes,
1219  const wordList& actualPatchTypes,
1220  const HashTable<word>& fieldSourceTypes
1221 )
1222 {
1223  const bool cacheTmp = gf.db().cacheTemporaryObject(newName);
1224 
1226  (
1228  (
1229  IOobject
1230  (
1231  newName,
1232  gf.instance(),
1233  gf.local(),
1234  gf.db(),
1237  cacheTmp
1238  ),
1239  gf,
1240  patchFieldTypes,
1241  actualPatchTypes,
1242  fieldSourceTypes
1243  ),
1244  cacheTmp
1245  );
1246 }
1247 
1248 
1249 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1252 (
1253  const word& newName,
1255  const wordList& patchFieldTypes,
1256  const wordList& actualPatchTypes,
1257  const HashTable<word>& fieldSourceTypes
1258 )
1259 {
1260  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
1261 
1263  (
1265  (
1266  IOobject
1267  (
1268  newName,
1269  tgf().instance(),
1270  tgf().local(),
1271  tgf().db(),
1274  cacheTmp
1275  ),
1276  tgf,
1277  patchFieldTypes,
1278  actualPatchTypes,
1279  fieldSourceTypes
1280  ),
1281  cacheTmp
1282  );
1283 }
1284 
1285 
1286 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1287 template<template<class> class PrimitiveField2>
1290 (
1291  const word& newName,
1293  const wordList& patchFieldTypes,
1294  const wordList& actualPatchTypes,
1295  const HashTable<word>& fieldSourceTypes
1296 )
1297 {
1298  const bool cacheTmp = df.db().cacheTemporaryObject(newName);
1299 
1301  (
1303  (
1304  IOobject
1305  (
1306  newName,
1307  df.instance(),
1308  df.local(),
1309  df.db(),
1312  cacheTmp
1313  ),
1314  df,
1315  patchFieldTypes,
1316  actualPatchTypes,
1317  fieldSourceTypes
1318  ),
1319  cacheTmp
1320  );
1321 }
1322 
1323 
1324 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1327 (
1328  const word& newName,
1329  const tmp<Internal>& tdf,
1330  const wordList& patchFieldTypes,
1331  const wordList& actualPatchTypes,
1332  const HashTable<word>& fieldSourceTypes
1333 )
1334 {
1335  const bool cacheTmp = tdf().db().cacheTemporaryObject(newName);
1336 
1338  (
1340  (
1341  IOobject
1342  (
1343  newName,
1344  tdf().instance(),
1345  tdf().local(),
1346  tdf().db(),
1349  cacheTmp
1350  ),
1351  tdf,
1352  patchFieldTypes,
1353  actualPatchTypes,
1354  fieldSourceTypes
1355  ),
1356  cacheTmp
1357  );
1358 }
1359 
1360 
1361 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
1362 
1363 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1365 {
1366  this->db().cacheTemporaryObject(*this);
1367 
1368  clearPrevIter();
1369 }
1370 
1371 
1372 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1373 
1374 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1377 {
1378  this->setUpToDate();
1379  storeOldTimes();
1380  return *this;
1381 }
1382 
1383 
1384 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1385 typename
1388 {
1389  this->setUpToDate();
1390  storeOldTimes();
1391  return *this;
1392 }
1393 
1394 
1395 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1396 typename
1399 {
1400  this->setUpToDate();
1401  storeOldTimes();
1402  return boundaryField_;
1403 }
1404 
1405 
1406 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1407 typename
1411 {
1412  this->setUpToDate();
1413  return boundaryField_;
1414 }
1415 
1416 
1417 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1418 typename
1421 {
1422  this->setUpToDate();
1423  storeOldTimes();
1424  return sources_;
1425 }
1426 
1427 
1428 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1430 {
1431  if (!fieldPrevIterPtr_)
1432  {
1433  if (debug)
1434  {
1436  << "Allocating previous iteration field" << endl
1437  << this->info() << endl;
1438  }
1439 
1440  fieldPrevIterPtr_ =
1442  (
1443  this->name() + "PrevIter",
1444  *this
1445  );
1446  }
1447  else
1448  {
1449  *fieldPrevIterPtr_ == *this;
1450  }
1451 }
1452 
1453 
1454 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1456 {
1457  deleteDemandDrivenData(fieldPrevIterPtr_);
1458 }
1459 
1460 
1461 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1464 {
1465  if (!fieldPrevIterPtr_)
1466  {
1468  << "previous iteration field" << endl << this->info() << endl
1469  << " not stored."
1470  << " Use field.storePrevIter() at start of iteration."
1471  << abort(FatalError);
1472  }
1473 
1474  return *fieldPrevIterPtr_;
1475 }
1476 
1477 
1478 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1481 {
1482  this->setUpToDate();
1483  storeOldTimes();
1484  boundaryField_.evaluate();
1485 }
1486 
1487 
1488 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1489 template<template<class> class PrimitiveField2>
1491 (
1493 )
1494 {
1495  Internal::reset(gf);
1496 
1497  boundaryField_.reset(gf.boundaryField());
1498  sources_.reset(*this, gf.sources());
1499 }
1500 
1501 
1502 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1504 (
1506 )
1507 {
1509 
1510  checkFieldAssignment(*this, gf);
1511 
1512  this->dimensions() = gf.dimensions();
1513 
1514  if (tgf.isTmp())
1515  {
1516  PrimitiveField<Type>::transfer(tgf.ref());
1517  }
1518  else
1519  {
1520  PrimitiveField<Type>::operator=(gf.primitiveField());
1521  }
1522 
1523  boundaryField_.reset(gf.boundaryField());
1524  sources_.reset(*this, gf.sources());
1525 
1526  tgf.clear();
1527 }
1528 
1529 
1530 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1531 template<template<class> class PrimitiveField2>
1533 (
1535 )
1536 {
1537  reset(tgf());
1538 
1539  tgf.clear();
1540 }
1541 
1542 
1543 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1545 {
1546  // Search all boundary conditions, if any are
1547  // fixed-value or mixed (Robin) do not set reference level for solution.
1548 
1549  bool needRef = true;
1550 
1551  forAll(boundaryField_, patchi)
1552  {
1553  if (boundaryField_[patchi].fixesValue())
1554  {
1555  needRef = false;
1556  break;
1557  }
1558  }
1559 
1560  reduce(needRef, andOp<bool>());
1561 
1562  return needRef;
1563 }
1564 
1565 
1566 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1568 (
1569  const scalar alpha
1570 )
1571 {
1572  if (alpha < 1)
1573  {
1574  if (debug)
1575  {
1577  << "Relaxing" << endl << this->info()
1578  << " by " << alpha << endl;
1579  }
1580 
1581  operator==(prevIter() + alpha*(*this - prevIter()));
1582  }
1583 }
1584 
1585 
1586 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1587 Foam::scalar
1589 {
1590  if
1591  (
1593  && this->mesh().solution().relaxField(this->name() + "Final")
1594  )
1595  {
1596  return this->mesh().solution().fieldRelaxationFactor
1597  (
1598  this->name() + "Final"
1599  );
1600  }
1601  else if (this->mesh().solution().relaxField(this->name()))
1602  {
1603  return this->mesh().solution().fieldRelaxationFactor(this->name());
1604  }
1605  else
1606  {
1607  return 1;
1608  }
1609 }
1610 
1611 
1612 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1614 {
1615  relax(relaxationFactor());
1616 }
1617 
1618 
1619 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1620 template<template<class> class PrimitiveField2>
1622 (
1624  const scalar alpha
1625 )
1626 {
1627  if (alpha < 1)
1628  {
1629  if (debug)
1630  {
1632  << "Relaxing" << endl << this->info()
1633  << " by " << alpha << endl;
1634  }
1635 
1636  operator==(*this + alpha*(tgf - *this));
1637  }
1638  else
1639  {
1640  operator==(tgf);
1641  }
1642 }
1643 
1644 
1645 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1646 template<template<class> class PrimitiveField2>
1648 (
1650 )
1651 {
1652  relax(tgf, relaxationFactor());
1653 }
1654 
1655 
1656 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1658 (
1659  bool final
1660 ) const
1661 {
1662  if (final)
1663  {
1664  return this->name() + "Final";
1665  }
1666  else
1667  {
1668  return this->name();
1669  }
1670 }
1671 
1672 
1673 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1675 (
1676  Ostream& os
1677 ) const
1678 {
1679  os << "min/max(" << this->name() << ") = "
1680  << Foam::min(*this).value() << ", "
1681  << Foam::max(*this).value()
1682  << endl;
1683 }
1684 
1685 
1686 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1688 (
1689  Ostream& os
1690 ) const
1691 {
1692  os << *this;
1693  return os.good();
1694 }
1695 
1696 
1697 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1698 
1699 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1702 {
1704  (
1706  (
1707  this->name() + ".T()",
1708  this->mesh(),
1709  this->dimensions()
1710  )
1711  );
1712 
1713  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1714  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1715 
1716  return result;
1717 }
1718 
1719 
1720 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1721 Foam::tmp
1722 <
1724  <
1726  GeoMesh,
1727  Foam::Field
1728  >
1729 >
1731 (
1732  const direction d
1733 ) const
1734 {
1736  (
1738  (
1739  this->name() + ".component(" + Foam::name(d) + ')',
1740  this->mesh(),
1741  this->dimensions()
1742  )
1743  );
1744 
1745  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1746  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1747 
1748  return Component;
1749 }
1750 
1751 
1752 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1753 template<template<class> class PrimitiveField2>
1755 (
1756  const direction d,
1757  const GeometricField
1758  <
1760  GeoMesh,
1761  PrimitiveField2
1762  >& gcf
1763 )
1764 {
1765  primitiveFieldRef().replace(d, gcf.primitiveField());
1766  boundaryFieldRef().replace(d, gcf.boundaryField());
1767 }
1768 
1769 
1770 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1771 template<template<class> class PrimitiveField2>
1773 (
1774  const direction d,
1775  const tmp
1776  <
1778  <
1780  GeoMesh,
1781  PrimitiveField2
1782  >
1783  >& gcf
1784 )
1785 {
1786  replace(d, gcf());
1787  gcf.clear();
1788 }
1789 
1790 
1791 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1793 (
1794  const direction d,
1795  const dimensioned<cmptType>& ds
1796 )
1797 {
1798  primitiveFieldRef().replace(d, ds.value());
1799  boundaryFieldRef().replace(d, ds.value());
1800 }
1801 
1802 
1803 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1805 (
1806  const dimensioned<Type>& dt
1807 )
1808 {
1809  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1810  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1811 }
1812 
1813 
1814 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1816 (
1817  const dimensioned<Type>& dt
1818 )
1819 {
1820  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1821  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1822 }
1823 
1824 
1825 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1827 (
1828  const dimensioned<Type>& minDt,
1829  const dimensioned<Type>& maxDt
1830 )
1831 {
1832  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1833  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1834  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1835  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1836 }
1837 
1838 
1839 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1841 {
1842  primitiveFieldRef().negate();
1843  boundaryFieldRef().negate();
1844 }
1845 
1846 
1847 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1848 
1849 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1851 (
1853 )
1854 {
1855  checkFieldAssignment(*this, gf);
1856  checkFieldOperation(*this, gf, "=");
1857 
1858  internalFieldRef() = gf.internalField();
1859  boundaryFieldRef() = gf.boundaryField();
1860 }
1861 
1862 
1863 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1865 (
1867 )
1868 {
1869  checkFieldAssignment(*this, gf);
1870  checkFieldOperation(*this, gf, "=");
1871 
1872  internalFieldRef() = move(gf.internalField());
1873  boundaryFieldRef() = move(gf.boundaryField());
1874 }
1875 
1876 
1877 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1878 template<template<class> class PrimitiveField2>
1880 (
1882 )
1883 {
1884  checkFieldOperation(*this, gf, "=");
1885 
1886  internalFieldRef() = gf.internalField();
1887  boundaryFieldRef() = gf.boundaryField();
1888 }
1889 
1890 
1891 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1893 (
1895 )
1896 {
1898 
1899  checkFieldAssignment(*this, gf);
1900  checkFieldOperation(*this, gf, "=");
1901 
1902  this->dimensions() = gf.dimensions();
1903 
1904  if (tgf.isTmp())
1905  {
1906  primitiveFieldRef().transfer(tgf.ref());
1907  }
1908  else
1909  {
1911  }
1912 
1913  boundaryFieldRef() = gf.boundaryField();
1914 
1915  tgf.clear();
1916 }
1917 
1918 
1919 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1920 template<template<class> class PrimitiveField2>
1922 (
1924 )
1925 {
1927 
1928  checkFieldOperation(*this, gf, "=");
1929 
1930  internalFieldRef() = gf.internalField();
1931  boundaryFieldRef() = gf.boundaryField();
1932 
1933  tgf.clear();
1934 }
1935 
1936 
1937 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1939 (
1940  const dimensioned<Type>& dt
1941 )
1942 {
1943  internalFieldRef() = dt;
1944  boundaryFieldRef() = dt.value();
1945 }
1946 
1947 
1948 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1950 (
1951  const zero&
1952 )
1953 {
1954  internalFieldRef() = Zero;
1955  boundaryFieldRef() = Zero;
1956 }
1957 
1958 
1959 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1960 template<template<class> class PrimitiveField2>
1962 (
1964 )
1965 {
1966  checkFieldOperation(*this, gf, "==");
1967 
1968  internalFieldRef() = gf.internalField();
1969  boundaryFieldRef() == gf.boundaryField();
1970 }
1971 
1972 
1973 template<class Type, class GeoMesh, template<class> class PrimitiveField>
1975 (
1977 )
1978 {
1980 
1981  checkFieldOperation(*this, gf, "==");
1982 
1983  this->dimensions() = gf.dimensions();
1984 
1985  if (tgf.isTmp())
1986  {
1987  primitiveFieldRef().transfer(tgf.ref());
1988  }
1989  else
1990  {
1992  }
1993 
1994  boundaryFieldRef() == gf.boundaryField();
1995 
1996  tgf.clear();
1997 }
1998 
1999 
2000 template<class Type, class GeoMesh, template<class> class PrimitiveField>
2001 template<template<class> class PrimitiveField2>
2003 (
2005 )
2006 {
2008 
2009  checkFieldOperation(*this, gf, "=");
2010 
2011  internalFieldRef() = gf.internalField();
2012  boundaryFieldRef() == gf.boundaryField();
2013 
2014  tgf.clear();
2015 }
2016 
2017 
2018 template<class Type, class GeoMesh, template<class> class PrimitiveField>
2020 (
2021  const dimensioned<Type>& dt
2022 )
2023 {
2024  internalFieldRef() = dt;
2025  boundaryFieldRef() == dt.value();
2026 }
2027 
2028 
2029 template<class Type, class GeoMesh, template<class> class PrimitiveField>
2031 (
2032  const zero&
2033 )
2034 {
2035  internalFieldRef() = Zero;
2036  boundaryFieldRef() == Zero;
2037 }
2038 
2039 
2040 #define COMPUTED_ASSIGNMENT(TYPE, op) \
2041  \
2042 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
2043 template<template<class> class PrimitiveField2> \
2044 void Foam::GeometricField<Type, GeoMesh, PrimitiveField>::operator op \
2045 ( \
2046  const GeometricField<TYPE, GeoMesh, PrimitiveField2>& gf \
2047 ) \
2048 { \
2049  checkFieldOperation(*this, gf, #op); \
2050  \
2051  internalFieldRef() op gf.internalField(); \
2052  boundaryFieldRef() op gf.boundaryField(); \
2053 } \
2054  \
2055 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
2056 template<template<class> class PrimitiveField2> \
2057 void Foam::GeometricField<Type, GeoMesh, PrimitiveField>::operator op \
2058 ( \
2059  const tmp<GeometricField<TYPE, GeoMesh, PrimitiveField2>>& tgf \
2060 ) \
2061 { \
2062  operator op(tgf()); \
2063  tgf.clear(); \
2064 } \
2065  \
2066 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
2067 void Foam::GeometricField<Type, GeoMesh, PrimitiveField>::operator op \
2068 ( \
2069  const dimensioned<TYPE>& dt \
2070 ) \
2071 { \
2072  internalFieldRef() op dt; \
2073  boundaryFieldRef() op dt.value(); \
2074 }
2075 
2080 
2081 #undef COMPUTED_ASSIGNMENT
2082 
2083 
2084 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2085 
2086 template<class Type, class GeoMesh, template<class> class PrimitiveField>
2087 Foam::Ostream& Foam::operator<<
2088 (
2089  Ostream& os,
2091 )
2092 {
2093  gf().writeData(os, "internalField");
2094  os << nl;
2095  gf.boundaryField().writeEntry("boundaryField", os);
2096 
2097  if (!gf.sources_.empty())
2098  {
2099  os << nl;
2100  gf.sources().writeEntry("sources", os);
2101  }
2102 
2103  // Check state of IOstream
2104  os.check
2105  (
2106  "Ostream& operator<<(Ostream&, "
2107  "const GeometricField<Type, GeoMesh, PrimitiveField>&)"
2108  );
2109 
2110  return (os);
2111 }
2112 
2113 
2114 template<class Type, class GeoMesh, template<class> class PrimitiveField>
2115 Foam::Ostream& Foam::operator<<
2116 (
2117  Ostream& os,
2119 )
2120 {
2121  os << tgf();
2122  tgf.clear();
2123  return os;
2124 }
2125 
2126 
2127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2128 
2129 #undef checkFieldAssignment
2130 #undef checkFieldOperation
2131 
2132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2133 
2134 #include "GeometricFieldFunctions.C"
2135 
2136 // ************************************************************************* //
EaEqn relax()
#define checkFieldOperation(gf1, gf2, op)
#define COMPUTED_ASSIGNMENT(TYPE, op)
#define checkFieldAssignment(gf1, gf2)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const dimensionSet & dimensions() const
Return dimensions.
PrimitiveField< Type > FieldType
Type of the field from which this DimensionedField is derived.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricBoundaryField class.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
Part of a geometric field used for setting the values associated with optional sources.
void writeEntry(const word &keyword, Ostream &os) const
Write sources as dictionary entry.
Generic GeometricField class.
void max(const dimensioned< Type > &)
tmp< GeometricField< Type, GeoMesh, Field > > T() const
Return transpose (only if it is a tensor field)
Sources & sourcesRef()
Return a reference to the sources.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
void relax()
Relax current field with respect to the cached previous iteration.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
const GeometricField< Type, GeoMesh, Field > & prevIter() const
Return previous iteration field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
Boundary & boundaryFieldRefNoStoreOldTimes()
Return a reference to the boundary field without storing old times.
tmp< GeometricField< cmptType, GeoMesh, Field > > component(const direction) const
Return a component of the field.
const Sources & sources() const
Return const-reference to the sources.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
void min(const dimensioned< Type > &)
friend class GeometricField
Declare friendship with other geometric fields.
Internal & internalFieldRef()
Return a reference to the dimensioned internal field.
Field< Type >::cmptType cmptType
Component type of the elements of the field.
void replace(const direction, const GeometricField< cmptType, GeoMesh, PrimitiveField2 > &)
Replace a component field of the field.
tmp< GeometricField< Type, GeoMesh, PrimitiveField > > clone() const
Clone.
bool needReference() const
Does the field need a reference level for solution.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
virtual ~GeometricField()
Destructor.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
void reset(const GeometricField< Type, GeoMesh, PrimitiveField2 > &)
Reset the field contents to the given field.
scalar relaxationFactor() const
Return the field relaxation factor read from fvSolution.
void clearPrevIter()
Delete the previous iteration field.
void storePrevIter() const
Store the field as the previous iteration value.
word select(bool final) const
Select the final iteration parameters if `final' is true.
void correctBoundaryConditions()
Correct boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
tmp< GeometricField< Type, GeoMesh, PrimitiveField > > cloneUnSliced() const
Clone un-sliced.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
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
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:485
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Class to add into field types to provide old-time storage and retrieval.
Definition: OldTimeField.H:113
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:198
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:426
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
const Time & time() const
Return time.
bool cacheTemporaryObject(const word &name) const
Return true if given name is in the cacheTemporaryObjects set.
static bool finalIteration(const objectRegistry &registry)
Lookup solutionControl from the objectRegistry and return finalIter.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
scalar fieldRelaxationFactor(const word &name) const
Return the relaxation factor for the given field.
Definition: solution.C:271
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:253
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)
Template functions to aid in the implementation of demand driven data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
#define InfoInFunction
Report an information message using Foam::Info.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word patchFieldType(const PatchField &pf)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
void operator+=(fvMatrix< Type > &fvEqn, const CarrierEqn< Type > &cEqn)
Add to a finite-volume equation.
Definition: CarrierEqn.C:71
static const char nl
Definition: Ostream.H:267
uint8_t direction
Definition: direction.H:45
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
dictionary dict
conserve primitiveFieldRef()+