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-2023 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"
31 #include "demandDrivenData.H"
32 #include "dictionary.H"
34 #include "data.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 // Check that both fields use the same mesh
39 #undef checkField
40 #define checkField(fld1, fld2, op) \
41 if (&(fld1).mesh() != &(fld2).mesh()) \
42 { \
43  FatalErrorInFunction \
44  << "Different mesh for fields " \
45  << (fld1).name() << " and " << (fld2).name() \
46  << " during operation " << op \
47  << abort(FatalError); \
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
52 
53 template<class Type, template<class> class PatchField, class GeoMesh>
55 (
56  const dictionary& dict
57 )
58 {
59  Internal::readField(dict, "internalField");
60 
61  boundaryField_.readField(*this, dict.subDict("boundaryField"));
62 
63  Type refLevel;
64 
65  if (dict.readIfPresent("referenceLevel", refLevel))
66  {
67  Field<Type>::operator+=(refLevel);
68 
69  forAll(boundaryField_, patchi)
70  {
71  boundaryField_[patchi] == boundaryField_[patchi] + refLevel;
72  }
73  }
74 }
75 
76 
77 template<class Type, template<class> class PatchField, class GeoMesh>
79 {
80  const localIOdictionary dict
81  (
82  IOobject
83  (
84  this->name(),
85  this->instance(),
86  this->local(),
87  this->db(),
88  IOobject::MUST_READ,
89  IOobject::NO_WRITE,
90  IOobject::NO_REGISTER
91  ),
92  typeName
93  );
94 
95  this->close();
96 
98 }
99 
100 
101 template<class Type, template<class> class PatchField, class GeoMesh>
103 {
104  if (this->isReadRequired())
105  {
107  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
108  << " suggests that a read constructor for field " << this->name()
109  << " would be more appropriate." << endl;
110  }
111  else if
112  (
113  this->isReadOptional()
114  && this->template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
115  (
116  true
117  )
118  )
119  {
120  readFields();
121 
122  // Check compatibility between field and mesh
123  if (this->size() != GeoMesh::size(this->mesh()))
124  {
125  FatalIOErrorInFunction(this->readStream(typeName))
126  << " number of field elements = " << this->size()
127  << " number of mesh elements = "
128  << GeoMesh::size(this->mesh())
129  << exit(FatalIOError);
130  }
131 
132  readOldTimeIfPresent();
133 
134  return true;
135  }
136 
137  return false;
138 }
139 
140 
141 template<class Type, template<class> class PatchField, class GeoMesh>
143 {
144  // Read the old time field if present
145  IOobject field0
146  (
147  this->name() + "_0",
148  this->time().timeName(),
149  this->db(),
150  IOobject::READ_IF_PRESENT,
151  IOobject::AUTO_WRITE,
152  this->registerObject()
153  );
154 
155  if
156  (
157  field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
158  (
159  true
160  )
161  )
162  {
164  << "Reading old time level for field" << nl << this->info() << endl;
165 
166  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
167  (
168  field0,
169  this->mesh()
170  );
171 
172  // Ensure the old time field oriented flag is set to the parent's state
173  // Note: required for backwards compatibility in case of restarting from
174  // an old run where the oriented state may not have been set
175  field0Ptr_->oriented() = this->oriented();
176 
177  field0Ptr_->timeIndex_ = timeIndex_ - 1;
178 
179  if (!field0Ptr_->readOldTimeIfPresent())
180  {
181  field0Ptr_->oldTime();
182  }
183 
184  return true;
185  }
186 
187  return false;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
192 
193 template<class Type, template<class> class PatchField, class GeoMesh>
195 (
196  const IOobject& io,
197  const Mesh& mesh,
198  const dimensionSet& dims,
199  const word& patchFieldType
200 )
201 :
202  Internal(io, mesh, dims, false),
203  timeIndex_(this->time().timeIndex()),
204  field0Ptr_(nullptr),
205  fieldPrevIterPtr_(nullptr),
206  boundaryField_(mesh.boundary(), *this, patchFieldType)
207 {
209  << "Creating" << nl << this->info() << endl;
211  readIfPresent();
212 }
213 
214 
215 template<class Type, template<class> class PatchField, class GeoMesh>
217 (
218  const IOobject& io,
219  const Mesh& mesh,
220  const dimensionSet& dims,
221  const wordList& patchFieldTypes,
222  const wordList& actualPatchTypes
223 )
224 :
225  Internal(io, mesh, dims, false),
226  timeIndex_(this->time().timeIndex()),
227  field0Ptr_(nullptr),
228  fieldPrevIterPtr_(nullptr),
229  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
230 {
232  << "Creating" << nl << this->info() << endl;
233 
235 }
236 
237 
238 
239 template<class Type, template<class> class PatchField, class GeoMesh>
241 (
242  const IOobject& io,
243  const Mesh& mesh,
244  const Type& value,
245  const dimensionSet& dims,
246  const word& patchFieldType
247 )
248 :
249  Internal(io, mesh, value, dims, false),
250  timeIndex_(this->time().timeIndex()),
251  field0Ptr_(nullptr),
252  fieldPrevIterPtr_(nullptr),
253  boundaryField_(mesh.boundary(), *this, patchFieldType)
254 {
256  << "Creating" << nl << this->info() << endl;
257 
258  boundaryField_ == value;
260  readIfPresent();
261 }
262 
263 
264 template<class Type, template<class> class PatchField, class GeoMesh>
266 (
267  const IOobject& io,
268  const Mesh& mesh,
269  const Type& value,
270  const dimensionSet& dims,
271  const wordList& patchFieldTypes,
272  const wordList& actualPatchTypes
273 )
274 :
275  Internal(io, mesh, value, dims, false),
276  timeIndex_(this->time().timeIndex()),
277  field0Ptr_(nullptr),
278  fieldPrevIterPtr_(nullptr),
279  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
280 {
282  << "Creating" << nl << this->info() << endl;
283 
284  boundaryField_ == value;
286  readIfPresent();
287 }
288 
289 
290 template<class Type, template<class> class PatchField, class GeoMesh>
292 (
293  const IOobject& io,
294  const Mesh& mesh,
295  const dimensioned<Type>& dt,
296  const word& patchFieldType
297 )
298 :
299  GeometricField<Type, PatchField, GeoMesh>
300  (
301  io,
302  mesh,
303  dt.value(),
304  dt.dimensions(),
305  patchFieldType
306  )
307 {}
308 
309 
310 template<class Type, template<class> class PatchField, class GeoMesh>
312 (
313  const IOobject& io,
314  const Mesh& mesh,
315  const dimensioned<Type>& dt,
316  const wordList& patchFieldTypes,
317  const wordList& actualPatchTypes
318 )
319 :
320  GeometricField<Type, PatchField, GeoMesh>
321  (
322  io,
323  mesh,
324  dt.value(),
325  dt.dimensions(),
326  patchFieldTypes,
327  actualPatchTypes
328  )
329 {}
330 
331 
332 template<class Type, template<class> class PatchField, class GeoMesh>
334 (
335  const IOobject& io,
336  const Internal& diField,
337  const PtrList<PatchField<Type>>& ptfl
338 )
339 :
340  Internal(io, diField),
341  timeIndex_(this->time().timeIndex()),
342  field0Ptr_(nullptr),
343  fieldPrevIterPtr_(nullptr),
344  boundaryField_(this->mesh().boundary(), *this, ptfl)
345 {
347  << "Copy construct from components" << nl << this->info() << endl;
349  readIfPresent();
350 }
351 
352 
353 template<class Type, template<class> class PatchField, class GeoMesh>
355 (
356  const IOobject& io,
357  Internal&& diField,
358  const PtrList<PatchField<Type>>& ptfl
359 )
360 :
361  Internal(io, std::move(diField)),
362  timeIndex_(this->time().timeIndex()),
363  field0Ptr_(nullptr),
364  fieldPrevIterPtr_(nullptr),
365  boundaryField_(this->mesh().boundary(), *this, ptfl)
366 {
368  << "Move construct from components" << nl << this->info() << endl;
370  readIfPresent();
371 }
372 
373 
374 template<class Type, template<class> class PatchField, class GeoMesh>
376 (
377  const IOobject& io,
378  const tmp<Internal>& tfield,
379  const PtrList<PatchField<Type>>& ptfl
380 )
381 :
382  Internal(io, tfield),
383  timeIndex_(this->time().timeIndex()),
384  field0Ptr_(nullptr),
385  fieldPrevIterPtr_(nullptr),
386  boundaryField_(this->mesh().boundary(), *this, ptfl)
387 {
389  << "Construct from tmp internalField" << nl << this->info() << endl;
391  readIfPresent();
392 }
393 
394 
395 template<class Type, template<class> class PatchField, class GeoMesh>
397 (
398  const Internal& diField,
399  const PtrList<PatchField<Type>>& ptfl
400 )
401 :
402  Internal(diField),
403  timeIndex_(this->time().timeIndex()),
404  field0Ptr_(nullptr),
405  fieldPrevIterPtr_(nullptr),
406  boundaryField_(this->mesh().boundary(), *this, ptfl)
407 {
409  << "Copy construct from components" << nl << this->info() << endl;
411  readIfPresent();
412 }
413 
414 
415 template<class Type, template<class> class PatchField, class GeoMesh>
417 (
418  Internal&& diField,
419  const PtrList<PatchField<Type>>& ptfl
420 )
421 :
422  Internal(std::move(diField)),
423  timeIndex_(this->time().timeIndex()),
424  field0Ptr_(nullptr),
425  fieldPrevIterPtr_(nullptr),
426  boundaryField_(this->mesh().boundary(), *this, ptfl)
427 {
429  << "Move construct from components" << nl << this->info() << endl;
431  readIfPresent();
432 }
433 
434 
435 template<class Type, template<class> class PatchField, class GeoMesh>
437 (
438  const IOobject& io,
439  const Mesh& mesh,
440  const dimensionSet& dims,
441  const Field<Type>& iField,
442  const word& patchFieldType
443 )
444 :
445  Internal(io, mesh, dims, iField),
446  timeIndex_(this->time().timeIndex()),
447  field0Ptr_(nullptr),
448  fieldPrevIterPtr_(nullptr),
449  boundaryField_(mesh.boundary(), *this, patchFieldType)
450 {
452  << "Copy construct from internal field" << nl << this->info() << endl;
454  readIfPresent();
455 }
456 
457 
458 template<class Type, template<class> class PatchField, class GeoMesh>
460 (
461  const IOobject& io,
462  const Mesh& mesh,
463  const dimensionSet& dims,
464  Field<Type>&& iField,
465  const word& patchFieldType
466 )
467 :
468  Internal(io, mesh, dims, std::move(iField)),
469  timeIndex_(this->time().timeIndex()),
470  field0Ptr_(nullptr),
471  fieldPrevIterPtr_(nullptr),
472  boundaryField_(mesh.boundary(), *this, patchFieldType)
473 {
475  << "Move construct from internal field" << nl << this->info() << endl;
477  readIfPresent();
478 }
479 
480 
481 template<class Type, template<class> class PatchField, class GeoMesh>
483 (
484  const IOobject& io,
485  const Mesh& mesh,
486  const dimensionSet& dims,
487  const Field<Type>& iField,
488  const PtrList<PatchField<Type>>& ptfl
489 )
490 :
491  Internal(io, mesh, dims, iField),
492  timeIndex_(this->time().timeIndex()),
493  field0Ptr_(nullptr),
494  fieldPrevIterPtr_(nullptr),
495  boundaryField_(mesh.boundary(), *this, ptfl)
496 {
498  << "Copy construct from components" << nl << this->info() << endl;
500  readIfPresent();
501 }
502 
503 
504 template<class Type, template<class> class PatchField, class GeoMesh>
506 (
507  const IOobject& io,
508  const Mesh& mesh,
509  const dimensionSet& dims,
510  Field<Type>&& iField,
511  const PtrList<PatchField<Type>>& ptfl
512 )
513 :
514  Internal(io, mesh, dims, std::move(iField)),
515  timeIndex_(this->time().timeIndex()),
516  field0Ptr_(nullptr),
517  fieldPrevIterPtr_(nullptr),
518  boundaryField_(mesh.boundary(), *this, ptfl)
519 {
521  << "Move construct from components" << nl << this->info() << endl;
523  readIfPresent();
524 }
525 
526 
527 template<class Type, template<class> class PatchField, class GeoMesh>
529 (
530  const IOobject& io,
531  const Mesh& mesh,
532  const dimensionSet& dims,
533  const tmp<Field<Type>>& tfield,
534  const PtrList<PatchField<Type>>& ptfl
535 )
536 :
537  Internal(io, mesh, dims, tfield),
538  timeIndex_(this->time().timeIndex()),
539  field0Ptr_(nullptr),
540  fieldPrevIterPtr_(nullptr),
541  boundaryField_(mesh.boundary(), *this, ptfl)
542 {
544  << "Construct from tmp internalField" << nl << this->info() << endl;
546  readIfPresent();
547 }
548 
549 
550 template<class Type, template<class> class PatchField, class GeoMesh>
552 (
553  const IOobject& io,
554  const Mesh& mesh,
555  const bool readOldTime
556 )
557 :
558  Internal(io, mesh, dimless, false),
559  timeIndex_(this->time().timeIndex()),
560  field0Ptr_(nullptr),
561  fieldPrevIterPtr_(nullptr),
562  boundaryField_(mesh.boundary())
563 {
564  readFields();
565 
566  // Check compatibility between field and mesh
567 
568  if (this->size() != GeoMesh::size(this->mesh()))
569  {
570  FatalIOErrorInFunction(this->readStream(typeName))
571  << " number of field elements = " << this->size()
572  << " number of mesh elements = " << GeoMesh::size(this->mesh())
573  << exit(FatalIOError);
574  }
575 
576  if (readOldTime)
577  {
578  readOldTimeIfPresent();
579  }
580 
582  << "Finishing read-construction" << nl << this->info() << endl;
583 }
584 
585 
586 template<class Type, template<class> class PatchField, class GeoMesh>
588 (
589  const IOobject& io,
590  const Mesh& mesh,
591  const dictionary& dict
592 )
593 :
594  Internal(io, mesh, dimless, false),
595  timeIndex_(this->time().timeIndex()),
596  field0Ptr_(nullptr),
597  fieldPrevIterPtr_(nullptr),
598  boundaryField_(mesh.boundary())
599 {
600  readFields(dict);
601 
602  // Check compatibility between field and mesh
603 
604  if (this->size() != GeoMesh::size(this->mesh()))
605  {
607  << " number of field elements = " << this->size()
608  << " number of mesh elements = " << GeoMesh::size(this->mesh())
609  << exit(FatalIOError);
610  }
611 
613  << "Finishing dictionary-construct" << nl << this->info() << endl;
614 }
615 
616 
617 template<class Type, template<class> class PatchField, class GeoMesh>
619 (
621 )
622 :
623  Internal(gf),
624  timeIndex_(gf.timeIndex()),
625  field0Ptr_(nullptr),
626  fieldPrevIterPtr_(nullptr),
627  boundaryField_(*this, gf.boundaryField_)
628 {
630  << "Copy construct" << nl << this->info() << endl;
631 
632  if (gf.field0Ptr_)
633  {
635  (
636  *gf.field0Ptr_
637  );
638  }
640  this->writeOpt(IOobject::NO_WRITE);
641 }
642 
643 
644 template<class Type, template<class> class PatchField, class GeoMesh>
646 (
648 )
649 :
650  Internal(tgf.constCast(), tgf.movable()),
651  timeIndex_(tgf().timeIndex()),
652  field0Ptr_(nullptr),
653  fieldPrevIterPtr_(nullptr),
654  boundaryField_(*this, tgf().boundaryField_)
655 {
657  << "Constructing from tmp" << nl << this->info() << endl;
658 
659  this->writeOpt(IOobject::NO_WRITE);
661  tgf.clear();
662 }
663 
664 
665 template<class Type, template<class> class PatchField, class GeoMesh>
667 (
668  const IOobject& io,
670 )
671 :
672  Internal(io, gf),
673  timeIndex_(gf.timeIndex()),
674  field0Ptr_(nullptr),
675  fieldPrevIterPtr_(nullptr),
676  boundaryField_(*this, gf.boundaryField_)
677 {
679  << "Copy construct, resetting IO params" << nl
680  << this->info() << endl;
681 
682  if (!readIfPresent() && gf.field0Ptr_)
683  {
685  (
686  io.name() + "_0",
687  *gf.field0Ptr_
688  );
689  }
690 }
691 
692 
693 template<class Type, template<class> class PatchField, class GeoMesh>
695 (
696  const IOobject& io,
697  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
698 )
699 :
700  Internal(io, tgf.constCast(), tgf.movable()),
701  timeIndex_(tgf().timeIndex()),
702  field0Ptr_(nullptr),
703  fieldPrevIterPtr_(nullptr),
704  boundaryField_(*this, tgf().boundaryField_)
705 {
707  << "Constructing from tmp resetting IO params" << nl
708  << this->info() << endl;
709 
710  tgf.clear();
712  readIfPresent();
713 }
714 
715 
716 template<class Type, template<class> class PatchField, class GeoMesh>
718 (
719  const word& newName,
721 )
722 :
723  Internal(newName, gf),
724  timeIndex_(gf.timeIndex()),
725  field0Ptr_(nullptr),
726  fieldPrevIterPtr_(nullptr),
727  boundaryField_(*this, gf.boundaryField_)
728 {
730  << "Copy construct, resetting name" << nl
731  << this->info() << endl;
732 
733  if (!readIfPresent() && gf.field0Ptr_)
734  {
736  (
737  newName + "_0",
738  *gf.field0Ptr_
739  );
740  }
741 }
742 
743 
744 template<class Type, template<class> class PatchField, class GeoMesh>
746 (
747  const word& newName,
748  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
749 )
750 :
751  Internal(newName, tgf.constCast(), tgf.movable()),
752  timeIndex_(tgf().timeIndex()),
753  field0Ptr_(nullptr),
754  fieldPrevIterPtr_(nullptr),
755  boundaryField_(*this, tgf().boundaryField_)
756 {
758  << "Constructing from tmp resetting name" << nl
759  << this->info() << endl;
761  tgf.clear();
762 }
763 
764 
765 template<class Type, template<class> class PatchField, class GeoMesh>
767 (
768  const IOobject& io,
770  const word& patchFieldType
771 )
772 :
773  Internal(io, gf),
774  timeIndex_(gf.timeIndex()),
775  field0Ptr_(nullptr),
776  fieldPrevIterPtr_(nullptr),
777  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
778 {
780  << "Copy construct, resetting IO params" << nl
781  << this->info() << endl;
782 
783  boundaryField_ == gf.boundaryField_;
784 
785  if (!readIfPresent() && gf.field0Ptr_)
786  {
788  (
789  io.name() + "_0",
790  *gf.field0Ptr_
791  );
792  }
793 }
794 
795 
796 template<class Type, template<class> class PatchField, class GeoMesh>
798 (
799  const IOobject& io,
800  const GeometricField<Type, PatchField, GeoMesh>& gf,
801  const wordList& patchFieldTypes,
802  const wordList& actualPatchTypes
803 )
804 :
805  Internal(io, gf),
806  timeIndex_(gf.timeIndex()),
807  field0Ptr_(nullptr),
808  fieldPrevIterPtr_(nullptr),
809  boundaryField_
810  (
811  this->mesh().boundary(),
812  *this,
813  patchFieldTypes,
814  actualPatchTypes
815  )
816 {
818  << "Copy construct, resetting IO params and patch types" << nl
819  << this->info() << endl;
820 
821  boundaryField_ == gf.boundaryField_;
822 
823  if (!readIfPresent() && gf.field0Ptr_)
824  {
825  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
826  (
827  io.name() + "_0",
828  *gf.field0Ptr_
829  );
830  }
831 }
832 
833 
834 template<class Type, template<class> class PatchField, class GeoMesh>
836 (
837  const IOobject& io,
838  const GeometricField<Type, PatchField, GeoMesh>& gf,
839  const labelList& patchIDs,
840  const word& patchFieldType
841 )
842 :
843  Internal(io, gf),
844  timeIndex_(gf.timeIndex()),
845  field0Ptr_(nullptr),
846  fieldPrevIterPtr_(nullptr),
847  boundaryField_(*this, gf.boundaryField_, patchIDs, patchFieldType)
848 {
850  << "Copy construct, resetting IO params and setting patchFieldType "
851  << "for patchIDs" << nl
852  << this->info() << endl;
853 
854  if (!readIfPresent() && gf.field0Ptr_)
855  {
856  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
857  (
858  io.name() + "_0",
859  *gf.field0Ptr_
860  );
861  }
862 }
863 
864 
865 template<class Type, template<class> class PatchField, class GeoMesh>
867 (
868  const IOobject& io,
869  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf,
870  const wordList& patchFieldTypes,
871  const wordList& actualPatchTypes
872 )
873 :
874  Internal(io, tgf.constCast(), tgf.movable()),
875  timeIndex_(tgf().timeIndex()),
876  field0Ptr_(nullptr),
877  fieldPrevIterPtr_(nullptr),
878  boundaryField_
879  (
880  this->mesh().boundary(),
881  *this,
882  patchFieldTypes,
883  actualPatchTypes
884  )
885 {
887  << "Constructing from tmp resetting IO params and patch types" << nl
888  << this->info() << endl;
889 
890  boundaryField_ == tgf().boundaryField_;
892  tgf.clear();
893 }
894 
895 
896 template<class Type, template<class> class PatchField, class GeoMesh>
899 {
901 }
902 
903 
904 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
905 
906 template<class Type, template<class> class PatchField, class GeoMesh>
908 {
909  deleteDemandDrivenData(field0Ptr_);
910  deleteDemandDrivenData(fieldPrevIterPtr_);
911 }
912 
914 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
915 
916 template<class Type, template<class> class PatchField, class GeoMesh>
917 typename
920 (
921  const bool updateAccessTime
922 )
923 {
924  if (updateAccessTime)
925  {
926  this->setUpToDate();
927  storeOldTimes();
928  }
929  return *this;
930 }
931 
932 
933 template<class Type, template<class> class PatchField, class GeoMesh>
934 typename
937 (
938  const bool updateAccessTime
939 )
940 {
941  if (updateAccessTime)
942  {
943  this->setUpToDate();
944  storeOldTimes();
945  }
946  return *this;
947 }
948 
949 
950 template<class Type, template<class> class PatchField, class GeoMesh>
951 typename
954 (
955  const bool updateAccessTime
956 )
957 {
958  if (updateAccessTime)
959  {
960  this->setUpToDate();
961  storeOldTimes();
962  }
963  return boundaryField_;
964 }
965 
966 
967 template<class Type, template<class> class PatchField, class GeoMesh>
969 {
970  if
971  (
972  field0Ptr_
973  && timeIndex_ != this->time().timeIndex()
974  && !this->name().ends_with("_0")
975  )
976  {
977  storeOldTime();
978  timeIndex_ = this->time().timeIndex();
979  }
981  // Correct time index
982  //timeIndex_ = this->time().timeIndex();
983 }
984 
985 
986 template<class Type, template<class> class PatchField, class GeoMesh>
988 {
989  if (field0Ptr_)
990  {
991  field0Ptr_->storeOldTime();
992 
994  << "Storing old time field for field" << nl << this->info() << endl;
995 
996  *field0Ptr_ == *this;
997  field0Ptr_->timeIndex_ = timeIndex_;
998 
999  if (field0Ptr_->field0Ptr_)
1000  {
1001  field0Ptr_->writeOpt(this->writeOpt());
1002  }
1003  }
1004 }
1005 
1006 
1007 template<class Type, template<class> class PatchField, class GeoMesh>
1009 {
1010  if (field0Ptr_)
1011  {
1012  return field0Ptr_->nOldTimes() + 1;
1013  }
1015  return 0;
1016 }
1017 
1018 
1019 template<class Type, template<class> class PatchField, class GeoMesh>
1022 {
1023  if (!field0Ptr_)
1024  {
1026  (
1027  IOobject
1028  (
1029  this->name() + "_0",
1030  this->time().timeName(),
1031  this->db(),
1032  IOobject::NO_READ,
1033  IOobject::NO_WRITE,
1034  this->registerObject()
1035  ),
1036  *this
1037  );
1038 
1039  if (debug)
1040  {
1042  << "created old time field " << field0Ptr_->info() << endl;
1043 
1044  if (debug&2)
1045  {
1046  error::printStack(Info);
1047  }
1048  }
1049  }
1050  else
1051  {
1052  storeOldTimes();
1053  }
1055  return *field0Ptr_;
1056 }
1057 
1058 
1059 template<class Type, template<class> class PatchField, class GeoMesh>
1062 {
1063  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1065 
1066  return *field0Ptr_;
1067 }
1068 
1069 
1070 template<class Type, template<class> class PatchField, class GeoMesh>
1072 {
1073  if (!fieldPrevIterPtr_)
1074  {
1076  << "Allocating previous iteration field" << nl
1077  << this->info() << endl;
1078 
1079  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1080  (
1081  this->name() + "PrevIter",
1082  *this
1083  );
1084  }
1085  else
1086  {
1087  *fieldPrevIterPtr_ == *this;
1088  }
1089 }
1090 
1091 
1092 template<class Type, template<class> class PatchField, class GeoMesh>
1095 {
1096  if (!fieldPrevIterPtr_)
1097  {
1099  << "previous iteration field" << endl << this->info() << endl
1100  << " not stored."
1101  << " Use field.storePrevIter() at start of iteration."
1102  << abort(FatalError);
1103  }
1105  return *fieldPrevIterPtr_;
1106 }
1107 
1108 
1109 template<class Type, template<class> class PatchField, class GeoMesh>
1112 {
1113  this->setUpToDate();
1114  storeOldTimes();
1115  boundaryField_.evaluate();
1116 }
1117 
1118 
1119 template<class Type, template<class> class PatchField, class GeoMesh>
1121 {
1122  // Search all boundary conditions, if any are
1123  // fixed-value or mixed (Robin) do not set reference level for solution.
1124 
1125  bool needRef = true;
1126 
1127  for (const auto& pf : boundaryField_)
1128  {
1129  if (pf.fixesValue())
1130  {
1131  needRef = false;
1132  break;
1133  }
1134  }
1135 
1136  return returnReduceAnd(needRef);
1137 }
1138 
1139 
1140 template<class Type, template<class> class PatchField, class GeoMesh>
1142 {
1144  << "Relaxing" << nl << this->info() << " by " << alpha << endl;
1145 
1146  operator==(prevIter() + alpha*(*this - prevIter()));
1147 }
1148 
1149 
1150 template<class Type, template<class> class PatchField, class GeoMesh>
1152 {
1153  word name = this->name();
1154 
1155  if
1156  (
1157  this->mesh().data::template getOrDefault<bool>
1158  (
1159  "finalIteration",
1160  false
1161  )
1162  )
1163  {
1164  name += "Final";
1165  }
1166 
1167  if (this->mesh().relaxField(name))
1168  {
1169  relax(this->mesh().fieldRelaxationFactor(name));
1170  }
1171 }
1172 
1173 
1174 template<class Type, template<class> class PatchField, class GeoMesh>
1176 (
1177  bool final
1178 ) const
1179 {
1180  if (final)
1181  {
1182  return this->name() + "Final";
1183  }
1185  return this->name();
1186 }
1187 
1188 
1189 template<class Type, template<class> class PatchField, class GeoMesh>
1191 (
1192  Ostream& os
1193 ) const
1194 {
1195  MinMax<Type> range = Foam::minMax(*this).value();
1196 
1197  os << "min/max(" << this->name() << ") = "
1198  << range.min() << ", " << range.max() << endl;
1199 }
1200 
1201 
1202 template<class Type, template<class> class PatchField, class GeoMesh>
1204 writeData(Ostream& os) const
1205 {
1206  os << *this;
1207  return os.good();
1209 
1210 
1211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1212 
1213 template<class Type, template<class> class PatchField, class GeoMesh>
1216 {
1218  (
1219  this->name() + ".T()",
1220  this->mesh(),
1221  this->dimensions()
1222  );
1223 
1224  Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1225  Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1226 
1227  return tresult;
1228 }
1229 
1230 
1231 template<class Type, template<class> class PatchField, class GeoMesh>
1232 Foam::tmp
1233 <
1235  <
1237  PatchField,
1238  GeoMesh
1239  >
1240 >
1242 (
1243  const direction d
1244 ) const
1245 {
1247  (
1248  this->name() + ".component(" + Foam::name(d) + ')',
1249  this->mesh(),
1250  this->dimensions()
1251  );
1252 
1253  Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1254  Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1255 
1256  return tresult;
1257 }
1258 
1259 
1260 template<class Type, template<class> class PatchField, class GeoMesh>
1262 (
1263  const direction d,
1264  const GeometricField
1265  <
1266  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1267  PatchField,
1268  GeoMesh
1269  >& gcf
1270 )
1271 {
1272  primitiveFieldRef().replace(d, gcf.primitiveField());
1273  boundaryFieldRef().replace(d, gcf.boundaryField());
1274 }
1275 
1276 
1277 template<class Type, template<class> class PatchField, class GeoMesh>
1279 (
1280  const direction d,
1281  const dimensioned<cmptType>& ds
1282 )
1283 {
1284  primitiveFieldRef().replace(d, ds.value());
1285  boundaryFieldRef().replace(d, ds.value());
1286 }
1287 
1288 
1289 template<class Type, template<class> class PatchField, class GeoMesh>
1291 (
1292  const Type& lower
1293 )
1294 {
1295  primitiveFieldRef().clamp_min(lower);
1296  boundaryFieldRef().clamp_min(lower);
1297 }
1298 
1299 
1300 template<class Type, template<class> class PatchField, class GeoMesh>
1302 (
1303  const Type& upper
1304 )
1305 {
1306  primitiveFieldRef().clamp_max(upper);
1307  boundaryFieldRef().clamp_max(upper);
1308 }
1309 
1310 
1311 template<class Type, template<class> class PatchField, class GeoMesh>
1313 (
1314  const dimensioned<Type>& lower
1315 )
1317  this->clamp_min(lower.value());
1318 }
1319 
1320 
1321 template<class Type, template<class> class PatchField, class GeoMesh>
1323 (
1324  const dimensioned<Type>& upper
1325 )
1327  this->clamp_max(upper.value());
1328 }
1329 
1330 
1331 template<class Type, template<class> class PatchField, class GeoMesh>
1333 (
1334  const Type& lower,
1335  const Type& upper
1336 )
1337 {
1338  primitiveFieldRef().clamp_range(lower, upper);
1339  boundaryFieldRef().clamp_range(lower, upper);
1340 }
1341 
1342 
1343 template<class Type, template<class> class PatchField, class GeoMesh>
1345 (
1346  const MinMax<Type>& range
1347 )
1348 {
1349  primitiveFieldRef().clamp_range(range.min(), range.max());
1350  boundaryFieldRef().clamp_range(range.min(), range.max());
1351 }
1352 
1353 
1354 template<class Type, template<class> class PatchField, class GeoMesh>
1356 (
1357  const dimensioned<Type>& lower,
1358  const dimensioned<Type>& upper
1359 )
1361  this->clamp_range(lower.value(), upper.value());
1362 }
1363 
1364 
1365 template<class Type, template<class> class PatchField, class GeoMesh>
1367 (
1370 {
1371  this->clamp_range(range.value());
1372 }
1373 
1374 
1375 template<class Type, template<class> class PatchField, class GeoMesh>
1378  primitiveFieldRef().negate();
1379  boundaryFieldRef().negate();
1380 }
1381 
1382 
1383 template<class Type, template<class> class PatchField, class GeoMesh>
1385 {
1386  primitiveFieldRef().normalise();
1387  boundaryFieldRef().normalise();
1389 
1390 
1391 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1392 
1393 template<class Type, template<class> class PatchField, class GeoMesh>
1395 (
1397 )
1398 {
1399  if (this == &gf)
1400  {
1401  return; // Self-assignment is a no-op
1402  }
1403 
1404  checkField(*this, gf, "=");
1405 
1406  // Only assign field contents not ID
1407 
1408  internalFieldRef() = gf.internalField();
1410 }
1411 
1412 
1413 template<class Type, template<class> class PatchField, class GeoMesh>
1415 (
1417 )
1418 {
1419  const auto& gf = tgf();
1420 
1421  if (this == &gf)
1422  {
1423  return; // Self-assignment is a no-op
1424  }
1425 
1426  checkField(*this, gf, "=");
1427 
1428  // Only assign field contents not ID
1429 
1430  this->dimensions() = gf.dimensions();
1431  this->oriented() = gf.oriented();
1432 
1433  if (tgf.movable())
1434  {
1435  // Transfer storage from the tmp
1436  primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1437  }
1438  else
1439  {
1441  }
1442 
1445  tgf.clear();
1446 }
1447 
1448 
1449 template<class Type, template<class> class PatchField, class GeoMesh>
1451 (
1452  const dimensioned<Type>& dt
1453 )
1454 {
1455  internalFieldRef() = dt;
1456  boundaryFieldRef() = dt.value();
1457 }
1458 
1459 
1460 template<class Type, template<class> class PatchField, class GeoMesh>
1462 (
1464 )
1465 {
1466  const auto& gf = tgf();
1467 
1468  checkField(*this, gf, "==");
1469 
1470  // Only assign field contents not ID
1471 
1472  internalFieldRef() = gf.internalField();
1473  boundaryFieldRef() == gf.boundaryField();
1475  tgf.clear();
1476 }
1477 
1478 
1479 template<class Type, template<class> class PatchField, class GeoMesh>
1481 (
1482  const dimensioned<Type>& dt
1484 {
1485  internalFieldRef() = dt;
1486  boundaryFieldRef() == dt.value();
1487 }
1488 
1489 
1490 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1491  \
1492 template<class Type, template<class> class PatchField, class GeoMesh> \
1493 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1494 ( \
1495  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1496 ) \
1497 { \
1498  checkField(*this, gf, #op); \
1499  \
1500  internalFieldRef() op gf.internalField(); \
1501  boundaryFieldRef() op gf.boundaryField(); \
1502 } \
1503  \
1504 template<class Type, template<class> class PatchField, class GeoMesh> \
1505 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1506 ( \
1507  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1508 ) \
1509 { \
1510  operator op(tgf()); \
1511  tgf.clear(); \
1512 } \
1513  \
1514 template<class Type, template<class> class PatchField, class GeoMesh> \
1515 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1516 ( \
1517  const dimensioned<TYPE>& dt \
1518 ) \
1519 { \
1520  internalFieldRef() op dt; \
1521  boundaryFieldRef() op dt.value(); \
1522 }
1523 
1524 COMPUTED_ASSIGNMENT(Type, +=)
1525 COMPUTED_ASSIGNMENT(Type, -=)
1526 COMPUTED_ASSIGNMENT(scalar, *=)
1527 COMPUTED_ASSIGNMENT(scalar, /=)
1528 
1529 #undef COMPUTED_ASSIGNMENT
1530 
1531 
1532 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1533 
1534 template<class Type, template<class> class PatchField, class GeoMesh>
1535 Foam::Ostream& Foam::operator<<
1536 (
1537  Ostream& os,
1539 )
1540 {
1541  gf.internalField().writeData(os, "internalField");
1542  os << nl;
1543  gf.boundaryField().writeEntry("boundaryField", os);
1544 
1545  os.check(FUNCTION_NAME);
1546  return os;
1547 }
1548 
1549 
1550 template<class Type, template<class> class PatchField, class GeoMesh>
1551 Foam::Ostream& Foam::operator<<
1552 (
1553  Ostream& os,
1555 )
1556 {
1557  os << tgf();
1558  tgf.clear();
1559  return os;
1560 }
1561 
1562 
1563 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1564 
1565 #undef checkField
1566 
1567 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1568 
1569 #include "GeometricFieldNew.C"
1570 #include "GeometricFieldFunctions.C"
1571 
1572 // ************************************************************************* //
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place)
faceListList boundary
void clamp_range(const dimensioned< MinMax< Type >> &range)
Clamp field values (in-place) to the specified range.
#define COMPUTED_ASSIGNMENT(TYPE, op)
const Type & value() const noexcept
Return const reference to value.
dictionary dict
const labelList patchIDs(pbm.patchSet(polyPatchNames, false, true).sortedToc())
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
uint8_t direction
Definition: direction.H:48
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
Y [inertIndex] clamp_min(0)
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Field< Type >::cmptType cmptType
Component type of the field elements.
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
GeometricField(const IOobject &io, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Construct given IOobject, mesh, dimensions and patch type.
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
orientedType oriented() const noexcept
Return oriented type.
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1186
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.
word select(bool final) const
Select the final iteration parameters if final is true by returning the field name + "Final" otherwis...
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:414
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:107
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
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:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
void storeOldTime() const
Store the old-time field.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
bool local
Definition: EEqn.H:20
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
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
void clamp_max(const Type &upper)
Impose upper (ceiling) clamp on the field values (in-place)
List< word > wordList
List of word.
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.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1170
#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
meshDefDict readIfPresent("polyMeshPatches", polyPatchNames)
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
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
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:171
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.
#define checkField(fld1, fld2, op)
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.