GeometricField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "GeometricField.H"
27 #include "Time.H"
28 #include "demandDrivenData.H"
29 #include "dictionary.H"
30 #include "localIOdictionary.H"
31 #include "solutionControl.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 #define checkField(gf1, gf2, op) \
36 if ((gf1).mesh() != (gf2).mesh()) \
37 { \
38  FatalErrorInFunction \
39  << "different mesh for fields " \
40  << (gf1).name() << " and " << (gf2).name() \
41  << " during operatrion " << op \
42  << abort(FatalError); \
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47 
48 template<class Type, template<class> class PatchField, class GeoMesh>
50 (
51  const dictionary& dict
52 )
53 {
54  Internal::readField(dict, "internalField");
55 
56  boundaryField_.readField(*this, dict.subDict("boundaryField"));
57 
58  // Don't use subOrEmptyDict here, or line numbers will be lost from any IO
59  // error messages. Use an actual sub-dict reference here if possible.
60  if (dict.found("sources"))
61  {
62  sources_.readField(*this, dict.subDict("sources"));
63  }
64  else
65  {
66  sources_.readField(*this, dictionary(dict, dictionary()));
67  }
68 
69  if (dict.found("referenceLevel"))
70  {
71  Type fieldAverage(pTraits<Type>(dict.lookup("referenceLevel")));
72 
73  Field<Type>::operator+=(fieldAverage);
74 
75  forAll(boundaryField_, patchi)
76  {
77  boundaryField_[patchi] == boundaryField_[patchi] + fieldAverage;
78  }
79  }
80 }
81 
82 
83 template<class Type, template<class> class PatchField, class GeoMesh>
85 {
86  const localIOdictionary dict
87  (
88  IOobject
89  (
90  this->name(),
91  this->instance(),
92  this->local(),
93  this->db(),
94  IOobject::MUST_READ,
95  IOobject::NO_WRITE,
96  false
97  ),
98  typeName
99  );
100 
101  this->close();
102 
103  readFields(dict);
104 }
105 
106 
107 template<class Type, template<class> class PatchField, class GeoMesh>
109 {
110  if
111  (
112  this->readOpt() == IOobject::MUST_READ
113  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
114  )
115  {
117  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
118  << " suggests that a read constructor for field " << this->name()
119  << " would be more appropriate." << endl;
120  }
121  else if
122  (
123  this->readOpt() == IOobject::READ_IF_PRESENT
124  && this->headerOk()
125  )
126  {
127  readFields();
128 
129  // Check compatibility between field and mesh
130  if (this->size() != GeoMesh::size(this->mesh()))
131  {
132  FatalIOErrorInFunction(this->readStream(typeName))
133  << " number of field elements = " << this->size()
134  << " number of mesh elements = "
135  << GeoMesh::size(this->mesh())
136  << exit(FatalIOError);
137  }
138 
139  readOldTimeIfPresent();
140 
141  return true;
142  }
143 
144  return false;
145 }
146 
147 
148 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
149 
150 template<class Type, template<class> class PatchField, class GeoMesh>
152 (
153  const IOobject& io,
154  const Mesh& mesh,
155  const dimensionSet& ds,
156  const word& patchFieldType
157 )
158 :
159  Internal(io, mesh, ds, false),
160  OldTimeField<GeometricField>(this->time().timeIndex()),
161  fieldPrevIterPtr_(nullptr),
162  boundaryField_(mesh.boundary(), *this, patchFieldType),
163  sources_()
164 {
165  if (debug)
166  {
167  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
168  }
169 
170  readIfPresent();
171 }
172 
173 
174 template<class Type, template<class> class PatchField, class GeoMesh>
176 (
177  const IOobject& io,
178  const Mesh& mesh,
179  const dimensionSet& ds,
180  const wordList& patchFieldTypes,
181  const wordList& actualPatchTypes,
182  const HashTable<word>& fieldSourceTypes
183 )
184 :
185  Internal(io, mesh, ds, false),
186  OldTimeField<GeometricField>(this->time().timeIndex()),
187  fieldPrevIterPtr_(nullptr),
188  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes),
189  sources_(*this, fieldSourceTypes)
190 {
191  if (debug)
192  {
193  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
194  }
195 
196  readIfPresent();
197 }
198 
199 
200 template<class Type, template<class> class PatchField, class GeoMesh>
202 (
203  const IOobject& io,
204  const Mesh& mesh,
205  const dimensioned<Type>& dt,
206  const word& patchFieldType
207 )
208 :
209  Internal(io, mesh, dt, false),
211  fieldPrevIterPtr_(nullptr),
212  boundaryField_(mesh.boundary(), *this, patchFieldType),
213  sources_()
214 {
215  if (debug)
216  {
217  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
218  }
219 
220  boundaryField_ == dt.value();
221 
222  readIfPresent();
223 }
224 
225 
226 template<class Type, template<class> class PatchField, class GeoMesh>
228 (
229  const IOobject& io,
230  const Mesh& mesh,
231  const dimensioned<Type>& dt,
232  const wordList& patchFieldTypes,
233  const wordList& actualPatchTypes,
234  const HashTable<word>& fieldSourceTypes
235 )
236 :
237  Internal(io, mesh, dt, false),
238  OldTimeField<GeometricField>(this->time().timeIndex()),
239  fieldPrevIterPtr_(nullptr),
240  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes),
241  sources_(*this, fieldSourceTypes)
242 {
243  if (debug)
244  {
245  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
246  }
247 
248  boundaryField_ == dt.value();
249 
250  readIfPresent();
251 }
252 
253 
254 template<class Type, template<class> class PatchField, class GeoMesh>
256 (
257  const IOobject& io,
258  const Internal& diField,
259  const PtrList<PatchField<Type>>& ptfl,
260  const HashPtrTable<Source>& stft
261 )
262 :
263  Internal(io, diField),
264  OldTimeField<GeometricField>(this->time().timeIndex()),
265  fieldPrevIterPtr_(nullptr),
266  boundaryField_(this->mesh().boundary(), *this, ptfl),
267  sources_(*this, stft)
268 {
269  if (debug)
270  {
272  << "Constructing from components" << endl << this->info() << endl;
273  }
274 
275  readIfPresent();
276 }
277 
278 
279 template<class Type, template<class> class PatchField, class GeoMesh>
281 (
282  const IOobject& io,
283  const Mesh& mesh,
284  const dimensionSet& ds,
285  const Field<Type>& iField,
286  const PtrList<PatchField<Type>>& ptfl,
287  const HashPtrTable<Source>& stft
288 )
289 :
290  Internal(io, mesh, ds, iField),
292  fieldPrevIterPtr_(nullptr),
293  boundaryField_(mesh.boundary(), *this, ptfl),
294  sources_(*this, stft)
295 {
296  if (debug)
297  {
299  << "Constructing from components" << endl << this->info() << endl;
300  }
301 
302  readIfPresent();
303 }
304 
305 
306 template<class Type, template<class> class PatchField, class GeoMesh>
308 (
309  const IOobject& io,
310  const Mesh& mesh
311 )
312 :
313  Internal(io, mesh, dimless, false),
314  OldTimeField<GeometricField>(this->time().timeIndex()),
315  fieldPrevIterPtr_(nullptr),
316  boundaryField_(mesh.boundary()),
317  sources_()
318 {
319  readFields();
320 
321  // Check compatibility between field and mesh
322 
323  if (this->size() != GeoMesh::size(this->mesh()))
324  {
325  FatalIOErrorInFunction(this->readStream(typeName))
326  << " number of field elements = " << this->size()
327  << " number of mesh elements = " << GeoMesh::size(this->mesh())
329  }
330 
331  readOldTimeIfPresent();
332 
333  if (debug)
334  {
336  << "Finishing read-construction of" << endl << this->info() << endl;
337  }
338 }
339 
340 
341 template<class Type, template<class> class PatchField, class GeoMesh>
343 (
344  const IOobject& io,
345  const Mesh& mesh,
346  const dictionary& dict
347 )
348 :
349  Internal(io, mesh, dimless, false),
350  OldTimeField<GeometricField>(this->time().timeIndex()),
351  fieldPrevIterPtr_(nullptr),
352  boundaryField_(mesh.boundary()),
353  sources_()
354 {
355  readFields(dict);
356 
357  // Check compatibility between field and mesh
358 
359  if (this->size() != GeoMesh::size(this->mesh()))
360  {
362  << " number of field elements = " << this->size()
363  << " number of mesh elements = " << GeoMesh::size(this->mesh())
365  }
366 
367  if (debug)
368  {
370  << "Finishing dictionary-construct of "
371  << endl << this->info() << endl;
372  }
373 }
374 
375 
376 template<class Type, template<class> class PatchField, class GeoMesh>
378 (
380 )
381 :
382  Internal(gf),
384  fieldPrevIterPtr_(nullptr),
385  boundaryField_(*this, gf.boundaryField_),
386  sources_(*this, gf.sources_)
387 {
388  if (debug)
389  {
391  << "Constructing as copy" << endl << this->info() << endl;
392  }
393 
394  this->writeOpt() = IOobject::NO_WRITE;
395 }
396 
397 
398 template<class Type, template<class> class PatchField, class GeoMesh>
400 (
402 )
403 :
404  Internal(move(gf)),
405  OldTimeField<GeometricField>(move(gf)),
406  fieldPrevIterPtr_(nullptr),
407  boundaryField_(*this, gf.boundaryField_),
408  sources_(*this, gf.sources_)
409 {
410  if (debug)
411  {
413  << "Constructing by moving" << endl << this->info() << endl;
414  }
415 
416  this->writeOpt() = IOobject::NO_WRITE;
417 }
418 
419 
420 template<class Type, template<class> class PatchField, class GeoMesh>
422 (
424 )
425 :
426  Internal
427  (
428  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
429  tgf.isTmp()
430  ),
432  fieldPrevIterPtr_(nullptr),
433  boundaryField_(*this, tgf().boundaryField_),
434  sources_(*this, tgf().sources_)
435 {
436  if (debug)
437  {
439  << "Constructing from tmp" << endl << this->info() << endl;
440  }
441 
442  this->writeOpt() = IOobject::NO_WRITE;
443 
444  tgf.clear();
445 }
446 
447 
448 template<class Type, template<class> class PatchField, class GeoMesh>
450 (
451  const IOobject& io,
453 )
454 :
455  Internal(io, gf, false),
457  fieldPrevIterPtr_(nullptr),
458  boundaryField_(*this, gf.boundaryField_),
459  sources_(*this, gf.sources_)
460 {
461  if (debug)
462  {
464  << "Constructing as copy resetting IO params"
465  << endl << this->info() << endl;
466  }
467 
468  if (!readIfPresent())
469  {
470  copyOldTimes(io, gf);
471  }
472 }
473 
474 
475 template<class Type, template<class> class PatchField, class GeoMesh>
477 (
478  const IOobject& io,
480 )
481 :
482  Internal
483  (
484  io,
485  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
486  tgf.isTmp(),
487  false
488  ),
490  fieldPrevIterPtr_(nullptr),
491  boundaryField_(*this, tgf().boundaryField_),
492  sources_(*this, tgf().sources_)
493 {
494  if (debug)
495  {
497  << "Constructing from tmp resetting IO params"
498  << endl << this->info() << endl;
499  }
500 
501  tgf.clear();
502 
503  readIfPresent();
504 }
505 
506 
507 template<class Type, template<class> class PatchField, class GeoMesh>
509 (
510  const word& newName,
512 )
513 :
514  Internal(newName, gf),
516  fieldPrevIterPtr_(nullptr),
517  boundaryField_(*this, gf.boundaryField_),
518  sources_(*this, gf.sources_)
519 {
520  if (debug)
521  {
523  << "Constructing as copy resetting name"
524  << endl << this->info() << endl;
525  }
526 
527  copyOldTimes(newName, gf);
528 }
529 
530 
531 template<class Type, template<class> class PatchField, class GeoMesh>
533 (
534  const word& newName,
536 )
537 :
538  Internal
539  (
540  newName,
541  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
542  tgf.isTmp()
543  ),
545  fieldPrevIterPtr_(nullptr),
546  boundaryField_(*this, tgf().boundaryField_),
547  sources_(*this, tgf().sources_)
548 {
549  if (debug)
550  {
552  << "Constructing from tmp resetting name"
553  << endl << this->info() << endl;
554  }
555 
556  tgf.clear();
557 }
558 
559 
560 template<class Type, template<class> class PatchField, class GeoMesh>
562 (
563  const IOobject& io,
565  const word& patchFieldType
566 )
567 :
568  Internal(io, gf, false),
570  fieldPrevIterPtr_(nullptr),
571  boundaryField_(this->mesh().boundary(), *this, patchFieldType),
572  sources_(*this, gf.sources_)
573 {
574  if (debug)
575  {
577  << "Constructing as copy resetting IO params"
578  << endl << this->info() << endl;
579  }
580 
581  boundaryField_ == gf.boundaryField_;
582 
583  if (!readIfPresent())
584  {
585  copyOldTimes(io, gf);
586  }
587 }
588 
589 
590 template<class Type, template<class> class PatchField, class GeoMesh>
592 (
593  const IOobject& io,
595  const wordList& patchFieldTypes,
596  const wordList& actualPatchTypes,
597  const HashTable<word>& fieldSourceTypes
598 )
599 :
600  Internal(io, gf, false),
602  fieldPrevIterPtr_(nullptr),
603  boundaryField_
604  (
605  this->mesh().boundary(),
606  *this,
607  patchFieldTypes,
608  actualPatchTypes
609  ),
610  sources_(*this, fieldSourceTypes)
611 {
612  if (debug)
613  {
615  << "Constructing as copy resetting IO params and patch types"
616  << endl << this->info() << endl;
617  }
618 
619  boundaryField_ == gf.boundaryField_;
620 
621  if (!readIfPresent())
622  {
623  copyOldTimes(io, gf);
624  }
625 }
626 
627 
628 template<class Type, template<class> class PatchField, class GeoMesh>
630 (
631  const IOobject& io,
633  const wordList& patchFieldTypes,
634  const wordList& actualPatchTypes,
635  const HashTable<word>& fieldSourceTypes
636 )
637 :
638  Internal
639  (
640  io,
641  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
642  tgf.isTmp(),
643  false
644  ),
646  fieldPrevIterPtr_(nullptr),
647  boundaryField_
648  (
649  this->mesh().boundary(),
650  *this,
651  patchFieldTypes,
652  actualPatchTypes
653  ),
654  sources_(*this, fieldSourceTypes)
655 {
656  if (debug)
657  {
659  << "Constructing from tmp resetting IO params and patch types"
660  << endl << this->info() << endl;
661  }
662 
663  boundaryField_ == tgf().boundaryField_;
664 
665  tgf.clear();
666 
667  readIfPresent();
668 }
669 
670 
671 template<class Type, template<class> class PatchField, class GeoMesh>
674 {
676  (
678  );
679 }
680 
681 
682 template<class Type, template<class> class PatchField, class GeoMesh>
685 {
687  (
689  (
690  IOobject
691  (
692  this->name(),
693  this->mesh().thisDb().time().name(),
694  this->mesh().thisDb(),
697  false
698  ),
699  *this,
701  )
702  );
703 }
704 
705 
706 template<class Type, template<class> class PatchField, class GeoMesh>
709 (
710  const word& name,
711  const Internal& diField,
712  const PtrList<PatchField<Type>>& ptfl,
713  const HashPtrTable<Source>& stft
714 )
715 {
716  const bool cacheTmp = diField.mesh().thisDb().cacheTemporaryObject(name);
717 
719  (
721  (
722  IOobject
723  (
724  name,
725  diField.mesh().thisDb().time().name(),
726  diField.mesh().thisDb(),
729  cacheTmp
730  ),
731  diField,
732  ptfl,
733  stft
734  ),
735  cacheTmp
736  );
737 }
738 
739 
740 template<class Type, template<class> class PatchField, class GeoMesh>
743 (
744  const word& name,
745  const Mesh& mesh,
746  const dimensionSet& ds,
747  const word& patchFieldType
748 )
749 {
750  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
751 
753  (
755  (
756  IOobject
757  (
758  name,
759  mesh.thisDb().time().name(),
760  mesh.thisDb(),
763  cacheTmp
764  ),
765  mesh,
766  ds,
768  ),
769  cacheTmp
770  );
771 }
772 
773 
774 template<class Type, template<class> class PatchField, class GeoMesh>
777 (
778  const word& name,
779  const Mesh& mesh,
780  const dimensioned<Type>& dt,
781  const word& patchFieldType
782 )
783 {
784  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
785 
787  (
789  (
790  IOobject
791  (
792  name,
793  mesh.thisDb().time().name(),
794  mesh.thisDb(),
797  cacheTmp
798  ),
799  mesh,
800  dt,
802  ),
803  cacheTmp
804  );
805 }
806 
807 
808 
809 template<class Type, template<class> class PatchField, class GeoMesh>
812 (
813  const word& name,
814  const Mesh& mesh,
815  const dimensioned<Type>& dt,
816  const wordList& patchFieldTypes,
817  const wordList& actualPatchTypes,
818  const HashTable<word>& fieldSourceTypes
819 )
820 {
821  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
822 
824  (
826  (
827  IOobject
828  (
829  name,
830  mesh.thisDb().time().name(),
831  mesh.thisDb(),
834  cacheTmp
835  ),
836  mesh,
837  dt,
838  patchFieldTypes,
839  actualPatchTypes,
840  fieldSourceTypes
841  ),
842  cacheTmp
843  );
844 }
845 
846 
847 template<class Type, template<class> class PatchField, class GeoMesh>
850 (
851  const word& newName,
853 )
854 {
855  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
856 
858  (
860  (
861  IOobject
862  (
863  newName,
864  tgf().instance(),
865  tgf().local(),
866  tgf().db(),
869  cacheTmp
870  ),
871  tgf
872  ),
873  cacheTmp
874  );
875 }
876 
877 
878 template<class Type, template<class> class PatchField, class GeoMesh>
881 (
882  const word& newName,
884  const word& patchFieldType
885 )
886 {
887  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
888 
890  (
892  (
893  IOobject
894  (
895  newName,
896  tgf().instance(),
897  tgf().local(),
898  tgf().db(),
901  cacheTmp
902  ),
903  tgf,
905  ),
906  cacheTmp
907  );
908 }
909 
910 
911 template<class Type, template<class> class PatchField, class GeoMesh>
914 (
915  const word& newName,
917  const wordList& patchFieldTypes,
918  const wordList& actualPatchTypes,
919  const HashTable<word>& fieldSourceTypes
920 )
921 {
922  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
923 
925  (
927  (
928  IOobject
929  (
930  newName,
931  tgf().instance(),
932  tgf().local(),
933  tgf().db(),
936  cacheTmp
937  ),
938  tgf,
939  patchFieldTypes,
940  actualPatchTypes,
941  fieldSourceTypes
942  ),
943  cacheTmp
944  );
945 }
946 
947 
948 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
949 
950 template<class Type, template<class> class PatchField, class GeoMesh>
952 {
953  this->db().cacheTemporaryObject(*this);
954 
955  clearPrevIter();
956 }
957 
958 
959 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
960 
961 template<class Type, template<class> class PatchField, class GeoMesh>
962 typename
965 {
966  this->setUpToDate();
967  storeOldTimes();
968  return *this;
969 }
970 
971 
972 template<class Type, template<class> class PatchField, class GeoMesh>
973 typename
976 {
977  this->setUpToDate();
978  storeOldTimes();
979  return *this;
980 }
981 
982 
983 template<class Type, template<class> class PatchField, class GeoMesh>
984 typename
987 {
988  this->setUpToDate();
989  storeOldTimes();
990  return boundaryField_;
991 }
992 
993 
994 template<class Type, template<class> class PatchField, class GeoMesh>
996 {
997  if (!fieldPrevIterPtr_)
998  {
999  if (debug)
1000  {
1002  << "Allocating previous iteration field" << endl
1003  << this->info() << endl;
1004  }
1005 
1006  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1007  (
1008  this->name() + "PrevIter",
1009  *this
1010  );
1011  }
1012  else
1013  {
1014  *fieldPrevIterPtr_ == *this;
1015  }
1016 }
1017 
1018 
1019 template<class Type, template<class> class PatchField, class GeoMesh>
1021 {
1022  deleteDemandDrivenData(fieldPrevIterPtr_);
1023 }
1024 
1025 
1026 template<class Type, template<class> class PatchField, class GeoMesh>
1029 {
1030  if (!fieldPrevIterPtr_)
1031  {
1033  << "previous iteration field" << endl << this->info() << endl
1034  << " not stored."
1035  << " Use field.storePrevIter() at start of iteration."
1036  << abort(FatalError);
1037  }
1038 
1039  return *fieldPrevIterPtr_;
1040 }
1041 
1042 
1043 template<class Type, template<class> class PatchField, class GeoMesh>
1046 {
1047  this->setUpToDate();
1048  storeOldTimes();
1049  boundaryField_.evaluate();
1050 }
1051 
1052 
1053 template<class Type, template<class> class PatchField, class GeoMesh>
1055 (
1057 )
1058 {
1060 
1061  Internal::reset(gf);
1062  boundaryField_.reset(gf.boundaryField());
1063  sources_.reset(gf.sources());
1064 
1065  tgf.clear();
1066 }
1067 
1068 
1069 template<class Type, template<class> class PatchField, class GeoMesh>
1071 {
1072  // Search all boundary conditions, if any are
1073  // fixed-value or mixed (Robin) do not set reference level for solution.
1074 
1075  bool needRef = true;
1076 
1077  forAll(boundaryField_, patchi)
1078  {
1079  if (boundaryField_[patchi].fixesValue())
1080  {
1081  needRef = false;
1082  break;
1083  }
1084  }
1085 
1086  reduce(needRef, andOp<bool>());
1087 
1088  return needRef;
1089 }
1090 
1091 
1092 template<class Type, template<class> class PatchField, class GeoMesh>
1094 {
1095  if (alpha < 1)
1096  {
1097  if (debug)
1098  {
1100  << "Relaxing" << endl << this->info()
1101  << " by " << alpha << endl;
1102  }
1103 
1104  operator==(prevIter() + alpha*(*this - prevIter()));
1105  }
1106 }
1107 
1108 
1109 template<class Type, template<class> class PatchField, class GeoMesh>
1110 Foam::scalar
1112 {
1113  if
1114  (
1115  solutionControl::finalIteration(this->mesh())
1116  && this->mesh().solution().relaxField(this->name() + "Final")
1117  )
1118  {
1119  return this->mesh().solution().fieldRelaxationFactor
1120  (
1121  this->name() + "Final"
1122  );
1123  }
1124  else if (this->mesh().solution().relaxField(this->name()))
1125  {
1126  return this->mesh().solution().fieldRelaxationFactor(this->name());
1127  }
1128  else
1129  {
1130  return 1;
1131  }
1132 }
1133 
1134 
1135 template<class Type, template<class> class PatchField, class GeoMesh>
1137 {
1138  relax(relaxationFactor());
1139 }
1140 
1141 
1142 template<class Type, template<class> class PatchField, class GeoMesh>
1144 (
1146  const scalar alpha
1147 )
1148 {
1149  if (alpha < 1)
1150  {
1151  if (debug)
1152  {
1154  << "Relaxing" << endl << this->info()
1155  << " by " << alpha << endl;
1156  }
1157 
1158  operator==(*this + alpha*(tgf - *this));
1159  }
1160  else
1161  {
1162  operator==(tgf);
1163  }
1164 }
1165 
1166 
1167 template<class Type, template<class> class PatchField, class GeoMesh>
1169 (
1171 )
1172 {
1173  relax(tgf, relaxationFactor());
1174 }
1175 
1176 
1177 template<class Type, template<class> class PatchField, class GeoMesh>
1179 (
1180  bool final
1181 ) const
1182 {
1183  if (final)
1184  {
1185  return this->name() + "Final";
1186  }
1187  else
1188  {
1189  return this->name();
1190  }
1191 }
1192 
1193 
1194 template<class Type, template<class> class PatchField, class GeoMesh>
1196 (
1197  Ostream& os
1198 ) const
1199 {
1200  os << "min/max(" << this->name() << ") = "
1201  << Foam::min(*this).value() << ", "
1202  << Foam::max(*this).value()
1203  << endl;
1204 }
1205 
1206 
1207 template<class Type, template<class> class PatchField, class GeoMesh>
1209 writeData(Ostream& os) const
1210 {
1211  os << *this;
1212  return os.good();
1213 }
1214 
1215 
1216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1217 
1218 template<class Type, template<class> class PatchField, class GeoMesh>
1221 {
1223  (
1225  (
1226  this->name() + ".T()",
1227  this->mesh(),
1228  this->dimensions()
1229  )
1230  );
1231 
1232  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1233  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1234 
1235  return result;
1236 }
1237 
1238 
1239 template<class Type, template<class> class PatchField, class GeoMesh>
1240 Foam::tmp
1241 <
1243  <
1245  PatchField,
1246  GeoMesh
1247  >
1248 >
1250 (
1251  const direction d
1252 ) const
1253 {
1255  (
1257  (
1258  this->name() + ".component(" + Foam::name(d) + ')',
1259  this->mesh(),
1260  this->dimensions()
1261  )
1262  );
1263 
1264  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1265  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1266 
1267  return Component;
1268 }
1269 
1270 
1271 template<class Type, template<class> class PatchField, class GeoMesh>
1273 (
1274  const direction d,
1275  const GeometricField
1276  <
1278  PatchField,
1279  GeoMesh
1280  >& gcf
1281 )
1282 {
1283  primitiveFieldRef().replace(d, gcf.primitiveField());
1284  boundaryFieldRef().replace(d, gcf.boundaryField());
1285 }
1286 
1287 
1288 template<class Type, template<class> class PatchField, class GeoMesh>
1290 (
1291  const direction d,
1292  const dimensioned<cmptType>& ds
1293 )
1294 {
1295  primitiveFieldRef().replace(d, ds.value());
1296  boundaryFieldRef().replace(d, ds.value());
1297 }
1298 
1299 
1300 template<class Type, template<class> class PatchField, class GeoMesh>
1302 (
1303  const dimensioned<Type>& dt
1304 )
1305 {
1306  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1307  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1308 }
1309 
1310 
1311 template<class Type, template<class> class PatchField, class GeoMesh>
1313 (
1314  const dimensioned<Type>& dt
1315 )
1316 {
1317  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1318  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1319 }
1320 
1321 
1322 template<class Type, template<class> class PatchField, class GeoMesh>
1324 (
1325  const dimensioned<Type>& minDt,
1326  const dimensioned<Type>& maxDt
1327 )
1328 {
1329  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1330  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1331  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1332  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1333 }
1334 
1335 
1336 template<class Type, template<class> class PatchField, class GeoMesh>
1338 {
1339  primitiveFieldRef().negate();
1340  boundaryFieldRef().negate();
1341 }
1342 
1343 
1344 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1345 
1346 template<class Type, template<class> class PatchField, class GeoMesh>
1348 (
1350 )
1351 {
1352  if (this == &gf)
1353  {
1355  << "attempted assignment to self"
1356  << abort(FatalError);
1357  }
1358 
1359  checkField(*this, gf, "=");
1360 
1361  // Only assign field contents not ID
1362 
1363  internalFieldRef() = gf.internalField();
1364  boundaryFieldRef() = gf.boundaryField();
1365 }
1366 
1367 
1368 template<class Type, template<class> class PatchField, class GeoMesh>
1370 (
1372 )
1373 {
1374  if (this == &gf)
1375  {
1377  << "attempted assignment to self"
1378  << abort(FatalError);
1379  }
1380 
1381  checkField(*this, gf, "=");
1382 
1383  // Only assign field contents not ID
1384 
1385  internalFieldRef() = move(gf.internalField());
1386  boundaryFieldRef() = move(gf.boundaryField());
1387 }
1388 
1389 
1390 template<class Type, template<class> class PatchField, class GeoMesh>
1392 (
1394 )
1395 {
1396  if (this == &(tgf()))
1397  {
1399  << "attempted assignment to self"
1400  << abort(FatalError);
1401  }
1402 
1404 
1405  checkField(*this, gf, "=");
1406 
1407  // Only assign field contents not ID
1408 
1409  this->dimensions() = gf.dimensions();
1410 
1411  if (tgf.isTmp())
1412  {
1413  primitiveFieldRef().transfer(tgf.ref());
1414  }
1415  else
1416  {
1418  }
1419 
1420  boundaryFieldRef() = gf.boundaryField();
1421 
1422  tgf.clear();
1423 }
1424 
1425 
1426 template<class Type, template<class> class PatchField, class GeoMesh>
1428 (
1429  const dimensioned<Type>& dt
1430 )
1431 {
1432  internalFieldRef() = dt;
1433  boundaryFieldRef() = dt.value();
1434 }
1435 
1436 
1437 template<class Type, template<class> class PatchField, class GeoMesh>
1439 (
1440  const zero&
1441 )
1442 {
1443  internalFieldRef() = Zero;
1444  boundaryFieldRef() = Zero;
1445 }
1446 
1447 
1448 template<class Type, template<class> class PatchField, class GeoMesh>
1450 (
1452 )
1453 {
1455 
1456  checkField(*this, gf, "==");
1457 
1458  // Only assign field contents not ID
1459 
1460  this->dimensions() = gf.dimensions();
1461 
1462  if (tgf.isTmp())
1463  {
1464  primitiveFieldRef().transfer(tgf.ref());
1465  }
1466  else
1467  {
1469  }
1470 
1471  boundaryFieldRef() == gf.boundaryField();
1472 
1473  tgf.clear();
1474 }
1475 
1476 
1477 template<class Type, template<class> class PatchField, class GeoMesh>
1479 (
1480  const dimensioned<Type>& dt
1481 )
1482 {
1483  internalFieldRef() = dt;
1484  boundaryFieldRef() == dt.value();
1485 }
1486 
1487 
1488 template<class Type, template<class> class PatchField, class GeoMesh>
1490 (
1491  const zero&
1492 )
1493 {
1494  internalFieldRef() = Zero;
1495  boundaryFieldRef() == Zero;
1496 }
1497 
1498 
1499 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1500  \
1501 template<class Type, template<class> class PatchField, class GeoMesh> \
1502 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1503 ( \
1504  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1505 ) \
1506 { \
1507  checkField(*this, gf, #op); \
1508  \
1509  internalFieldRef() op gf.internalField(); \
1510  boundaryFieldRef() op gf.boundaryField(); \
1511 } \
1512  \
1513 template<class Type, template<class> class PatchField, class GeoMesh> \
1514 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1515 ( \
1516  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1517 ) \
1518 { \
1519  operator op(tgf()); \
1520  tgf.clear(); \
1521 } \
1522  \
1523 template<class Type, template<class> class PatchField, class GeoMesh> \
1524 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1525 ( \
1526  const dimensioned<TYPE>& dt \
1527 ) \
1528 { \
1529  internalFieldRef() op dt; \
1530  boundaryFieldRef() op dt.value(); \
1531 }
1532 
1537 
1538 #undef COMPUTED_ASSIGNMENT
1539 
1540 
1541 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1542 
1543 template<class Type, template<class> class PatchField, class GeoMesh>
1544 Foam::Ostream& Foam::operator<<
1545 (
1546  Ostream& os,
1548 )
1549 {
1550  gf().writeData(os, "internalField");
1551  os << nl;
1552  gf.boundaryField().writeEntry("boundaryField", os);
1553 
1554  if (!gf.sources_.empty())
1555  {
1556  os << nl;
1557  gf.sources().writeEntry("sources", os);
1558  }
1559 
1560  // Check state of IOstream
1561  os.check
1562  (
1563  "Ostream& operator<<(Ostream&, "
1564  "const GeometricField<Type, PatchField, GeoMesh>&)"
1565  );
1566 
1567  return (os);
1568 }
1569 
1570 
1571 template<class Type, template<class> class PatchField, class GeoMesh>
1572 Foam::Ostream& Foam::operator<<
1573 (
1574  Ostream& os,
1576 )
1577 {
1578  os << tgf();
1579  tgf.clear();
1580  return os;
1581 }
1582 
1583 
1584 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1585 
1586 #undef checkField
1587 
1588 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1589 
1590 #include "GeometricFieldFunctions.C"
1591 
1592 // ************************************************************************* //
EaEqn relax()
#define checkField(gf1, gf2, op)
#define COMPUTED_ASSIGNMENT(TYPE, op)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricBoundaryField class.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
void writeEntry(const word &keyword, Ostream &os) const
Write sources as dictionary entry.
Generic GeometricField class.
void max(const dimensioned< Type > &)
Internal & internalFieldRef()
Return a reference to the dimensioned internal field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
const Sources & sources() const
Return const-reference to the sources.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
void relax()
Relax current field with respect to the cached previous iteration.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
void min(const dimensioned< Type > &)
Field< Type >::cmptType cmptType
Component type of the elements of the field.
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch field type.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
bool needReference() const
Does the field need a reference level for solution.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
virtual ~GeometricField()
Destructor.
scalar relaxationFactor() const
Return the field relaxation factor read from fvSolution.
void clearPrevIter()
Delete the previous iteration field.
void storePrevIter() const
Store the field as the previous iteration value.
word select(bool final) const
Select the final iteration parameters if `final' is true.
void correctBoundaryConditions()
Correct boundary field.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:486
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
Class to add into field types to provide old-time storage and retrieval.
Definition: OldTimeField.H:93
void copyOldTimes(const IOobject &io, const OldTimeField< GeometricField< Type, PatchField, GeoMesh > > &otf)
Copy the old times from the given field.
Definition: OldTimeField.C:161
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:198
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const Type & value() const
Return const reference to value.
static bool finalIteration(const objectRegistry &registry)
Lookup solutionControl from the objectRegistry and return finalIter.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
Template functions to aid in the implementation of demand driven data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
#define InfoInFunction
Report an information message using Foam::Info.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word patchFieldType(const PatchField &pf)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
uint8_t direction
Definition: direction.H:45
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
dictionary dict
conserve primitiveFieldRef()+