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-2020 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_(*this, 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 Internal& diField,
752  const PtrList<PatchField<Type>>& ptfl
753 )
754 {
755  const bool cacheTmp = diField.mesh().thisDb().cacheTemporaryObject(name);
756 
758  (
760  (
761  IOobject
762  (
763  name,
764  diField.mesh().thisDb().time().timeName(),
765  diField.mesh().thisDb(),
766  IOobject::NO_READ,
767  IOobject::NO_WRITE,
768  cacheTmp
769  ),
770  diField,
771  ptfl
772  ),
773  cacheTmp
774  );
775 }
776 
777 
778 template<class Type, template<class> class PatchField, class GeoMesh>
781 (
782  const word& name,
783  const Mesh& mesh,
784  const dimensionSet& ds,
785  const word& patchFieldType
786 )
787 {
788  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
789 
791  (
793  (
794  IOobject
795  (
796  name,
797  mesh.thisDb().time().timeName(),
798  mesh.thisDb(),
799  IOobject::NO_READ,
800  IOobject::NO_WRITE,
801  cacheTmp
802  ),
803  mesh,
804  ds,
805  patchFieldType
806  ),
807  cacheTmp
808  );
809 }
810 
811 
812 template<class Type, template<class> class PatchField, class GeoMesh>
815 (
816  const word& name,
817  const Mesh& mesh,
818  const dimensioned<Type>& dt,
819  const word& patchFieldType
820 )
821 {
822  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
823 
825  (
827  (
828  IOobject
829  (
830  name,
831  mesh.thisDb().time().timeName(),
832  mesh.thisDb(),
833  IOobject::NO_READ,
834  IOobject::NO_WRITE,
835  cacheTmp
836  ),
837  mesh,
838  dt,
839  patchFieldType
840  ),
841  cacheTmp
842  );
843 }
844 
845 
846 
847 template<class Type, template<class> class PatchField, class GeoMesh>
850 (
851  const word& name,
852  const Mesh& mesh,
853  const dimensioned<Type>& dt,
854  const wordList& patchFieldTypes,
855  const wordList& actualPatchTypes
856 )
857 {
858  const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
859 
861  (
863  (
864  IOobject
865  (
866  name,
867  mesh.thisDb().time().timeName(),
868  mesh.thisDb(),
869  IOobject::NO_READ,
870  IOobject::NO_WRITE,
871  cacheTmp
872  ),
873  mesh,
874  dt,
875  patchFieldTypes,
876  actualPatchTypes
877  ),
878  cacheTmp
879  );
880 }
881 
882 
883 template<class Type, template<class> class PatchField, class GeoMesh>
886 (
887  const word& newName,
889 )
890 {
891  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
892 
894  (
896  (
897  IOobject
898  (
899  newName,
900  tgf().instance(),
901  tgf().local(),
902  tgf().db(),
903  IOobject::NO_READ,
904  IOobject::NO_WRITE,
905  cacheTmp
906  ),
907  tgf
908  ),
909  cacheTmp
910  );
911 }
912 
913 
914 template<class Type, template<class> class PatchField, class GeoMesh>
917 (
918  const word& newName,
920  const word& patchFieldType
921 )
922 {
923  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
924 
926  (
928  (
929  IOobject
930  (
931  newName,
932  tgf().instance(),
933  tgf().local(),
934  tgf().db(),
935  IOobject::NO_READ,
936  IOobject::NO_WRITE,
937  cacheTmp
938  ),
939  tgf,
941  ),
942  cacheTmp
943  );
944 }
945 
946 
947 template<class Type, template<class> class PatchField, class GeoMesh>
950 (
951  const word& newName,
953  const wordList& patchFieldTypes,
954  const wordList& actualPatchTypes
955 )
956 {
957  const bool cacheTmp = tgf().db().cacheTemporaryObject(newName);
958 
960  (
962  (
963  IOobject
964  (
965  newName,
966  tgf().instance(),
967  tgf().local(),
968  tgf().db(),
969  IOobject::NO_READ,
970  IOobject::NO_WRITE,
971  cacheTmp
972  ),
973  tgf,
974  patchFieldTypes,
975  actualPatchTypes
976  ),
977  cacheTmp
978  );
979 }
980 
981 
982 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
983 
984 template<class Type, template<class> class PatchField, class GeoMesh>
986 {
987  this->db().cacheTemporaryObject(*this);
988 
989  deleteDemandDrivenData(field0Ptr_);
990  deleteDemandDrivenData(fieldPrevIterPtr_);
991 }
992 
993 
994 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
995 
996 template<class Type, template<class> class PatchField, class GeoMesh>
997 typename
1000 {
1001  this->setUpToDate();
1002  storeOldTimes();
1003  return *this;
1004 }
1005 
1006 
1007 template<class Type, template<class> class PatchField, class GeoMesh>
1008 typename
1011 {
1012  this->setUpToDate();
1013  storeOldTimes();
1014  return *this;
1015 }
1016 
1017 
1018 template<class Type, template<class> class PatchField, class GeoMesh>
1019 typename
1022 {
1023  this->setUpToDate();
1024  storeOldTimes();
1025  return boundaryField_;
1026 }
1027 
1028 
1029 template<class Type, template<class> class PatchField, class GeoMesh>
1031 {
1032  if
1033  (
1034  field0Ptr_
1035  && timeIndex_ != this->time().timeIndex()
1036  && !(
1037  this->name().size() > 2
1038  && this->name()(this->name().size()-2, 2) == "_0"
1039  )
1040  )
1041  {
1042  storeOldTime();
1043  }
1044 
1045  // Correct time index
1046  timeIndex_ = this->time().timeIndex();
1047 }
1048 
1049 
1050 template<class Type, template<class> class PatchField, class GeoMesh>
1052 {
1053  if (field0Ptr_)
1054  {
1055  field0Ptr_->storeOldTime();
1056 
1057  if (debug)
1058  {
1060  << "Storing old time field for field" << endl
1061  << this->info() << endl;
1062  }
1063 
1064  *field0Ptr_ == *this;
1065  field0Ptr_->timeIndex_ = timeIndex_;
1066 
1067  if (field0Ptr_->field0Ptr_)
1068  {
1069  field0Ptr_->writeOpt() = this->writeOpt();
1070  }
1071  }
1072 }
1073 
1074 
1075 template<class Type, template<class> class PatchField, class GeoMesh>
1077 {
1078  if (field0Ptr_)
1079  {
1080  return field0Ptr_->nOldTimes() + 1;
1081  }
1082  else
1083  {
1084  return 0;
1085  }
1086 }
1087 
1088 
1089 template<class Type, template<class> class PatchField, class GeoMesh>
1092 {
1093  if (!field0Ptr_)
1094  {
1096  (
1097  IOobject
1098  (
1099  this->name() + "_0",
1100  this->time().timeName(),
1101  this->db(),
1102  IOobject::NO_READ,
1103  IOobject::NO_WRITE,
1104  this->registerObject()
1105  ),
1106  *this
1107  );
1108  }
1109  else
1110  {
1111  storeOldTimes();
1112  }
1113 
1114  return *field0Ptr_;
1115 }
1116 
1117 
1118 template<class Type, template<class> class PatchField, class GeoMesh>
1121 {
1122  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1123  .oldTime();
1124 
1125  return *field0Ptr_;
1126 }
1127 
1128 
1129 template<class Type, template<class> class PatchField, class GeoMesh>
1131 {
1132  if (!fieldPrevIterPtr_)
1133  {
1134  if (debug)
1135  {
1137  << "Allocating previous iteration field" << endl
1138  << this->info() << endl;
1139  }
1140 
1141  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1142  (
1143  this->name() + "PrevIter",
1144  *this
1145  );
1146  }
1147  else
1148  {
1149  *fieldPrevIterPtr_ == *this;
1150  }
1151 }
1152 
1153 
1154 template<class Type, template<class> class PatchField, class GeoMesh>
1157 {
1158  if (!fieldPrevIterPtr_)
1159  {
1161  << "previous iteration field" << endl << this->info() << endl
1162  << " not stored."
1163  << " Use field.storePrevIter() at start of iteration."
1164  << abort(FatalError);
1165  }
1166 
1167  return *fieldPrevIterPtr_;
1168 }
1169 
1170 
1171 template<class Type, template<class> class PatchField, class GeoMesh>
1174 {
1175  this->setUpToDate();
1176  storeOldTimes();
1177  boundaryField_.evaluate();
1178 }
1179 
1180 
1181 template<class Type, template<class> class PatchField, class GeoMesh>
1183 {
1184  // Search all boundary conditions, if any are
1185  // fixed-value or mixed (Robin) do not set reference level for solution.
1186 
1187  bool needRef = true;
1188 
1189  forAll(boundaryField_, patchi)
1190  {
1191  if (boundaryField_[patchi].fixesValue())
1192  {
1193  needRef = false;
1194  break;
1195  }
1196  }
1197 
1198  reduce(needRef, andOp<bool>());
1199 
1200  return needRef;
1201 }
1202 
1203 
1204 template<class Type, template<class> class PatchField, class GeoMesh>
1206 {
1207  if (debug)
1208  {
1210  << "Relaxing" << endl << this->info() << " by " << alpha << endl;
1211  }
1212 
1213  operator==(prevIter() + alpha*(*this - prevIter()));
1214 }
1215 
1216 
1217 template<class Type, template<class> class PatchField, class GeoMesh>
1219 {
1220  if
1221  (
1222  this->mesh().data::template lookupOrDefault<bool>
1223  (
1224  "finalIteration",
1225  false
1226  )
1227  && this->mesh().relaxField(this->name() + "Final")
1228  )
1229  {
1230  relax(this->mesh().fieldRelaxationFactor(this->name() + "Final"));
1231  }
1232  else if (this->mesh().relaxField(this->name()))
1233  {
1234  relax(this->mesh().fieldRelaxationFactor(this->name()));
1235  }
1236 }
1237 
1238 
1239 template<class Type, template<class> class PatchField, class GeoMesh>
1242  bool final
1243 ) const
1244 {
1245  if (final)
1246  {
1247  return this->name() + "Final";
1248  }
1249  else
1250  {
1251  return this->name();
1252  }
1253 }
1254 
1255 
1256 template<class Type, template<class> class PatchField, class GeoMesh>
1259  Ostream& os
1260 ) const
1261 {
1262  os << "min/max(" << this->name() << ") = "
1263  << Foam::min(*this).value() << ", "
1264  << Foam::max(*this).value()
1265  << endl;
1266 }
1267 
1268 
1269 template<class Type, template<class> class PatchField, class GeoMesh>
1271 writeData(Ostream& os) const
1272 {
1273  os << *this;
1274  return os.good();
1275 }
1276 
1277 
1278 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1279 
1280 template<class Type, template<class> class PatchField, class GeoMesh>
1283 {
1285  (
1287  (
1288  this->name() + ".T()",
1289  this->mesh(),
1290  this->dimensions()
1291  )
1292  );
1293 
1294  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1295  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1296 
1297  return result;
1298 }
1299 
1300 
1301 template<class Type, template<class> class PatchField, class GeoMesh>
1302 Foam::tmp
1303 <
1305  <
1306  typename Foam::GeometricField<Type, PatchField, GeoMesh>::cmptType,
1307  PatchField,
1308  GeoMesh
1309  >
1310 >
1312 (
1313  const direction d
1314 ) const
1315 {
1317  (
1319  (
1320  this->name() + ".component(" + Foam::name(d) + ')',
1321  this->mesh(),
1322  this->dimensions()
1323  )
1324  );
1325 
1326  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1327  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1328 
1329  return Component;
1330 }
1331 
1332 
1333 template<class Type, template<class> class PatchField, class GeoMesh>
1335 (
1336  const direction d,
1337  const GeometricField
1338  <
1340  PatchField,
1341  GeoMesh
1342  >& gcf
1343 )
1344 {
1345  primitiveFieldRef().replace(d, gcf.primitiveField());
1346  boundaryFieldRef().replace(d, gcf.boundaryField());
1347 }
1348 
1349 
1350 template<class Type, template<class> class PatchField, class GeoMesh>
1353  const direction d,
1354  const dimensioned<cmptType>& ds
1355 )
1356 {
1357  primitiveFieldRef().replace(d, ds.value());
1358  boundaryFieldRef().replace(d, ds.value());
1359 }
1360 
1361 
1362 template<class Type, template<class> class PatchField, class GeoMesh>
1365  const dimensioned<Type>& dt
1366 )
1367 {
1368  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1369  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1370 }
1371 
1372 
1373 template<class Type, template<class> class PatchField, class GeoMesh>
1376  const dimensioned<Type>& dt
1377 )
1378 {
1379  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1380  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1381 }
1382 
1383 
1384 template<class Type, template<class> class PatchField, class GeoMesh>
1387  const dimensioned<Type>& minDt,
1388  const dimensioned<Type>& maxDt
1389 )
1390 {
1391  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1392  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1393  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1394  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1395 }
1396 
1397 
1398 template<class Type, template<class> class PatchField, class GeoMesh>
1400 {
1401  primitiveFieldRef().negate();
1402  boundaryFieldRef().negate();
1403 }
1404 
1405 
1406 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1407 
1408 template<class Type, template<class> class PatchField, class GeoMesh>
1409 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1412 )
1413 {
1414  if (this == &gf)
1415  {
1417  << "attempted assignment to self"
1418  << abort(FatalError);
1419  }
1420 
1421  checkField(*this, gf, "=");
1422 
1423  // Only assign field contents not ID
1424 
1425  ref() = gf();
1426  boundaryFieldRef() = gf.boundaryField();
1427 }
1428 
1429 
1430 template<class Type, template<class> class PatchField, class GeoMesh>
1431 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1434 )
1435 {
1436  if (this == &gf)
1437  {
1439  << "attempted assignment to self"
1440  << abort(FatalError);
1441  }
1442 
1443  checkField(*this, gf, "=");
1444 
1445  // Only assign field contents not ID
1446 
1447  ref() = move(gf());
1448  boundaryFieldRef() = move(gf.boundaryField());
1449 }
1450 
1451 
1452 template<class Type, template<class> class PatchField, class GeoMesh>
1453 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1456 )
1457 {
1458  if (this == &(tgf()))
1459  {
1461  << "attempted assignment to self"
1462  << abort(FatalError);
1463  }
1464 
1466 
1467  checkField(*this, gf, "=");
1468 
1469  // Only assign field contents not ID
1470 
1471  this->dimensions() = gf.dimensions();
1472 
1473  if (tgf.isTmp())
1474  {
1475  // Transfer the storage from the tmp
1476  primitiveFieldRef().transfer
1477  (
1478  const_cast<Field<Type>&>(gf.primitiveField())
1479  );
1480  }
1481  else
1482  {
1484  }
1485 
1486  boundaryFieldRef() = gf.boundaryField();
1487 
1488  tgf.clear();
1489 }
1490 
1491 
1492 template<class Type, template<class> class PatchField, class GeoMesh>
1493 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1495  const dimensioned<Type>& dt
1496 )
1497 {
1498  ref() = dt;
1499  boundaryFieldRef() = dt.value();
1500 }
1501 
1502 
1503 template<class Type, template<class> class PatchField, class GeoMesh>
1504 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1506  const zero&
1507 )
1508 {
1509  ref() = Zero;
1510  boundaryFieldRef() = Zero;
1511 }
1512 
1513 
1514 template<class Type, template<class> class PatchField, class GeoMesh>
1515 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1518 )
1519 {
1521 
1522  checkField(*this, gf, "==");
1523 
1524  // Only assign field contents not ID
1525 
1526  ref() = gf();
1527  boundaryFieldRef() == gf.boundaryField();
1528 
1529  tgf.clear();
1530 }
1531 
1532 
1533 template<class Type, template<class> class PatchField, class GeoMesh>
1534 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1536  const dimensioned<Type>& dt
1537 )
1538 {
1539  ref() = dt;
1540  boundaryFieldRef() == dt.value();
1541 }
1542 
1543 
1544 template<class Type, template<class> class PatchField, class GeoMesh>
1545 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1547  const zero&
1548 )
1549 {
1550  ref() = Zero;
1551  boundaryFieldRef() == Zero;
1552 }
1553 
1554 
1555 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1556  \
1557 template<class Type, template<class> class PatchField, class GeoMesh> \
1558 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1559 ( \
1560  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1561 ) \
1562 { \
1563  checkField(*this, gf, #op); \
1564  \
1565  ref() op gf(); \
1566  boundaryFieldRef() op gf.boundaryField(); \
1567 } \
1568  \
1569 template<class Type, template<class> class PatchField, class GeoMesh> \
1570 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1571 ( \
1572  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1573 ) \
1574 { \
1575  operator op(tgf()); \
1576  tgf.clear(); \
1577 } \
1578  \
1579 template<class Type, template<class> class PatchField, class GeoMesh> \
1580 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1581 ( \
1582  const dimensioned<TYPE>& dt \
1583 ) \
1584 { \
1585  ref() op dt; \
1586  boundaryFieldRef() op dt.value(); \
1587 }
1588 
1593 
1594 #undef COMPUTED_ASSIGNMENT
1595 
1596 
1597 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1598 
1599 template<class Type, template<class> class PatchField, class GeoMesh>
1600 Foam::Ostream& Foam::operator<<
1601 (
1602  Ostream& os,
1604 )
1605 {
1606  gf().writeData(os, "internalField");
1607  os << nl;
1608  gf.boundaryField().writeEntry("boundaryField", os);
1609 
1610  // Check state of IOstream
1611  os.check
1612  (
1613  "Ostream& operator<<(Ostream&, "
1614  "const GeometricField<Type, PatchField, GeoMesh>&)"
1615  );
1616 
1617  return (os);
1618 }
1619 
1620 
1621 template<class Type, template<class> class PatchField, class GeoMesh>
1622 Foam::Ostream& Foam::operator<<
1623 (
1624  Ostream& os,
1626 )
1627 {
1628  os << tgf();
1629  tgf.clear();
1630  return os;
1631 }
1632 
1633 
1634 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1635 
1636 #undef checkField
1637 
1638 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1639 
1640 #include "GeometricBoundaryField.C"
1641 #include "GeometricFieldFunctions.C"
1642 
1643 // ************************************************************************* //
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:303
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
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:181
uint8_t direction
Definition: direction.H:45
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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.
const dimensionSet dimless
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:330
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
conserve primitiveFieldRef()+
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
EaEqn relax()
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch field type.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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:260
const Mesh & mesh() const
Return mesh.
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
#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:335
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.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
rDeltaT ref()
void deleteDemandDrivenData(DataPtr &dataPtr)
word patchFieldType(const PatchField &pf)
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.