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-2022 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 
92  readFields(dict);
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().timeName(),
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_)
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_)
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_)
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 :
527  Internal
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_)
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_)
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_)
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().timeName(),
751  this->mesh().thisDb(),
752  IOobject::NO_READ,
753  IOobject::NO_WRITE,
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().timeName(),
782  diField.mesh().thisDb(),
783  IOobject::NO_READ,
784  IOobject::NO_WRITE,
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().timeName(),
815  mesh.thisDb(),
816  IOobject::NO_READ,
817  IOobject::NO_WRITE,
818  cacheTmp
819  ),
820  mesh,
821  ds,
822  patchFieldType
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().timeName(),
849  mesh.thisDb(),
850  IOobject::NO_READ,
851  IOobject::NO_WRITE,
852  cacheTmp
853  ),
854  mesh,
855  dt,
856  patchFieldType
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().timeName(),
885  mesh.thisDb(),
886  IOobject::NO_READ,
887  IOobject::NO_WRITE,
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(),
920  IOobject::NO_READ,
921  IOobject::NO_WRITE,
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(),
952  IOobject::NO_READ,
953  IOobject::NO_WRITE,
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(),
986  IOobject::NO_READ,
987  IOobject::NO_WRITE,
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  field0Ptr_->storeOldTime();
1073 
1074  if (debug)
1075  {
1077  << "Storing old time field for field" << endl
1078  << this->info() << endl;
1079  }
1080 
1081  *field0Ptr_ == *this;
1082  field0Ptr_->timeIndex_ = timeIndex_;
1083 
1084  if (field0Ptr_->field0Ptr_)
1085  {
1086  field0Ptr_->writeOpt() = this->writeOpt();
1087  }
1088  }
1089 }
1090 
1091 
1092 template<class Type, template<class> class PatchField, class GeoMesh>
1094 {
1095  if (field0Ptr_)
1096  {
1097  return field0Ptr_->nOldTimes() + 1;
1098  }
1099  else
1100  {
1101  return 0;
1102  }
1103 }
1104 
1105 
1106 template<class Type, template<class> class PatchField, class GeoMesh>
1109 {
1110  if (!field0Ptr_)
1111  {
1113  (
1114  IOobject
1115  (
1116  this->name() + "_0",
1117  this->time().timeName(),
1118  this->db(),
1119  IOobject::NO_READ,
1120  IOobject::NO_WRITE,
1121  this->registerObject()
1122  ),
1123  *this
1124  );
1125  }
1126  else
1127  {
1128  storeOldTimes();
1129  }
1130 
1131  return *field0Ptr_;
1132 }
1133 
1134 
1135 template<class Type, template<class> class PatchField, class GeoMesh>
1138 {
1139  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1140  .oldTime();
1141 
1142  return *field0Ptr_;
1143 }
1144 
1145 
1146 template<class Type, template<class> class PatchField, class GeoMesh>
1149 {
1150  if (n == 0)
1151  {
1152  return *this;
1153  }
1154  else
1155  {
1156  return oldTime().oldTime(n - 1);
1157  }
1158 }
1159 
1160 
1161 template<class Type, template<class> class PatchField, class GeoMesh>
1164 {
1165  if (n == 0)
1166  {
1167  return *this;
1168  }
1169  else
1170  {
1171  return oldTime().oldTime(n - 1);
1172  }
1173 }
1174 
1175 
1176 template<class Type, template<class> class PatchField, class GeoMesh>
1178 {
1179  if (!fieldPrevIterPtr_)
1180  {
1181  if (debug)
1182  {
1184  << "Allocating previous iteration field" << endl
1185  << this->info() << endl;
1186  }
1187 
1188  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1189  (
1190  this->name() + "PrevIter",
1191  *this
1192  );
1193  }
1194  else
1195  {
1196  *fieldPrevIterPtr_ == *this;
1197  }
1198 }
1199 
1200 
1201 template<class Type, template<class> class PatchField, class GeoMesh>
1204 {
1205  if (!fieldPrevIterPtr_)
1206  {
1208  << "previous iteration field" << endl << this->info() << endl
1209  << " not stored."
1210  << " Use field.storePrevIter() at start of iteration."
1211  << abort(FatalError);
1212  }
1213 
1214  return *fieldPrevIterPtr_;
1215 }
1216 
1217 
1218 template<class Type, template<class> class PatchField, class GeoMesh>
1220 {
1221  deleteDemandDrivenData(field0Ptr_);
1222  deleteDemandDrivenData(fieldPrevIterPtr_);
1223 }
1224 
1225 
1226 template<class Type, template<class> class PatchField, class GeoMesh>
1229 {
1230  this->setUpToDate();
1231  storeOldTimes();
1232  boundaryField_.evaluate();
1233 }
1234 
1235 
1236 template<class Type, template<class> class PatchField, class GeoMesh>
1240 )
1241 {
1243 
1244  Internal::reset(gf);
1245  boundaryField_.reset(*this, gf.boundaryField());
1246 
1247  tgf.clear();
1248 }
1249 
1250 
1251 template<class Type, template<class> class PatchField, class GeoMesh>
1253 {
1254  // Search all boundary conditions, if any are
1255  // fixed-value or mixed (Robin) do not set reference level for solution.
1256 
1257  bool needRef = true;
1258 
1259  forAll(boundaryField_, patchi)
1260  {
1261  if (boundaryField_[patchi].fixesValue())
1262  {
1263  needRef = false;
1264  break;
1265  }
1266  }
1267 
1268  reduce(needRef, andOp<bool>());
1269 
1270  return needRef;
1271 }
1272 
1273 
1274 template<class Type, template<class> class PatchField, class GeoMesh>
1276 {
1277  if (debug)
1278  {
1280  << "Relaxing" << endl << this->info() << " by " << alpha << endl;
1281  }
1282 
1283  operator==(prevIter() + alpha*(*this - prevIter()));
1284 }
1285 
1286 
1287 template<class Type, template<class> class PatchField, class GeoMesh>
1289 {
1290  if
1291  (
1292  this->mesh().data::template lookupOrDefault<bool>
1293  (
1294  "finalIteration",
1295  false
1296  )
1297  && this->mesh().solution().relaxField(this->name() + "Final")
1298  )
1299  {
1300  relax
1301  (
1302  this->mesh().solution().fieldRelaxationFactor
1303  (
1304  this->name() + "Final"
1305  )
1306  );
1307  }
1308  else if (this->mesh().solution().relaxField(this->name()))
1309  {
1310  relax(this->mesh().solution().fieldRelaxationFactor(this->name()));
1311  }
1312 }
1313 
1314 
1315 template<class Type, template<class> class PatchField, class GeoMesh>
1318  bool final
1319 ) const
1320 {
1321  if (final)
1322  {
1323  return this->name() + "Final";
1324  }
1325  else
1326  {
1327  return this->name();
1328  }
1329 }
1330 
1331 
1332 template<class Type, template<class> class PatchField, class GeoMesh>
1335  Ostream& os
1336 ) const
1337 {
1338  os << "min/max(" << this->name() << ") = "
1339  << Foam::min(*this).value() << ", "
1340  << Foam::max(*this).value()
1341  << endl;
1342 }
1343 
1344 
1345 template<class Type, template<class> class PatchField, class GeoMesh>
1347 writeData(Ostream& os) const
1348 {
1349  os << *this;
1350  return os.good();
1351 }
1352 
1353 
1354 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1355 
1356 template<class Type, template<class> class PatchField, class GeoMesh>
1359 {
1361  (
1363  (
1364  this->name() + ".T()",
1365  this->mesh(),
1366  this->dimensions()
1367  )
1368  );
1369 
1370  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1371  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1372 
1373  return result;
1374 }
1375 
1376 
1377 template<class Type, template<class> class PatchField, class GeoMesh>
1378 Foam::tmp
1379 <
1381  <
1382  typename Foam::GeometricField<Type, PatchField, GeoMesh>::cmptType,
1383  PatchField,
1384  GeoMesh
1385  >
1386 >
1388 (
1389  const direction d
1390 ) const
1391 {
1393  (
1395  (
1396  this->name() + ".component(" + Foam::name(d) + ')',
1397  this->mesh(),
1398  this->dimensions()
1399  )
1400  );
1401 
1402  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1403  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1404 
1405  return Component;
1406 }
1407 
1408 
1409 template<class Type, template<class> class PatchField, class GeoMesh>
1411 (
1412  const direction d,
1413  const GeometricField
1414  <
1416  PatchField,
1417  GeoMesh
1418  >& gcf
1419 )
1420 {
1421  primitiveFieldRef().replace(d, gcf.primitiveField());
1422  boundaryFieldRef().replace(d, gcf.boundaryField());
1423 }
1424 
1425 
1426 template<class Type, template<class> class PatchField, class GeoMesh>
1429  const direction d,
1430  const dimensioned<cmptType>& ds
1431 )
1432 {
1433  primitiveFieldRef().replace(d, ds.value());
1434  boundaryFieldRef().replace(d, ds.value());
1435 }
1436 
1437 
1438 template<class Type, template<class> class PatchField, class GeoMesh>
1441  const dimensioned<Type>& dt
1442 )
1443 {
1444  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1445  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1446 }
1447 
1448 
1449 template<class Type, template<class> class PatchField, class GeoMesh>
1452  const dimensioned<Type>& dt
1453 )
1454 {
1455  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1456  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1457 }
1458 
1459 
1460 template<class Type, template<class> class PatchField, class GeoMesh>
1463  const dimensioned<Type>& minDt,
1464  const dimensioned<Type>& maxDt
1465 )
1466 {
1467  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1468  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1469  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1470  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1471 }
1472 
1473 
1474 template<class Type, template<class> class PatchField, class GeoMesh>
1476 {
1477  primitiveFieldRef().negate();
1478  boundaryFieldRef().negate();
1479 }
1480 
1481 
1482 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1483 
1484 template<class Type, template<class> class PatchField, class GeoMesh>
1485 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1488 )
1489 {
1490  if (this == &gf)
1491  {
1493  << "attempted assignment to self"
1494  << abort(FatalError);
1495  }
1496 
1497  checkField(*this, gf, "=");
1498 
1499  // Only assign field contents not ID
1500 
1501  ref() = gf();
1502  boundaryFieldRef() = gf.boundaryField();
1503 }
1504 
1505 
1506 template<class Type, template<class> class PatchField, class GeoMesh>
1507 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1510 )
1511 {
1512  if (this == &gf)
1513  {
1515  << "attempted assignment to self"
1516  << abort(FatalError);
1517  }
1518 
1519  checkField(*this, gf, "=");
1520 
1521  // Only assign field contents not ID
1522 
1523  ref() = move(gf());
1524  boundaryFieldRef() = move(gf.boundaryField());
1525 }
1526 
1527 
1528 template<class Type, template<class> class PatchField, class GeoMesh>
1529 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1532 )
1533 {
1534  if (this == &(tgf()))
1535  {
1537  << "attempted assignment to self"
1538  << abort(FatalError);
1539  }
1540 
1542 
1543  checkField(*this, gf, "=");
1544 
1545  // Only assign field contents not ID
1546 
1547  this->dimensions() = gf.dimensions();
1548 
1549  if (tgf.isTmp())
1550  {
1551  // Transfer the storage from the tmp
1552  primitiveFieldRef().transfer
1553  (
1554  const_cast<Field<Type>&>(gf.primitiveField())
1555  );
1556  }
1557  else
1558  {
1560  }
1561 
1562  boundaryFieldRef() = gf.boundaryField();
1563 
1564  tgf.clear();
1565 }
1566 
1567 
1568 template<class Type, template<class> class PatchField, class GeoMesh>
1569 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1571  const dimensioned<Type>& dt
1572 )
1573 {
1574  ref() = dt;
1575  boundaryFieldRef() = dt.value();
1576 }
1577 
1578 
1579 template<class Type, template<class> class PatchField, class GeoMesh>
1580 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1582  const zero&
1583 )
1584 {
1585  ref() = Zero;
1586  boundaryFieldRef() = Zero;
1587 }
1588 
1589 
1590 template<class Type, template<class> class PatchField, class GeoMesh>
1591 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1594 )
1595 {
1597 
1598  checkField(*this, gf, "==");
1599 
1600  // Only assign field contents not ID
1601 
1602  ref() = gf();
1603  boundaryFieldRef() == gf.boundaryField();
1604 
1605  tgf.clear();
1606 }
1607 
1608 
1609 template<class Type, template<class> class PatchField, class GeoMesh>
1610 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1612  const dimensioned<Type>& dt
1613 )
1614 {
1615  ref() = dt;
1616  boundaryFieldRef() == dt.value();
1617 }
1618 
1619 
1620 template<class Type, template<class> class PatchField, class GeoMesh>
1621 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1623  const zero&
1624 )
1625 {
1626  ref() = Zero;
1627  boundaryFieldRef() == Zero;
1628 }
1629 
1630 
1631 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1632  \
1633 template<class Type, template<class> class PatchField, class GeoMesh> \
1634 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1635 ( \
1636  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1637 ) \
1638 { \
1639  checkField(*this, gf, #op); \
1640  \
1641  ref() op gf(); \
1642  boundaryFieldRef() op gf.boundaryField(); \
1643 } \
1644  \
1645 template<class Type, template<class> class PatchField, class GeoMesh> \
1646 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1647 ( \
1648  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1649 ) \
1650 { \
1651  operator op(tgf()); \
1652  tgf.clear(); \
1653 } \
1654  \
1655 template<class Type, template<class> class PatchField, class GeoMesh> \
1656 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1657 ( \
1658  const dimensioned<TYPE>& dt \
1659 ) \
1660 { \
1661  ref() op dt; \
1662  boundaryFieldRef() op dt.value(); \
1663 }
1664 
1669 
1670 #undef COMPUTED_ASSIGNMENT
1671 
1672 
1673 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1674 
1675 template<class Type, template<class> class PatchField, class GeoMesh>
1676 Foam::Ostream& Foam::operator<<
1677 (
1678  Ostream& os,
1680 )
1681 {
1682  gf().writeData(os, "internalField");
1683  os << nl;
1684  gf.boundaryField().writeEntry("boundaryField", os);
1685 
1686  // Check state of IOstream
1687  os.check
1688  (
1689  "Ostream& operator<<(Ostream&, "
1690  "const GeometricField<Type, PatchField, GeoMesh>&)"
1691  );
1692 
1693  return (os);
1694 }
1695 
1696 
1697 template<class Type, template<class> class PatchField, class GeoMesh>
1698 Foam::Ostream& Foam::operator<<
1699 (
1700  Ostream& os,
1702 )
1703 {
1704  os << tgf();
1705  tgf.clear();
1706  return os;
1707 }
1708 
1709 
1710 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1711 
1712 #undef checkField
1713 
1714 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1715 
1716 #include "GeometricFieldFunctions.C"
1717 
1718 // ************************************************************************* //
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
#define COMPUTED_ASSIGNMENT(TYPE, op)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
#define checkField(gf1, gf2, op)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
Definition: IOobject.H:315
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
void clearOldTimes()
Delete old time and previous iteration fields.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Field< Type >::cmptType cmptType
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
uint8_t direction
Definition: direction.H:45
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
const dimensionSet dimless
Generic dimensioned Type class.
label nOldTimes() const
Return the number of old time fields stored.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
fvMesh & mesh
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
conserve primitiveFieldRef()+
bool isOldTime() const
Return whether or not this is an old-time field.
bool needReference() const
Does the field need a reference level for solution.
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
Definition: dimensionSet.H:121
void storeOldTimes() const
Store the old-time fields.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
Pre-declare SubField and related Field type.
Definition: Field.H:56
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
void storeOldTime() const
Store the old-time field.
const Type & value() const
Return const reference to value.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
word timeName
Definition: getTimeIndex.H:3
void min(const dimensioned< Type > &)
static const zero Zero
Definition: zero.H:97
faceListList boundary(nPatches)
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch field type.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual ~GeometricField()
Destructor.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Generic GeometricBoundaryField class.
Definition: fvMesh.H:80
label timeIndex() const
Return the time index of the field.
UEqn relax()
Internal & ref()
Return a reference to the dimensioned internal field.
rho oldTime()
Template functions to aid in the implementation of demand driven data.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void max(const dimensioned< Type > &)
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void correctBoundaryConditions()
Correct boundary field.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
label n
void storePrevIter() const
Store the field as the previous iteration value.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
rDeltaT ref()
void deleteDemandDrivenData(DataPtr &dataPtr)
word patchFieldType(const PatchField &pf)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
void relax()
Relax field (for steady-state solution).
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.