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-2024 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 "dictionary.H"
33 #include "meshState.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 // Check that both fields use the same mesh
38 #undef checkField
39 #define checkField(fld1, fld2, op) \
40 if (&(fld1).mesh() != &(fld2).mesh()) \
41 { \
42  FatalErrorInFunction \
43  << "Different mesh for fields " \
44  << (fld1).name() << " and " << (fld2).name() \
45  << " during operation " << op \
46  << abort(FatalError); \
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
51 
52 template<class Type, template<class> class PatchField, class GeoMesh>
54 (
55  const dictionary& dict
56 )
57 {
58  Internal::readField(dict, "internalField"); // Includes size check
59 
60  boundaryField_.readField(*this, dict.subDict("boundaryField"));
61 
62  Type refLevel;
63 
64  if (dict.readIfPresent("referenceLevel", refLevel))
65  {
66  Field<Type>::operator+=(refLevel);
67 
68  forAll(boundaryField_, patchi)
69  {
70  boundaryField_[patchi] == boundaryField_[patchi] + refLevel;
71  }
72  }
73 }
74 
75 
76 template<class Type, template<class> class PatchField, class GeoMesh>
78 {
80  (
81  localIOdictionary::readContents
82  (
83  IOobject
84  (
85  this->name(),
86  this->instance(),
87  this->local(),
88  this->db(),
89  IOobjectOption::MUST_READ,
90  IOobjectOption::NO_WRITE,
91  IOobjectOption::NO_REGISTER
92  ),
93  typeName
94  )
95  );
96 
97  this->close();
98 
100 }
101 
102 
103 template<class Type, template<class> class PatchField, class GeoMesh>
105 {
106  if (this->isReadRequired())
107  {
109  << "The readOption MUST_READ or READ_MODIFIED"
110  << " suggests that a read constructor for field " << this->name()
111  << " would be more appropriate." << endl;
112  }
113  else if
114  (
115  this->isReadOptional()
116  && this->template typeHeaderOk<this_type>(true) // checkType = true
117  )
118  {
119  readFields();
120  readOldTimeIfPresent();
121 
122  return true;
123  }
124 
125  return false;
126 }
127 
128 
129 template<class Type, template<class> class PatchField, class GeoMesh>
131 {
132  // Read the old time field if present
133  IOobject field0
134  (
135  this->name() + "_0",
136  this->time().timeName(),
137  this->db(),
138  IOobjectOption::LAZY_READ,
139  IOobjectOption::AUTO_WRITE,
140  this->registerObject()
141  );
142 
143  if
144  (
145  field0.template typeHeaderOk<this_type>(true) // checkType = true
146  )
147  {
149  << "Reading old time level for field" << nl << this->info() << endl;
150 
151  field0Ptr_ = std::make_unique<this_type>(field0, this->mesh());
152 
153  // Ensure the old time field oriented flag is set to the parent's state
154  // Note: required for backwards compatibility in case of restarting from
155  // an old run where the oriented state may not have been set
156  field0Ptr_->oriented() = this->oriented();
157 
158  field0Ptr_->timeIndex_ = timeIndex_ - 1;
159 
160  if (!field0Ptr_->readOldTimeIfPresent())
161  {
162  field0Ptr_->oldTime();
163  }
164 
165  return true;
166  }
167 
168  return false;
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
173 
174 template<class Type, template<class> class PatchField, class GeoMesh>
176 (
177  const IOobject& io,
178  const Mesh& mesh,
179  const dimensionSet& dims,
180  const word& patchFieldType
181 )
182 :
183  Internal(io, mesh, dims, false),
184  timeIndex_(this->time().timeIndex()),
185  boundaryField_(mesh.boundary(), *this, patchFieldType)
186 {
188  << "Creating" << nl << this->info() << endl;
190  readIfPresent();
191 }
192 
193 
194 template<class Type, template<class> class PatchField, class GeoMesh>
196 (
197  const IOobject& io,
198  const Mesh& mesh,
199  const dimensionSet& dims,
200  const wordList& patchFieldTypes,
201  const wordList& actualPatchTypes
202 )
203 :
204  Internal(io, mesh, dims, false),
205  timeIndex_(this->time().timeIndex()),
206  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
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 Type& value,
221  const dimensionSet& dims,
222  const word& patchFieldType
223 )
224 :
225  Internal(io, mesh, value, dims, false),
226  timeIndex_(this->time().timeIndex()),
227  boundaryField_(mesh.boundary(), *this, patchFieldType)
228 {
230  << "Creating" << nl << this->info() << endl;
231 
232  boundaryField_ == value;
234  readIfPresent();
235 }
236 
237 
238 template<class Type, template<class> class PatchField, class GeoMesh>
240 (
241  const IOobject& io,
242  const Mesh& mesh,
243  const Type& value,
244  const dimensionSet& dims,
245  const wordList& patchFieldTypes,
246  const wordList& actualPatchTypes
247 )
248 :
249  Internal(io, mesh, value, dims, false),
250  timeIndex_(this->time().timeIndex()),
251  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
252 {
254  << "Creating" << nl << this->info() << endl;
255 
256  boundaryField_ == value;
258  readIfPresent();
259 }
260 
261 
262 template<class Type, template<class> class PatchField, class GeoMesh>
264 (
265  const IOobject& io,
266  const Mesh& mesh,
267  const dimensioned<Type>& dt,
268  const word& patchFieldType
269 )
270 :
271  GeometricField<Type, PatchField, GeoMesh>
272  (
273  io,
274  mesh,
275  dt.value(),
276  dt.dimensions(),
277  patchFieldType
278  )
279 {}
280 
281 
282 template<class Type, template<class> class PatchField, class GeoMesh>
284 (
285  const IOobject& io,
286  const Mesh& mesh,
287  const dimensioned<Type>& dt,
288  const wordList& patchFieldTypes,
289  const wordList& actualPatchTypes
290 )
291 :
292  GeometricField<Type, PatchField, GeoMesh>
293  (
294  io,
295  mesh,
296  dt.value(),
297  dt.dimensions(),
298  patchFieldTypes,
299  actualPatchTypes
300  )
301 {}
302 
303 
304 template<class Type, template<class> class PatchField, class GeoMesh>
306 (
307  const IOobject& io,
308  const Internal& diField,
309  const PtrList<PatchField<Type>>& ptfl
310 )
311 :
312  Internal(io, diField),
313  timeIndex_(this->time().timeIndex()),
314  boundaryField_(this->mesh().boundary(), *this, ptfl)
315 {
317  << "Copy construct from components" << nl << this->info() << endl;
319  readIfPresent();
320 }
321 
322 
323 template<class Type, template<class> class PatchField, class GeoMesh>
325 (
326  const IOobject& io,
327  Internal&& diField,
328  const PtrList<PatchField<Type>>& ptfl
329 )
330 :
331  Internal(io, std::move(diField)),
332  timeIndex_(this->time().timeIndex()),
333  boundaryField_(this->mesh().boundary(), *this, ptfl)
334 {
336  << "Move construct from components" << nl << this->info() << endl;
338  readIfPresent();
339 }
340 
341 
342 template<class Type, template<class> class PatchField, class GeoMesh>
344 (
345  const IOobject& io,
346  const tmp<Internal>& tfield,
347  const PtrList<PatchField<Type>>& ptfl
348 )
349 :
350  Internal(io, tfield),
351  timeIndex_(this->time().timeIndex()),
352  boundaryField_(this->mesh().boundary(), *this, ptfl)
353 {
355  << "Construct from tmp internalField" << nl << this->info() << endl;
357  readIfPresent();
358 }
359 
360 
361 template<class Type, template<class> class PatchField, class GeoMesh>
363 (
364  const Internal& diField,
365  const PtrList<PatchField<Type>>& ptfl
366 )
367 :
368  Internal(diField),
369  timeIndex_(this->time().timeIndex()),
370  boundaryField_(this->mesh().boundary(), *this, ptfl)
371 {
373  << "Copy construct from components" << nl << this->info() << endl;
375  readIfPresent();
376 }
377 
378 
379 template<class Type, template<class> class PatchField, class GeoMesh>
381 (
382  Internal&& diField,
383  const PtrList<PatchField<Type>>& ptfl
384 )
385 :
386  Internal(std::move(diField)),
387  timeIndex_(this->time().timeIndex()),
388  boundaryField_(this->mesh().boundary(), *this, ptfl)
389 {
391  << "Move construct from components" << nl << this->info() << endl;
393  readIfPresent();
394 }
395 
396 
397 template<class Type, template<class> class PatchField, class GeoMesh>
399 (
400  const IOobject& io,
401  const Mesh& mesh,
402  const dimensionSet& dims,
403  const Field<Type>& iField,
404  const word& patchFieldType
405 )
406 :
407  Internal(io, mesh, dims, iField),
408  timeIndex_(this->time().timeIndex()),
409  boundaryField_(mesh.boundary(), *this, patchFieldType)
410 {
412  << "Copy construct from internal field" << nl << this->info() << endl;
414  readIfPresent();
415 }
416 
417 
418 template<class Type, template<class> class PatchField, class GeoMesh>
420 (
421  const IOobject& io,
422  const Mesh& mesh,
423  const dimensionSet& dims,
424  Field<Type>&& iField,
425  const word& patchFieldType
426 )
427 :
428  Internal(io, mesh, dims, std::move(iField)),
429  timeIndex_(this->time().timeIndex()),
430  boundaryField_(mesh.boundary(), *this, patchFieldType)
431 {
433  << "Move construct from internal field" << nl << this->info() << endl;
435  readIfPresent();
436 }
437 
438 
439 template<class Type, template<class> class PatchField, class GeoMesh>
441 (
442  const IOobject& io,
443  const Mesh& mesh,
444  const dimensionSet& dims,
445  const Field<Type>& iField,
446  const PtrList<PatchField<Type>>& ptfl
447 )
448 :
449  Internal(io, mesh, dims, iField),
450  timeIndex_(this->time().timeIndex()),
451  boundaryField_(mesh.boundary(), *this, ptfl)
452 {
454  << "Copy construct from components" << nl << this->info() << endl;
456  readIfPresent();
457 }
458 
459 
460 template<class Type, template<class> class PatchField, class GeoMesh>
462 (
463  const IOobject& io,
464  const Mesh& mesh,
465  const dimensionSet& dims,
466  Field<Type>&& iField,
467  const PtrList<PatchField<Type>>& ptfl
468 )
469 :
470  Internal(io, mesh, dims, std::move(iField)),
471  timeIndex_(this->time().timeIndex()),
472  boundaryField_(mesh.boundary(), *this, ptfl)
473 {
475  << "Move construct from components" << 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 tmp<Field<Type>>& tfield,
488  const PtrList<PatchField<Type>>& ptfl
489 )
490 :
491  Internal(io, mesh, dims, tfield),
492  timeIndex_(this->time().timeIndex()),
493  boundaryField_(mesh.boundary(), *this, ptfl)
494 {
496  << "Construct from tmp internalField" << nl << this->info() << endl;
498  readIfPresent();
499 }
500 
501 
502 template<class Type, template<class> class PatchField, class GeoMesh>
504 (
505  const IOobject& io,
506  const Mesh& mesh,
507  const bool readOldTime
508 )
509 :
510  Internal(io, mesh, dimless, false),
511  timeIndex_(this->time().timeIndex()),
512  boundaryField_(mesh.boundary())
513 {
515  << "Read construct" << nl << this->info() << endl;
516 
517  if (!this->isAnyRead())
518  {
519  // Do not warn about LAZY_READ since we may have already checked
520  // the IOobject before calling.
522  << "Had readOption NO_READ for field "
523  << this->name() << ", but constructor always reads field!"
524  << endl;
525  }
526 
527  readFields();
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  boundaryField_(mesh.boundary())
550 {
551  readFields(dict);
552 
554  << "Finishing dictionary-construct" << nl << this->info() << endl;
555 }
556 
557 
558 template<class Type, template<class> class PatchField, class GeoMesh>
560 (
562 )
563 :
564  Internal(gf),
565  timeIndex_(gf.timeIndex()),
566  boundaryField_(*this, gf.boundaryField_)
567 {
569  << "Copy construct" << nl << this->info() << endl;
570 
571  if (gf.field0Ptr_)
572  {
573  field0Ptr_ = std::make_unique<this_type>(*gf.field0Ptr_);
574  }
576  this->writeOpt(IOobjectOption::NO_WRITE);
577 }
578 
579 
580 template<class Type, template<class> class PatchField, class GeoMesh>
582 (
584 )
585 :
586  Internal(tgf.constCast(), tgf.movable()),
587  timeIndex_(tgf().timeIndex()),
588  boundaryField_(*this, tgf().boundaryField_)
589 {
591  << "Constructing from tmp" << nl << this->info() << endl;
592 
593  this->writeOpt(IOobjectOption::NO_WRITE);
595  tgf.clear();
596 }
597 
598 
599 template<class Type, template<class> class PatchField, class GeoMesh>
601 (
602  const IOobject& io,
604 )
605 :
606  Internal(io, gf),
607  timeIndex_(gf.timeIndex()),
608  boundaryField_(*this, gf.boundaryField_)
609 {
611  << "Copy construct, resetting IO params" << nl
612  << this->info() << endl;
613 
614  if (!readIfPresent() && gf.field0Ptr_)
615  {
616  field0Ptr_ = std::make_unique<this_type>
617  (
618  io.name() + "_0",
619  *gf.field0Ptr_
620  );
621  }
622 }
623 
624 
625 template<class Type, template<class> class PatchField, class GeoMesh>
627 (
628  const IOobject& io,
629  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
630 )
631 :
632  Internal(io, tgf.constCast(), tgf.movable()),
633  timeIndex_(tgf().timeIndex()),
634  boundaryField_(*this, tgf().boundaryField_)
635 {
637  << "Constructing from tmp resetting IO params" << nl
638  << this->info() << endl;
639 
640  tgf.clear();
642  readIfPresent();
643 }
644 
645 
646 template<class Type, template<class> class PatchField, class GeoMesh>
648 (
649  const word& newName,
651 )
652 :
653  Internal(newName, gf),
654  timeIndex_(gf.timeIndex()),
655  boundaryField_(*this, gf.boundaryField_)
656 {
658  << "Copy construct, resetting name" << nl
659  << this->info() << endl;
660 
661  if (!readIfPresent() && gf.field0Ptr_)
662  {
663  field0Ptr_ = std::make_unique<this_type>
664  (
665  newName + "_0",
666  *gf.field0Ptr_
667  );
668  }
669 }
670 
671 
672 template<class Type, template<class> class PatchField, class GeoMesh>
674 (
675  const word& newName,
676  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
677 )
678 :
679  Internal(newName, tgf.constCast(), tgf.movable()),
680  timeIndex_(tgf().timeIndex()),
681  boundaryField_(*this, tgf().boundaryField_)
682 {
684  << "Constructing from tmp resetting name" << nl
685  << this->info() << endl;
687  tgf.clear();
688 }
689 
690 
691 template<class Type, template<class> class PatchField, class GeoMesh>
693 (
694  const IOobject& io,
696  const word& patchFieldType
697 )
698 :
699  Internal(io, gf),
700  timeIndex_(gf.timeIndex()),
701  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
702 {
704  << "Copy construct, resetting IO params" << nl
705  << this->info() << endl;
706 
707  boundaryField_ == gf.boundaryField_;
708 
709  if (!readIfPresent() && gf.field0Ptr_)
710  {
711  field0Ptr_ = std::make_unique<this_type>
712  (
713  io.name() + "_0",
714  *gf.field0Ptr_
715  );
716  }
717 }
718 
719 
720 template<class Type, template<class> class PatchField, class GeoMesh>
722 (
723  const IOobject& io,
724  const GeometricField<Type, PatchField, GeoMesh>& gf,
725  const wordList& patchFieldTypes,
726  const wordList& actualPatchTypes
727 )
728 :
729  Internal(io, gf),
730  timeIndex_(gf.timeIndex()),
731  boundaryField_
732  (
733  this->mesh().boundary(),
734  *this,
735  patchFieldTypes,
736  actualPatchTypes
737  )
738 {
740  << "Copy construct, resetting IO params and patch types" << nl
741  << this->info() << endl;
742 
743  boundaryField_ == gf.boundaryField_;
744 
745  if (!readIfPresent() && gf.field0Ptr_)
746  {
747  field0Ptr_ = std::make_unique<this_type>
748  (
749  io.name() + "_0",
750  *gf.field0Ptr_
751  );
752  }
753 }
754 
755 
756 template<class Type, template<class> class PatchField, class GeoMesh>
758 (
759  const IOobject& io,
760  const GeometricField<Type, PatchField, GeoMesh>& gf,
761  const labelList& patchIDs,
762  const word& patchFieldType
763 )
764 :
765  Internal(io, gf),
766  timeIndex_(gf.timeIndex()),
767  boundaryField_(*this, gf.boundaryField_, patchIDs, patchFieldType)
768 {
770  << "Copy construct, resetting IO params and setting patchFieldType "
771  << "for patchIDs" << nl
772  << this->info() << endl;
773 
774  if (!readIfPresent() && gf.field0Ptr_)
775  {
776  field0Ptr_ = std::make_unique<this_type>
777  (
778  io.name() + "_0",
779  *gf.field0Ptr_
780  );
781  }
782 }
783 
784 
785 template<class Type, template<class> class PatchField, class GeoMesh>
787 (
788  const IOobject& io,
789  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf,
790  const wordList& patchFieldTypes,
791  const wordList& actualPatchTypes
792 )
793 :
794  Internal(io, tgf.constCast(), tgf.movable()),
795  timeIndex_(tgf().timeIndex()),
796  boundaryField_
797  (
798  this->mesh().boundary(),
799  *this,
800  patchFieldTypes,
801  actualPatchTypes
802  )
803 {
805  << "Constructing from tmp resetting IO params and patch types" << nl
806  << this->info() << endl;
807 
808  boundaryField_ == tgf().boundaryField_;
810  tgf.clear();
811 }
812 
813 
814 template<class Type, template<class> class PatchField, class GeoMesh>
817 {
819 }
820 
821 
822 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
823 
824 template<class Type, template<class> class PatchField, class GeoMesh>
826 {
827  /*
828  if (debug)
829  {
830  // Problem: temporary fields might have their internal field
831  // already stolen so boundary fields will not be able to access the
832  // internal field anymore
833  boundaryField_.check();
834  }
835  */
836 
837  // FUTURE: register cache field info
838  // // this->db().cacheTemporaryObject(*this);
839 }
840 
842 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
843 
844 template<class Type, template<class> class PatchField, class GeoMesh>
845 typename
848 (
849  const bool updateAccessTime
850 )
851 {
852  if (updateAccessTime)
853  {
854  this->setUpToDate();
855  storeOldTimes();
856  }
857  return *this;
858 }
859 
860 
861 template<class Type, template<class> class PatchField, class GeoMesh>
862 typename
865 (
866  const bool updateAccessTime
867 )
868 {
869  if (updateAccessTime)
870  {
871  this->setUpToDate();
872  storeOldTimes();
873  }
874  return *this;
875 }
876 
877 
878 template<class Type, template<class> class PatchField, class GeoMesh>
879 typename
882 (
883  const bool updateAccessTime
884 )
885 {
886  if (updateAccessTime)
887  {
888  this->setUpToDate();
889  storeOldTimes();
890  }
891  return boundaryField_;
892 }
893 
894 
895 template<class Type, template<class> class PatchField, class GeoMesh>
896 Foam::label
898 {
899  return (field0Ptr_ ? (field0Ptr_->nOldTimes() + 1) : 0);
900 }
901 
902 
903 template<class Type, template<class> class PatchField, class GeoMesh>
905 {
906  if
907  (
908  field0Ptr_
909  && timeIndex_ != this->time().timeIndex()
910  && !this->name().ends_with("_0")
911  )
912  {
913  storeOldTime();
914  timeIndex_ = this->time().timeIndex();
915  }
917  // Correct time index
918  //timeIndex_ = this->time().timeIndex();
919 }
920 
921 
922 template<class Type, template<class> class PatchField, class GeoMesh>
924 {
925  if (field0Ptr_)
926  {
927  field0Ptr_->storeOldTime();
928 
930  << "Storing old time field for field" << nl << this->info() << endl;
931 
932  *field0Ptr_ == *this;
933  field0Ptr_->timeIndex_ = timeIndex_;
934 
935  if (field0Ptr_->field0Ptr_)
936  {
937  field0Ptr_->writeOpt(this->writeOpt());
938  }
939  }
940 }
941 
942 
943 template<class Type, template<class> class PatchField, class GeoMesh>
946 {
947  if (!field0Ptr_)
948  {
949  field0Ptr_ = std::make_unique<this_type>
950  (
951  IOobject
952  (
953  this->name() + "_0",
954  this->time().timeName(),
955  this->db(),
956  IOobjectOption::NO_READ,
957  IOobjectOption::NO_WRITE,
958  this->registerObject()
959  ),
960  *this
961  );
962 
963  if (debug)
964  {
966  << "created old time field " << field0Ptr_->info() << endl;
967 
968  if (debug&2)
969  {
970  error::printStack(Info);
971  }
972  }
973  }
974  else
975  {
976  storeOldTimes();
977  }
979  return *field0Ptr_;
980 }
981 
982 
983 template<class Type, template<class> class PatchField, class GeoMesh>
986 {
987  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
989 
990  return *field0Ptr_;
991 }
992 
993 
994 template<class Type, template<class> class PatchField, class GeoMesh>
996 {
997  if (!fieldPrevIterPtr_)
998  {
1000  << "Allocating previous iteration field" << nl
1001  << this->info() << endl;
1002 
1003  fieldPrevIterPtr_ = std::make_unique<this_type>
1004  (
1005  this->name() + "PrevIter",
1006  *this
1007  );
1008  }
1009  else
1010  {
1011  *fieldPrevIterPtr_ == *this;
1012  }
1013 }
1014 
1015 
1016 template<class Type, template<class> class PatchField, class GeoMesh>
1019 {
1020  if (!fieldPrevIterPtr_)
1021  {
1023  << "previous iteration field" << endl << this->info() << endl
1024  << " not stored."
1025  << " Use field.storePrevIter() at start of iteration."
1026  << abort(FatalError);
1027  }
1028 
1029  return *fieldPrevIterPtr_;
1030 }
1031 
1032 
1033 template<class Type, template<class> class PatchField, class GeoMesh>
1035 {
1036  field0Ptr_.reset(nullptr);
1037  fieldPrevIterPtr_.reset(nullptr);
1038 }
1039 
1040 
1041 template<class Type, template<class> class PatchField, class GeoMesh>
1044 {
1045  // updateAccessTime
1046  {
1047  this->setUpToDate();
1048  storeOldTimes();
1049  }
1050  boundaryField_.evaluate();
1051 }
1052 
1053 
1054 template<class Type, template<class> class PatchField, class GeoMesh>
1057 {
1058  // updateAccessTime
1059  {
1060  this->setUpToDate();
1061  storeOldTimes();
1062  }
1063  boundaryField_.evaluateLocal();
1064 }
1065 
1066 
1067 template<class Type, template<class> class PatchField, class GeoMesh>
1069 {
1070  // Search all boundary conditions, if any are
1071  // fixed-value or mixed (Robin) do not set reference level for solution.
1072 
1073  bool needRef = true;
1074 
1075  for (const auto& pf : boundaryField_)
1076  {
1077  if (pf.fixesValue())
1078  {
1079  needRef = false;
1080  break;
1081  }
1082  }
1083 
1084  return returnReduceAnd(needRef);
1085 }
1086 
1087 
1088 template<class Type, template<class> class PatchField, class GeoMesh>
1090 {
1092  << "Relaxing" << nl << this->info() << " by " << alpha << endl;
1093 
1094  operator==(prevIter() + alpha*(*this - prevIter()));
1095 }
1096 
1097 
1098 template<class Type, template<class> class PatchField, class GeoMesh>
1100 {
1101  word name = this->name();
1102 
1103  if (this->mesh().data().isFinalIteration())
1104  {
1105  name += "Final";
1106  }
1107 
1108  scalar relaxCoeff = 1;
1109 
1110  if (this->mesh().relaxField(name, relaxCoeff))
1111  {
1112  relax(relaxCoeff);
1113  }
1114 }
1115 
1116 
1117 template<class Type, template<class> class PatchField, class GeoMesh>
1119 (
1120  bool final
1121 ) const
1122 {
1123  if (final)
1124  {
1125  return this->name() + "Final";
1126  }
1128  return this->name();
1129 }
1130 
1131 
1132 template<class Type, template<class> class PatchField, class GeoMesh>
1134 (
1135  Ostream& os
1136 ) const
1137 {
1138  MinMax<Type> range = Foam::minMax(*this).value();
1139 
1140  os << "min/max(" << this->name() << ") = "
1141  << range.min() << ", " << range.max() << endl;
1142 }
1143 
1144 
1145 template<class Type, template<class> class PatchField, class GeoMesh>
1147 writeData(Ostream& os) const
1148 {
1149  this->internalField().writeData(os, "internalField");
1150  os << nl;
1151  this->boundaryField().writeEntry("boundaryField", os);
1152 
1153  os.check(FUNCTION_NAME);
1154  return os.good();
1156 
1157 
1158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1159 
1160 template<class Type, template<class> class PatchField, class GeoMesh>
1163 {
1165  (
1166  this->name() + ".T()",
1167  this->mesh(),
1168  this->dimensions()
1169  );
1170 
1171  Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1172  Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1173 
1174  return tresult;
1175 }
1176 
1177 
1178 template<class Type, template<class> class PatchField, class GeoMesh>
1179 Foam::tmp
1180 <
1182  <
1184  PatchField,
1185  GeoMesh
1186  >
1187 >
1189 (
1190  const direction d
1191 ) const
1192 {
1194  (
1195  this->name() + ".component(" + Foam::name(d) + ')',
1196  this->mesh(),
1197  this->dimensions()
1198  );
1199 
1200  Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1201  Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1202 
1203  return tresult;
1204 }
1205 
1206 
1207 template<class Type, template<class> class PatchField, class GeoMesh>
1209 (
1210  const direction d,
1211  const GeometricField
1212  <
1213  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1214  PatchField,
1215  GeoMesh
1216  >& gcf
1217 )
1218 {
1219  primitiveFieldRef().replace(d, gcf.primitiveField());
1220  boundaryFieldRef().replace(d, gcf.boundaryField());
1221 }
1222 
1223 
1224 template<class Type, template<class> class PatchField, class GeoMesh>
1226 (
1227  const direction d,
1228  const dimensioned<cmptType>& ds
1229 )
1230 {
1231  primitiveFieldRef().replace(d, ds.value());
1232  boundaryFieldRef().replace(d, ds.value());
1233 }
1234 
1235 
1236 template<class Type, template<class> class PatchField, class GeoMesh>
1238 (
1239  const Type& lower
1240 )
1241 {
1242  primitiveFieldRef().clamp_min(lower);
1243  boundaryFieldRef().clamp_min(lower);
1244 }
1245 
1246 
1247 template<class Type, template<class> class PatchField, class GeoMesh>
1249 (
1250  const Type& upper
1251 )
1252 {
1253  primitiveFieldRef().clamp_max(upper);
1254  boundaryFieldRef().clamp_max(upper);
1255 }
1256 
1257 
1258 template<class Type, template<class> class PatchField, class GeoMesh>
1260 (
1261  const dimensioned<Type>& lower
1262 )
1264  this->clamp_min(lower.value());
1265 }
1266 
1267 
1268 template<class Type, template<class> class PatchField, class GeoMesh>
1270 (
1271  const dimensioned<Type>& upper
1272 )
1274  this->clamp_max(upper.value());
1275 }
1276 
1277 
1278 template<class Type, template<class> class PatchField, class GeoMesh>
1280 (
1281  const Type& lower,
1282  const Type& upper
1283 )
1284 {
1285  primitiveFieldRef().clamp_range(lower, upper);
1286  boundaryFieldRef().clamp_range(lower, upper);
1287 }
1288 
1289 
1290 template<class Type, template<class> class PatchField, class GeoMesh>
1292 (
1293  const MinMax<Type>& range
1294 )
1295 {
1296  primitiveFieldRef().clamp_range(range.min(), range.max());
1297  boundaryFieldRef().clamp_range(range.min(), range.max());
1298 }
1299 
1300 
1301 template<class Type, template<class> class PatchField, class GeoMesh>
1303 (
1304  const dimensioned<Type>& lower,
1305  const dimensioned<Type>& upper
1306 )
1308  this->clamp_range(lower.value(), upper.value());
1309 }
1310 
1311 
1312 template<class Type, template<class> class PatchField, class GeoMesh>
1314 (
1317 {
1318  this->clamp_range(range.value());
1319 }
1320 
1321 
1322 template<class Type, template<class> class PatchField, class GeoMesh>
1325  primitiveFieldRef().negate();
1326  boundaryFieldRef().negate();
1327 }
1328 
1329 
1330 template<class Type, template<class> class PatchField, class GeoMesh>
1332 {
1333  primitiveFieldRef().normalise();
1334  boundaryFieldRef().normalise();
1336 
1337 
1338 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1339 
1340 template<class Type, template<class> class PatchField, class GeoMesh>
1342 (
1344 )
1345 {
1346  if (this == &gf)
1347  {
1348  return; // Self-assignment is a no-op
1349  }
1350 
1351  checkField(*this, gf, "=");
1352 
1353  // Only assign field contents not ID
1354 
1355  internalFieldRef() = gf.internalField();
1357 }
1358 
1359 
1360 template<class Type, template<class> class PatchField, class GeoMesh>
1362 (
1364 )
1365 {
1366  const auto& gf = tgf();
1367 
1368  if (this == &gf)
1369  {
1370  return; // Self-assignment is a no-op
1371  }
1372 
1373  checkField(*this, gf, "=");
1374 
1375  // Only assign field contents not ID
1376 
1377  this->dimensions() = gf.dimensions();
1378  this->oriented() = gf.oriented();
1379 
1380  if (tgf.movable())
1381  {
1382  // Transfer storage from the tmp
1383  primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1384  }
1385  else
1386  {
1388  }
1389 
1392  tgf.clear();
1393 }
1394 
1395 
1396 template<class Type, template<class> class PatchField, class GeoMesh>
1398 (
1399  const dimensioned<Type>& dt
1400 )
1401 {
1402  internalFieldRef() = dt;
1403  boundaryFieldRef() = dt.value();
1404 }
1405 
1406 
1407 template<class Type, template<class> class PatchField, class GeoMesh>
1409 (
1411 )
1412 {
1413  const auto& gf = tgf();
1414 
1415  checkField(*this, gf, "==");
1416 
1417  // Only assign field contents not ID
1418 
1419  internalFieldRef() = gf.internalField();
1420  boundaryFieldRef() == gf.boundaryField();
1422  tgf.clear();
1423 }
1424 
1425 
1426 template<class Type, template<class> class PatchField, class GeoMesh>
1428 (
1429  const dimensioned<Type>& dt
1431 {
1432  internalFieldRef() = dt;
1433  boundaryFieldRef() == dt.value();
1434 }
1435 
1436 
1437 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1438  \
1439 template<class Type, template<class> class PatchField, class GeoMesh> \
1440 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1441 ( \
1442  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1443 ) \
1444 { \
1445  checkField(*this, gf, #op); \
1446  \
1447  internalFieldRef() op gf.internalField(); \
1448  boundaryFieldRef() op gf.boundaryField(); \
1449 } \
1450  \
1451 template<class Type, template<class> class PatchField, class GeoMesh> \
1452 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1453 ( \
1454  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1455 ) \
1456 { \
1457  operator op(tgf()); \
1458  tgf.clear(); \
1459 } \
1460  \
1461 template<class Type, template<class> class PatchField, class GeoMesh> \
1462 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1463 ( \
1464  const dimensioned<TYPE>& dt \
1465 ) \
1466 { \
1467  internalFieldRef() op dt; \
1468  boundaryFieldRef() op dt.value(); \
1469 }
1470 
1471 COMPUTED_ASSIGNMENT(Type, +=)
1472 COMPUTED_ASSIGNMENT(Type, -=)
1473 COMPUTED_ASSIGNMENT(scalar, *=)
1474 COMPUTED_ASSIGNMENT(scalar, /=)
1475 
1476 #undef COMPUTED_ASSIGNMENT
1477 
1478 
1479 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1480 
1481 template<class Type, template<class> class PatchField, class GeoMesh>
1482 Foam::Ostream& Foam::operator<<
1483 (
1484  Ostream& os,
1486 )
1487 {
1488  fld.writeData(os);
1489  return os;
1490 }
1491 
1492 
1493 template<class Type, template<class> class PatchField, class GeoMesh>
1494 Foam::Ostream& Foam::operator<<
1495 (
1496  Ostream& os,
1498 )
1499 {
1500  tfld().writeData(os);
1501  tfld.clear();
1502 
1503  return os;
1504 }
1505 
1506 
1507 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1508 
1509 #undef checkField
1510 
1511 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1512 
1513 #include "GeometricFieldNew.C"
1514 #include "GeometricFieldFunctions.C"
1515 
1516 // ************************************************************************* //
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:608
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.
Generic dimensioned Type class.
const dimensionSet dimless
Dimensionless.
bool writeData(Ostream &os) const
The writeData function (required by regIOobject)
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.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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
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...
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: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
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.
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.