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-2019 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 wordList& patchFieldTypes,
669  const wordList& actualPatchTypes
670 )
671 :
672  Internal(io, tgf.constCast(), tgf.movable()),
673  timeIndex_(tgf().timeIndex()),
674  field0Ptr_(nullptr),
675  fieldPrevIterPtr_(nullptr),
676  boundaryField_
677  (
678  this->mesh().boundary(),
679  *this,
680  patchFieldTypes,
681  actualPatchTypes
682  )
683 {
685  << "Constructing from tmp resetting IO params and patch types" << nl
686  << this->info() << endl;
687 
688  boundaryField_ == tgf().boundaryField_;
689 
690  tgf.clear();
691 }
692 
693 
694 template<class Type, template<class> class PatchField, class GeoMesh>
697 {
699 }
700 
701 
702 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
703 
704 template<class Type, template<class> class PatchField, class GeoMesh>
706 {
707  deleteDemandDrivenData(field0Ptr_);
708  deleteDemandDrivenData(fieldPrevIterPtr_);
709 }
710 
711 
712 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
713 
714 template<class Type, template<class> class PatchField, class GeoMesh>
715 typename
718 (
719  const bool updateAccessTime
720 )
721 {
722  if (updateAccessTime)
723  {
724  this->setUpToDate();
725  storeOldTimes();
726  }
727  return *this;
728 }
729 
730 
731 template<class Type, template<class> class PatchField, class GeoMesh>
732 typename
735 (
736  const bool updateAccessTime
737 )
738 {
739  if (updateAccessTime)
740  {
741  this->setUpToDate();
742  storeOldTimes();
743  }
744  return *this;
745 }
746 
747 
748 template<class Type, template<class> class PatchField, class GeoMesh>
749 typename
752 (
753  const bool updateAccessTime
754 )
755 {
756  if (updateAccessTime)
757  {
758  this->setUpToDate();
759  storeOldTimes();
760  }
761  return boundaryField_;
762 }
763 
764 
765 template<class Type, template<class> class PatchField, class GeoMesh>
767 {
768  if
769  (
770  field0Ptr_
771  && timeIndex_ != this->time().timeIndex()
772  && !this->name().ends_with("_0")
773  )
774  {
775  storeOldTime();
776  timeIndex_ = this->time().timeIndex();
777  }
778 
779  // Correct time index
780  //timeIndex_ = this->time().timeIndex();
781 }
782 
783 
784 template<class Type, template<class> class PatchField, class GeoMesh>
786 {
787  if (field0Ptr_)
788  {
789  field0Ptr_->storeOldTime();
790 
792  << "Storing old time field for field" << nl << this->info() << endl;
793 
794  *field0Ptr_ == *this;
795  field0Ptr_->timeIndex_ = timeIndex_;
796 
797  if (field0Ptr_->field0Ptr_)
798  {
799  field0Ptr_->writeOpt() = this->writeOpt();
800  }
801  }
802 }
803 
804 
805 template<class Type, template<class> class PatchField, class GeoMesh>
807 {
808  if (field0Ptr_)
809  {
810  return field0Ptr_->nOldTimes() + 1;
811  }
812 
813  return 0;
814 }
815 
816 
817 template<class Type, template<class> class PatchField, class GeoMesh>
820 {
821  if (!field0Ptr_)
822  {
824  (
825  IOobject
826  (
827  this->name() + "_0",
828  this->time().timeName(),
829  this->db(),
830  IOobject::NO_READ,
831  IOobject::NO_WRITE,
832  this->registerObject()
833  ),
834  *this
835  );
836 
837  if (debug)
838  {
840  << "created old time field " << field0Ptr_->info() << endl;
841 
842  if (debug&2)
843  {
844  error::printStack(Info);
845  }
846  }
847  }
848  else
849  {
850  storeOldTimes();
851  }
852 
853  return *field0Ptr_;
854 }
855 
856 
857 template<class Type, template<class> class PatchField, class GeoMesh>
860 {
861  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
862  .oldTime();
863 
864  return *field0Ptr_;
865 }
866 
867 
868 template<class Type, template<class> class PatchField, class GeoMesh>
870 {
871  if (!fieldPrevIterPtr_)
872  {
874  << "Allocating previous iteration field" << nl
875  << this->info() << endl;
876 
877  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
878  (
879  this->name() + "PrevIter",
880  *this
881  );
882  }
883  else
884  {
885  *fieldPrevIterPtr_ == *this;
886  }
887 }
888 
889 
890 template<class Type, template<class> class PatchField, class GeoMesh>
893 {
894  if (!fieldPrevIterPtr_)
895  {
897  << "previous iteration field" << endl << this->info() << endl
898  << " not stored."
899  << " Use field.storePrevIter() at start of iteration."
900  << abort(FatalError);
901  }
902 
903  return *fieldPrevIterPtr_;
904 }
905 
906 
907 template<class Type, template<class> class PatchField, class GeoMesh>
910 {
911  this->setUpToDate();
912  storeOldTimes();
913  boundaryField_.evaluate();
914 }
915 
916 
917 template<class Type, template<class> class PatchField, class GeoMesh>
919 {
920  // Search all boundary conditions, if any are
921  // fixed-value or mixed (Robin) do not set reference level for solution.
922 
923  bool needRef = true;
924 
925  forAll(boundaryField_, patchi)
926  {
927  if (boundaryField_[patchi].fixesValue())
928  {
929  needRef = false;
930  break;
931  }
932  }
933 
934  reduce(needRef, andOp<bool>());
935 
936  return needRef;
937 }
938 
939 
940 template<class Type, template<class> class PatchField, class GeoMesh>
942 {
944  << "Relaxing" << nl << this->info() << " by " << alpha << endl;
945 
946  operator==(prevIter() + alpha*(*this - prevIter()));
947 }
948 
949 
950 template<class Type, template<class> class PatchField, class GeoMesh>
952 {
953  word name = this->name();
954 
955  if
956  (
957  this->mesh().data::template lookupOrDefault<bool>
958  (
959  "finalIteration",
960  false
961  )
962  )
963  {
964  name += "Final";
965  }
966 
967  if (this->mesh().relaxField(name))
968  {
969  relax(this->mesh().fieldRelaxationFactor(name));
970  }
971 }
972 
973 
974 template<class Type, template<class> class PatchField, class GeoMesh>
976 (
977  bool final
978 ) const
979 {
980  if (final)
981  {
982  return this->name() + "Final";
983  }
984 
985  return this->name();
986 }
987 
988 
989 template<class Type, template<class> class PatchField, class GeoMesh>
991 (
992  Ostream& os
993 ) const
994 {
995  MinMax<Type> range = Foam::minMax(*this).value();
996 
997  os << "min/max(" << this->name() << ") = "
998  << range.min() << ", " << range.max() << endl;
999 }
1000 
1001 
1002 template<class Type, template<class> class PatchField, class GeoMesh>
1004 writeData(Ostream& os) const
1005 {
1006  os << *this;
1007  return os.good();
1008 }
1009 
1010 
1011 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1012 
1013 template<class Type, template<class> class PatchField, class GeoMesh>
1016 {
1018  (
1019  IOobject
1020  (
1021  this->name() + ".T()",
1022  this->instance(),
1023  this->db()
1024  ),
1025  this->mesh(),
1026  this->dimensions()
1027  );
1028 
1029  Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1030  Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1031 
1032  return tresult;
1033 }
1034 
1035 
1036 template<class Type, template<class> class PatchField, class GeoMesh>
1037 Foam::tmp
1038 <
1040  <
1042  PatchField,
1043  GeoMesh
1044  >
1045 >
1047 (
1048  const direction d
1049 ) const
1050 {
1052  (
1053  IOobject
1054  (
1055  this->name() + ".component(" + Foam::name(d) + ')',
1056  this->instance(),
1057  this->db()
1058  ),
1059  this->mesh(),
1060  this->dimensions()
1061  );
1062 
1063  Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1064  Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1065 
1066  return tresult;
1067 }
1068 
1069 
1070 template<class Type, template<class> class PatchField, class GeoMesh>
1072 (
1073  const direction d,
1074  const GeometricField
1075  <
1076  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1077  PatchField,
1078  GeoMesh
1079  >& gcf
1080 )
1081 {
1082  primitiveFieldRef().replace(d, gcf.primitiveField());
1083  boundaryFieldRef().replace(d, gcf.boundaryField());
1084 }
1085 
1086 
1087 template<class Type, template<class> class PatchField, class GeoMesh>
1090  const direction d,
1091  const dimensioned<cmptType>& ds
1092 )
1093 {
1094  primitiveFieldRef().replace(d, ds.value());
1095  boundaryFieldRef().replace(d, ds.value());
1096 }
1097 
1098 
1099 template<class Type, template<class> class PatchField, class GeoMesh>
1102  const dimensioned<Type>& dt
1103 )
1104 {
1105  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1106  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1107 }
1108 
1109 
1110 template<class Type, template<class> class PatchField, class GeoMesh>
1113  const dimensioned<Type>& dt
1114 )
1115 {
1116  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1117  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1118 }
1119 
1120 
1121 template<class Type, template<class> class PatchField, class GeoMesh>
1125 )
1126 {
1127  Foam::clip(primitiveFieldRef(), primitiveField(), range.value());
1128  Foam::clip(boundaryFieldRef(), boundaryField(), range.value());
1129 }
1130 
1131 
1132 template<class Type, template<class> class PatchField, class GeoMesh>
1135  const dimensioned<Type>& minVal,
1136  const dimensioned<Type>& maxVal
1137 )
1138 {
1139  MinMax<Type> range(minVal.value(), maxVal.value());
1140 
1141  Foam::clip(primitiveFieldRef(), primitiveField(), range);
1142  Foam::clip(boundaryFieldRef(), boundaryField(), range);
1143 }
1144 
1145 
1146 template<class Type, template<class> class PatchField, class GeoMesh>
1149  const dimensioned<Type>& minVal,
1150  const dimensioned<Type>& maxVal
1151 )
1152 {
1153  this->clip(minVal, maxVal);
1154 }
1155 
1156 
1157 template<class Type, template<class> class PatchField, class GeoMesh>
1159 {
1160  primitiveFieldRef().negate();
1161  boundaryFieldRef().negate();
1162 }
1163 
1164 
1165 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1166 
1167 template<class Type, template<class> class PatchField, class GeoMesh>
1171 )
1172 {
1173  if (this == &gf)
1174  {
1175  return; // Self-assignment is a no-op
1176  }
1177 
1178  checkField(*this, gf, "=");
1179 
1180  // Only assign field contents not ID
1181 
1182  ref() = gf();
1183  boundaryFieldRef() = gf.boundaryField();
1184 }
1185 
1186 
1187 template<class Type, template<class> class PatchField, class GeoMesh>
1191 )
1192 {
1193  const auto& gf = tgf();
1194 
1195  if (this == &gf)
1196  {
1197  return; // Self-assignment is a no-op
1198  }
1199 
1200  checkField(*this, gf, "=");
1201 
1202  // Only assign field contents not ID
1203 
1204  this->dimensions() = gf.dimensions();
1205  this->oriented() = gf.oriented();
1206 
1207  if (tgf.movable())
1208  {
1209  // Transfer storage from the tmp
1210  primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1211  }
1212  else
1213  {
1215  }
1216 
1217  boundaryFieldRef() = gf.boundaryField();
1218 
1219  tgf.clear();
1220 }
1221 
1222 
1223 template<class Type, template<class> class PatchField, class GeoMesh>
1226  const dimensioned<Type>& dt
1227 )
1228 {
1229  ref() = dt;
1230  boundaryFieldRef() = dt.value();
1231 }
1232 
1233 
1234 template<class Type, template<class> class PatchField, class GeoMesh>
1238 )
1239 {
1240  const auto& gf = tgf();
1241 
1242  checkField(*this, gf, "==");
1243 
1244  // Only assign field contents not ID
1245 
1246  ref() = gf();
1247  boundaryFieldRef() == gf.boundaryField();
1248 
1249  tgf.clear();
1250 }
1251 
1252 
1253 template<class Type, template<class> class PatchField, class GeoMesh>
1256  const dimensioned<Type>& dt
1257 )
1258 {
1259  ref() = dt;
1260  boundaryFieldRef() == dt.value();
1261 }
1262 
1263 
1264 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1265  \
1266 template<class Type, template<class> class PatchField, class GeoMesh> \
1267 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1268 ( \
1269  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1270 ) \
1271 { \
1272  checkField(*this, gf, #op); \
1273  \
1274  ref() op gf(); \
1275  boundaryFieldRef() op gf.boundaryField(); \
1276 } \
1277  \
1278 template<class Type, template<class> class PatchField, class GeoMesh> \
1279 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1280 ( \
1281  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1282 ) \
1283 { \
1284  operator op(tgf()); \
1285  tgf.clear(); \
1286 } \
1287  \
1288 template<class Type, template<class> class PatchField, class GeoMesh> \
1289 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1290 ( \
1291  const dimensioned<TYPE>& dt \
1292 ) \
1293 { \
1294  ref() op dt; \
1295  boundaryFieldRef() op dt.value(); \
1296 }
1297 
1302 
1303 #undef COMPUTED_ASSIGNMENT
1304 
1305 
1306 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1307 
1308 template<class Type, template<class> class PatchField, class GeoMesh>
1309 Foam::Ostream& Foam::operator<<
1310 (
1311  Ostream& os,
1313 )
1314 {
1315  gf().writeData(os, "internalField");
1316  os << nl;
1317  gf.boundaryField().writeEntry("boundaryField", os);
1318 
1319  os.check(FUNCTION_NAME);
1320  return os;
1321 }
1322 
1323 
1324 template<class Type, template<class> class PatchField, class GeoMesh>
1325 Foam::Ostream& Foam::operator<<
1326 (
1327  Ostream& os,
1329 )
1330 {
1331  os << tgf();
1332  tgf.clear();
1333  return os;
1334 }
1335 
1336 
1337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1338 
1339 #undef checkField
1340 
1341 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1342 
1343 #include "GeometricFieldNew.C"
1344 #include "GeometricBoundaryField.C"
1345 #include "GeometricFieldFunctions.C"
1346 
1347 // ************************************************************************* //
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:104
ref
rDeltaT ref()
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:316
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::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
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:62
Foam::GeometricField::writeData
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Definition: GeometricField.C:1004
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
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:1264
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:819
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:337
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::GeometricField::select
word select(bool final) const
Select the final iteration parameters if `final' is true.
Definition: GeometricField.C:976
primitiveFieldRef
conserve primitiveFieldRef()+
Foam::GeometricField::cmptType
Field< Type >::cmptType cmptType
Definition: GeometricField.H:285
Foam::GeometricField::min
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
Definition: GeometricField.C:1101
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:918
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::GeometricField::nOldTimes
label nOldTimes() const
Return the number of old time fields stored.
Definition: GeometricField.C:806
Foam::GeometricField::maxMin
void maxMin(const dimensioned< Type > &minVal, const dimensioned< Type > &maxVal) FOAM_DEPRECATED_FOR(2019-01
Deprecated(2019-01) identical to clip()
Definition: GeometricField.C:1148
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::GeometricField::storeOldTimes
void storeOldTimes() const
Store the old-time fields.
Definition: GeometricField.C:766
Foam::GeometricField< Type, fvPatchField, volMesh >::Mesh
volMesh ::Mesh Mesh
Type of mesh on which this GeometricField is instantiated.
Definition: GeometricField.H:100
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::GeometricField::negate
void negate()
Negate the field inplace.
Definition: GeometricField.C:1158
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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 (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:356
Foam::GeometricField::writeMinMax
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
Definition: GeometricField.C:991
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:892
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:65
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:785
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
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:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
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:137
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:735
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:696
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:909
Foam::GeometricField::T
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Definition: GeometricField.C:1015
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
Time.H
Foam::GeometricField::storePrevIter
void storePrevIter() const
Store the field as the previous iteration value.
Definition: GeometricField.C:869
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::GeometricField::Boundary
Definition: GeometricField.H:114
Foam::clip
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:276
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::GeometricField::clip
void clip() method")
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:718
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::List< word >
Foam::GeometricField::relax
void relax()
Relax field (for steady-state solution).
Definition: GeometricField.C:951
dictionary.H
Foam::direction
uint8_t direction
Definition: direction.H:47
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:261
Foam::GeometricField::~GeometricField
virtual ~GeometricField()
Destructor.
Definition: GeometricField.C:705
localIOdictionary.H
Foam::GeometricField::max
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Definition: GeometricField.C:1112
data.H
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
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
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:216
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:294
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
relax
UEqn relax()