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-2019 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 "data.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  if (dict.found("referenceLevel"))
59  {
60  Type fieldAverage(pTraits<Type>(dict.lookup("referenceLevel")));
61 
62  Field<Type>::operator+=(fieldAverage);
63 
64  forAll(boundaryField_, patchi)
65  {
66  boundaryField_[patchi] == boundaryField_[patchi] + fieldAverage;
67  }
68  }
69 }
70 
71 
72 template<class Type, template<class> class PatchField, class GeoMesh>
74 {
75  const localIOdictionary dict
76  (
77  IOobject
78  (
79  this->name(),
80  this->instance(),
81  this->local(),
82  this->db(),
83  IOobject::MUST_READ,
84  IOobject::NO_WRITE,
85  false
86  ),
87  typeName
88  );
89 
90  this->close();
91 
92  readFields(dict);
93 }
94 
95 
96 template<class Type, template<class> class PatchField, class GeoMesh>
98 {
99  if
100  (
101  this->readOpt() == IOobject::MUST_READ
102  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
103  )
104  {
106  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
107  << " suggests that a read constructor for field " << this->name()
108  << " would be more appropriate." << endl;
109  }
110  else if
111  (
112  this->readOpt() == IOobject::READ_IF_PRESENT
113  && this->template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
114  (
115  true
116  )
117  )
118  {
119  readFields();
120 
121  // Check compatibility between field and mesh
122  if (this->size() != GeoMesh::size(this->mesh()))
123  {
124  FatalIOErrorInFunction(this->readStream(typeName))
125  << " number of field elements = " << this->size()
126  << " number of mesh elements = "
127  << GeoMesh::size(this->mesh())
128  << exit(FatalIOError);
129  }
130 
131  readOldTimeIfPresent();
132 
133  return true;
134  }
135 
136  return false;
137 }
138 
139 
140 template<class Type, template<class> class PatchField, class GeoMesh>
142 {
143  // Read the old time field if present
144  IOobject field0
145  (
146  this->name() + "_0",
147  this->time().timeName(),
148  this->db(),
149  IOobject::READ_IF_PRESENT,
150  IOobject::AUTO_WRITE,
151  this->registerObject()
152  );
153 
154  if
155  (
156  field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
157  (
158  true
159  )
160  )
161  {
162  if (debug)
163  {
164  InfoInFunction << "Reading old time level for field"
165  << endl << this->info() << endl;
166  }
167 
168  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
169  (
170  field0,
171  this->mesh()
172  );
173 
174  field0Ptr_->timeIndex_ = timeIndex_ - 1;
175 
176  if (!field0Ptr_->readOldTimeIfPresent())
177  {
178  field0Ptr_->oldTime();
179  }
180 
181  return true;
182  }
183 
184  return false;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
189 
190 template<class Type, template<class> class PatchField, class GeoMesh>
192 (
193  const IOobject& io,
194  const Mesh& mesh,
195  const dimensionSet& ds,
196  const word& patchFieldType
197 )
198 :
199  Internal(io, mesh, ds, false),
200  timeIndex_(this->time().timeIndex()),
201  field0Ptr_(nullptr),
202  fieldPrevIterPtr_(nullptr),
203  boundaryField_(mesh.boundary(), *this, patchFieldType)
204 {
205  if (debug)
206  {
207  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
208  }
209 
210  readIfPresent();
211 }
212 
213 
214 template<class Type, template<class> class PatchField, class GeoMesh>
216 (
217  const IOobject& io,
218  const Mesh& mesh,
219  const dimensionSet& ds,
220  const wordList& patchFieldTypes,
221  const wordList& actualPatchTypes
222 )
223 :
224  Internal(io, mesh, ds, false),
225  timeIndex_(this->time().timeIndex()),
226  field0Ptr_(nullptr),
227  fieldPrevIterPtr_(nullptr),
228  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
229 {
230  if (debug)
231  {
232  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
233  }
234 
235  readIfPresent();
236 }
237 
238 
239 template<class Type, template<class> class PatchField, class GeoMesh>
241 (
242  const IOobject& io,
243  const Mesh& mesh,
244  const dimensioned<Type>& dt,
245  const word& patchFieldType
246 )
247 :
248  Internal(io, mesh, dt, false),
249  timeIndex_(this->time().timeIndex()),
250  field0Ptr_(nullptr),
251  fieldPrevIterPtr_(nullptr),
252  boundaryField_(mesh.boundary(), *this, patchFieldType)
253 {
254  if (debug)
255  {
256  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
257  }
258 
259  boundaryField_ == dt.value();
260 
261  readIfPresent();
262 }
263 
264 
265 template<class Type, template<class> class PatchField, class GeoMesh>
267 (
268  const IOobject& io,
269  const Mesh& mesh,
270  const dimensioned<Type>& dt,
271  const wordList& patchFieldTypes,
272  const wordList& actualPatchTypes
273 )
274 :
275  Internal(io, mesh, dt, false),
276  timeIndex_(this->time().timeIndex()),
277  field0Ptr_(nullptr),
278  fieldPrevIterPtr_(nullptr),
279  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
280 {
281  if (debug)
282  {
283  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
284  }
285 
286  boundaryField_ == dt.value();
287 
288  readIfPresent();
289 }
290 
291 
292 template<class Type, template<class> class PatchField, class GeoMesh>
294 (
295  const IOobject& io,
296  const Internal& diField,
297  const PtrList<PatchField<Type>>& ptfl
298 )
299 :
300  Internal(io, diField),
301  timeIndex_(this->time().timeIndex()),
302  field0Ptr_(nullptr),
303  fieldPrevIterPtr_(nullptr),
304  boundaryField_(this->mesh().boundary(), *this, ptfl)
305 {
306  if (debug)
307  {
309  << "Constructing from components" << endl << this->info() << endl;
310  }
311 
312  readIfPresent();
313 }
314 
315 
316 template<class Type, template<class> class PatchField, class GeoMesh>
318 (
319  const IOobject& io,
320  const Mesh& mesh,
321  const dimensionSet& ds,
322  const Field<Type>& iField,
323  const PtrList<PatchField<Type>>& ptfl
324 )
325 :
326  Internal(io, mesh, ds, iField),
327  timeIndex_(this->time().timeIndex()),
328  field0Ptr_(nullptr),
329  fieldPrevIterPtr_(nullptr),
330  boundaryField_(mesh.boundary(), *this, ptfl)
331 {
332  if (debug)
333  {
335  << "Constructing from components" << endl << this->info() << endl;
336  }
337 
338  readIfPresent();
339 }
340 
341 
342 template<class Type, template<class> class PatchField, class GeoMesh>
344 (
345  const IOobject& io,
346  const Mesh& mesh
347 )
348 :
349  Internal(io, mesh, dimless, false),
350  timeIndex_(this->time().timeIndex()),
351  field0Ptr_(nullptr),
352  fieldPrevIterPtr_(nullptr),
353  boundaryField_(mesh.boundary())
354 {
355  readFields();
356 
357  // Check compatibility between field and mesh
358 
359  if (this->size() != GeoMesh::size(this->mesh()))
360  {
361  FatalIOErrorInFunction(this->readStream(typeName))
362  << " number of field elements = " << this->size()
363  << " number of mesh elements = " << GeoMesh::size(this->mesh())
364  << exit(FatalIOError);
365  }
366 
367  readOldTimeIfPresent();
368 
369  if (debug)
370  {
372  << "Finishing read-construction of" << endl << this->info() << endl;
373  }
374 }
375 
376 
377 template<class Type, template<class> class PatchField, class GeoMesh>
379 (
380  const IOobject& io,
381  const Mesh& mesh,
382  const dictionary& dict
383 )
384 :
385  Internal(io, mesh, dimless, false),
386  timeIndex_(this->time().timeIndex()),
387  field0Ptr_(nullptr),
388  fieldPrevIterPtr_(nullptr),
389  boundaryField_(mesh.boundary())
390 {
391  readFields(dict);
392 
393  // Check compatibility between field and mesh
394 
395  if (this->size() != GeoMesh::size(this->mesh()))
396  {
398  << " number of field elements = " << this->size()
399  << " number of mesh elements = " << GeoMesh::size(this->mesh())
400  << exit(FatalIOError);
401  }
402 
403  if (debug)
404  {
406  << "Finishing dictionary-construct of "
407  << endl << this->info() << endl;
408  }
409 }
410 
411 
412 template<class Type, template<class> class PatchField, class GeoMesh>
414 (
416 )
417 :
418  Internal(gf),
419  timeIndex_(gf.timeIndex()),
420  field0Ptr_(nullptr),
421  fieldPrevIterPtr_(nullptr),
422  boundaryField_(*this, gf.boundaryField_)
423 {
424  if (debug)
425  {
427  << "Constructing as copy" << endl << this->info() << endl;
428  }
429 
430  if (gf.field0Ptr_)
431  {
433  (
434  *gf.field0Ptr_
435  );
436  }
437 
438  this->writeOpt() = IOobject::NO_WRITE;
439 }
440 
441 
442 template<class Type, template<class> class PatchField, class GeoMesh>
444 (
446 )
447 :
448  Internal(move(gf)),
449  timeIndex_(gf.timeIndex()),
450  field0Ptr_(nullptr),
451  fieldPrevIterPtr_(nullptr),
452  boundaryField_(move(gf.boundaryField_))
453 {
454  if (debug)
455  {
457  << "Constructing by moving" << endl << this->info() << endl;
458  }
459 
460  if (gf.field0Ptr_)
461  {
462  field0Ptr_ = gf.field0Ptr_;
463  gf.field0Ptr_ = nullptr;
464  }
465 
466  this->writeOpt() = IOobject::NO_WRITE;
467 }
468 
469 
470 template<class Type, template<class> class PatchField, class GeoMesh>
472 (
474 )
475 :
476  Internal
477  (
478  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
479  tgf.isTmp()
480  ),
481  timeIndex_(tgf().timeIndex()),
482  field0Ptr_(nullptr),
483  fieldPrevIterPtr_(nullptr),
484  boundaryField_(*this, tgf().boundaryField_)
485 {
486  if (debug)
487  {
489  << "Constructing from tmp" << endl << this->info() << endl;
490  }
491 
492  this->writeOpt() = IOobject::NO_WRITE;
493 
494  tgf.clear();
495 }
496 
497 
498 template<class Type, template<class> class PatchField, class GeoMesh>
500 (
501  const IOobject& io,
503 )
504 :
505  Internal(io, gf),
506  timeIndex_(gf.timeIndex()),
507  field0Ptr_(nullptr),
508  fieldPrevIterPtr_(nullptr),
509  boundaryField_(*this, gf.boundaryField_)
510 {
511  if (debug)
512  {
514  << "Constructing as copy resetting IO params"
515  << endl << this->info() << endl;
516  }
517 
518  if (!readIfPresent() && gf.field0Ptr_)
519  {
521  (
522  io.name() + "_0",
523  *gf.field0Ptr_
524  );
525  }
526 }
527 
528 
529 template<class Type, template<class> class PatchField, class GeoMesh>
531 (
532  const IOobject& io,
534 )
535 :
536  Internal
537  (
538  io,
539  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
540  tgf.isTmp()
541  ),
542  timeIndex_(tgf().timeIndex()),
543  field0Ptr_(nullptr),
544  fieldPrevIterPtr_(nullptr),
545  boundaryField_(*this, tgf().boundaryField_)
546 {
547  if (debug)
548  {
550  << "Constructing from tmp resetting IO params"
551  << endl << this->info() << endl;
552  }
553 
554  tgf.clear();
555 
556  readIfPresent();
557 }
558 
559 
560 template<class Type, template<class> class PatchField, class GeoMesh>
562 (
563  const word& newName,
565 )
566 :
567  Internal(newName, gf),
568  timeIndex_(gf.timeIndex()),
569  field0Ptr_(nullptr),
570  fieldPrevIterPtr_(nullptr),
571  boundaryField_(*this, gf.boundaryField_)
572 {
573  if (debug)
574  {
576  << "Constructing as copy resetting name"
577  << endl << this->info() << endl;
578  }
579 
580  if (!readIfPresent() && gf.field0Ptr_)
581  {
583  (
584  newName + "_0",
585  *gf.field0Ptr_
586  );
587  }
588 }
589 
590 
591 template<class Type, template<class> class PatchField, class GeoMesh>
593 (
594  const word& newName,
596 )
597 :
598  Internal
599  (
600  newName,
601  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
602  tgf.isTmp()
603  ),
604  timeIndex_(tgf().timeIndex()),
605  field0Ptr_(nullptr),
606  fieldPrevIterPtr_(nullptr),
607  boundaryField_(*this, tgf().boundaryField_)
608 {
609  if (debug)
610  {
612  << "Constructing from tmp resetting name"
613  << endl << this->info() << endl;
614  }
615 
616  tgf.clear();
617 }
618 
619 
620 template<class Type, template<class> class PatchField, class GeoMesh>
622 (
623  const IOobject& io,
625  const word& patchFieldType
626 )
627 :
628  Internal(io, gf),
629  timeIndex_(gf.timeIndex()),
630  field0Ptr_(nullptr),
631  fieldPrevIterPtr_(nullptr),
632  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
633 {
634  if (debug)
635  {
637  << "Constructing as copy resetting IO params"
638  << endl << this->info() << endl;
639  }
640 
641  boundaryField_ == gf.boundaryField_;
642 
643  if (!readIfPresent() && gf.field0Ptr_)
644  {
646  (
647  io.name() + "_0",
648  *gf.field0Ptr_
649  );
650  }
651 }
652 
653 
654 template<class Type, template<class> class PatchField, class GeoMesh>
656 (
657  const IOobject& io,
659  const wordList& patchFieldTypes,
660  const wordList& actualPatchTypes
661 
662 )
663 :
664  Internal(io, gf),
665  timeIndex_(gf.timeIndex()),
666  field0Ptr_(nullptr),
667  fieldPrevIterPtr_(nullptr),
668  boundaryField_
669  (
670  this->mesh().boundary(),
671  *this,
672  patchFieldTypes,
673  actualPatchTypes
674  )
675 {
676  if (debug)
677  {
679  << "Constructing as copy resetting IO params and patch types"
680  << endl << this->info() << endl;
681  }
682 
683  boundaryField_ == gf.boundaryField_;
684 
685  if (!readIfPresent() && gf.field0Ptr_)
686  {
688  (
689  io.name() + "_0",
690  *gf.field0Ptr_
691  );
692  }
693 }
694 
695 
696 template<class Type, template<class> class PatchField, class GeoMesh>
698 (
699  const IOobject& io,
701  const wordList& patchFieldTypes,
702  const wordList& actualPatchTypes
703 )
704 :
705  Internal
706  (
707  io,
708  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
709  tgf.isTmp()
710  ),
711  timeIndex_(tgf().timeIndex()),
712  field0Ptr_(nullptr),
713  fieldPrevIterPtr_(nullptr),
714  boundaryField_
715  (
716  this->mesh().boundary(),
717  *this,
718  patchFieldTypes,
719  actualPatchTypes
720  )
721 {
722  if (debug)
723  {
725  << "Constructing from tmp resetting IO params and patch types"
726  << endl << this->info() << endl;
727  }
728 
729  boundaryField_ == tgf().boundaryField_;
730 
731  tgf.clear();
732 }
733 
734 
735 template<class Type, template<class> class PatchField, class GeoMesh>
738 {
740  (
742  );
743 }
744 
745 
746 template<class Type, template<class> class PatchField, class GeoMesh>
749 (
750  const word& name,
751  const Mesh& mesh,
752  const dimensionSet& ds,
753  const word& patchFieldType
754 )
755 {
757  (
759  (
760  IOobject
761  (
762  name,
763  mesh.thisDb().time().timeName(),
764  mesh.thisDb(),
765  IOobject::NO_READ,
766  IOobject::NO_WRITE,
767  false
768  ),
769  mesh,
770  ds,
771  patchFieldType
772  )
773  );
774 }
775 
776 
777 template<class Type, template<class> class PatchField, class GeoMesh>
780 (
781  const word& name,
782  const Mesh& mesh,
783  const dimensioned<Type>& dt,
784  const word& patchFieldType
785 )
786 {
788  (
790  (
791  IOobject
792  (
793  name,
794  mesh.thisDb().time().timeName(),
795  mesh.thisDb(),
796  IOobject::NO_READ,
797  IOobject::NO_WRITE,
798  false
799  ),
800  mesh,
801  dt,
802  patchFieldType
803  )
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 )
819 {
821  (
823  (
824  IOobject
825  (
826  name,
827  mesh.thisDb().time().timeName(),
828  mesh.thisDb(),
829  IOobject::NO_READ,
830  IOobject::NO_WRITE,
831  false
832  ),
833  mesh,
834  dt,
835  patchFieldTypes,
836  actualPatchTypes
837  )
838  );
839 }
840 
841 
842 template<class Type, template<class> class PatchField, class GeoMesh>
845 (
846  const word& newName,
848 )
849 {
851  (
853  (
854  IOobject
855  (
856  newName,
857  tgf().instance(),
858  tgf().local(),
859  tgf().db(),
860  IOobject::NO_READ,
861  IOobject::NO_WRITE,
862  false
863  ),
864  tgf
865  )
866  );
867 }
868 
869 
870 template<class Type, template<class> class PatchField, class GeoMesh>
873 (
874  const word& newName,
876  const wordList& patchFieldTypes,
877  const wordList& actualPatchTypes
878 )
879 {
881  (
883  (
884  IOobject
885  (
886  newName,
887  tgf().instance(),
888  tgf().local(),
889  tgf().db(),
890  IOobject::NO_READ,
891  IOobject::NO_WRITE,
892  false
893  ),
894  tgf,
895  patchFieldTypes,
896  actualPatchTypes
897  )
898  );
899 }
900 
901 
902 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
903 
904 template<class Type, template<class> class PatchField, class GeoMesh>
906 {
907  deleteDemandDrivenData(field0Ptr_);
908  deleteDemandDrivenData(fieldPrevIterPtr_);
909 }
910 
911 
912 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
913 
914 template<class Type, template<class> class PatchField, class GeoMesh>
915 typename
918 {
919  this->setUpToDate();
920  storeOldTimes();
921  return *this;
922 }
923 
924 
925 template<class Type, template<class> class PatchField, class GeoMesh>
926 typename
929 {
930  this->setUpToDate();
931  storeOldTimes();
932  return *this;
933 }
934 
935 
936 template<class Type, template<class> class PatchField, class GeoMesh>
937 typename
940 {
941  this->setUpToDate();
942  storeOldTimes();
943  return boundaryField_;
944 }
945 
946 
947 template<class Type, template<class> class PatchField, class GeoMesh>
949 {
950  if
951  (
952  field0Ptr_
953  && timeIndex_ != this->time().timeIndex()
954  && !(
955  this->name().size() > 2
956  && this->name()(this->name().size()-2, 2) == "_0"
957  )
958  )
959  {
960  storeOldTime();
961  }
962 
963  // Correct time index
964  timeIndex_ = this->time().timeIndex();
965 }
966 
967 
968 template<class Type, template<class> class PatchField, class GeoMesh>
970 {
971  if (field0Ptr_)
972  {
973  field0Ptr_->storeOldTime();
974 
975  if (debug)
976  {
978  << "Storing old time field for field" << endl
979  << this->info() << endl;
980  }
981 
982  *field0Ptr_ == *this;
983  field0Ptr_->timeIndex_ = timeIndex_;
984 
985  if (field0Ptr_->field0Ptr_)
986  {
987  field0Ptr_->writeOpt() = this->writeOpt();
988  }
989  }
990 }
991 
992 
993 template<class Type, template<class> class PatchField, class GeoMesh>
995 {
996  if (field0Ptr_)
997  {
998  return field0Ptr_->nOldTimes() + 1;
999  }
1000  else
1001  {
1002  return 0;
1003  }
1004 }
1005 
1006 
1007 template<class Type, template<class> class PatchField, class GeoMesh>
1010 {
1011  if (!field0Ptr_)
1012  {
1014  (
1015  IOobject
1016  (
1017  this->name() + "_0",
1018  this->time().timeName(),
1019  this->db(),
1020  IOobject::NO_READ,
1021  IOobject::NO_WRITE,
1022  this->registerObject()
1023  ),
1024  *this
1025  );
1026  }
1027  else
1028  {
1029  storeOldTimes();
1030  }
1031 
1032  return *field0Ptr_;
1033 }
1034 
1035 
1036 template<class Type, template<class> class PatchField, class GeoMesh>
1039 {
1040  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1041  .oldTime();
1042 
1043  return *field0Ptr_;
1044 }
1045 
1046 
1047 template<class Type, template<class> class PatchField, class GeoMesh>
1049 {
1050  if (!fieldPrevIterPtr_)
1051  {
1052  if (debug)
1053  {
1055  << "Allocating previous iteration field" << endl
1056  << this->info() << endl;
1057  }
1058 
1059  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1060  (
1061  this->name() + "PrevIter",
1062  *this
1063  );
1064  }
1065  else
1066  {
1067  *fieldPrevIterPtr_ == *this;
1068  }
1069 }
1070 
1071 
1072 template<class Type, template<class> class PatchField, class GeoMesh>
1075 {
1076  if (!fieldPrevIterPtr_)
1077  {
1079  << "previous iteration field" << endl << this->info() << endl
1080  << " not stored."
1081  << " Use field.storePrevIter() at start of iteration."
1082  << abort(FatalError);
1083  }
1084 
1085  return *fieldPrevIterPtr_;
1086 }
1087 
1088 
1089 template<class Type, template<class> class PatchField, class GeoMesh>
1092 {
1093  this->setUpToDate();
1094  storeOldTimes();
1095  boundaryField_.evaluate();
1096 }
1097 
1098 
1099 template<class Type, template<class> class PatchField, class GeoMesh>
1101 {
1102  // Search all boundary conditions, if any are
1103  // fixed-value or mixed (Robin) do not set reference level for solution.
1104 
1105  bool needRef = true;
1106 
1107  forAll(boundaryField_, patchi)
1108  {
1109  if (boundaryField_[patchi].fixesValue())
1110  {
1111  needRef = false;
1112  break;
1113  }
1114  }
1115 
1116  reduce(needRef, andOp<bool>());
1117 
1118  return needRef;
1119 }
1120 
1121 
1122 template<class Type, template<class> class PatchField, class GeoMesh>
1124 {
1125  if (debug)
1126  {
1128  << "Relaxing" << endl << this->info() << " by " << alpha << endl;
1129  }
1130 
1131  operator==(prevIter() + alpha*(*this - prevIter()));
1132 }
1133 
1134 
1135 template<class Type, template<class> class PatchField, class GeoMesh>
1137 {
1138  if
1139  (
1140  this->mesh().data::template lookupOrDefault<bool>
1141  (
1142  "finalIteration",
1143  false
1144  )
1145  && this->mesh().relaxField(this->name() + "Final")
1146  )
1147  {
1148  relax(this->mesh().fieldRelaxationFactor(this->name() + "Final"));
1149  }
1150  else if (this->mesh().relaxField(this->name()))
1151  {
1152  relax(this->mesh().fieldRelaxationFactor(this->name()));
1153  }
1154 }
1155 
1156 
1157 template<class Type, template<class> class PatchField, class GeoMesh>
1160  bool final
1161 ) const
1162 {
1163  if (final)
1164  {
1165  return this->name() + "Final";
1166  }
1167  else
1168  {
1169  return this->name();
1170  }
1171 }
1172 
1173 
1174 template<class Type, template<class> class PatchField, class GeoMesh>
1177  Ostream& os
1178 ) const
1179 {
1180  os << "min/max(" << this->name() << ") = "
1181  << Foam::min(*this).value() << ", "
1182  << Foam::max(*this).value()
1183  << endl;
1184 }
1185 
1186 
1187 template<class Type, template<class> class PatchField, class GeoMesh>
1189 writeData(Ostream& os) const
1190 {
1191  os << *this;
1192  return os.good();
1193 }
1194 
1195 
1196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1197 
1198 template<class Type, template<class> class PatchField, class GeoMesh>
1201 {
1203  (
1205  (
1206  this->name() + ".T()",
1207  this->mesh(),
1208  this->dimensions()
1209  )
1210  );
1211 
1212  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1213  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1214 
1215  return result;
1216 }
1217 
1218 
1219 template<class Type, template<class> class PatchField, class GeoMesh>
1220 Foam::tmp
1221 <
1223  <
1224  typename Foam::GeometricField<Type, PatchField, GeoMesh>::cmptType,
1225  PatchField,
1226  GeoMesh
1227  >
1228 >
1230 (
1231  const direction d
1232 ) const
1233 {
1235  (
1237  (
1238  this->name() + ".component(" + Foam::name(d) + ')',
1239  this->mesh(),
1240  this->dimensions()
1241  )
1242  );
1243 
1244  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1245  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1246 
1247  return Component;
1248 }
1249 
1250 
1251 template<class Type, template<class> class PatchField, class GeoMesh>
1253 (
1254  const direction d,
1255  const GeometricField
1256  <
1258  PatchField,
1259  GeoMesh
1260  >& gcf
1261 )
1262 {
1263  primitiveFieldRef().replace(d, gcf.primitiveField());
1264  boundaryFieldRef().replace(d, gcf.boundaryField());
1265 }
1266 
1267 
1268 template<class Type, template<class> class PatchField, class GeoMesh>
1271  const direction d,
1272  const dimensioned<cmptType>& ds
1273 )
1274 {
1275  primitiveFieldRef().replace(d, ds.value());
1276  boundaryFieldRef().replace(d, ds.value());
1277 }
1278 
1279 
1280 template<class Type, template<class> class PatchField, class GeoMesh>
1283  const dimensioned<Type>& dt
1284 )
1285 {
1286  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1287  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1288 }
1289 
1290 
1291 template<class Type, template<class> class PatchField, class GeoMesh>
1294  const dimensioned<Type>& dt
1295 )
1296 {
1297  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1298  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1299 }
1300 
1301 
1302 template<class Type, template<class> class PatchField, class GeoMesh>
1305  const dimensioned<Type>& minDt,
1306  const dimensioned<Type>& maxDt
1307 )
1308 {
1309  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1310  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1311  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1312  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1313 }
1314 
1315 
1316 template<class Type, template<class> class PatchField, class GeoMesh>
1318 {
1319  primitiveFieldRef().negate();
1320  boundaryFieldRef().negate();
1321 }
1322 
1323 
1324 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1325 
1326 template<class Type, template<class> class PatchField, class GeoMesh>
1327 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1330 )
1331 {
1332  if (this == &gf)
1333  {
1335  << "attempted assignment to self"
1336  << abort(FatalError);
1337  }
1338 
1339  checkField(*this, gf, "=");
1340 
1341  // Only assign field contents not ID
1342 
1343  ref() = gf();
1344  boundaryFieldRef() = gf.boundaryField();
1345 }
1346 
1347 
1348 template<class Type, template<class> class PatchField, class GeoMesh>
1349 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1352 )
1353 {
1354  if (this == &gf)
1355  {
1357  << "attempted assignment to self"
1358  << abort(FatalError);
1359  }
1360 
1361  checkField(*this, gf, "=");
1362 
1363  // Only assign field contents not ID
1364 
1365  ref() = move(gf());
1366  boundaryFieldRef() = move(gf.boundaryField());
1367 }
1368 
1369 
1370 template<class Type, template<class> class PatchField, class GeoMesh>
1371 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1374 )
1375 {
1376  if (this == &(tgf()))
1377  {
1379  << "attempted assignment to self"
1380  << abort(FatalError);
1381  }
1382 
1384 
1385  checkField(*this, gf, "=");
1386 
1387  // Only assign field contents not ID
1388 
1389  this->dimensions() = gf.dimensions();
1390 
1391  if (tgf.isTmp())
1392  {
1393  // Transfer the storage from the tmp
1394  primitiveFieldRef().transfer
1395  (
1396  const_cast<Field<Type>&>(gf.primitiveField())
1397  );
1398  }
1399  else
1400  {
1402  }
1403 
1404  boundaryFieldRef() = gf.boundaryField();
1405 
1406  tgf.clear();
1407 }
1408 
1409 
1410 template<class Type, template<class> class PatchField, class GeoMesh>
1411 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1413  const dimensioned<Type>& dt
1414 )
1415 {
1416  ref() = dt;
1417  boundaryFieldRef() = dt.value();
1418 }
1419 
1420 
1421 template<class Type, template<class> class PatchField, class GeoMesh>
1422 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1424  const zero&
1425 )
1426 {
1427  ref() = Zero;
1428  boundaryFieldRef() = Zero;
1429 }
1430 
1431 
1432 template<class Type, template<class> class PatchField, class GeoMesh>
1433 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1436 )
1437 {
1439 
1440  checkField(*this, gf, "==");
1441 
1442  // Only assign field contents not ID
1443 
1444  ref() = gf();
1445  boundaryFieldRef() == gf.boundaryField();
1446 
1447  tgf.clear();
1448 }
1449 
1450 
1451 template<class Type, template<class> class PatchField, class GeoMesh>
1452 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1454  const dimensioned<Type>& dt
1455 )
1456 {
1457  ref() = dt;
1458  boundaryFieldRef() == dt.value();
1459 }
1460 
1461 
1462 template<class Type, template<class> class PatchField, class GeoMesh>
1463 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1465  const zero&
1466 )
1467 {
1468  ref() = Zero;
1469  boundaryFieldRef() == Zero;
1470 }
1471 
1472 
1473 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1474  \
1475 template<class Type, template<class> class PatchField, class GeoMesh> \
1476 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1477 ( \
1478  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1479 ) \
1480 { \
1481  checkField(*this, gf, #op); \
1482  \
1483  ref() op gf(); \
1484  boundaryFieldRef() op gf.boundaryField(); \
1485 } \
1486  \
1487 template<class Type, template<class> class PatchField, class GeoMesh> \
1488 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1489 ( \
1490  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1491 ) \
1492 { \
1493  operator op(tgf()); \
1494  tgf.clear(); \
1495 } \
1496  \
1497 template<class Type, template<class> class PatchField, class GeoMesh> \
1498 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1499 ( \
1500  const dimensioned<TYPE>& dt \
1501 ) \
1502 { \
1503  ref() op dt; \
1504  boundaryFieldRef() op dt.value(); \
1505 }
1506 
1511 
1512 #undef COMPUTED_ASSIGNMENT
1513 
1514 
1515 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1516 
1517 template<class Type, template<class> class PatchField, class GeoMesh>
1518 Foam::Ostream& Foam::operator<<
1519 (
1520  Ostream& os,
1522 )
1523 {
1524  gf().writeData(os, "internalField");
1525  os << nl;
1526  gf.boundaryField().writeEntry("boundaryField", os);
1527 
1528  // Check state of IOstream
1529  os.check
1530  (
1531  "Ostream& operator<<(Ostream&, "
1532  "const GeometricField<Type, PatchField, GeoMesh>&)"
1533  );
1534 
1535  return (os);
1536 }
1537 
1538 
1539 template<class Type, template<class> class PatchField, class GeoMesh>
1540 Foam::Ostream& Foam::operator<<
1541 (
1542  Ostream& os,
1544 )
1545 {
1546  os << tgf();
1547  tgf.clear();
1548  return os;
1549 }
1550 
1551 
1552 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1553 
1554 #undef checkField
1555 
1556 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1557 
1558 #include "GeometricBoundaryField.C"
1559 #include "GeometricFieldFunctions.C"
1560 
1561 // ************************************************************************* //
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
#define COMPUTED_ASSIGNMENT(TYPE, op)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
#define checkField(gf1, gf2, op)
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:59
const word & name() const
Return name.
Definition: IOobject.H:295
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
UEqn relax()
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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Field< Type >::cmptType cmptType
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
Generic dimensioned Type class.
label nOldTimes() const
Return the number of old time fields stored.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
conserve primitiveFieldRef()+
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
bool needReference() const
Does the field need a reference level for solution.
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
void storeOldTimes() const
Store the old-time fields.
Pre-declare SubField and related Field type.
Definition: Field.H:56
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
void storeOldTime() const
Store the old-time field.
const Type & value() const
Return const reference to value.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
word timeName
Definition: getTimeIndex.H:3
void min(const dimensioned< Type > &)
static const zero Zero
Definition: zero.H:97
faceListList boundary(nPatches)
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch type.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:265
rDeltaT ref()
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual ~GeometricField()
Destructor.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label timeIndex() const
Return the time index of the field.
Internal & ref()
Return a reference to the dimensioned internal field.
Template functions to aid in the implementation of demand driven data.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void max(const dimensioned< Type > &)
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void correctBoundaryConditions()
Correct boundary field.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
void storePrevIter() const
Store the field as the previous iteration value.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
void deleteDemandDrivenData(DataPtr &dataPtr)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
void relax()
Relax field (for steady-state solution).
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.