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-2021 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"
33 #include "localIOdictionary.H"
34 #include "data.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 #define checkField(gf1, gf2, op) \
39 if ((gf1).mesh() != (gf2).mesh()) \
40 { \
41  FatalErrorInFunction \
42  << "different mesh for fields " \
43  << (gf1).name() << " and " << (gf2).name() \
44  << " during operation " << op \
45  << abort(FatalError); \
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50 
51 template<class Type, template<class> class PatchField, class GeoMesh>
53 (
54  const dictionary& dict
55 )
56 {
57  Internal::readField(dict, "internalField");
58 
59  boundaryField_.readField(*this, dict.subDict("boundaryField"));
60 
61  Type refLevel;
62 
63  if (dict.readIfPresent("referenceLevel", refLevel))
64  {
65  Field<Type>::operator+=(refLevel);
66 
67  forAll(boundaryField_, patchi)
68  {
69  boundaryField_[patchi] == boundaryField_[patchi] + refLevel;
70  }
71  }
72 }
73 
74 
75 template<class Type, template<class> class PatchField, class GeoMesh>
77 {
78  const localIOdictionary dict
79  (
80  IOobject
81  (
82  this->name(),
83  this->instance(),
84  this->local(),
85  this->db(),
86  IOobject::MUST_READ,
87  IOobject::NO_WRITE,
88  false
89  ),
90  typeName
91  );
92 
93  this->close();
94 
96 }
97 
98 
99 template<class Type, template<class> class PatchField, class GeoMesh>
101 {
102  if
103  (
104  this->readOpt() == IOobject::MUST_READ
105  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
106  )
107  {
109  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
110  << " suggests that a read constructor for field " << this->name()
111  << " would be more appropriate." << endl;
112  }
113  else if
114  (
115  this->readOpt() == IOobject::READ_IF_PRESENT
116  && this->template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
117  (
118  true
119  )
120  )
121  {
122  readFields();
123 
124  // Check compatibility between field and mesh
125  if (this->size() != GeoMesh::size(this->mesh()))
126  {
127  FatalIOErrorInFunction(this->readStream(typeName))
128  << " number of field elements = " << this->size()
129  << " number of mesh elements = "
130  << GeoMesh::size(this->mesh())
131  << exit(FatalIOError);
132  }
133 
134  readOldTimeIfPresent();
135 
136  return true;
137  }
138 
139  return false;
140 }
141 
142 
143 template<class Type, template<class> class PatchField, class GeoMesh>
145 {
146  // Read the old time field if present
147  IOobject field0
148  (
149  this->name() + "_0",
150  this->time().timeName(),
151  this->db(),
152  IOobject::READ_IF_PRESENT,
153  IOobject::AUTO_WRITE,
154  this->registerObject()
155  );
156 
157  if
158  (
159  field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
160  (
161  true
162  )
163  )
164  {
166  << "Reading old time level for field" << nl << this->info() << endl;
167 
168  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
169  (
170  field0,
171  this->mesh()
172  );
173 
174  // Ensure the old time field oriented flag is set to the parent's state
175  // Note: required for backwards compatibility in case of restarting from
176  // an old run where the oriented state may not have been set
177  field0Ptr_->oriented() = this->oriented();
178 
179  field0Ptr_->timeIndex_ = timeIndex_ - 1;
180 
181  if (!field0Ptr_->readOldTimeIfPresent())
182  {
183  field0Ptr_->oldTime();
184  }
185 
186  return true;
187  }
188 
189  return false;
190 }
191 
192 
193 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
194 
195 template<class Type, template<class> class PatchField, class GeoMesh>
197 (
198  const IOobject& io,
199  const Mesh& mesh,
200  const dimensionSet& ds,
201  const word& patchFieldType
202 )
203 :
204  Internal(io, mesh, ds, false),
205  timeIndex_(this->time().timeIndex()),
206  field0Ptr_(nullptr),
207  fieldPrevIterPtr_(nullptr),
208  boundaryField_(mesh.boundary(), *this, patchFieldType)
209 {
211  << "Creating temporary" << nl << this->info() << endl;
212 
213  readIfPresent();
214 }
215 
216 
217 template<class Type, template<class> class PatchField, class GeoMesh>
219 (
220  const IOobject& io,
221  const Mesh& mesh,
222  const dimensionSet& ds,
223  const wordList& patchFieldTypes,
224  const wordList& actualPatchTypes
225 )
226 :
227  Internal(io, mesh, ds, false),
228  timeIndex_(this->time().timeIndex()),
229  field0Ptr_(nullptr),
230  fieldPrevIterPtr_(nullptr),
231  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
232 {
234  << "Creating temporary" << nl << this->info() << endl;
235 
236  readIfPresent();
237 }
238 
239 
240 template<class Type, template<class> class PatchField, class GeoMesh>
242 (
243  const IOobject& io,
244  const Mesh& mesh,
245  const dimensioned<Type>& dt,
246  const word& patchFieldType
247 )
248 :
249  Internal(io, mesh, dt, false),
250  timeIndex_(this->time().timeIndex()),
251  field0Ptr_(nullptr),
252  fieldPrevIterPtr_(nullptr),
253  boundaryField_(mesh.boundary(), *this, patchFieldType)
254 {
256  << "Creating temporary" << nl << this->info() << endl;
257 
258  boundaryField_ == dt.value();
259 
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 dimensioned<Type>& dt,
270  const wordList& patchFieldTypes,
271  const wordList& actualPatchTypes
272 )
273 :
274  Internal(io, mesh, dt, false),
275  timeIndex_(this->time().timeIndex()),
276  field0Ptr_(nullptr),
277  fieldPrevIterPtr_(nullptr),
278  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
279 {
281  << "Creating temporary" << nl << this->info() << endl;
282 
283  boundaryField_ == dt.value();
284 
285  readIfPresent();
286 }
287 
288 
289 template<class Type, template<class> class PatchField, class GeoMesh>
291 (
292  const IOobject& io,
293  const Internal& diField,
294  const PtrList<PatchField<Type>>& ptfl
295 )
296 :
297  Internal(io, diField),
298  timeIndex_(this->time().timeIndex()),
299  field0Ptr_(nullptr),
300  fieldPrevIterPtr_(nullptr),
301  boundaryField_(this->mesh().boundary(), *this, ptfl)
302 {
304  << "Copy construct from components" << nl << this->info() << endl;
305 
306  readIfPresent();
307 }
308 
309 
310 template<class Type, template<class> class PatchField, class GeoMesh>
312 (
313  const IOobject& io,
314  const Mesh& mesh,
315  const dimensionSet& ds,
316  const Field<Type>& iField,
317  const word& patchFieldType
318 )
319 :
320  Internal(io, mesh, ds, iField),
321  timeIndex_(this->time().timeIndex()),
322  field0Ptr_(nullptr),
323  fieldPrevIterPtr_(nullptr),
324  boundaryField_(mesh.boundary(), *this, patchFieldType)
325 {
327  << "Copy construct from internal field" << nl << this->info() << endl;
328 
329  readIfPresent();
330 }
331 
332 
333 template<class Type, template<class> class PatchField, class GeoMesh>
335 (
336  const IOobject& io,
337  const Mesh& mesh,
338  const dimensionSet& ds,
339  Field<Type>&& iField,
340  const word& patchFieldType
341 )
342 :
343  Internal(io, mesh, ds, std::move(iField)),
344  timeIndex_(this->time().timeIndex()),
345  field0Ptr_(nullptr),
346  fieldPrevIterPtr_(nullptr),
347  boundaryField_(mesh.boundary(), *this, patchFieldType)
348 {
350  << "Move construct from internal field" << nl << this->info() << endl;
351 
352  readIfPresent();
353 }
354 
355 
356 template<class Type, template<class> class PatchField, class GeoMesh>
358 (
359  const IOobject& io,
360  const Mesh& mesh,
361  const dimensionSet& ds,
362  const Field<Type>& iField,
363  const PtrList<PatchField<Type>>& ptfl
364 )
365 :
366  Internal(io, mesh, ds, iField),
367  timeIndex_(this->time().timeIndex()),
368  field0Ptr_(nullptr),
369  fieldPrevIterPtr_(nullptr),
370  boundaryField_(mesh.boundary(), *this, ptfl)
371 {
373  << "Copy construct from components" << nl << this->info() << endl;
374 
375  readIfPresent();
376 }
377 
378 
379 template<class Type, template<class> class PatchField, class GeoMesh>
381 (
382  const IOobject& io,
383  const Mesh& mesh,
384  const bool readOldTime
385 )
386 :
387  Internal(io, mesh, dimless, false),
388  timeIndex_(this->time().timeIndex()),
389  field0Ptr_(nullptr),
390  fieldPrevIterPtr_(nullptr),
391  boundaryField_(mesh.boundary())
392 {
393  readFields();
394 
395  // Check compatibility between field and mesh
396 
397  if (this->size() != GeoMesh::size(this->mesh()))
398  {
399  FatalIOErrorInFunction(this->readStream(typeName))
400  << " number of field elements = " << this->size()
401  << " number of mesh elements = " << GeoMesh::size(this->mesh())
402  << exit(FatalIOError);
403  }
404 
405  if (readOldTime)
406  {
407  readOldTimeIfPresent();
408  }
409 
411  << "Finishing read-construction" << nl << this->info() << endl;
412 }
413 
414 
415 template<class Type, template<class> class PatchField, class GeoMesh>
417 (
418  const IOobject& io,
419  const Mesh& mesh,
420  const dictionary& dict
421 )
422 :
423  Internal(io, mesh, dimless, false),
424  timeIndex_(this->time().timeIndex()),
425  field0Ptr_(nullptr),
426  fieldPrevIterPtr_(nullptr),
427  boundaryField_(mesh.boundary())
428 {
429  readFields(dict);
430 
431  // Check compatibility between field and mesh
432 
433  if (this->size() != GeoMesh::size(this->mesh()))
434  {
436  << " number of field elements = " << this->size()
437  << " number of mesh elements = " << GeoMesh::size(this->mesh())
438  << exit(FatalIOError);
439  }
440 
442  << "Finishing dictionary-construct" << nl << this->info() << endl;
443 }
444 
445 
446 template<class Type, template<class> class PatchField, class GeoMesh>
448 (
450 )
451 :
452  Internal(gf),
453  timeIndex_(gf.timeIndex()),
454  field0Ptr_(nullptr),
455  fieldPrevIterPtr_(nullptr),
456  boundaryField_(*this, gf.boundaryField_)
457 {
459  << "Copy construct" << nl << this->info() << endl;
460 
461  if (gf.field0Ptr_)
462  {
464  (
465  *gf.field0Ptr_
466  );
467  }
468 
469  this->writeOpt(IOobject::NO_WRITE);
470 }
471 
472 
473 template<class Type, template<class> class PatchField, class GeoMesh>
475 (
477 )
478 :
479  Internal(tgf.constCast(), tgf.movable()),
480  timeIndex_(tgf().timeIndex()),
481  field0Ptr_(nullptr),
482  fieldPrevIterPtr_(nullptr),
483  boundaryField_(*this, tgf().boundaryField_)
484 {
486  << "Constructing from tmp" << nl << this->info() << endl;
487 
488  this->writeOpt(IOobject::NO_WRITE);
489 
490  tgf.clear();
491 }
492 
493 
494 template<class Type, template<class> class PatchField, class GeoMesh>
496 (
497  const IOobject& io,
499 )
500 :
501  Internal(io, gf),
502  timeIndex_(gf.timeIndex()),
503  field0Ptr_(nullptr),
504  fieldPrevIterPtr_(nullptr),
505  boundaryField_(*this, gf.boundaryField_)
506 {
508  << "Copy construct, resetting IO params" << nl
509  << this->info() << endl;
510 
511  if (!readIfPresent() && gf.field0Ptr_)
512  {
514  (
515  io.name() + "_0",
516  *gf.field0Ptr_
517  );
518  }
519 }
520 
521 
522 template<class Type, template<class> class PatchField, class GeoMesh>
524 (
525  const IOobject& io,
527 )
528 :
529  Internal(io, tgf.constCast(), tgf.movable()),
530  timeIndex_(tgf().timeIndex()),
531  field0Ptr_(nullptr),
532  fieldPrevIterPtr_(nullptr),
533  boundaryField_(*this, tgf().boundaryField_)
534 {
536  << "Constructing from tmp resetting IO params" << nl
537  << this->info() << endl;
538 
539  tgf.clear();
540 
541  readIfPresent();
542 }
543 
544 
545 template<class Type, template<class> class PatchField, class GeoMesh>
547 (
548  const word& newName,
550 )
551 :
552  Internal(newName, gf),
553  timeIndex_(gf.timeIndex()),
554  field0Ptr_(nullptr),
555  fieldPrevIterPtr_(nullptr),
556  boundaryField_(*this, gf.boundaryField_)
557 {
559  << "Copy construct, resetting name" << nl
560  << this->info() << endl;
561 
562  if (!readIfPresent() && gf.field0Ptr_)
563  {
565  (
566  newName + "_0",
567  *gf.field0Ptr_
568  );
569  }
570 }
571 
572 
573 template<class Type, template<class> class PatchField, class GeoMesh>
575 (
576  const word& newName,
578 )
579 :
580  Internal(newName, tgf.constCast(), tgf.movable()),
581  timeIndex_(tgf().timeIndex()),
582  field0Ptr_(nullptr),
583  fieldPrevIterPtr_(nullptr),
584  boundaryField_(*this, tgf().boundaryField_)
585 {
587  << "Constructing from tmp resetting name" << nl
588  << this->info() << endl;
589 
590  tgf.clear();
591 }
592 
593 
594 template<class Type, template<class> class PatchField, class GeoMesh>
596 (
597  const IOobject& io,
599  const word& patchFieldType
600 )
601 :
602  Internal(io, gf),
603  timeIndex_(gf.timeIndex()),
604  field0Ptr_(nullptr),
605  fieldPrevIterPtr_(nullptr),
606  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
607 {
609  << "Copy construct, resetting IO params" << nl
610  << this->info() << endl;
611 
612  boundaryField_ == gf.boundaryField_;
613 
614  if (!readIfPresent() && gf.field0Ptr_)
615  {
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,
630  const wordList& patchFieldTypes,
631  const wordList& actualPatchTypes
632 )
633 :
634  Internal(io, gf),
635  timeIndex_(gf.timeIndex()),
636  field0Ptr_(nullptr),
637  fieldPrevIterPtr_(nullptr),
638  boundaryField_
639  (
640  this->mesh().boundary(),
641  *this,
642  patchFieldTypes,
643  actualPatchTypes
644  )
645 {
647  << "Copy construct, resetting IO params and patch types" << nl
648  << this->info() << endl;
649 
650  boundaryField_ == gf.boundaryField_;
651 
652  if (!readIfPresent() && gf.field0Ptr_)
653  {
655  (
656  io.name() + "_0",
657  *gf.field0Ptr_
658  );
659  }
660 }
661 
662 
663 template<class Type, template<class> class PatchField, class GeoMesh>
665 (
666  const IOobject& io,
668  const labelList& patchIDs,
669  const word& patchFieldType
670 )
671 :
672  Internal(io, gf),
673  timeIndex_(gf.timeIndex()),
674  field0Ptr_(nullptr),
675  fieldPrevIterPtr_(nullptr),
676  boundaryField_(*this, gf.boundaryField_, patchIDs, patchFieldType)
677 {
679  << "Copy construct, resetting IO params and setting patchFieldType "
680  << "for patchIDs" << nl
681  << this->info() << endl;
682 
683  if (!readIfPresent() && gf.field0Ptr_)
684  {
686  (
687  io.name() + "_0",
688  *gf.field0Ptr_
689  );
690  }
691 }
692 
693 
694 template<class Type, template<class> class PatchField, class GeoMesh>
696 (
697  const IOobject& io,
699  const wordList& patchFieldTypes,
700  const wordList& actualPatchTypes
701 )
702 :
703  Internal(io, tgf.constCast(), tgf.movable()),
704  timeIndex_(tgf().timeIndex()),
705  field0Ptr_(nullptr),
706  fieldPrevIterPtr_(nullptr),
707  boundaryField_
708  (
709  this->mesh().boundary(),
710  *this,
711  patchFieldTypes,
712  actualPatchTypes
713  )
714 {
716  << "Constructing from tmp resetting IO params and patch types" << nl
717  << this->info() << endl;
718 
719  boundaryField_ == tgf().boundaryField_;
720 
721  tgf.clear();
722 }
723 
724 
725 template<class Type, template<class> class PatchField, class GeoMesh>
728 {
730 }
731 
732 
733 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
734 
735 template<class Type, template<class> class PatchField, class GeoMesh>
737 {
738  deleteDemandDrivenData(field0Ptr_);
739  deleteDemandDrivenData(fieldPrevIterPtr_);
740 }
741 
742 
743 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
744 
745 template<class Type, template<class> class PatchField, class GeoMesh>
746 typename
749 (
750  const bool updateAccessTime
751 )
752 {
753  if (updateAccessTime)
754  {
755  this->setUpToDate();
756  storeOldTimes();
757  }
758  return *this;
759 }
760 
761 
762 template<class Type, template<class> class PatchField, class GeoMesh>
763 typename
766 (
767  const bool updateAccessTime
768 )
769 {
770  if (updateAccessTime)
771  {
772  this->setUpToDate();
773  storeOldTimes();
774  }
775  return *this;
776 }
777 
778 
779 template<class Type, template<class> class PatchField, class GeoMesh>
780 typename
783 (
784  const bool updateAccessTime
785 )
786 {
787  if (updateAccessTime)
788  {
789  this->setUpToDate();
790  storeOldTimes();
791  }
792  return boundaryField_;
793 }
794 
795 
796 template<class Type, template<class> class PatchField, class GeoMesh>
798 {
799  if
800  (
801  field0Ptr_
802  && timeIndex_ != this->time().timeIndex()
803  && !this->name().ends_with("_0")
804  )
805  {
806  storeOldTime();
807  timeIndex_ = this->time().timeIndex();
808  }
809 
810  // Correct time index
811  //timeIndex_ = this->time().timeIndex();
812 }
813 
814 
815 template<class Type, template<class> class PatchField, class GeoMesh>
817 {
818  if (field0Ptr_)
819  {
820  field0Ptr_->storeOldTime();
821 
823  << "Storing old time field for field" << nl << this->info() << endl;
824 
825  *field0Ptr_ == *this;
826  field0Ptr_->timeIndex_ = timeIndex_;
827 
828  if (field0Ptr_->field0Ptr_)
829  {
830  field0Ptr_->writeOpt(this->writeOpt());
831  }
832  }
833 }
834 
835 
836 template<class Type, template<class> class PatchField, class GeoMesh>
838 {
839  if (field0Ptr_)
840  {
841  return field0Ptr_->nOldTimes() + 1;
842  }
843 
844  return 0;
845 }
846 
847 
848 template<class Type, template<class> class PatchField, class GeoMesh>
851 {
852  if (!field0Ptr_)
853  {
855  (
856  IOobject
857  (
858  this->name() + "_0",
859  this->time().timeName(),
860  this->db(),
861  IOobject::NO_READ,
862  IOobject::NO_WRITE,
863  this->registerObject()
864  ),
865  *this
866  );
867 
868  if (debug)
869  {
871  << "created old time field " << field0Ptr_->info() << endl;
872 
873  if (debug&2)
874  {
875  error::printStack(Info);
876  }
877  }
878  }
879  else
880  {
881  storeOldTimes();
882  }
883 
884  return *field0Ptr_;
885 }
886 
887 
888 template<class Type, template<class> class PatchField, class GeoMesh>
891 {
892  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
893  .oldTime();
894 
895  return *field0Ptr_;
896 }
897 
898 
899 template<class Type, template<class> class PatchField, class GeoMesh>
901 {
902  if (!fieldPrevIterPtr_)
903  {
905  << "Allocating previous iteration field" << nl
906  << this->info() << endl;
907 
908  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
909  (
910  this->name() + "PrevIter",
911  *this
912  );
913  }
914  else
915  {
916  *fieldPrevIterPtr_ == *this;
917  }
918 }
919 
920 
921 template<class Type, template<class> class PatchField, class GeoMesh>
924 {
925  if (!fieldPrevIterPtr_)
926  {
928  << "previous iteration field" << endl << this->info() << endl
929  << " not stored."
930  << " Use field.storePrevIter() at start of iteration."
931  << abort(FatalError);
932  }
933 
934  return *fieldPrevIterPtr_;
935 }
936 
937 
938 template<class Type, template<class> class PatchField, class GeoMesh>
941 {
942  this->setUpToDate();
943  storeOldTimes();
944  boundaryField_.evaluate();
945 }
946 
947 
948 template<class Type, template<class> class PatchField, class GeoMesh>
950 {
951  // Search all boundary conditions, if any are
952  // fixed-value or mixed (Robin) do not set reference level for solution.
953 
954  bool needRef = true;
955 
956  forAll(boundaryField_, patchi)
957  {
958  if (boundaryField_[patchi].fixesValue())
959  {
960  needRef = false;
961  break;
962  }
963  }
964 
965  reduce(needRef, andOp<bool>());
966 
967  return needRef;
968 }
969 
970 
971 template<class Type, template<class> class PatchField, class GeoMesh>
973 {
975  << "Relaxing" << nl << this->info() << " by " << alpha << endl;
976 
977  operator==(prevIter() + alpha*(*this - prevIter()));
978 }
979 
980 
981 template<class Type, template<class> class PatchField, class GeoMesh>
983 {
984  word name = this->name();
985 
986  if
987  (
988  this->mesh().data::template getOrDefault<bool>
989  (
990  "finalIteration",
991  false
992  )
993  )
994  {
995  name += "Final";
996  }
997 
998  if (this->mesh().relaxField(name))
999  {
1000  relax(this->mesh().fieldRelaxationFactor(name));
1001  }
1002 }
1003 
1004 
1005 template<class Type, template<class> class PatchField, class GeoMesh>
1008  bool final
1009 ) const
1010 {
1011  if (final)
1012  {
1013  return this->name() + "Final";
1014  }
1015 
1016  return this->name();
1017 }
1018 
1019 
1020 template<class Type, template<class> class PatchField, class GeoMesh>
1023  Ostream& os
1024 ) const
1025 {
1026  MinMax<Type> range = Foam::minMax(*this).value();
1027 
1028  os << "min/max(" << this->name() << ") = "
1029  << range.min() << ", " << range.max() << endl;
1030 }
1031 
1032 
1033 template<class Type, template<class> class PatchField, class GeoMesh>
1036 {
1037  os << *this;
1038  return os.good();
1039 }
1040 
1041 
1042 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1043 
1044 template<class Type, template<class> class PatchField, class GeoMesh>
1047 {
1049  (
1050  IOobject
1051  (
1052  this->name() + ".T()",
1053  this->instance(),
1054  this->db()
1055  ),
1056  this->mesh(),
1057  this->dimensions()
1058  );
1059 
1060  Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1061  Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1062 
1063  return tresult;
1064 }
1065 
1066 
1067 template<class Type, template<class> class PatchField, class GeoMesh>
1068 Foam::tmp
1069 <
1071  <
1073  PatchField,
1074  GeoMesh
1075  >
1076 >
1078 (
1079  const direction d
1080 ) const
1081 {
1083  (
1084  IOobject
1085  (
1086  this->name() + ".component(" + Foam::name(d) + ')',
1087  this->instance(),
1088  this->db()
1089  ),
1090  this->mesh(),
1091  this->dimensions()
1092  );
1093 
1094  Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1095  Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1096 
1097  return tresult;
1098 }
1099 
1100 
1101 template<class Type, template<class> class PatchField, class GeoMesh>
1103 (
1104  const direction d,
1105  const GeometricField
1106  <
1107  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1108  PatchField,
1109  GeoMesh
1110  >& gcf
1111 )
1112 {
1113  primitiveFieldRef().replace(d, gcf.primitiveField());
1114  boundaryFieldRef().replace(d, gcf.boundaryField());
1115 }
1116 
1117 
1118 template<class Type, template<class> class PatchField, class GeoMesh>
1121  const direction d,
1122  const dimensioned<cmptType>& ds
1123 )
1124 {
1125  primitiveFieldRef().replace(d, ds.value());
1126  boundaryFieldRef().replace(d, ds.value());
1127 }
1128 
1129 
1130 template<class Type, template<class> class PatchField, class GeoMesh>
1133  const dimensioned<Type>& dt
1134 )
1135 {
1136  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1137  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1138 }
1139 
1140 
1141 template<class Type, template<class> class PatchField, class GeoMesh>
1144  const dimensioned<Type>& dt
1145 )
1146 {
1147  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1148  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1149 }
1150 
1151 
1152 template<class Type, template<class> class PatchField, class GeoMesh>
1156 )
1157 {
1158  Foam::clip(primitiveFieldRef(), primitiveField(), range.value());
1159  Foam::clip(boundaryFieldRef(), boundaryField(), range.value());
1160 }
1161 
1162 
1163 template<class Type, template<class> class PatchField, class GeoMesh>
1166  const dimensioned<Type>& minVal,
1167  const dimensioned<Type>& maxVal
1168 )
1169 {
1170  MinMax<Type> range(minVal.value(), maxVal.value());
1171 
1172  Foam::clip(primitiveFieldRef(), primitiveField(), range);
1173  Foam::clip(boundaryFieldRef(), boundaryField(), range);
1174 }
1175 
1176 
1177 template<class Type, template<class> class PatchField, class GeoMesh>
1180  const dimensioned<Type>& minVal,
1181  const dimensioned<Type>& maxVal
1182 )
1183 {
1184  this->clip(minVal, maxVal);
1185 }
1186 
1187 
1188 template<class Type, template<class> class PatchField, class GeoMesh>
1190 {
1191  primitiveFieldRef().negate();
1192  boundaryFieldRef().negate();
1193 }
1194 
1195 
1196 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1197 
1198 template<class Type, template<class> class PatchField, class GeoMesh>
1202 )
1203 {
1204  if (this == &gf)
1205  {
1206  return; // Self-assignment is a no-op
1207  }
1208 
1209  checkField(*this, gf, "=");
1210 
1211  // Only assign field contents not ID
1212 
1213  ref() = gf();
1214  boundaryFieldRef() = gf.boundaryField();
1215 }
1216 
1217 
1218 template<class Type, template<class> class PatchField, class GeoMesh>
1222 )
1223 {
1224  const auto& gf = tgf();
1225 
1226  if (this == &gf)
1227  {
1228  return; // Self-assignment is a no-op
1229  }
1230 
1231  checkField(*this, gf, "=");
1232 
1233  // Only assign field contents not ID
1234 
1235  this->dimensions() = gf.dimensions();
1236  this->oriented() = gf.oriented();
1237 
1238  if (tgf.movable())
1239  {
1240  // Transfer storage from the tmp
1241  primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1242  }
1243  else
1244  {
1246  }
1247 
1248  boundaryFieldRef() = gf.boundaryField();
1249 
1250  tgf.clear();
1251 }
1252 
1253 
1254 template<class Type, template<class> class PatchField, class GeoMesh>
1257  const dimensioned<Type>& dt
1258 )
1259 {
1260  ref() = dt;
1261  boundaryFieldRef() = dt.value();
1262 }
1263 
1264 
1265 template<class Type, template<class> class PatchField, class GeoMesh>
1269 )
1270 {
1271  const auto& gf = tgf();
1272 
1273  checkField(*this, gf, "==");
1274 
1275  // Only assign field contents not ID
1276 
1277  ref() = gf();
1278  boundaryFieldRef() == gf.boundaryField();
1279 
1280  tgf.clear();
1281 }
1282 
1283 
1284 template<class Type, template<class> class PatchField, class GeoMesh>
1287  const dimensioned<Type>& dt
1288 )
1289 {
1290  ref() = dt;
1291  boundaryFieldRef() == dt.value();
1292 }
1293 
1294 
1295 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1296  \
1297 template<class Type, template<class> class PatchField, class GeoMesh> \
1298 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1299 ( \
1300  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1301 ) \
1302 { \
1303  checkField(*this, gf, #op); \
1304  \
1305  ref() op gf(); \
1306  boundaryFieldRef() op gf.boundaryField(); \
1307 } \
1308  \
1309 template<class Type, template<class> class PatchField, class GeoMesh> \
1310 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1311 ( \
1312  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1313 ) \
1314 { \
1315  operator op(tgf()); \
1316  tgf.clear(); \
1317 } \
1318  \
1319 template<class Type, template<class> class PatchField, class GeoMesh> \
1320 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1321 ( \
1322  const dimensioned<TYPE>& dt \
1323 ) \
1324 { \
1325  ref() op dt; \
1326  boundaryFieldRef() op dt.value(); \
1327 }
1328 
1333 
1334 #undef COMPUTED_ASSIGNMENT
1335 
1336 
1337 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1338 
1339 template<class Type, template<class> class PatchField, class GeoMesh>
1340 Foam::Ostream& Foam::operator<<
1341 (
1342  Ostream& os,
1344 )
1345 {
1346  gf().writeData(os, "internalField");
1347  os << nl;
1348  gf.boundaryField().writeEntry("boundaryField", os);
1349 
1350  os.check(FUNCTION_NAME);
1351  return os;
1352 }
1353 
1354 
1355 template<class Type, template<class> class PatchField, class GeoMesh>
1356 Foam::Ostream& Foam::operator<<
1357 (
1358  Ostream& os,
1360 )
1361 {
1362  os << tgf();
1363  tgf.clear();
1364  return os;
1365 }
1366 
1367 
1368 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1369 
1370 #undef checkField
1371 
1372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1373 
1374 #include "GeometricFieldNew.C"
1375 #include "GeometricBoundaryField.C"
1376 #include "GeometricFieldFunctions.C"
1377 
1378 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::GeometricField::replace
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::GeometricField::writeData
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Definition: GeometricField.C:1035
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::MinMax::min
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
COMPUTED_ASSIGNMENT
#define COMPUTED_ASSIGNMENT(TYPE, op)
Definition: GeometricField.C:1295
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:850
checkField
#define checkField(gf1, gf2, op)
Definition: GeometricField.C:38
Foam::GeometricField::GeometricField
GeometricField(const IOobject &io, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=PatchField< Type >::calculatedType())
Construct given IOobject, mesh, dimensions and patch type.
Definition: GeometricField.C:197
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
ref
rDeltaT ref()
Foam::GeometricField::select
word select(bool final) const
Select the final iteration parameters if `final' is true.
Definition: GeometricField.C:1007
primitiveFieldRef
conserve primitiveFieldRef()+
Foam::GeometricField::clip
void clip(const dimensioned< MinMax< Type >> &range)
Clip the field to be bounded within the specified range.
Definition: GeometricField.C:1154
Foam::GeometricField::cmptType
Field< Type >::cmptType cmptType
Definition: GeometricField.H:290
Foam::GeometricField::min
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
Definition: GeometricField.C:1132
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::GeometricField::needReference
bool needReference() const
Does the field need a reference level for solution.
Definition: GeometricField.C:949
Foam::GeometricField::maxMin
void maxMin(const dimensioned< Type > &minVal, const dimensioned< Type > &maxVal)
Deprecated(2019-01) identical to clip()
Definition: GeometricField.C:1179
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::GeometricField::nOldTimes
label nOldTimes() const
Return the number of old time fields stored.
Definition: GeometricField.C:837
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::GeometricField::storeOldTimes
void storeOldTimes() const
Store the old-time fields.
Definition: GeometricField.C:797
Foam::GeometricField< Type, fvPatchField, volMesh >::Mesh
volMesh ::Mesh Mesh
Type of mesh on which this GeometricField is instantiated.
Definition: GeometricField.H:100
Foam::GeometricField::negate
void negate()
Negate the field inplace.
Definition: GeometricField.C:1189
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::GeometricField::writeMinMax
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
Definition: GeometricField.C:1022
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::GeometricField::prevIter
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
Definition: GeometricField.C:923
Foam::andOp
Definition: ops.H:233
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::GeometricField::storeOldTime
void storeOldTime() const
Store the old-time field.
Definition: GeometricField.C:816
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
reduce
reduce(hasMovingMesh, orOp< bool >())
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:48
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::GeometricField::timeIndex
label timeIndex() const
Return the time index of the field.
Definition: GeometricFieldI.H:70
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::GeometricField::clone
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
Definition: GeometricField.C:727
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::GeometricField::T
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Definition: GeometricField.C:1046
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
GeometricField.H
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Time.H
Foam::GeometricField::storePrevIter
void storePrevIter() const
Store the field as the previous iteration value.
Definition: GeometricField.C:900
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::GeometricField::Boundary
The boundary fields.
Definition: GeometricField.H:115
Foam::clip
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:278
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::List< word >
Foam::GeometricField::relax
void relax()
Relax field (for steady-state solution).
Definition: GeometricField.C:982
dictionary.H
Foam::direction
uint8_t direction
Definition: direction.H:52
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::GeometricField::~GeometricField
virtual ~GeometricField()
Destructor.
Definition: GeometricField.C:736
localIOdictionary.H
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::GeometricField::max
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Definition: GeometricField.C:1143
data.H
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
Foam::minMax
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:61
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
GeometricFieldNew.C
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
GeometricBoundaryField.C
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::MinMax
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition: HashSet.H:76
GeometricFieldFunctions.C
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
relax
UEqn relax()