GeometricField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "GeometricField.H"
30 #include "Time.H"
32 #include "dictionary.H"
33 #include "localIOdictionary.H"
34 #include "data.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 #define checkField(gf1, gf2, op) \
39 if ((gf1).mesh() != (gf2).mesh()) \
40 { \
41  FatalErrorInFunction \
42  << "different mesh for fields " \
43  << (gf1).name() << " and " << (gf2).name() \
44  << " during operation " << op \
45  << abort(FatalError); \
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50 
51 template<class Type, template<class> class PatchField, class GeoMesh>
53 (
54  const dictionary& dict
55 )
56 {
57  Internal::readField(dict, "internalField");
58 
59  boundaryField_.readField(*this, dict.subDict("boundaryField"));
60 
61  Type refLevel;
62 
63  if (dict.readIfPresent("referenceLevel", refLevel))
64  {
65  Field<Type>::operator+=(refLevel);
66 
67  forAll(boundaryField_, patchi)
68  {
69  boundaryField_[patchi] == boundaryField_[patchi] + refLevel;
70  }
71  }
72 }
73 
74 
75 template<class Type, template<class> class PatchField, class GeoMesh>
77 {
78  const localIOdictionary dict
79  (
80  IOobject
81  (
82  this->name(),
83  this->instance(),
84  this->local(),
85  this->db(),
86  IOobject::MUST_READ,
87  IOobject::NO_WRITE,
88  false
89  ),
90  typeName
91  );
92 
93  this->close();
94 
96 }
97 
98 
99 template<class Type, template<class> class PatchField, class GeoMesh>
101 {
102  if (this->isReadRequired())
103  {
105  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
106  << " suggests that a read constructor for field " << this->name()
107  << " would be more appropriate." << endl;
108  }
109  else if
110  (
111  this->isReadOptional()
112  && this->template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
113  (
114  true
115  )
116  )
117  {
118  readFields();
119 
120  // Check compatibility between field and mesh
121  if (this->size() != GeoMesh::size(this->mesh()))
122  {
123  FatalIOErrorInFunction(this->readStream(typeName))
124  << " number of field elements = " << this->size()
125  << " number of mesh elements = "
126  << GeoMesh::size(this->mesh())
127  << exit(FatalIOError);
128  }
129 
130  readOldTimeIfPresent();
131 
132  return true;
133  }
134 
135  return false;
136 }
137 
138 
139 template<class Type, template<class> class PatchField, class GeoMesh>
141 {
142  // Read the old time field if present
143  IOobject field0
144  (
145  this->name() + "_0",
146  this->time().timeName(),
147  this->db(),
148  IOobject::READ_IF_PRESENT,
149  IOobject::AUTO_WRITE,
150  this->registerObject()
151  );
152 
153  if
154  (
155  field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
156  (
157  true
158  )
159  )
160  {
162  << "Reading old time level for field" << nl << this->info() << endl;
163 
164  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
165  (
166  field0,
167  this->mesh()
168  );
169 
170  // Ensure the old time field oriented flag is set to the parent's state
171  // Note: required for backwards compatibility in case of restarting from
172  // an old run where the oriented state may not have been set
173  field0Ptr_->oriented() = this->oriented();
174 
175  field0Ptr_->timeIndex_ = timeIndex_ - 1;
176 
177  if (!field0Ptr_->readOldTimeIfPresent())
178  {
179  field0Ptr_->oldTime();
180  }
181 
182  return true;
183  }
184 
185  return false;
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
190 
191 template<class Type, template<class> class PatchField, class GeoMesh>
193 (
194  const IOobject& io,
195  const Mesh& mesh,
196  const dimensionSet& ds,
197  const word& patchFieldType
198 )
199 :
200  Internal(io, mesh, ds, false),
201  timeIndex_(this->time().timeIndex()),
202  field0Ptr_(nullptr),
203  fieldPrevIterPtr_(nullptr),
204  boundaryField_(mesh.boundary(), *this, patchFieldType)
205 {
207  << "Creating temporary" << nl << this->info() << endl;
209  readIfPresent();
210 }
211 
212 
213 template<class Type, template<class> class PatchField, class GeoMesh>
215 (
216  const IOobject& io,
217  const Mesh& mesh,
218  const dimensionSet& ds,
219  const wordList& patchFieldTypes,
220  const wordList& actualPatchTypes
221 )
222 :
223  Internal(io, mesh, ds, false),
224  timeIndex_(this->time().timeIndex()),
225  field0Ptr_(nullptr),
226  fieldPrevIterPtr_(nullptr),
227  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
228 {
230  << "Creating temporary" << nl << this->info() << endl;
232  readIfPresent();
233 }
234 
235 
236 template<class Type, template<class> class PatchField, class GeoMesh>
238 (
239  const IOobject& io,
240  const Mesh& mesh,
241  const dimensioned<Type>& dt,
242  const word& patchFieldType
243 )
244 :
245  Internal(io, mesh, dt, false),
246  timeIndex_(this->time().timeIndex()),
247  field0Ptr_(nullptr),
248  fieldPrevIterPtr_(nullptr),
249  boundaryField_(mesh.boundary(), *this, patchFieldType)
250 {
252  << "Creating temporary" << nl << this->info() << endl;
253 
254  boundaryField_ == dt.value();
256  readIfPresent();
257 }
258 
259 
260 template<class Type, template<class> class PatchField, class GeoMesh>
262 (
263  const IOobject& io,
264  const Mesh& mesh,
265  const dimensioned<Type>& dt,
266  const wordList& patchFieldTypes,
267  const wordList& actualPatchTypes
268 )
269 :
270  Internal(io, mesh, dt, false),
271  timeIndex_(this->time().timeIndex()),
272  field0Ptr_(nullptr),
273  fieldPrevIterPtr_(nullptr),
274  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
275 {
277  << "Creating temporary" << nl << this->info() << endl;
278 
279  boundaryField_ == dt.value();
281  readIfPresent();
282 }
283 
284 
285 template<class Type, template<class> class PatchField, class GeoMesh>
287 (
288  const IOobject& io,
289  const Internal& diField,
290  const PtrList<PatchField<Type>>& ptfl
291 )
292 :
293  Internal(io, diField),
294  timeIndex_(this->time().timeIndex()),
295  field0Ptr_(nullptr),
296  fieldPrevIterPtr_(nullptr),
297  boundaryField_(this->mesh().boundary(), *this, ptfl)
298 {
300  << "Copy construct from components" << nl << this->info() << endl;
302  readIfPresent();
303 }
304 
305 
306 template<class Type, template<class> class PatchField, class GeoMesh>
308 (
309  const IOobject& io,
310  Internal&& diField,
311  const PtrList<PatchField<Type>>& ptfl
312 )
313 :
314  Internal(io, std::move(diField)),
315  timeIndex_(this->time().timeIndex()),
316  field0Ptr_(nullptr),
317  fieldPrevIterPtr_(nullptr),
318  boundaryField_(this->mesh().boundary(), *this, ptfl)
319 {
321  << "Move construct from components" << nl << this->info() << endl;
323  readIfPresent();
324 }
325 
326 
327 template<class Type, template<class> class PatchField, class GeoMesh>
329 (
330  const IOobject& io,
331  const tmp<Internal>& tfield,
332  const PtrList<PatchField<Type>>& ptfl
333 )
334 :
335  Internal(io, tfield),
336  timeIndex_(this->time().timeIndex()),
337  field0Ptr_(nullptr),
338  fieldPrevIterPtr_(nullptr),
339  boundaryField_(this->mesh().boundary(), *this, ptfl)
340 {
342  << "Construct from tmp internalField" << nl << this->info() << endl;
344  readIfPresent();
345 }
346 
347 
348 template<class Type, template<class> class PatchField, class GeoMesh>
350 (
351  const Internal& diField,
352  const PtrList<PatchField<Type>>& ptfl
353 )
354 :
355  Internal(diField),
356  timeIndex_(this->time().timeIndex()),
357  field0Ptr_(nullptr),
358  fieldPrevIterPtr_(nullptr),
359  boundaryField_(this->mesh().boundary(), *this, ptfl)
360 {
362  << "Copy construct from components" << nl << this->info() << endl;
364  readIfPresent();
365 }
366 
367 
368 template<class Type, template<class> class PatchField, class GeoMesh>
370 (
371  Internal&& diField,
372  const PtrList<PatchField<Type>>& ptfl
373 )
374 :
375  Internal(std::move(diField)),
376  timeIndex_(this->time().timeIndex()),
377  field0Ptr_(nullptr),
378  fieldPrevIterPtr_(nullptr),
379  boundaryField_(this->mesh().boundary(), *this, ptfl)
380 {
382  << "Move construct from components" << nl << this->info() << endl;
384  readIfPresent();
385 }
386 
387 
388 template<class Type, template<class> class PatchField, class GeoMesh>
390 (
391  const IOobject& io,
392  const Mesh& mesh,
393  const dimensionSet& ds,
394  const Field<Type>& iField,
395  const word& patchFieldType
396 )
397 :
398  Internal(io, mesh, ds, iField),
399  timeIndex_(this->time().timeIndex()),
400  field0Ptr_(nullptr),
401  fieldPrevIterPtr_(nullptr),
402  boundaryField_(mesh.boundary(), *this, patchFieldType)
403 {
405  << "Copy construct from internal field" << nl << this->info() << endl;
407  readIfPresent();
408 }
409 
410 
411 template<class Type, template<class> class PatchField, class GeoMesh>
413 (
414  const IOobject& io,
415  const Mesh& mesh,
416  const dimensionSet& ds,
417  Field<Type>&& iField,
418  const word& patchFieldType
419 )
420 :
421  Internal(io, mesh, ds, std::move(iField)),
422  timeIndex_(this->time().timeIndex()),
423  field0Ptr_(nullptr),
424  fieldPrevIterPtr_(nullptr),
425  boundaryField_(mesh.boundary(), *this, patchFieldType)
426 {
428  << "Move construct from internal field" << nl << this->info() << endl;
430  readIfPresent();
431 }
432 
433 
434 template<class Type, template<class> class PatchField, class GeoMesh>
436 (
437  const IOobject& io,
438  const Mesh& mesh,
439  const dimensionSet& ds,
440  const Field<Type>& iField,
441  const PtrList<PatchField<Type>>& ptfl
442 )
443 :
444  Internal(io, mesh, ds, iField),
445  timeIndex_(this->time().timeIndex()),
446  field0Ptr_(nullptr),
447  fieldPrevIterPtr_(nullptr),
448  boundaryField_(mesh.boundary(), *this, ptfl)
449 {
451  << "Copy construct from components" << nl << this->info() << endl;
453  readIfPresent();
454 }
455 
456 
457 template<class Type, template<class> class PatchField, class GeoMesh>
459 (
460  const IOobject& io,
461  const Mesh& mesh,
462  const dimensionSet& ds,
463  Field<Type>&& iField,
464  const PtrList<PatchField<Type>>& ptfl
465 )
466 :
467  Internal(io, mesh, ds, std::move(iField)),
468  timeIndex_(this->time().timeIndex()),
469  field0Ptr_(nullptr),
470  fieldPrevIterPtr_(nullptr),
471  boundaryField_(mesh.boundary(), *this, ptfl)
472 {
474  << "Move construct from components" << nl << this->info() << endl;
476  readIfPresent();
477 }
478 
479 
480 template<class Type, template<class> class PatchField, class GeoMesh>
482 (
483  const IOobject& io,
484  const Mesh& mesh,
485  const dimensionSet& ds,
486  const tmp<Field<Type>>& tfield,
487  const PtrList<PatchField<Type>>& ptfl
488 )
489 :
490  Internal(io, mesh, ds, tfield),
491  timeIndex_(this->time().timeIndex()),
492  field0Ptr_(nullptr),
493  fieldPrevIterPtr_(nullptr),
494  boundaryField_(mesh.boundary(), *this, ptfl)
495 {
497  << "Construct from tmp internalField" << nl << this->info() << endl;
499  readIfPresent();
500 }
501 
502 
503 template<class Type, template<class> class PatchField, class GeoMesh>
505 (
506  const IOobject& io,
507  const Mesh& mesh,
508  const bool readOldTime
509 )
510 :
511  Internal(io, mesh, dimless, false),
512  timeIndex_(this->time().timeIndex()),
513  field0Ptr_(nullptr),
514  fieldPrevIterPtr_(nullptr),
515  boundaryField_(mesh.boundary())
516 {
517  readFields();
518 
519  // Check compatibility between field and mesh
520 
521  if (this->size() != GeoMesh::size(this->mesh()))
522  {
523  FatalIOErrorInFunction(this->readStream(typeName))
524  << " number of field elements = " << this->size()
525  << " number of mesh elements = " << GeoMesh::size(this->mesh())
526  << exit(FatalIOError);
527  }
528 
529  if (readOldTime)
530  {
531  readOldTimeIfPresent();
532  }
533 
535  << "Finishing read-construction" << nl << this->info() << endl;
536 }
537 
538 
539 template<class Type, template<class> class PatchField, class GeoMesh>
541 (
542  const IOobject& io,
543  const Mesh& mesh,
544  const dictionary& dict
545 )
546 :
547  Internal(io, mesh, dimless, false),
548  timeIndex_(this->time().timeIndex()),
549  field0Ptr_(nullptr),
550  fieldPrevIterPtr_(nullptr),
551  boundaryField_(mesh.boundary())
552 {
553  readFields(dict);
554 
555  // Check compatibility between field and mesh
556 
557  if (this->size() != GeoMesh::size(this->mesh()))
558  {
560  << " number of field elements = " << this->size()
561  << " number of mesh elements = " << GeoMesh::size(this->mesh())
562  << exit(FatalIOError);
563  }
564 
566  << "Finishing dictionary-construct" << nl << this->info() << endl;
567 }
568 
569 
570 template<class Type, template<class> class PatchField, class GeoMesh>
572 (
574 )
575 :
576  Internal(gf),
577  timeIndex_(gf.timeIndex()),
578  field0Ptr_(nullptr),
579  fieldPrevIterPtr_(nullptr),
580  boundaryField_(*this, gf.boundaryField_)
581 {
583  << "Copy construct" << nl << this->info() << endl;
584 
585  if (gf.field0Ptr_)
586  {
588  (
589  *gf.field0Ptr_
590  );
591  }
593  this->writeOpt(IOobject::NO_WRITE);
594 }
595 
596 
597 template<class Type, template<class> class PatchField, class GeoMesh>
599 (
601 )
602 :
603  Internal(tgf.constCast(), tgf.movable()),
604  timeIndex_(tgf().timeIndex()),
605  field0Ptr_(nullptr),
606  fieldPrevIterPtr_(nullptr),
607  boundaryField_(*this, tgf().boundaryField_)
608 {
610  << "Constructing from tmp" << nl << this->info() << endl;
611 
612  this->writeOpt(IOobject::NO_WRITE);
614  tgf.clear();
615 }
616 
617 
618 template<class Type, template<class> class PatchField, class GeoMesh>
620 (
621  const IOobject& io,
623 )
624 :
625  Internal(io, gf),
626  timeIndex_(gf.timeIndex()),
627  field0Ptr_(nullptr),
628  fieldPrevIterPtr_(nullptr),
629  boundaryField_(*this, gf.boundaryField_)
630 {
632  << "Copy construct, resetting IO params" << nl
633  << this->info() << endl;
634 
635  if (!readIfPresent() && gf.field0Ptr_)
636  {
638  (
639  io.name() + "_0",
640  *gf.field0Ptr_
641  );
642  }
643 }
644 
645 
646 template<class Type, template<class> class PatchField, class GeoMesh>
648 (
649  const IOobject& io,
650  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
651 )
652 :
653  Internal(io, tgf.constCast(), tgf.movable()),
654  timeIndex_(tgf().timeIndex()),
655  field0Ptr_(nullptr),
656  fieldPrevIterPtr_(nullptr),
657  boundaryField_(*this, tgf().boundaryField_)
658 {
660  << "Constructing from tmp resetting IO params" << nl
661  << this->info() << endl;
662 
663  tgf.clear();
665  readIfPresent();
666 }
667 
668 
669 template<class Type, template<class> class PatchField, class GeoMesh>
671 (
672  const word& newName,
674 )
675 :
676  Internal(newName, gf),
677  timeIndex_(gf.timeIndex()),
678  field0Ptr_(nullptr),
679  fieldPrevIterPtr_(nullptr),
680  boundaryField_(*this, gf.boundaryField_)
681 {
683  << "Copy construct, resetting name" << nl
684  << this->info() << endl;
685 
686  if (!readIfPresent() && gf.field0Ptr_)
687  {
689  (
690  newName + "_0",
691  *gf.field0Ptr_
692  );
693  }
694 }
695 
696 
697 template<class Type, template<class> class PatchField, class GeoMesh>
699 (
700  const word& newName,
701  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
702 )
703 :
704  Internal(newName, tgf.constCast(), tgf.movable()),
705  timeIndex_(tgf().timeIndex()),
706  field0Ptr_(nullptr),
707  fieldPrevIterPtr_(nullptr),
708  boundaryField_(*this, tgf().boundaryField_)
709 {
711  << "Constructing from tmp resetting name" << nl
712  << this->info() << endl;
714  tgf.clear();
715 }
716 
717 
718 template<class Type, template<class> class PatchField, class GeoMesh>
720 (
721  const IOobject& io,
723  const word& patchFieldType
724 )
725 :
726  Internal(io, gf),
727  timeIndex_(gf.timeIndex()),
728  field0Ptr_(nullptr),
729  fieldPrevIterPtr_(nullptr),
730  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
731 {
733  << "Copy construct, resetting IO params" << nl
734  << this->info() << endl;
735 
736  boundaryField_ == gf.boundaryField_;
737 
738  if (!readIfPresent() && gf.field0Ptr_)
739  {
741  (
742  io.name() + "_0",
743  *gf.field0Ptr_
744  );
745  }
746 }
747 
748 
749 template<class Type, template<class> class PatchField, class GeoMesh>
751 (
752  const IOobject& io,
753  const GeometricField<Type, PatchField, GeoMesh>& gf,
754  const wordList& patchFieldTypes,
755  const wordList& actualPatchTypes
756 )
757 :
758  Internal(io, gf),
759  timeIndex_(gf.timeIndex()),
760  field0Ptr_(nullptr),
761  fieldPrevIterPtr_(nullptr),
762  boundaryField_
763  (
764  this->mesh().boundary(),
765  *this,
766  patchFieldTypes,
767  actualPatchTypes
768  )
769 {
771  << "Copy construct, resetting IO params and patch types" << nl
772  << this->info() << endl;
773 
774  boundaryField_ == gf.boundaryField_;
775 
776  if (!readIfPresent() && gf.field0Ptr_)
777  {
778  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
779  (
780  io.name() + "_0",
781  *gf.field0Ptr_
782  );
783  }
784 }
785 
786 
787 template<class Type, template<class> class PatchField, class GeoMesh>
789 (
790  const IOobject& io,
791  const GeometricField<Type, PatchField, GeoMesh>& gf,
792  const labelList& patchIDs,
793  const word& patchFieldType
794 )
795 :
796  Internal(io, gf),
797  timeIndex_(gf.timeIndex()),
798  field0Ptr_(nullptr),
799  fieldPrevIterPtr_(nullptr),
800  boundaryField_(*this, gf.boundaryField_, patchIDs, patchFieldType)
801 {
803  << "Copy construct, resetting IO params and setting patchFieldType "
804  << "for patchIDs" << nl
805  << this->info() << endl;
806 
807  if (!readIfPresent() && gf.field0Ptr_)
808  {
809  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
810  (
811  io.name() + "_0",
812  *gf.field0Ptr_
813  );
814  }
815 }
816 
817 
818 template<class Type, template<class> class PatchField, class GeoMesh>
820 (
821  const IOobject& io,
822  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf,
823  const wordList& patchFieldTypes,
824  const wordList& actualPatchTypes
825 )
826 :
827  Internal(io, tgf.constCast(), tgf.movable()),
828  timeIndex_(tgf().timeIndex()),
829  field0Ptr_(nullptr),
830  fieldPrevIterPtr_(nullptr),
831  boundaryField_
832  (
833  this->mesh().boundary(),
834  *this,
835  patchFieldTypes,
836  actualPatchTypes
837  )
838 {
840  << "Constructing from tmp resetting IO params and patch types" << nl
841  << this->info() << endl;
842 
843  boundaryField_ == tgf().boundaryField_;
845  tgf.clear();
846 }
847 
848 
849 template<class Type, template<class> class PatchField, class GeoMesh>
852 {
854 }
855 
856 
857 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
858 
859 template<class Type, template<class> class PatchField, class GeoMesh>
861 {
862  deleteDemandDrivenData(field0Ptr_);
863  deleteDemandDrivenData(fieldPrevIterPtr_);
864 }
865 
867 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
868 
869 template<class Type, template<class> class PatchField, class GeoMesh>
870 typename
873 (
874  const bool updateAccessTime
875 )
876 {
877  if (updateAccessTime)
878  {
879  this->setUpToDate();
880  storeOldTimes();
881  }
882  return *this;
883 }
884 
885 
886 template<class Type, template<class> class PatchField, class GeoMesh>
887 typename
890 (
891  const bool updateAccessTime
892 )
893 {
894  if (updateAccessTime)
895  {
896  this->setUpToDate();
897  storeOldTimes();
898  }
899  return *this;
900 }
901 
902 
903 template<class Type, template<class> class PatchField, class GeoMesh>
904 typename
907 (
908  const bool updateAccessTime
909 )
910 {
911  if (updateAccessTime)
912  {
913  this->setUpToDate();
914  storeOldTimes();
915  }
916  return boundaryField_;
917 }
918 
919 
920 template<class Type, template<class> class PatchField, class GeoMesh>
922 {
923  if
924  (
925  field0Ptr_
926  && timeIndex_ != this->time().timeIndex()
927  && !this->name().ends_with("_0")
928  )
929  {
930  storeOldTime();
931  timeIndex_ = this->time().timeIndex();
932  }
934  // Correct time index
935  //timeIndex_ = this->time().timeIndex();
936 }
937 
938 
939 template<class Type, template<class> class PatchField, class GeoMesh>
941 {
942  if (field0Ptr_)
943  {
944  field0Ptr_->storeOldTime();
945 
947  << "Storing old time field for field" << nl << this->info() << endl;
948 
949  *field0Ptr_ == *this;
950  field0Ptr_->timeIndex_ = timeIndex_;
951 
952  if (field0Ptr_->field0Ptr_)
953  {
954  field0Ptr_->writeOpt(this->writeOpt());
955  }
956  }
957 }
958 
959 
960 template<class Type, template<class> class PatchField, class GeoMesh>
962 {
963  if (field0Ptr_)
964  {
965  return field0Ptr_->nOldTimes() + 1;
966  }
968  return 0;
969 }
970 
971 
972 template<class Type, template<class> class PatchField, class GeoMesh>
975 {
976  if (!field0Ptr_)
977  {
979  (
980  IOobject
981  (
982  this->name() + "_0",
983  this->time().timeName(),
984  this->db(),
985  IOobject::NO_READ,
986  IOobject::NO_WRITE,
987  this->registerObject()
988  ),
989  *this
990  );
991 
992  if (debug)
993  {
995  << "created old time field " << field0Ptr_->info() << endl;
996 
997  if (debug&2)
998  {
999  error::printStack(Info);
1000  }
1001  }
1002  }
1003  else
1004  {
1005  storeOldTimes();
1006  }
1008  return *field0Ptr_;
1009 }
1010 
1011 
1012 template<class Type, template<class> class PatchField, class GeoMesh>
1015 {
1016  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1018 
1019  return *field0Ptr_;
1020 }
1021 
1022 
1023 template<class Type, template<class> class PatchField, class GeoMesh>
1025 {
1026  if (!fieldPrevIterPtr_)
1027  {
1029  << "Allocating previous iteration field" << nl
1030  << this->info() << endl;
1031 
1032  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1033  (
1034  this->name() + "PrevIter",
1035  *this
1036  );
1037  }
1038  else
1039  {
1040  *fieldPrevIterPtr_ == *this;
1041  }
1042 }
1043 
1044 
1045 template<class Type, template<class> class PatchField, class GeoMesh>
1048 {
1049  if (!fieldPrevIterPtr_)
1050  {
1052  << "previous iteration field" << endl << this->info() << endl
1053  << " not stored."
1054  << " Use field.storePrevIter() at start of iteration."
1055  << abort(FatalError);
1056  }
1058  return *fieldPrevIterPtr_;
1059 }
1060 
1061 
1062 template<class Type, template<class> class PatchField, class GeoMesh>
1065 {
1066  this->setUpToDate();
1067  storeOldTimes();
1068  boundaryField_.evaluate();
1069 }
1070 
1071 
1072 template<class Type, template<class> class PatchField, class GeoMesh>
1074 {
1075  // Search all boundary conditions, if any are
1076  // fixed-value or mixed (Robin) do not set reference level for solution.
1077 
1078  bool needRef = true;
1079 
1080  for (const auto& pf : boundaryField_)
1081  {
1082  if (pf.fixesValue())
1083  {
1084  needRef = false;
1085  break;
1086  }
1087  }
1088 
1089  return returnReduceAnd(needRef);
1090 }
1091 
1092 
1093 template<class Type, template<class> class PatchField, class GeoMesh>
1095 {
1097  << "Relaxing" << nl << this->info() << " by " << alpha << endl;
1098 
1099  operator==(prevIter() + alpha*(*this - prevIter()));
1100 }
1101 
1102 
1103 template<class Type, template<class> class PatchField, class GeoMesh>
1105 {
1106  word name = this->name();
1107 
1108  if
1109  (
1110  this->mesh().data::template getOrDefault<bool>
1111  (
1112  "finalIteration",
1113  false
1114  )
1115  )
1116  {
1117  name += "Final";
1118  }
1119 
1120  if (this->mesh().relaxField(name))
1121  {
1122  relax(this->mesh().fieldRelaxationFactor(name));
1123  }
1124 }
1125 
1126 
1127 template<class Type, template<class> class PatchField, class GeoMesh>
1129 (
1130  bool final
1131 ) const
1132 {
1133  if (final)
1134  {
1135  return this->name() + "Final";
1136  }
1138  return this->name();
1139 }
1140 
1141 
1142 template<class Type, template<class> class PatchField, class GeoMesh>
1144 (
1145  Ostream& os
1146 ) const
1147 {
1148  MinMax<Type> range = Foam::minMax(*this).value();
1149 
1150  os << "min/max(" << this->name() << ") = "
1151  << range.min() << ", " << range.max() << endl;
1152 }
1153 
1154 
1155 template<class Type, template<class> class PatchField, class GeoMesh>
1157 writeData(Ostream& os) const
1158 {
1159  os << *this;
1160  return os.good();
1162 
1163 
1164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1165 
1166 template<class Type, template<class> class PatchField, class GeoMesh>
1169 {
1171  (
1172  IOobject
1173  (
1174  this->name() + ".T()",
1175  this->instance(),
1176  this->db()
1177  ),
1178  this->mesh(),
1179  this->dimensions()
1180  );
1181 
1182  Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1183  Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1184 
1185  return tresult;
1186 }
1187 
1188 
1189 template<class Type, template<class> class PatchField, class GeoMesh>
1190 Foam::tmp
1191 <
1193  <
1195  PatchField,
1196  GeoMesh
1197  >
1198 >
1200 (
1201  const direction d
1202 ) const
1203 {
1204  auto tresult = tmp<GeometricField<cmptType, PatchField, GeoMesh>>::New
1205  (
1206  IOobject
1207  (
1208  this->name() + ".component(" + Foam::name(d) + ')',
1209  this->instance(),
1210  this->db()
1211  ),
1212  this->mesh(),
1213  this->dimensions()
1214  );
1215 
1216  Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1217  Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1218 
1219  return tresult;
1220 }
1221 
1222 
1223 template<class Type, template<class> class PatchField, class GeoMesh>
1225 (
1226  const direction d,
1227  const GeometricField
1228  <
1229  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1230  PatchField,
1231  GeoMesh
1232  >& gcf
1233 )
1234 {
1235  primitiveFieldRef().replace(d, gcf.primitiveField());
1236  boundaryFieldRef().replace(d, gcf.boundaryField());
1237 }
1238 
1239 
1240 template<class Type, template<class> class PatchField, class GeoMesh>
1242 (
1243  const direction d,
1244  const dimensioned<cmptType>& ds
1245 )
1246 {
1247  primitiveFieldRef().replace(d, ds.value());
1248  boundaryFieldRef().replace(d, ds.value());
1249 }
1250 
1251 
1252 template<class Type, template<class> class PatchField, class GeoMesh>
1254 (
1255  const dimensioned<Type>& dt
1256 )
1257 {
1258  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1259  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1260 }
1261 
1262 
1263 template<class Type, template<class> class PatchField, class GeoMesh>
1265 (
1266  const dimensioned<Type>& dt
1267 )
1268 {
1269  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1270  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1271 }
1272 
1273 
1274 template<class Type, template<class> class PatchField, class GeoMesh>
1276 (
1278 )
1279 {
1280  Foam::clip(primitiveFieldRef(), primitiveField(), range.value());
1281  Foam::clip(boundaryFieldRef(), boundaryField(), range.value());
1282 }
1283 
1284 
1285 template<class Type, template<class> class PatchField, class GeoMesh>
1287 (
1288  const dimensioned<Type>& minVal,
1289  const dimensioned<Type>& maxVal
1290 )
1291 {
1292  MinMax<Type> range(minVal.value(), maxVal.value());
1293 
1294  Foam::clip(primitiveFieldRef(), primitiveField(), range);
1295  Foam::clip(boundaryFieldRef(), boundaryField(), range);
1296 }
1297 
1298 
1299 template<class Type, template<class> class PatchField, class GeoMesh>
1301 (
1302  const dimensioned<Type>& minVal,
1303  const dimensioned<Type>& maxVal
1305 {
1306  this->clip(minVal, maxVal);
1307 }
1308 
1309 
1310 template<class Type, template<class> class PatchField, class GeoMesh>
1313  primitiveFieldRef().negate();
1314  boundaryFieldRef().negate();
1315 }
1316 
1317 
1318 template<class Type, template<class> class PatchField, class GeoMesh>
1320 {
1321  primitiveFieldRef().normalise();
1322  boundaryFieldRef().normalise();
1324 
1325 
1326 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1327 
1328 template<class Type, template<class> class PatchField, class GeoMesh>
1330 (
1332 )
1333 {
1334  if (this == &gf)
1335  {
1336  return; // Self-assignment is a no-op
1337  }
1338 
1339  checkField(*this, gf, "=");
1340 
1341  // Only assign field contents not ID
1342 
1343  internalFieldRef() = gf.internalField();
1344  boundaryFieldRef() = gf.boundaryField();
1345 }
1346 
1347 
1348 template<class Type, template<class> class PatchField, class GeoMesh>
1350 (
1352 )
1353 {
1354  const auto& gf = tgf();
1355 
1356  if (this == &gf)
1357  {
1358  return; // Self-assignment is a no-op
1359  }
1360 
1361  checkField(*this, gf, "=");
1362 
1363  // Only assign field contents not ID
1364 
1365  this->dimensions() = gf.dimensions();
1366  this->oriented() = gf.oriented();
1367 
1368  if (tgf.movable())
1369  {
1370  // Transfer storage from the tmp
1371  primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1372  }
1373  else
1374  {
1376  }
1377 
1378  boundaryFieldRef() = gf.boundaryField();
1380  tgf.clear();
1381 }
1382 
1383 
1384 template<class Type, template<class> class PatchField, class GeoMesh>
1386 (
1387  const dimensioned<Type>& dt
1388 )
1389 {
1390  internalFieldRef() = dt;
1391  boundaryFieldRef() = dt.value();
1392 }
1393 
1394 
1395 template<class Type, template<class> class PatchField, class GeoMesh>
1397 (
1399 )
1400 {
1401  const auto& gf = tgf();
1402 
1403  checkField(*this, gf, "==");
1404 
1405  // Only assign field contents not ID
1406 
1407  internalFieldRef() = gf.internalField();
1408  boundaryFieldRef() == gf.boundaryField();
1410  tgf.clear();
1411 }
1412 
1413 
1414 template<class Type, template<class> class PatchField, class GeoMesh>
1416 (
1417  const dimensioned<Type>& dt
1419 {
1420  internalFieldRef() = dt;
1421  boundaryFieldRef() == dt.value();
1422 }
1423 
1424 
1425 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1426  \
1427 template<class Type, template<class> class PatchField, class GeoMesh> \
1428 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1429 ( \
1430  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1431 ) \
1432 { \
1433  checkField(*this, gf, #op); \
1434  \
1435  internalFieldRef() op gf.internalField(); \
1436  boundaryFieldRef() op gf.boundaryField(); \
1437 } \
1438  \
1439 template<class Type, template<class> class PatchField, class GeoMesh> \
1440 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1441 ( \
1442  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1443 ) \
1444 { \
1445  operator op(tgf()); \
1446  tgf.clear(); \
1447 } \
1448  \
1449 template<class Type, template<class> class PatchField, class GeoMesh> \
1450 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1451 ( \
1452  const dimensioned<TYPE>& dt \
1453 ) \
1454 { \
1455  internalFieldRef() op dt; \
1456  boundaryFieldRef() op dt.value(); \
1457 }
1458 
1459 COMPUTED_ASSIGNMENT(Type, +=)
1460 COMPUTED_ASSIGNMENT(Type, -=)
1461 COMPUTED_ASSIGNMENT(scalar, *=)
1462 COMPUTED_ASSIGNMENT(scalar, /=)
1463 
1464 #undef COMPUTED_ASSIGNMENT
1465 
1466 
1467 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1468 
1469 template<class Type, template<class> class PatchField, class GeoMesh>
1470 Foam::Ostream& Foam::operator<<
1471 (
1472  Ostream& os,
1474 )
1475 {
1476  gf.internalField().writeData(os, "internalField");
1477  os << nl;
1478  gf.boundaryField().writeEntry("boundaryField", os);
1479 
1480  os.check(FUNCTION_NAME);
1481  return os;
1482 }
1483 
1484 
1485 template<class Type, template<class> class PatchField, class GeoMesh>
1486 Foam::Ostream& Foam::operator<<
1487 (
1488  Ostream& os,
1490 )
1491 {
1492  os << tgf();
1493  tgf.clear();
1494  return os;
1495 }
1496 
1497 
1498 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1499 
1500 #undef checkField
1501 
1502 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1503 
1504 #include "GeometricFieldNew.C"
1505 #include "GeometricFieldFunctions.C"
1506 
1507 // ************************************************************************* //
faceListList boundary
#define COMPUTED_ASSIGNMENT(TYPE, op)
const Type & value() const noexcept
Return const reference to value.
dictionary dict
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
#define checkField(gf1, gf2, op)
uint8_t direction
Definition: direction.H:46
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
UEqn relax()
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Field< Type >::cmptType cmptType
Component type of the field elements.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
A min/max value pair with additional methods. In addition to conveniently storing values...
Definition: HashSet.H:72
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
orientedType oriented() const noexcept
Return oriented type.
void clip(const dimensioned< MinMax< Type >> &range)
Clip the field to be bounded within the specified range.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Generic dimensioned Type class.
const dimensionSet dimless
Dimensionless.
label nOldTimes() const
Return the number of old time fields stored.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField(const IOobject &io, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=PatchField< Type >::calculatedType())
Construct given IOobject, mesh, dimensions and patch type.
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
conserve primitiveFieldRef()+
scalar range
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
word timeName
Definition: getTimeIndex.H:3
bool needReference() const
Does the field need a reference level for solution.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:100
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
dynamicFvMesh & mesh
void storeOldTimes() const
Store the old-time fields.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
Generic templated field type.
Definition: Field.H:61
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:271
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
void storeOldTime() const
Store the old-time field.
void maxMin(const dimensioned< Type > &minVal, const dimensioned< Type > &maxVal)
Deprecated(2019-01) identical to clip()
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
bool local
Definition: EEqn.H:20
propsDict readIfPresent("fields", acceptFields)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void negate()
Negate the field inplace. See notes in Field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
Internal & internalFieldRef(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
int debug
Static debugging option.
void normalise()
Normalise the field inplace. See notes in Field.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual ~GeometricField()
Destructor.
Generic GeometricBoundaryField class.
Definition: areaFieldsFwd.H:46
List< word > wordList
A List of words.
Definition: fileName.H:58
Template functions to aid in the implementation of demand driven data.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:607
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void storePrevIter() const
Store the field as the previous iteration value.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void deleteDemandDrivenData(DataPtr &dataPtr)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
void relax()
Relax field (for steady-state solution).
label timeIndex
Definition: getTimeIndex.H:24
const dimensionSet & dimensions() const noexcept
Return dimensions.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
#define InfoInFunction
Report an information message using Foam::Info.
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.