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-2023 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 "data.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 #define checkField(gf1, gf2, op) \
36 if ((gf1).mesh() != (gf2).mesh()) \
37 { \
38  FatalErrorInFunction \
39  << "different mesh for fields " \
40  << (gf1).name() << " and " << (gf2).name() \
41  << " during operatrion " << op \
42  << abort(FatalError); \
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47 
48 template<class Type, template<class> class PatchField, class GeoMesh>
50 (
51  const dictionary& dict
52 )
53 {
54  Internal::readField(dict, "internalField");
55 
56  boundaryField_.readField(*this, dict.subDict("boundaryField"));
57 
58  if (dict.found("referenceLevel"))
59  {
60  Type fieldAverage(pTraits<Type>(dict.lookup("referenceLevel")));
61 
62  Field<Type>::operator+=(fieldAverage);
63 
64  forAll(boundaryField_, patchi)
65  {
66  boundaryField_[patchi] == boundaryField_[patchi] + fieldAverage;
67  }
68  }
69 }
70 
71 
72 template<class Type, template<class> class PatchField, class GeoMesh>
74 {
75  const localIOdictionary dict
76  (
77  IOobject
78  (
79  this->name(),
80  this->instance(),
81  this->local(),
82  this->db(),
83  IOobject::MUST_READ,
84  IOobject::NO_WRITE,
85  false
86  ),
87  typeName
88  );
89 
90  this->close();
91 
93 }
94 
95 
96 template<class Type, template<class> class PatchField, class GeoMesh>
98 {
99  if
100  (
101  this->readOpt() == IOobject::MUST_READ
102  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
103  )
104  {
106  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
107  << " suggests that a read constructor for field " << this->name()
108  << " would be more appropriate." << endl;
109  }
110  else if
111  (
112  this->readOpt() == IOobject::READ_IF_PRESENT
113  && this->headerOk()
114  )
115  {
116  readFields();
117 
118  // Check compatibility between field and mesh
119  if (this->size() != GeoMesh::size(this->mesh()))
120  {
121  FatalIOErrorInFunction(this->readStream(typeName))
122  << " number of field elements = " << this->size()
123  << " number of mesh elements = "
124  << GeoMesh::size(this->mesh())
125  << exit(FatalIOError);
126  }
127 
128  readOldTimeIfPresent();
129 
130  return true;
131  }
132 
133  return false;
134 }
135 
136 
137 template<class Type, template<class> class PatchField, class GeoMesh>
139 {
140  // Read the old time field if present
141  typeIOobject<GeometricField<Type, PatchField, GeoMesh>> field0
142  (
143  this->name() + "_0",
144  this->time().name(),
145  this->db(),
146  IOobject::READ_IF_PRESENT,
147  IOobject::AUTO_WRITE,
148  this->registerObject()
149  );
150 
151  if (field0.headerOk())
152  {
153  if (debug)
154  {
155  InfoInFunction << "Reading old time level for field"
156  << endl << this->info() << endl;
157  }
158 
159  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
160  (
161  field0,
162  this->mesh()
163  );
164 
165  field0Ptr_->timeIndex_ = timeIndex_ - 1;
166 
167  if (!field0Ptr_->readOldTimeIfPresent())
168  {
169  field0Ptr_->oldTime();
170  }
171 
172  return true;
173  }
174 
175  return false;
176 }
177 
178 
179 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
180 
181 template<class Type, template<class> class PatchField, class GeoMesh>
183 (
184  const IOobject& io,
185  const Mesh& mesh,
186  const dimensionSet& ds,
187  const word& patchFieldType
188 )
189 :
190  Internal(io, mesh, ds, false),
191  timeIndex_(this->time().timeIndex()),
192  field0Ptr_(nullptr),
193  fieldPrevIterPtr_(nullptr),
194  boundaryField_(mesh.boundary(), *this, patchFieldType)
195 {
196  if (debug)
197  {
198  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
199  }
200 
201  readIfPresent();
202 }
203 
204 
205 template<class Type, template<class> class PatchField, class GeoMesh>
207 (
208  const IOobject& io,
209  const Mesh& mesh,
210  const dimensionSet& ds,
211  const wordList& patchFieldTypes,
212  const wordList& actualPatchTypes
213 )
214 :
215  Internal(io, mesh, ds, false),
216  timeIndex_(this->time().timeIndex()),
217  field0Ptr_(nullptr),
218  fieldPrevIterPtr_(nullptr),
219  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
220 {
221  if (debug)
222  {
223  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
224  }
225 
226  readIfPresent();
227 }
228 
229 
230 template<class Type, template<class> class PatchField, class GeoMesh>
232 (
233  const IOobject& io,
234  const Mesh& mesh,
235  const dimensioned<Type>& dt,
236  const word& patchFieldType
237 )
238 :
239  Internal(io, mesh, dt, false),
240  timeIndex_(this->time().timeIndex()),
241  field0Ptr_(nullptr),
242  fieldPrevIterPtr_(nullptr),
243  boundaryField_(mesh.boundary(), *this, patchFieldType)
244 {
245  if (debug)
246  {
247  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
248  }
249 
250  boundaryField_ == dt.value();
251 
252  readIfPresent();
253 }
254 
255 
256 template<class Type, template<class> class PatchField, class GeoMesh>
258 (
259  const IOobject& io,
260  const Mesh& mesh,
261  const dimensioned<Type>& dt,
262  const wordList& patchFieldTypes,
263  const wordList& actualPatchTypes
264 )
265 :
266  Internal(io, mesh, dt, false),
267  timeIndex_(this->time().timeIndex()),
268  field0Ptr_(nullptr),
269  fieldPrevIterPtr_(nullptr),
270  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
271 {
272  if (debug)
273  {
274  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
275  }
276 
277  boundaryField_ == dt.value();
278 
279  readIfPresent();
280 }
281 
282 
283 template<class Type, template<class> class PatchField, class GeoMesh>
285 (
286  const IOobject& io,
287  const Internal& diField,
288  const PtrList<PatchField<Type>>& ptfl
289 )
290 :
291  Internal(io, diField),
292  timeIndex_(this->time().timeIndex()),
293  field0Ptr_(nullptr),
294  fieldPrevIterPtr_(nullptr),
295  boundaryField_(this->mesh().boundary(), *this, ptfl)
296 {
297  if (debug)
298  {
300  << "Constructing from components" << endl << this->info() << endl;
301  }
302 
303  readIfPresent();
304 }
305 
306 
307 template<class Type, template<class> class PatchField, class GeoMesh>
309 (
310  const IOobject& io,
311  const Mesh& mesh,
312  const dimensionSet& ds,
313  const Field<Type>& iField,
314  const PtrList<PatchField<Type>>& ptfl
315 )
316 :
317  Internal(io, mesh, ds, iField),
318  timeIndex_(this->time().timeIndex()),
319  field0Ptr_(nullptr),
320  fieldPrevIterPtr_(nullptr),
321  boundaryField_(mesh.boundary(), *this, ptfl)
322 {
323  if (debug)
324  {
326  << "Constructing from components" << endl << this->info() << endl;
327  }
328 
329  readIfPresent();
330 }
331 
332 
333 template<class Type, template<class> class PatchField, class GeoMesh>
335 (
336  const IOobject& io,
337  const Mesh& mesh
338 )
339 :
340  Internal(io, mesh, dimless, false),
341  timeIndex_(this->time().timeIndex()),
342  field0Ptr_(nullptr),
343  fieldPrevIterPtr_(nullptr),
344  boundaryField_(mesh.boundary())
345 {
346  readFields();
347 
348  // Check compatibility between field and mesh
349 
350  if (this->size() != GeoMesh::size(this->mesh()))
351  {
352  FatalIOErrorInFunction(this->readStream(typeName))
353  << " number of field elements = " << this->size()
354  << " number of mesh elements = " << GeoMesh::size(this->mesh())
355  << exit(FatalIOError);
356  }
357 
358  readOldTimeIfPresent();
359 
360  if (debug)
361  {
363  << "Finishing read-construction of" << endl << this->info() << endl;
364  }
365 }
366 
367 
368 template<class Type, template<class> class PatchField, class GeoMesh>
370 (
371  const IOobject& io,
372  const Mesh& mesh,
373  const dictionary& dict
374 )
375 :
376  Internal(io, mesh, dimless, false),
377  timeIndex_(this->time().timeIndex()),
378  field0Ptr_(nullptr),
379  fieldPrevIterPtr_(nullptr),
380  boundaryField_(mesh.boundary())
381 {
382  readFields(dict);
383 
384  // Check compatibility between field and mesh
385 
386  if (this->size() != GeoMesh::size(this->mesh()))
387  {
389  << " number of field elements = " << this->size()
390  << " number of mesh elements = " << GeoMesh::size(this->mesh())
391  << exit(FatalIOError);
392  }
393 
394  if (debug)
395  {
397  << "Finishing dictionary-construct of "
398  << endl << this->info() << endl;
399  }
400 }
401 
402 
403 template<class Type, template<class> class PatchField, class GeoMesh>
405 (
407 )
408 :
409  Internal(gf),
410  timeIndex_(gf.timeIndex()),
411  field0Ptr_(nullptr),
412  fieldPrevIterPtr_(nullptr),
413  boundaryField_(*this, gf.boundaryField_)
414 {
415  if (debug)
416  {
418  << "Constructing as copy" << endl << this->info() << endl;
419  }
420 
421  if (gf.field0Ptr_ && notNull(gf.field0Ptr_))
422  {
424  (
425  *gf.field0Ptr_
426  );
427  }
428 
429  this->writeOpt() = IOobject::NO_WRITE;
430 }
431 
432 
433 template<class Type, template<class> class PatchField, class GeoMesh>
435 (
437 )
438 :
439  Internal(move(gf)),
440  timeIndex_(gf.timeIndex()),
441  field0Ptr_(nullptr),
442  fieldPrevIterPtr_(nullptr),
443  boundaryField_(*this, gf.boundaryField_)
444 {
445  if (debug)
446  {
448  << "Constructing by moving" << endl << this->info() << endl;
449  }
450 
451  if (gf.field0Ptr_ && notNull(gf.field0Ptr_))
452  {
453  field0Ptr_ = gf.field0Ptr_;
454  gf.field0Ptr_ = nullptr;
455  }
456 
457  this->writeOpt() = IOobject::NO_WRITE;
458 }
459 
460 
461 template<class Type, template<class> class PatchField, class GeoMesh>
463 (
465 )
466 :
467  Internal
468  (
469  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
470  tgf.isTmp()
471  ),
472  timeIndex_(tgf().timeIndex()),
473  field0Ptr_(nullptr),
474  fieldPrevIterPtr_(nullptr),
475  boundaryField_(*this, tgf().boundaryField_)
476 {
477  if (debug)
478  {
480  << "Constructing from tmp" << endl << this->info() << endl;
481  }
482 
483  this->writeOpt() = IOobject::NO_WRITE;
484 
485  tgf.clear();
486 }
487 
488 
489 template<class Type, template<class> class PatchField, class GeoMesh>
491 (
492  const IOobject& io,
494 )
495 :
496  Internal(io, gf),
497  timeIndex_(gf.timeIndex()),
498  field0Ptr_(nullptr),
499  fieldPrevIterPtr_(nullptr),
500  boundaryField_(*this, gf.boundaryField_)
501 {
502  if (debug)
503  {
505  << "Constructing as copy resetting IO params"
506  << endl << this->info() << endl;
507  }
508 
509  if (!readIfPresent() && gf.field0Ptr_ && notNull(gf.field0Ptr_))
510  {
512  (
513  io.name() + "_0",
514  *gf.field0Ptr_
515  );
516  }
517 }
518 
519 
520 template<class Type, template<class> class PatchField, class GeoMesh>
522 (
523  const IOobject& io,
525 )
526 :
528  (
529  io,
530  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
531  tgf.isTmp()
532  ),
533  timeIndex_(tgf().timeIndex()),
534  field0Ptr_(nullptr),
535  fieldPrevIterPtr_(nullptr),
536  boundaryField_(*this, tgf().boundaryField_)
537 {
538  if (debug)
539  {
541  << "Constructing from tmp resetting IO params"
542  << endl << this->info() << endl;
543  }
544 
545  tgf.clear();
546 
547  readIfPresent();
548 }
549 
550 
551 template<class Type, template<class> class PatchField, class GeoMesh>
553 (
554  const word& newName,
556 )
557 :
558  Internal(newName, gf),
559  timeIndex_(gf.timeIndex()),
560  field0Ptr_(nullptr),
561  fieldPrevIterPtr_(nullptr),
562  boundaryField_(*this, gf.boundaryField_)
563 {
564  if (debug)
565  {
567  << "Constructing as copy resetting name"
568  << endl << this->info() << endl;
569  }
570 
571  if (!readIfPresent() && gf.field0Ptr_ && notNull(gf.field0Ptr_))
572  {
574  (
575  newName + "_0",
576  *gf.field0Ptr_
577  );
578  }
579 }
580 
581 
582 template<class Type, template<class> class PatchField, class GeoMesh>
584 (
585  const word& newName,
587 )
588 :
589  Internal
590  (
591  newName,
592  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
593  tgf.isTmp()
594  ),
595  timeIndex_(tgf().timeIndex()),
596  field0Ptr_(nullptr),
597  fieldPrevIterPtr_(nullptr),
598  boundaryField_(*this, tgf().boundaryField_)
599 {
600  if (debug)
601  {
603  << "Constructing from tmp resetting name"
604  << endl << this->info() << endl;
605  }
606 
607  tgf.clear();
608 }
609 
610 
611 template<class Type, template<class> class PatchField, class GeoMesh>
613 (
614  const IOobject& io,
616  const word& patchFieldType
617 )
618 :
619  Internal(io, gf),
620  timeIndex_(gf.timeIndex()),
621  field0Ptr_(nullptr),
622  fieldPrevIterPtr_(nullptr),
623  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
624 {
625  if (debug)
626  {
628  << "Constructing as copy resetting IO params"
629  << endl << this->info() << endl;
630  }
631 
632  boundaryField_ == gf.boundaryField_;
633 
634  if (!readIfPresent() && gf.field0Ptr_ && notNull(gf.field0Ptr_))
635  {
637  (
638  io.name() + "_0",
639  *gf.field0Ptr_
640  );
641  }
642 }
643 
644 
645 template<class Type, template<class> class PatchField, class GeoMesh>
647 (
648  const IOobject& io,
650  const wordList& patchFieldTypes,
651  const wordList& actualPatchTypes
652 
653 )
654 :
655  Internal(io, gf),
656  timeIndex_(gf.timeIndex()),
657  field0Ptr_(nullptr),
658  fieldPrevIterPtr_(nullptr),
659  boundaryField_
660  (
661  this->mesh().boundary(),
662  *this,
663  patchFieldTypes,
664  actualPatchTypes
665  )
666 {
667  if (debug)
668  {
670  << "Constructing as copy resetting IO params and patch types"
671  << endl << this->info() << endl;
672  }
673 
674  boundaryField_ == gf.boundaryField_;
675 
676  if (!readIfPresent() && gf.field0Ptr_ && notNull(gf.field0Ptr_))
677  {
679  (
680  io.name() + "_0",
681  *gf.field0Ptr_
682  );
683  }
684 }
685 
686 
687 template<class Type, template<class> class PatchField, class GeoMesh>
689 (
690  const IOobject& io,
692  const wordList& patchFieldTypes,
693  const wordList& actualPatchTypes
694 )
695 :
696  Internal
697  (
698  io,
699  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
700  tgf.isTmp()
701  ),
702  timeIndex_(tgf().timeIndex()),
703  field0Ptr_(nullptr),
704  fieldPrevIterPtr_(nullptr),
705  boundaryField_
706  (
707  this->mesh().boundary(),
708  *this,
709  patchFieldTypes,
710  actualPatchTypes
711  )
712 {
713  if (debug)
714  {
716  << "Constructing from tmp resetting IO params and patch types"
717  << endl << this->info() << endl;
718  }
719 
720  boundaryField_ == tgf().boundaryField_;
721 
722  tgf.clear();
723 
724  readIfPresent();
725 }
726 
727 
728 template<class Type, template<class> class PatchField, class GeoMesh>
731 {
733  (
735  );
736 }
737 
738 
739 template<class Type, template<class> class PatchField, class GeoMesh>
742 {
744  (
746  (
747  IOobject
748  (
749  this->name(),
750  this->mesh().thisDb().time().name(),
751  this->mesh().thisDb(),
754  false
755  ),
756  *this,
758  )
759  );
760 }
761 
762 
763 template<class Type, template<class> class PatchField, class GeoMesh>
766 (
767  const word& name,
768  const Internal& diField,
769  const PtrList<PatchField<Type>>& ptfl
770 )
771 {
772  const bool cacheTmp = diField.mesh().thisDb().cacheTemporaryObject(name);
773 
775  (
777  (
778  IOobject
779  (
780  name,
781  diField.mesh().thisDb().time().name(),
782  diField.mesh().thisDb(),
785  cacheTmp
786  ),
787  diField,
788  ptfl
789  ),
790  cacheTmp
791  );
792 }
793 
794 
795 template<class Type, template<class> class PatchField, class GeoMesh>
798 (
799  const word& name,
800  const Mesh& mesh,
801  const dimensionSet& ds,
802  const word& patchFieldType
803 )
804 {
805  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
806 
808  (
810  (
811  IOobject
812  (
813  name,
814  mesh.thisDb().time().name(),
815  mesh.thisDb(),
818  cacheTmp
819  ),
820  mesh,
821  ds,
823  ),
824  cacheTmp
825  );
826 }
827 
828 
829 template<class Type, template<class> class PatchField, class GeoMesh>
832 (
833  const word& name,
834  const Mesh& mesh,
835  const dimensioned<Type>& dt,
836  const word& patchFieldType
837 )
838 {
839  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
840 
842  (
844  (
845  IOobject
846  (
847  name,
848  mesh.thisDb().time().name(),
849  mesh.thisDb(),
852  cacheTmp
853  ),
854  mesh,
855  dt,
857  ),
858  cacheTmp
859  );
860 }
861 
862 
863 
864 template<class Type, template<class> class PatchField, class GeoMesh>
867 (
868  const word& name,
869  const Mesh& mesh,
870  const dimensioned<Type>& dt,
871  const wordList& patchFieldTypes,
872  const wordList& actualPatchTypes
873 )
874 {
875  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
876 
878  (
880  (
881  IOobject
882  (
883  name,
884  mesh.thisDb().time().name(),
885  mesh.thisDb(),
888  cacheTmp
889  ),
890  mesh,
891  dt,
892  patchFieldTypes,
893  actualPatchTypes
894  ),
895  cacheTmp
896  );
897 }
898 
899 
900 template<class Type, template<class> class PatchField, class GeoMesh>
903 (
904  const word& newName,
906 )
907 {
908  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
909 
911  (
913  (
914  IOobject
915  (
916  newName,
917  tgf().instance(),
918  tgf().local(),
919  tgf().db(),
922  cacheTmp
923  ),
924  tgf
925  ),
926  cacheTmp
927  );
928 }
929 
930 
931 template<class Type, template<class> class PatchField, class GeoMesh>
934 (
935  const word& newName,
937  const word& patchFieldType
938 )
939 {
940  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
941 
943  (
945  (
946  IOobject
947  (
948  newName,
949  tgf().instance(),
950  tgf().local(),
951  tgf().db(),
954  cacheTmp
955  ),
956  tgf,
958  ),
959  cacheTmp
960  );
961 }
962 
963 
964 template<class Type, template<class> class PatchField, class GeoMesh>
967 (
968  const word& newName,
970  const wordList& patchFieldTypes,
971  const wordList& actualPatchTypes
972 )
973 {
974  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
975 
977  (
979  (
980  IOobject
981  (
982  newName,
983  tgf().instance(),
984  tgf().local(),
985  tgf().db(),
988  cacheTmp
989  ),
990  tgf,
991  patchFieldTypes,
992  actualPatchTypes
993  ),
994  cacheTmp
995  );
996 }
997 
998 
999 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
1000 
1001 template<class Type, template<class> class PatchField, class GeoMesh>
1003 {
1004  this->db().cacheTemporaryObject(*this);
1005 
1006  clearOldTimes();
1007 }
1008 
1009 
1010 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1011 
1012 template<class Type, template<class> class PatchField, class GeoMesh>
1013 typename
1016 {
1017  this->setUpToDate();
1018  storeOldTimes();
1019  return *this;
1020 }
1021 
1022 
1023 template<class Type, template<class> class PatchField, class GeoMesh>
1024 typename
1027 {
1028  this->setUpToDate();
1029  storeOldTimes();
1030  return *this;
1031 }
1032 
1033 
1034 template<class Type, template<class> class PatchField, class GeoMesh>
1035 typename
1038 {
1039  this->setUpToDate();
1040  storeOldTimes();
1041  return boundaryField_;
1042 }
1043 
1044 
1045 template<class Type, template<class> class PatchField, class GeoMesh>
1047 {
1048  return
1049  this->name().size() > 2
1050  && this->name()(this->name().size()-2, 2) == "_0";
1051 }
1052 
1053 
1054 template<class Type, template<class> class PatchField, class GeoMesh>
1056 {
1057  if (field0Ptr_ && timeIndex_ != this->time().timeIndex() && !isOldTime())
1058  {
1059  storeOldTime();
1060  }
1061 
1062  // Correct time index
1063  timeIndex_ = this->time().timeIndex();
1064 }
1065 
1066 
1067 template<class Type, template<class> class PatchField, class GeoMesh>
1069 {
1070  if (field0Ptr_)
1071  {
1072  if (notNull(field0Ptr_))
1073  {
1074  field0Ptr_->storeOldTime();
1075 
1076  if (debug)
1077  {
1079  << "Storing old time field for field" << endl
1080  << this->info() << endl;
1081  }
1082 
1083  *field0Ptr_ == *this;
1084  field0Ptr_->timeIndex_ = timeIndex_;
1085 
1086  if (field0Ptr_->field0Ptr_)
1087  {
1088  field0Ptr_->writeOpt() = this->writeOpt();
1089  }
1090  }
1091  else
1092  {
1093  // Reinstate old-time field
1094  oldTime();
1095  }
1096  }
1097 }
1098 
1099 
1100 template<class Type, template<class> class PatchField, class GeoMesh>
1102 {
1103  if (field0Ptr_)
1104  {
1105  if (isNull(field0Ptr_))
1106  {
1107  return 1;
1108  }
1109  else
1110  {
1111  return field0Ptr_->nOldTimes() + 1;
1112  }
1113  }
1114  else
1115  {
1116  return 0;
1117  }
1118 }
1119 
1120 
1121 template<class Type, template<class> class PatchField, class GeoMesh>
1124 {
1125  if (!field0Ptr_ || isNull(field0Ptr_))
1126  {
1127  // Set field0Ptr_ null to ensure the old-time field constructor
1128  // does not construct the old-old-time field
1129  field0Ptr_ = nullptr;
1130 
1132  (
1133  IOobject
1134  (
1135  this->name() + "_0",
1136  this->time().name(),
1137  this->db(),
1140  this->registerObject()
1141  ),
1142  *this
1143  );
1144  }
1145  else
1146  {
1147  storeOldTimes();
1148  }
1149 
1150  return *field0Ptr_;
1151 }
1152 
1153 
1154 template<class Type, template<class> class PatchField, class GeoMesh>
1157 {
1158  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1159  .oldTime();
1160 
1161  return *field0Ptr_;
1162 }
1163 
1164 
1165 template<class Type, template<class> class PatchField, class GeoMesh>
1168 {
1169  if (n == 0)
1170  {
1171  return *this;
1172  }
1173  else
1174  {
1175  return oldTime().oldTime(n - 1);
1176  }
1177 }
1178 
1179 
1180 template<class Type, template<class> class PatchField, class GeoMesh>
1183 {
1184  if (n == 0)
1185  {
1186  return *this;
1187  }
1188  else
1189  {
1190  return oldTime().oldTime(n - 1);
1191  }
1192 }
1193 
1194 
1195 template<class Type, template<class> class PatchField, class GeoMesh>
1197 {
1198  if (!fieldPrevIterPtr_)
1199  {
1200  if (debug)
1201  {
1203  << "Allocating previous iteration field" << endl
1204  << this->info() << endl;
1205  }
1206 
1207  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1208  (
1209  this->name() + "PrevIter",
1210  *this
1211  );
1212  }
1213  else
1214  {
1215  *fieldPrevIterPtr_ == *this;
1216  }
1217 }
1218 
1219 
1220 template<class Type, template<class> class PatchField, class GeoMesh>
1223 {
1224  if (!fieldPrevIterPtr_)
1225  {
1227  << "previous iteration field" << endl << this->info() << endl
1228  << " not stored."
1229  << " Use field.storePrevIter() at start of iteration."
1230  << abort(FatalError);
1231  }
1232 
1233  return *fieldPrevIterPtr_;
1234 }
1235 
1236 
1237 template<class Type, template<class> class PatchField, class GeoMesh>
1239 {
1240  if (field0Ptr_ && notNull(field0Ptr_))
1241  {
1242  deleteDemandDrivenData(field0Ptr_);
1243  }
1244 
1245  deleteDemandDrivenData(fieldPrevIterPtr_);
1246 }
1247 
1248 
1249 template<class Type, template<class> class PatchField, class GeoMesh>
1251 {
1252  // Check that the field is not an old-time field
1253  if (!isOldTime())
1254  {
1255  // Search for the oldest old-time field and set to nullObject
1256  nullOldestTimeFound();
1257  }
1258 }
1259 
1260 
1261 template<class Type, template<class> class PatchField, class GeoMesh>
1263 {
1264  if (field0Ptr_ && notNull(field0Ptr_))
1265  {
1266  if (field0Ptr_->field0Ptr_)
1267  {
1268  field0Ptr_->nullOldestTimeFound();
1269  }
1270  else
1271  {
1272  deleteDemandDrivenData(field0Ptr_);
1273  field0Ptr_ =
1274  NullObjectPtr<GeometricField<Type, PatchField, GeoMesh>>();
1275  }
1276  }
1277 }
1278 
1279 
1280 template<class Type, template<class> class PatchField, class GeoMesh>
1283 {
1284  this->setUpToDate();
1285  storeOldTimes();
1286  boundaryField_.evaluate();
1287 }
1288 
1289 
1290 template<class Type, template<class> class PatchField, class GeoMesh>
1292 (
1294 )
1295 {
1297 
1298  Internal::reset(gf);
1299  boundaryField_.reset(gf.boundaryField());
1300 
1301  tgf.clear();
1302 }
1303 
1304 
1305 template<class Type, template<class> class PatchField, class GeoMesh>
1307 {
1308  // Search all boundary conditions, if any are
1309  // fixed-value or mixed (Robin) do not set reference level for solution.
1310 
1311  bool needRef = true;
1312 
1313  forAll(boundaryField_, patchi)
1314  {
1315  if (boundaryField_[patchi].fixesValue())
1316  {
1317  needRef = false;
1318  break;
1319  }
1320  }
1321 
1322  reduce(needRef, andOp<bool>());
1323 
1324  return needRef;
1325 }
1326 
1327 
1328 template<class Type, template<class> class PatchField, class GeoMesh>
1330 {
1331  if (alpha < 1)
1332  {
1333  if (debug)
1334  {
1336  << "Relaxing" << endl << this->info()
1337  << " by " << alpha << endl;
1338  }
1339 
1340  operator==(prevIter() + alpha*(*this - prevIter()));
1341  }
1342 }
1343 
1344 
1345 template<class Type, template<class> class PatchField, class GeoMesh>
1346 Foam::scalar
1348 {
1349  if
1350  (
1351  this->mesh().data::template lookupOrDefault<bool>
1352  (
1353  "finalIteration",
1354  false
1355  )
1356  && this->mesh().solution().relaxField(this->name() + "Final")
1357  )
1358  {
1359  return this->mesh().solution().fieldRelaxationFactor
1360  (
1361  this->name() + "Final"
1362  );
1363  }
1364  else if (this->mesh().solution().relaxField(this->name()))
1365  {
1366  return this->mesh().solution().fieldRelaxationFactor(this->name());
1367  }
1368  else
1369  {
1370  return 1;
1371  }
1372 }
1373 
1374 
1375 template<class Type, template<class> class PatchField, class GeoMesh>
1377 {
1378  relax(relaxationFactor());
1379 }
1380 
1381 
1382 template<class Type, template<class> class PatchField, class GeoMesh>
1384 (
1386  const scalar alpha
1387 )
1388 {
1389  if (alpha < 1)
1390  {
1391  if (debug)
1392  {
1394  << "Relaxing" << endl << this->info()
1395  << " by " << alpha << endl;
1396  }
1397 
1398  operator==(*this + alpha*(tgf - *this));
1399  }
1400  else
1401  {
1402  operator==(tgf);
1403  }
1404 }
1405 
1406 
1407 template<class Type, template<class> class PatchField, class GeoMesh>
1409 (
1411 )
1412 {
1413  relax(tgf, relaxationFactor());
1414 }
1415 
1416 
1417 template<class Type, template<class> class PatchField, class GeoMesh>
1419 (
1420  bool final
1421 ) const
1422 {
1423  if (final)
1424  {
1425  return this->name() + "Final";
1426  }
1427  else
1428  {
1429  return this->name();
1430  }
1431 }
1432 
1433 
1434 template<class Type, template<class> class PatchField, class GeoMesh>
1436 (
1437  Ostream& os
1438 ) const
1439 {
1440  os << "min/max(" << this->name() << ") = "
1441  << Foam::min(*this).value() << ", "
1442  << Foam::max(*this).value()
1443  << endl;
1444 }
1445 
1446 
1447 template<class Type, template<class> class PatchField, class GeoMesh>
1449 writeData(Ostream& os) const
1450 {
1451  os << *this;
1452  return os.good();
1453 }
1454 
1455 
1456 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1457 
1458 template<class Type, template<class> class PatchField, class GeoMesh>
1461 {
1463  (
1465  (
1466  this->name() + ".T()",
1467  this->mesh(),
1468  this->dimensions()
1469  )
1470  );
1471 
1472  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1473  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1474 
1475  return result;
1476 }
1477 
1478 
1479 template<class Type, template<class> class PatchField, class GeoMesh>
1480 Foam::tmp
1481 <
1483  <
1485  PatchField,
1486  GeoMesh
1487  >
1488 >
1490 (
1491  const direction d
1492 ) const
1493 {
1495  (
1497  (
1498  this->name() + ".component(" + Foam::name(d) + ')',
1499  this->mesh(),
1500  this->dimensions()
1501  )
1502  );
1503 
1504  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1505  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1506 
1507  return Component;
1508 }
1509 
1510 
1511 template<class Type, template<class> class PatchField, class GeoMesh>
1513 (
1514  const direction d,
1515  const GeometricField
1516  <
1518  PatchField,
1519  GeoMesh
1520  >& gcf
1521 )
1522 {
1523  primitiveFieldRef().replace(d, gcf.primitiveField());
1524  boundaryFieldRef().replace(d, gcf.boundaryField());
1525 }
1526 
1527 
1528 template<class Type, template<class> class PatchField, class GeoMesh>
1530 (
1531  const direction d,
1532  const dimensioned<cmptType>& ds
1533 )
1534 {
1535  primitiveFieldRef().replace(d, ds.value());
1536  boundaryFieldRef().replace(d, ds.value());
1537 }
1538 
1539 
1540 template<class Type, template<class> class PatchField, class GeoMesh>
1542 (
1543  const dimensioned<Type>& dt
1544 )
1545 {
1546  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1547  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1548 }
1549 
1550 
1551 template<class Type, template<class> class PatchField, class GeoMesh>
1553 (
1554  const dimensioned<Type>& dt
1555 )
1556 {
1557  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1558  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1559 }
1560 
1561 
1562 template<class Type, template<class> class PatchField, class GeoMesh>
1564 (
1565  const dimensioned<Type>& minDt,
1566  const dimensioned<Type>& maxDt
1567 )
1568 {
1569  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1570  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1571  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1572  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1573 }
1574 
1575 
1576 template<class Type, template<class> class PatchField, class GeoMesh>
1578 {
1579  primitiveFieldRef().negate();
1580  boundaryFieldRef().negate();
1581 }
1582 
1583 
1584 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1585 
1586 template<class Type, template<class> class PatchField, class GeoMesh>
1588 (
1590 )
1591 {
1592  if (this == &gf)
1593  {
1595  << "attempted assignment to self"
1596  << abort(FatalError);
1597  }
1598 
1599  checkField(*this, gf, "=");
1600 
1601  // Only assign field contents not ID
1602 
1603  ref() = gf();
1604  boundaryFieldRef() = gf.boundaryField();
1605 }
1606 
1607 
1608 template<class Type, template<class> class PatchField, class GeoMesh>
1610 (
1612 )
1613 {
1614  if (this == &gf)
1615  {
1617  << "attempted assignment to self"
1618  << abort(FatalError);
1619  }
1620 
1621  checkField(*this, gf, "=");
1622 
1623  // Only assign field contents not ID
1624 
1625  ref() = move(gf());
1626  boundaryFieldRef() = move(gf.boundaryField());
1627 }
1628 
1629 
1630 template<class Type, template<class> class PatchField, class GeoMesh>
1632 (
1634 )
1635 {
1636  if (this == &(tgf()))
1637  {
1639  << "attempted assignment to self"
1640  << abort(FatalError);
1641  }
1642 
1644 
1645  checkField(*this, gf, "=");
1646 
1647  // Only assign field contents not ID
1648 
1649  this->dimensions() = gf.dimensions();
1650 
1651  if (tgf.isTmp())
1652  {
1653  // Transfer the storage from the tmp
1654  primitiveFieldRef().transfer
1655  (
1656  const_cast<Field<Type>&>(gf.primitiveField())
1657  );
1658  }
1659  else
1660  {
1662  }
1663 
1664  boundaryFieldRef() = gf.boundaryField();
1665 
1666  tgf.clear();
1667 }
1668 
1669 
1670 template<class Type, template<class> class PatchField, class GeoMesh>
1672 (
1673  const dimensioned<Type>& dt
1674 )
1675 {
1676  ref() = dt;
1677  boundaryFieldRef() = dt.value();
1678 }
1679 
1680 
1681 template<class Type, template<class> class PatchField, class GeoMesh>
1683 (
1684  const zero&
1685 )
1686 {
1687  ref() = Zero;
1688  boundaryFieldRef() = Zero;
1689 }
1690 
1691 
1692 template<class Type, template<class> class PatchField, class GeoMesh>
1694 (
1696 )
1697 {
1699 
1700  checkField(*this, gf, "==");
1701 
1702  // Only assign field contents not ID
1703 
1704  ref() = gf();
1705  boundaryFieldRef() == gf.boundaryField();
1706 
1707  tgf.clear();
1708 }
1709 
1710 
1711 template<class Type, template<class> class PatchField, class GeoMesh>
1713 (
1714  const dimensioned<Type>& dt
1715 )
1716 {
1717  ref() = dt;
1718  boundaryFieldRef() == dt.value();
1719 }
1720 
1721 
1722 template<class Type, template<class> class PatchField, class GeoMesh>
1724 (
1725  const zero&
1726 )
1727 {
1728  ref() = Zero;
1729  boundaryFieldRef() == Zero;
1730 }
1731 
1732 
1733 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1734  \
1735 template<class Type, template<class> class PatchField, class GeoMesh> \
1736 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1737 ( \
1738  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1739 ) \
1740 { \
1741  checkField(*this, gf, #op); \
1742  \
1743  ref() op gf(); \
1744  boundaryFieldRef() op gf.boundaryField(); \
1745 } \
1746  \
1747 template<class Type, template<class> class PatchField, class GeoMesh> \
1748 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1749 ( \
1750  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1751 ) \
1752 { \
1753  operator op(tgf()); \
1754  tgf.clear(); \
1755 } \
1756  \
1757 template<class Type, template<class> class PatchField, class GeoMesh> \
1758 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1759 ( \
1760  const dimensioned<TYPE>& dt \
1761 ) \
1762 { \
1763  ref() op dt; \
1764  boundaryFieldRef() op dt.value(); \
1765 }
1766 
1771 
1772 #undef COMPUTED_ASSIGNMENT
1773 
1774 
1775 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1776 
1777 template<class Type, template<class> class PatchField, class GeoMesh>
1778 Foam::Ostream& Foam::operator<<
1779 (
1780  Ostream& os,
1782 )
1783 {
1784  gf().writeData(os, "internalField");
1785  os << nl;
1786  gf.boundaryField().writeEntry("boundaryField", os);
1787 
1788  // Check state of IOstream
1789  os.check
1790  (
1791  "Ostream& operator<<(Ostream&, "
1792  "const GeometricField<Type, PatchField, GeoMesh>&)"
1793  );
1794 
1795  return (os);
1796 }
1797 
1798 
1799 template<class Type, template<class> class PatchField, class GeoMesh>
1800 Foam::Ostream& Foam::operator<<
1801 (
1802  Ostream& os,
1804 )
1805 {
1806  os << tgf();
1807  tgf.clear();
1808  return os;
1809 }
1810 
1811 
1812 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1813 
1814 #undef checkField
1815 
1816 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1817 
1818 #include "GeometricFieldFunctions.C"
1819 
1820 // ************************************************************************* //
EaEqn relax()
#define checkField(gf1, gf2, op)
#define COMPUTED_ASSIGNMENT(TYPE, op)
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:82
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.
Generic GeometricField class.
void max(const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
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.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
label timeIndex() const
Return the time index of the field.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
void min(const dimensioned< Type > &)
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Field< Type >::cmptType cmptType
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch field type.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
bool needReference() const
Does the field need a reference level for solution.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
virtual ~GeometricField()
Destructor.
scalar relaxationFactor() const
Return the field relaxation factor read from fvSolution.
bool isOldTime() const
Return whether or not this is an old-time field.
void storePrevIter() const
Store the field as the previous iteration value.
label nOldTimes() const
Return the number of old time fields stored.
word select(bool final) const
Select the final iteration parameters if `final' is true.
void storeOldTimes() const
Store the old-time fields.
void correctBoundaryConditions()
Correct boundary field.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
void nullOldestTime()
Set oldest time field pointer to nullObjectPtr.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void clearOldTimes()
Delete old time and previous iteration fields.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:486
const word & name() const
Return name.
Definition: IOobject.H:310
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
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:174
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Dimension set for the base types.
Definition: dimensionSet.H:122
Generic dimensioned Type class.
const Type & value() const
Return const reference to value.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
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:181
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
Template functions to aid in the implementation of demand driven data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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 component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void deleteDemandDrivenData(DataPtr &dataPtr)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:58
word patchFieldType(const PatchField &pf)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:52
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
uint8_t direction
Definition: direction.H:45
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
dictionary dict
trDeltaTf ref()
conserve primitiveFieldRef()+