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 "meshState.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"); // Includes size check
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  IOobjectOption::MUST_READ,
89  IOobjectOption::NO_WRITE,
90  IOobjectOption::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  << "The readOption MUST_READ or READ_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  readOldTimeIfPresent();
122 
123  return true;
124  }
125 
126  return false;
127 }
128 
129 
130 template<class Type, template<class> class PatchField, class GeoMesh>
132 {
133  // Read the old time field if present
134  IOobject field0
135  (
136  this->name() + "_0",
137  this->time().timeName(),
138  this->db(),
139  IOobjectOption::LAZY_READ,
140  IOobjectOption::AUTO_WRITE,
141  this->registerObject()
142  );
143 
144  if
145  (
146  field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
147  (
148  true
149  )
150  )
151  {
153  << "Reading old time level for field" << nl << this->info() << endl;
154 
155  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
156  (
157  field0,
158  this->mesh()
159  );
160 
161  // Ensure the old time field oriented flag is set to the parent's state
162  // Note: required for backwards compatibility in case of restarting from
163  // an old run where the oriented state may not have been set
164  field0Ptr_->oriented() = this->oriented();
165 
166  field0Ptr_->timeIndex_ = timeIndex_ - 1;
167 
168  if (!field0Ptr_->readOldTimeIfPresent())
169  {
170  field0Ptr_->oldTime();
171  }
172 
173  return true;
174  }
175 
176  return false;
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
181 
182 template<class Type, template<class> class PatchField, class GeoMesh>
184 (
185  const IOobject& io,
186  const Mesh& mesh,
187  const dimensionSet& dims,
188  const word& patchFieldType
189 )
190 :
191  Internal(io, mesh, dims, false),
192  timeIndex_(this->time().timeIndex()),
193  field0Ptr_(nullptr),
194  fieldPrevIterPtr_(nullptr),
195  boundaryField_(mesh.boundary(), *this, patchFieldType)
196 {
198  << "Creating" << nl << this->info() << endl;
200  readIfPresent();
201 }
202 
203 
204 template<class Type, template<class> class PatchField, class GeoMesh>
206 (
207  const IOobject& io,
208  const Mesh& mesh,
209  const dimensionSet& dims,
210  const wordList& patchFieldTypes,
211  const wordList& actualPatchTypes
212 )
213 :
214  Internal(io, mesh, dims, false),
215  timeIndex_(this->time().timeIndex()),
216  field0Ptr_(nullptr),
217  fieldPrevIterPtr_(nullptr),
218  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
219 {
221  << "Creating" << nl << this->info() << endl;
223  readIfPresent();
224 }
225 
226 
227 template<class Type, template<class> class PatchField, class GeoMesh>
229 (
230  const IOobject& io,
231  const Mesh& mesh,
232  const Type& value,
233  const dimensionSet& dims,
234  const word& patchFieldType
235 )
236 :
237  Internal(io, mesh, value, dims, false),
238  timeIndex_(this->time().timeIndex()),
239  field0Ptr_(nullptr),
240  fieldPrevIterPtr_(nullptr),
241  boundaryField_(mesh.boundary(), *this, patchFieldType)
242 {
244  << "Creating" << nl << this->info() << endl;
245 
246  boundaryField_ == value;
248  readIfPresent();
249 }
250 
251 
252 template<class Type, template<class> class PatchField, class GeoMesh>
254 (
255  const IOobject& io,
256  const Mesh& mesh,
257  const Type& value,
258  const dimensionSet& dims,
259  const wordList& patchFieldTypes,
260  const wordList& actualPatchTypes
261 )
262 :
263  Internal(io, mesh, value, dims, false),
264  timeIndex_(this->time().timeIndex()),
265  field0Ptr_(nullptr),
266  fieldPrevIterPtr_(nullptr),
267  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
268 {
270  << "Creating" << nl << this->info() << endl;
271 
272  boundaryField_ == value;
274  readIfPresent();
275 }
276 
277 
278 template<class Type, template<class> class PatchField, class GeoMesh>
280 (
281  const IOobject& io,
282  const Mesh& mesh,
283  const dimensioned<Type>& dt,
284  const word& patchFieldType
285 )
286 :
287  GeometricField<Type, PatchField, GeoMesh>
288  (
289  io,
290  mesh,
291  dt.value(),
292  dt.dimensions(),
293  patchFieldType
294  )
295 {}
296 
297 
298 template<class Type, template<class> class PatchField, class GeoMesh>
300 (
301  const IOobject& io,
302  const Mesh& mesh,
303  const dimensioned<Type>& dt,
304  const wordList& patchFieldTypes,
305  const wordList& actualPatchTypes
306 )
307 :
308  GeometricField<Type, PatchField, GeoMesh>
309  (
310  io,
311  mesh,
312  dt.value(),
313  dt.dimensions(),
314  patchFieldTypes,
315  actualPatchTypes
316  )
317 {}
318 
319 
320 template<class Type, template<class> class PatchField, class GeoMesh>
322 (
323  const IOobject& io,
324  const Internal& diField,
325  const PtrList<PatchField<Type>>& ptfl
326 )
327 :
328  Internal(io, diField),
329  timeIndex_(this->time().timeIndex()),
330  field0Ptr_(nullptr),
331  fieldPrevIterPtr_(nullptr),
332  boundaryField_(this->mesh().boundary(), *this, ptfl)
333 {
335  << "Copy construct from components" << nl << this->info() << endl;
337  readIfPresent();
338 }
339 
340 
341 template<class Type, template<class> class PatchField, class GeoMesh>
343 (
344  const IOobject& io,
345  Internal&& diField,
346  const PtrList<PatchField<Type>>& ptfl
347 )
348 :
349  Internal(io, std::move(diField)),
350  timeIndex_(this->time().timeIndex()),
351  field0Ptr_(nullptr),
352  fieldPrevIterPtr_(nullptr),
353  boundaryField_(this->mesh().boundary(), *this, ptfl)
354 {
356  << "Move construct from components" << nl << this->info() << endl;
358  readIfPresent();
359 }
360 
361 
362 template<class Type, template<class> class PatchField, class GeoMesh>
364 (
365  const IOobject& io,
366  const tmp<Internal>& tfield,
367  const PtrList<PatchField<Type>>& ptfl
368 )
369 :
370  Internal(io, tfield),
371  timeIndex_(this->time().timeIndex()),
372  field0Ptr_(nullptr),
373  fieldPrevIterPtr_(nullptr),
374  boundaryField_(this->mesh().boundary(), *this, ptfl)
375 {
377  << "Construct from tmp internalField" << nl << this->info() << endl;
379  readIfPresent();
380 }
381 
382 
383 template<class Type, template<class> class PatchField, class GeoMesh>
385 (
386  const Internal& diField,
387  const PtrList<PatchField<Type>>& ptfl
388 )
389 :
390  Internal(diField),
391  timeIndex_(this->time().timeIndex()),
392  field0Ptr_(nullptr),
393  fieldPrevIterPtr_(nullptr),
394  boundaryField_(this->mesh().boundary(), *this, ptfl)
395 {
397  << "Copy construct from components" << nl << this->info() << endl;
399  readIfPresent();
400 }
401 
402 
403 template<class Type, template<class> class PatchField, class GeoMesh>
405 (
406  Internal&& diField,
407  const PtrList<PatchField<Type>>& ptfl
408 )
409 :
410  Internal(std::move(diField)),
411  timeIndex_(this->time().timeIndex()),
412  field0Ptr_(nullptr),
413  fieldPrevIterPtr_(nullptr),
414  boundaryField_(this->mesh().boundary(), *this, ptfl)
415 {
417  << "Move construct from components" << nl << this->info() << endl;
419  readIfPresent();
420 }
421 
422 
423 template<class Type, template<class> class PatchField, class GeoMesh>
425 (
426  const IOobject& io,
427  const Mesh& mesh,
428  const dimensionSet& dims,
429  const Field<Type>& iField,
430  const word& patchFieldType
431 )
432 :
433  Internal(io, mesh, dims, iField),
434  timeIndex_(this->time().timeIndex()),
435  field0Ptr_(nullptr),
436  fieldPrevIterPtr_(nullptr),
437  boundaryField_(mesh.boundary(), *this, patchFieldType)
438 {
440  << "Copy construct from internal field" << nl << this->info() << endl;
442  readIfPresent();
443 }
444 
445 
446 template<class Type, template<class> class PatchField, class GeoMesh>
448 (
449  const IOobject& io,
450  const Mesh& mesh,
451  const dimensionSet& dims,
452  Field<Type>&& iField,
453  const word& patchFieldType
454 )
455 :
456  Internal(io, mesh, dims, std::move(iField)),
457  timeIndex_(this->time().timeIndex()),
458  field0Ptr_(nullptr),
459  fieldPrevIterPtr_(nullptr),
460  boundaryField_(mesh.boundary(), *this, patchFieldType)
461 {
463  << "Move construct from internal field" << nl << this->info() << endl;
465  readIfPresent();
466 }
467 
468 
469 template<class Type, template<class> class PatchField, class GeoMesh>
471 (
472  const IOobject& io,
473  const Mesh& mesh,
474  const dimensionSet& dims,
475  const Field<Type>& iField,
476  const PtrList<PatchField<Type>>& ptfl
477 )
478 :
479  Internal(io, mesh, dims, iField),
480  timeIndex_(this->time().timeIndex()),
481  field0Ptr_(nullptr),
482  fieldPrevIterPtr_(nullptr),
483  boundaryField_(mesh.boundary(), *this, ptfl)
484 {
486  << "Copy construct from components" << nl << this->info() << endl;
488  readIfPresent();
489 }
490 
491 
492 template<class Type, template<class> class PatchField, class GeoMesh>
494 (
495  const IOobject& io,
496  const Mesh& mesh,
497  const dimensionSet& dims,
498  Field<Type>&& iField,
499  const PtrList<PatchField<Type>>& ptfl
500 )
501 :
502  Internal(io, mesh, dims, std::move(iField)),
503  timeIndex_(this->time().timeIndex()),
504  field0Ptr_(nullptr),
505  fieldPrevIterPtr_(nullptr),
506  boundaryField_(mesh.boundary(), *this, ptfl)
507 {
509  << "Move construct from components" << nl << this->info() << endl;
511  readIfPresent();
512 }
513 
514 
515 template<class Type, template<class> class PatchField, class GeoMesh>
517 (
518  const IOobject& io,
519  const Mesh& mesh,
520  const dimensionSet& dims,
521  const tmp<Field<Type>>& tfield,
522  const PtrList<PatchField<Type>>& ptfl
523 )
524 :
525  Internal(io, mesh, dims, tfield),
526  timeIndex_(this->time().timeIndex()),
527  field0Ptr_(nullptr),
528  fieldPrevIterPtr_(nullptr),
529  boundaryField_(mesh.boundary(), *this, ptfl)
530 {
532  << "Construct from tmp internalField" << nl << this->info() << endl;
534  readIfPresent();
535 }
536 
537 
538 template<class Type, template<class> class PatchField, class GeoMesh>
540 (
541  const IOobject& io,
542  const Mesh& mesh,
543  const bool readOldTime
544 )
545 :
546  Internal(io, mesh, dimless, false),
547  timeIndex_(this->time().timeIndex()),
548  field0Ptr_(nullptr),
549  fieldPrevIterPtr_(nullptr),
550  boundaryField_(mesh.boundary())
551 {
553  << "Read construct" << nl << this->info() << endl;
554 
555  if (!this->isAnyRead())
556  {
557  // Do not warn about LAZY_READ since we may have already checked
558  // the IOobject before calling.
560  << "Had readOption NO_READ for field "
561  << this->name() << ", but constructor always reads field!"
562  << endl;
563  }
564 
565  readFields();
566 
567  if (readOldTime)
568  {
569  readOldTimeIfPresent();
570  }
571 
573  << "Finishing read-construction" << nl << this->info() << endl;
574 }
575 
576 
577 template<class Type, template<class> class PatchField, class GeoMesh>
579 (
580  const IOobject& io,
581  const Mesh& mesh,
582  const dictionary& dict
583 )
584 :
585  Internal(io, mesh, dimless, false),
586  timeIndex_(this->time().timeIndex()),
587  field0Ptr_(nullptr),
588  fieldPrevIterPtr_(nullptr),
589  boundaryField_(mesh.boundary())
590 {
591  readFields(dict);
592 
594  << "Finishing dictionary-construct" << nl << this->info() << endl;
595 }
596 
597 
598 template<class Type, template<class> class PatchField, class GeoMesh>
600 (
602 )
603 :
604  Internal(gf),
605  timeIndex_(gf.timeIndex()),
606  field0Ptr_(nullptr),
607  fieldPrevIterPtr_(nullptr),
608  boundaryField_(*this, gf.boundaryField_)
609 {
611  << "Copy construct" << nl << this->info() << endl;
612 
613  if (gf.field0Ptr_)
614  {
616  (
617  *gf.field0Ptr_
618  );
619  }
621  this->writeOpt(IOobjectOption::NO_WRITE);
622 }
623 
624 
625 template<class Type, template<class> class PatchField, class GeoMesh>
627 (
629 )
630 :
631  Internal(tgf.constCast(), tgf.movable()),
632  timeIndex_(tgf().timeIndex()),
633  field0Ptr_(nullptr),
634  fieldPrevIterPtr_(nullptr),
635  boundaryField_(*this, tgf().boundaryField_)
636 {
638  << "Constructing from tmp" << nl << this->info() << endl;
639 
640  this->writeOpt(IOobjectOption::NO_WRITE);
642  tgf.clear();
643 }
644 
645 
646 template<class Type, template<class> class PatchField, class GeoMesh>
648 (
649  const IOobject& io,
651 )
652 :
653  Internal(io, gf),
654  timeIndex_(gf.timeIndex()),
655  field0Ptr_(nullptr),
656  fieldPrevIterPtr_(nullptr),
657  boundaryField_(*this, gf.boundaryField_)
658 {
660  << "Copy construct, resetting IO params" << nl
661  << this->info() << endl;
662 
663  if (!readIfPresent() && gf.field0Ptr_)
664  {
666  (
667  io.name() + "_0",
668  *gf.field0Ptr_
669  );
670  }
671 }
672 
673 
674 template<class Type, template<class> class PatchField, class GeoMesh>
676 (
677  const IOobject& io,
678  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
679 )
680 :
681  Internal(io, tgf.constCast(), tgf.movable()),
682  timeIndex_(tgf().timeIndex()),
683  field0Ptr_(nullptr),
684  fieldPrevIterPtr_(nullptr),
685  boundaryField_(*this, tgf().boundaryField_)
686 {
688  << "Constructing from tmp resetting IO params" << nl
689  << this->info() << endl;
690 
691  tgf.clear();
693  readIfPresent();
694 }
695 
696 
697 template<class Type, template<class> class PatchField, class GeoMesh>
699 (
700  const word& newName,
702 )
703 :
704  Internal(newName, gf),
705  timeIndex_(gf.timeIndex()),
706  field0Ptr_(nullptr),
707  fieldPrevIterPtr_(nullptr),
708  boundaryField_(*this, gf.boundaryField_)
709 {
711  << "Copy construct, resetting name" << nl
712  << this->info() << endl;
713 
714  if (!readIfPresent() && gf.field0Ptr_)
715  {
717  (
718  newName + "_0",
719  *gf.field0Ptr_
720  );
721  }
722 }
723 
724 
725 template<class Type, template<class> class PatchField, class GeoMesh>
727 (
728  const word& newName,
729  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
730 )
731 :
732  Internal(newName, tgf.constCast(), tgf.movable()),
733  timeIndex_(tgf().timeIndex()),
734  field0Ptr_(nullptr),
735  fieldPrevIterPtr_(nullptr),
736  boundaryField_(*this, tgf().boundaryField_)
737 {
739  << "Constructing from tmp resetting name" << nl
740  << this->info() << endl;
742  tgf.clear();
743 }
744 
745 
746 template<class Type, template<class> class PatchField, class GeoMesh>
748 (
749  const IOobject& io,
751  const word& patchFieldType
752 )
753 :
754  Internal(io, gf),
755  timeIndex_(gf.timeIndex()),
756  field0Ptr_(nullptr),
757  fieldPrevIterPtr_(nullptr),
758  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
759 {
761  << "Copy construct, resetting IO params" << nl
762  << this->info() << endl;
763 
764  boundaryField_ == gf.boundaryField_;
765 
766  if (!readIfPresent() && gf.field0Ptr_)
767  {
769  (
770  io.name() + "_0",
771  *gf.field0Ptr_
772  );
773  }
774 }
775 
776 
777 template<class Type, template<class> class PatchField, class GeoMesh>
779 (
780  const IOobject& io,
781  const GeometricField<Type, PatchField, GeoMesh>& gf,
782  const wordList& patchFieldTypes,
783  const wordList& actualPatchTypes
784 )
785 :
786  Internal(io, gf),
787  timeIndex_(gf.timeIndex()),
788  field0Ptr_(nullptr),
789  fieldPrevIterPtr_(nullptr),
790  boundaryField_
791  (
792  this->mesh().boundary(),
793  *this,
794  patchFieldTypes,
795  actualPatchTypes
796  )
797 {
799  << "Copy construct, resetting IO params and patch types" << nl
800  << this->info() << endl;
801 
802  boundaryField_ == gf.boundaryField_;
803 
804  if (!readIfPresent() && gf.field0Ptr_)
805  {
806  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
807  (
808  io.name() + "_0",
809  *gf.field0Ptr_
810  );
811  }
812 }
813 
814 
815 template<class Type, template<class> class PatchField, class GeoMesh>
817 (
818  const IOobject& io,
819  const GeometricField<Type, PatchField, GeoMesh>& gf,
820  const labelList& patchIDs,
821  const word& patchFieldType
822 )
823 :
824  Internal(io, gf),
825  timeIndex_(gf.timeIndex()),
826  field0Ptr_(nullptr),
827  fieldPrevIterPtr_(nullptr),
828  boundaryField_(*this, gf.boundaryField_, patchIDs, patchFieldType)
829 {
831  << "Copy construct, resetting IO params and setting patchFieldType "
832  << "for patchIDs" << nl
833  << this->info() << endl;
834 
835  if (!readIfPresent() && gf.field0Ptr_)
836  {
837  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
838  (
839  io.name() + "_0",
840  *gf.field0Ptr_
841  );
842  }
843 }
844 
845 
846 template<class Type, template<class> class PatchField, class GeoMesh>
848 (
849  const IOobject& io,
850  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf,
851  const wordList& patchFieldTypes,
852  const wordList& actualPatchTypes
853 )
854 :
855  Internal(io, tgf.constCast(), tgf.movable()),
856  timeIndex_(tgf().timeIndex()),
857  field0Ptr_(nullptr),
858  fieldPrevIterPtr_(nullptr),
859  boundaryField_
860  (
861  this->mesh().boundary(),
862  *this,
863  patchFieldTypes,
864  actualPatchTypes
865  )
866 {
868  << "Constructing from tmp resetting IO params and patch types" << nl
869  << this->info() << endl;
870 
871  boundaryField_ == tgf().boundaryField_;
873  tgf.clear();
874 }
875 
876 
877 template<class Type, template<class> class PatchField, class GeoMesh>
880 {
882 }
883 
884 
885 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
886 
887 template<class Type, template<class> class PatchField, class GeoMesh>
889 {
890  /*
891  if (debug)
892  {
893  // Problem: temporary fields might have their internal field
894  // already stolen so boundary fields will not be able to access the
895  // internal field anymore
896  boundaryField_.check();
897  }
898  */
899 
900  deleteDemandDrivenData(field0Ptr_);
901  deleteDemandDrivenData(fieldPrevIterPtr_);
902 
903  // FUTURE: register cache field info
904  // // this->db().cacheTemporaryObject(*this);
905 
906  clearOldTimes();
907 }
908 
910 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
911 
912 template<class Type, template<class> class PatchField, class GeoMesh>
913 typename
916 (
917  const bool updateAccessTime
918 )
919 {
920  if (updateAccessTime)
921  {
922  this->setUpToDate();
923  storeOldTimes();
924  }
925  return *this;
926 }
927 
928 
929 template<class Type, template<class> class PatchField, class GeoMesh>
930 typename
933 (
934  const bool updateAccessTime
935 )
936 {
937  if (updateAccessTime)
938  {
939  this->setUpToDate();
940  storeOldTimes();
941  }
942  return *this;
943 }
944 
945 
946 template<class Type, template<class> class PatchField, class GeoMesh>
947 typename
950 (
951  const bool updateAccessTime
952 )
953 {
954  if (updateAccessTime)
955  {
956  this->setUpToDate();
957  storeOldTimes();
958  }
959  return boundaryField_;
960 }
961 
962 
963 template<class Type, template<class> class PatchField, class GeoMesh>
964 Foam::label
966 {
967  return (field0Ptr_ ? (field0Ptr_->nOldTimes() + 1) : 0);
968 }
969 
970 
971 template<class Type, template<class> class PatchField, class GeoMesh>
973 {
974  if
975  (
976  field0Ptr_
977  && timeIndex_ != this->time().timeIndex()
978  && !this->name().ends_with("_0")
979  )
980  {
981  storeOldTime();
982  timeIndex_ = this->time().timeIndex();
983  }
985  // Correct time index
986  //timeIndex_ = this->time().timeIndex();
987 }
988 
989 
990 template<class Type, template<class> class PatchField, class GeoMesh>
992 {
993  if (field0Ptr_)
994  {
995  field0Ptr_->storeOldTime();
996 
998  << "Storing old time field for field" << nl << this->info() << endl;
999 
1000  *field0Ptr_ == *this;
1001  field0Ptr_->timeIndex_ = timeIndex_;
1002 
1003  if (field0Ptr_->field0Ptr_)
1004  {
1005  field0Ptr_->writeOpt(this->writeOpt());
1006  }
1007  }
1008 }
1009 
1010 
1011 template<class Type, template<class> class PatchField, class GeoMesh>
1014 {
1015  if (!field0Ptr_)
1016  {
1017  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
1018  (
1019  IOobject
1020  (
1021  this->name() + "_0",
1022  this->time().timeName(),
1023  this->db(),
1024  IOobjectOption::NO_READ,
1025  IOobjectOption::NO_WRITE,
1026  this->registerObject()
1027  ),
1028  *this
1029  );
1030 
1031  if (debug)
1032  {
1034  << "created old time field " << field0Ptr_->info() << endl;
1035 
1036  if (debug&2)
1037  {
1038  error::printStack(Info);
1039  }
1040  }
1041  }
1042  else
1043  {
1044  storeOldTimes();
1045  }
1047  return *field0Ptr_;
1048 }
1049 
1050 
1051 template<class Type, template<class> class PatchField, class GeoMesh>
1054 {
1055  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1057 
1058  return *field0Ptr_;
1059 }
1060 
1061 
1062 template<class Type, template<class> class PatchField, class GeoMesh>
1064 {
1065  if (!fieldPrevIterPtr_)
1066  {
1068  << "Allocating previous iteration field" << nl
1069  << this->info() << endl;
1070 
1071  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1072  (
1073  this->name() + "PrevIter",
1074  *this
1075  );
1076  }
1077  else
1078  {
1079  *fieldPrevIterPtr_ == *this;
1080  }
1081 }
1082 
1083 
1084 template<class Type, template<class> class PatchField, class GeoMesh>
1087 {
1088  if (!fieldPrevIterPtr_)
1089  {
1091  << "previous iteration field" << endl << this->info() << endl
1092  << " not stored."
1093  << " Use field.storePrevIter() at start of iteration."
1094  << abort(FatalError);
1095  }
1096 
1097  return *fieldPrevIterPtr_;
1098 }
1099 
1100 
1101 template<class Type, template<class> class PatchField, class GeoMesh>
1103 {
1105  deleteDemandDrivenData(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>
1122 {
1123  this->setUpToDate();
1124  storeOldTimes();
1125  boundaryField_.evaluateLocal();
1126 }
1127 
1128 
1129 template<class Type, template<class> class PatchField, class GeoMesh>
1131 {
1132  // Search all boundary conditions, if any are
1133  // fixed-value or mixed (Robin) do not set reference level for solution.
1134 
1135  bool needRef = true;
1136 
1137  for (const auto& pf : boundaryField_)
1138  {
1139  if (pf.fixesValue())
1140  {
1141  needRef = false;
1142  break;
1143  }
1144  }
1145 
1146  return returnReduceAnd(needRef);
1147 }
1148 
1149 
1150 template<class Type, template<class> class PatchField, class GeoMesh>
1152 {
1154  << "Relaxing" << nl << this->info() << " by " << alpha << endl;
1155 
1156  operator==(prevIter() + alpha*(*this - prevIter()));
1157 }
1158 
1159 
1160 template<class Type, template<class> class PatchField, class GeoMesh>
1162 {
1163  word name = this->name();
1164 
1165  if (this->mesh().data().isFinalIteration())
1166  {
1167  name += "Final";
1168  }
1169 
1170  scalar relaxCoeff = 1;
1171 
1172  if (this->mesh().relaxField(name, relaxCoeff))
1173  {
1174  relax(relaxCoeff);
1175  }
1176 }
1177 
1178 
1179 template<class Type, template<class> class PatchField, class GeoMesh>
1181 (
1182  bool final
1183 ) const
1184 {
1185  if (final)
1186  {
1187  return this->name() + "Final";
1188  }
1190  return this->name();
1191 }
1192 
1193 
1194 template<class Type, template<class> class PatchField, class GeoMesh>
1196 (
1197  Ostream& os
1198 ) const
1199 {
1200  MinMax<Type> range = Foam::minMax(*this).value();
1201 
1202  os << "min/max(" << this->name() << ") = "
1203  << range.min() << ", " << range.max() << endl;
1204 }
1205 
1206 
1207 template<class Type, template<class> class PatchField, class GeoMesh>
1209 writeData(Ostream& os) const
1210 {
1211  os << *this;
1212  return os.good();
1214 
1215 
1216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1217 
1218 template<class Type, template<class> class PatchField, class GeoMesh>
1221 {
1223  (
1224  this->name() + ".T()",
1225  this->mesh(),
1226  this->dimensions()
1227  );
1228 
1229  Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1230  Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1231 
1232  return tresult;
1233 }
1234 
1235 
1236 template<class Type, template<class> class PatchField, class GeoMesh>
1237 Foam::tmp
1238 <
1240  <
1242  PatchField,
1243  GeoMesh
1244  >
1245 >
1247 (
1248  const direction d
1249 ) const
1250 {
1252  (
1253  this->name() + ".component(" + Foam::name(d) + ')',
1254  this->mesh(),
1255  this->dimensions()
1256  );
1257 
1258  Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1259  Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1260 
1261  return tresult;
1262 }
1263 
1264 
1265 template<class Type, template<class> class PatchField, class GeoMesh>
1267 (
1268  const direction d,
1269  const GeometricField
1270  <
1271  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1272  PatchField,
1273  GeoMesh
1274  >& gcf
1275 )
1276 {
1277  primitiveFieldRef().replace(d, gcf.primitiveField());
1278  boundaryFieldRef().replace(d, gcf.boundaryField());
1279 }
1280 
1281 
1282 template<class Type, template<class> class PatchField, class GeoMesh>
1284 (
1285  const direction d,
1286  const dimensioned<cmptType>& ds
1287 )
1288 {
1289  primitiveFieldRef().replace(d, ds.value());
1290  boundaryFieldRef().replace(d, ds.value());
1291 }
1292 
1293 
1294 template<class Type, template<class> class PatchField, class GeoMesh>
1296 (
1297  const Type& lower
1298 )
1299 {
1300  primitiveFieldRef().clamp_min(lower);
1301  boundaryFieldRef().clamp_min(lower);
1302 }
1303 
1304 
1305 template<class Type, template<class> class PatchField, class GeoMesh>
1307 (
1308  const Type& upper
1309 )
1310 {
1311  primitiveFieldRef().clamp_max(upper);
1312  boundaryFieldRef().clamp_max(upper);
1313 }
1314 
1315 
1316 template<class Type, template<class> class PatchField, class GeoMesh>
1318 (
1319  const dimensioned<Type>& lower
1320 )
1322  this->clamp_min(lower.value());
1323 }
1324 
1325 
1326 template<class Type, template<class> class PatchField, class GeoMesh>
1328 (
1329  const dimensioned<Type>& upper
1330 )
1332  this->clamp_max(upper.value());
1333 }
1334 
1335 
1336 template<class Type, template<class> class PatchField, class GeoMesh>
1338 (
1339  const Type& lower,
1340  const Type& upper
1341 )
1342 {
1343  primitiveFieldRef().clamp_range(lower, upper);
1344  boundaryFieldRef().clamp_range(lower, upper);
1345 }
1346 
1347 
1348 template<class Type, template<class> class PatchField, class GeoMesh>
1350 (
1351  const MinMax<Type>& range
1352 )
1353 {
1354  primitiveFieldRef().clamp_range(range.min(), range.max());
1355  boundaryFieldRef().clamp_range(range.min(), range.max());
1356 }
1357 
1358 
1359 template<class Type, template<class> class PatchField, class GeoMesh>
1361 (
1362  const dimensioned<Type>& lower,
1363  const dimensioned<Type>& upper
1364 )
1366  this->clamp_range(lower.value(), upper.value());
1367 }
1368 
1369 
1370 template<class Type, template<class> class PatchField, class GeoMesh>
1372 (
1375 {
1376  this->clamp_range(range.value());
1377 }
1378 
1379 
1380 template<class Type, template<class> class PatchField, class GeoMesh>
1383  primitiveFieldRef().negate();
1384  boundaryFieldRef().negate();
1385 }
1386 
1387 
1388 template<class Type, template<class> class PatchField, class GeoMesh>
1390 {
1391  primitiveFieldRef().normalise();
1392  boundaryFieldRef().normalise();
1394 
1395 
1396 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1397 
1398 template<class Type, template<class> class PatchField, class GeoMesh>
1400 (
1402 )
1403 {
1404  if (this == &gf)
1405  {
1406  return; // Self-assignment is a no-op
1407  }
1408 
1409  checkField(*this, gf, "=");
1410 
1411  // Only assign field contents not ID
1412 
1413  internalFieldRef() = gf.internalField();
1415 }
1416 
1417 
1418 template<class Type, template<class> class PatchField, class GeoMesh>
1420 (
1422 )
1423 {
1424  const auto& gf = tgf();
1425 
1426  if (this == &gf)
1427  {
1428  return; // Self-assignment is a no-op
1429  }
1430 
1431  checkField(*this, gf, "=");
1432 
1433  // Only assign field contents not ID
1434 
1435  this->dimensions() = gf.dimensions();
1436  this->oriented() = gf.oriented();
1437 
1438  if (tgf.movable())
1439  {
1440  // Transfer storage from the tmp
1441  primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1442  }
1443  else
1444  {
1446  }
1447 
1450  tgf.clear();
1451 }
1452 
1453 
1454 template<class Type, template<class> class PatchField, class GeoMesh>
1456 (
1457  const dimensioned<Type>& dt
1458 )
1459 {
1460  internalFieldRef() = dt;
1461  boundaryFieldRef() = dt.value();
1462 }
1463 
1464 
1465 template<class Type, template<class> class PatchField, class GeoMesh>
1467 (
1469 )
1470 {
1471  const auto& gf = tgf();
1472 
1473  checkField(*this, gf, "==");
1474 
1475  // Only assign field contents not ID
1476 
1477  internalFieldRef() = gf.internalField();
1478  boundaryFieldRef() == gf.boundaryField();
1480  tgf.clear();
1481 }
1482 
1483 
1484 template<class Type, template<class> class PatchField, class GeoMesh>
1486 (
1487  const dimensioned<Type>& dt
1489 {
1490  internalFieldRef() = dt;
1491  boundaryFieldRef() == dt.value();
1492 }
1493 
1494 
1495 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1496  \
1497 template<class Type, template<class> class PatchField, class GeoMesh> \
1498 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1499 ( \
1500  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1501 ) \
1502 { \
1503  checkField(*this, gf, #op); \
1504  \
1505  internalFieldRef() op gf.internalField(); \
1506  boundaryFieldRef() op gf.boundaryField(); \
1507 } \
1508  \
1509 template<class Type, template<class> class PatchField, class GeoMesh> \
1510 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1511 ( \
1512  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1513 ) \
1514 { \
1515  operator op(tgf()); \
1516  tgf.clear(); \
1517 } \
1518  \
1519 template<class Type, template<class> class PatchField, class GeoMesh> \
1520 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1521 ( \
1522  const dimensioned<TYPE>& dt \
1523 ) \
1524 { \
1525  internalFieldRef() op dt; \
1526  boundaryFieldRef() op dt.value(); \
1527 }
1528 
1529 COMPUTED_ASSIGNMENT(Type, +=)
1530 COMPUTED_ASSIGNMENT(Type, -=)
1531 COMPUTED_ASSIGNMENT(scalar, *=)
1532 COMPUTED_ASSIGNMENT(scalar, /=)
1533 
1534 #undef COMPUTED_ASSIGNMENT
1535 
1536 
1537 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1538 
1539 template<class Type, template<class> class PatchField, class GeoMesh>
1540 Foam::Ostream& Foam::operator<<
1541 (
1542  Ostream& os,
1544 )
1545 {
1546  gf.internalField().writeData(os, "internalField");
1547  os << nl;
1548  gf.boundaryField().writeEntry("boundaryField", os);
1549 
1550  os.check(FUNCTION_NAME);
1551  return os;
1552 }
1553 
1554 
1555 template<class Type, template<class> class PatchField, class GeoMesh>
1556 Foam::Ostream& Foam::operator<<
1557 (
1558  Ostream& os,
1560 )
1561 {
1562  os << tgf();
1563  tgf.clear();
1564  return os;
1565 }
1566 
1567 
1568 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1569 
1570 #undef checkField
1571 
1572 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1573 
1574 #include "GeometricFieldNew.C"
1575 #include "GeometricFieldFunctions.C"
1576 
1577 // ************************************************************************* //
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place)
faceListList boundary
const labelList patchIDs(pbm.indices(polyPatchNames, true))
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 GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
uint8_t direction
Definition: direction.H:46
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Y [inertIndex] clamp_min(0)
void clearOldTimes()
Remove old-time and prev-iter fields.
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:50
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:531
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
label nOldTimes() const noexcept
The number of old time fields stored.
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:1187
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Generic dimensioned Type class.
const dimensionSet dimless
Dimensionless.
bool writeData(Ostream &os) const
The writeData function (required by regIOobject), calls operator<<.
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:421
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Generic templated field type.
Definition: Field.H:62
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
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:56
const direction noexcept
Definition: Scalar.H:258
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:59
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:1171
#define WarningInFunction
Report a warning using Foam::Warning.
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.
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:172
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).
void correctLocalBoundaryConditions()
Correct boundary conditions after a purely local operation. Is.
label timeIndex
Definition: getTimeIndex.H:24
const dimensionSet & dimensions() const noexcept
Return dimensions.
#define checkField(fld1, fld2, op)
#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.