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_(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 Internal& diField,
752  const PtrList<PatchField<Type>>& ptfl
753 )
754 {
756  (
758  (
759  IOobject
760  (
761  name,
762  diField.mesh().thisDb().time().timeName(),
763  diField.mesh().thisDb(),
764  IOobject::NO_READ,
765  IOobject::NO_WRITE,
766  diField.mesh().thisDb().cacheTemporaryObject(name)
767  ),
768  diField,
769  ptfl
770  )
771  );
772 }
773 
774 
775 template<class Type, template<class> class PatchField, class GeoMesh>
778 (
779  const word& name,
780  const Mesh& mesh,
781  const dimensionSet& ds,
782  const word& patchFieldType
783 )
784 {
786  (
788  (
789  IOobject
790  (
791  name,
792  mesh.thisDb().time().timeName(),
793  mesh.thisDb(),
794  IOobject::NO_READ,
795  IOobject::NO_WRITE,
796  mesh.thisDb().cacheTemporaryObject(name)
797  ),
798  mesh,
799  ds,
800  patchFieldType
801  )
802  );
803 }
804 
805 
806 template<class Type, template<class> class PatchField, class GeoMesh>
809 (
810  const word& name,
811  const Mesh& mesh,
812  const dimensioned<Type>& dt,
813  const word& patchFieldType
814 )
815 {
817  (
819  (
820  IOobject
821  (
822  name,
823  mesh.thisDb().time().timeName(),
824  mesh.thisDb(),
825  IOobject::NO_READ,
826  IOobject::NO_WRITE,
827  mesh.thisDb().cacheTemporaryObject(name)
828  ),
829  mesh,
830  dt,
831  patchFieldType
832  )
833  );
834 }
835 
836 
837 
838 template<class Type, template<class> class PatchField, class GeoMesh>
841 (
842  const word& name,
843  const Mesh& mesh,
844  const dimensioned<Type>& dt,
845  const wordList& patchFieldTypes,
846  const wordList& actualPatchTypes
847 )
848 {
850  (
852  (
853  IOobject
854  (
855  name,
856  mesh.thisDb().time().timeName(),
857  mesh.thisDb(),
858  IOobject::NO_READ,
859  IOobject::NO_WRITE,
860  mesh.thisDb().cacheTemporaryObject(name)
861  ),
862  mesh,
863  dt,
864  patchFieldTypes,
865  actualPatchTypes
866  )
867  );
868 }
869 
870 
871 template<class Type, template<class> class PatchField, class GeoMesh>
874 (
875  const word& newName,
877 )
878 {
880  (
882  (
883  IOobject
884  (
885  newName,
886  tgf().instance(),
887  tgf().local(),
888  tgf().db(),
889  IOobject::NO_READ,
890  IOobject::NO_WRITE,
891  tgf().db().cacheTemporaryObject(newName)
892  ),
893  tgf
894  )
895  );
896 }
897 
898 
899 template<class Type, template<class> class PatchField, class GeoMesh>
902 (
903  const word& newName,
905  const word& patchFieldType
906 )
907 {
909  (
911  (
912  IOobject
913  (
914  newName,
915  tgf().instance(),
916  tgf().local(),
917  tgf().db(),
918  IOobject::NO_READ,
919  IOobject::NO_WRITE,
920  tgf().db().cacheTemporaryObject(newName)
921  ),
922  tgf,
924  )
925  );
926 }
927 
928 
929 template<class Type, template<class> class PatchField, class GeoMesh>
932 (
933  const word& newName,
935  const wordList& patchFieldTypes,
936  const wordList& actualPatchTypes
937 )
938 {
940  (
942  (
943  IOobject
944  (
945  newName,
946  tgf().instance(),
947  tgf().local(),
948  tgf().db(),
949  IOobject::NO_READ,
950  IOobject::NO_WRITE,
951  tgf().db().cacheTemporaryObject(newName)
952  ),
953  tgf,
954  patchFieldTypes,
955  actualPatchTypes
956  )
957  );
958 }
959 
960 
961 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
962 
963 template<class Type, template<class> class PatchField, class GeoMesh>
965 {
966  this->db().cacheTemporaryObject(*this);
967 
968  deleteDemandDrivenData(field0Ptr_);
969  deleteDemandDrivenData(fieldPrevIterPtr_);
970 }
971 
972 
973 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
974 
975 template<class Type, template<class> class PatchField, class GeoMesh>
976 typename
979 {
980  this->setUpToDate();
981  storeOldTimes();
982  return *this;
983 }
984 
985 
986 template<class Type, template<class> class PatchField, class GeoMesh>
987 typename
990 {
991  this->setUpToDate();
992  storeOldTimes();
993  return *this;
994 }
995 
996 
997 template<class Type, template<class> class PatchField, class GeoMesh>
998 typename
1001 {
1002  this->setUpToDate();
1003  storeOldTimes();
1004  return boundaryField_;
1005 }
1006 
1007 
1008 template<class Type, template<class> class PatchField, class GeoMesh>
1010 {
1011  if
1012  (
1013  field0Ptr_
1014  && timeIndex_ != this->time().timeIndex()
1015  && !(
1016  this->name().size() > 2
1017  && this->name()(this->name().size()-2, 2) == "_0"
1018  )
1019  )
1020  {
1021  storeOldTime();
1022  }
1023 
1024  // Correct time index
1025  timeIndex_ = this->time().timeIndex();
1026 }
1027 
1028 
1029 template<class Type, template<class> class PatchField, class GeoMesh>
1031 {
1032  if (field0Ptr_)
1033  {
1034  field0Ptr_->storeOldTime();
1035 
1036  if (debug)
1037  {
1039  << "Storing old time field for field" << endl
1040  << this->info() << endl;
1041  }
1042 
1043  *field0Ptr_ == *this;
1044  field0Ptr_->timeIndex_ = timeIndex_;
1045 
1046  if (field0Ptr_->field0Ptr_)
1047  {
1048  field0Ptr_->writeOpt() = this->writeOpt();
1049  }
1050  }
1051 }
1052 
1053 
1054 template<class Type, template<class> class PatchField, class GeoMesh>
1056 {
1057  if (field0Ptr_)
1058  {
1059  return field0Ptr_->nOldTimes() + 1;
1060  }
1061  else
1062  {
1063  return 0;
1064  }
1065 }
1066 
1067 
1068 template<class Type, template<class> class PatchField, class GeoMesh>
1071 {
1072  if (!field0Ptr_)
1073  {
1075  (
1076  IOobject
1077  (
1078  this->name() + "_0",
1079  this->time().timeName(),
1080  this->db(),
1081  IOobject::NO_READ,
1082  IOobject::NO_WRITE,
1083  this->registerObject()
1084  ),
1085  *this
1086  );
1087  }
1088  else
1089  {
1090  storeOldTimes();
1091  }
1092 
1093  return *field0Ptr_;
1094 }
1095 
1096 
1097 template<class Type, template<class> class PatchField, class GeoMesh>
1100 {
1101  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1102  .oldTime();
1103 
1104  return *field0Ptr_;
1105 }
1106 
1107 
1108 template<class Type, template<class> class PatchField, class GeoMesh>
1110 {
1111  if (!fieldPrevIterPtr_)
1112  {
1113  if (debug)
1114  {
1116  << "Allocating previous iteration field" << endl
1117  << this->info() << endl;
1118  }
1119 
1120  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1121  (
1122  this->name() + "PrevIter",
1123  *this
1124  );
1125  }
1126  else
1127  {
1128  *fieldPrevIterPtr_ == *this;
1129  }
1130 }
1131 
1132 
1133 template<class Type, template<class> class PatchField, class GeoMesh>
1136 {
1137  if (!fieldPrevIterPtr_)
1138  {
1140  << "previous iteration field" << endl << this->info() << endl
1141  << " not stored."
1142  << " Use field.storePrevIter() at start of iteration."
1143  << abort(FatalError);
1144  }
1145 
1146  return *fieldPrevIterPtr_;
1147 }
1148 
1149 
1150 template<class Type, template<class> class PatchField, class GeoMesh>
1153 {
1154  this->setUpToDate();
1155  storeOldTimes();
1156  boundaryField_.evaluate();
1157 }
1158 
1159 
1160 template<class Type, template<class> class PatchField, class GeoMesh>
1162 {
1163  // Search all boundary conditions, if any are
1164  // fixed-value or mixed (Robin) do not set reference level for solution.
1165 
1166  bool needRef = true;
1167 
1168  forAll(boundaryField_, patchi)
1169  {
1170  if (boundaryField_[patchi].fixesValue())
1171  {
1172  needRef = false;
1173  break;
1174  }
1175  }
1176 
1177  reduce(needRef, andOp<bool>());
1178 
1179  return needRef;
1180 }
1181 
1182 
1183 template<class Type, template<class> class PatchField, class GeoMesh>
1185 {
1186  if (debug)
1187  {
1189  << "Relaxing" << endl << this->info() << " by " << alpha << endl;
1190  }
1191 
1192  operator==(prevIter() + alpha*(*this - prevIter()));
1193 }
1194 
1195 
1196 template<class Type, template<class> class PatchField, class GeoMesh>
1198 {
1199  if
1200  (
1201  this->mesh().data::template lookupOrDefault<bool>
1202  (
1203  "finalIteration",
1204  false
1205  )
1206  && this->mesh().relaxField(this->name() + "Final")
1207  )
1208  {
1209  relax(this->mesh().fieldRelaxationFactor(this->name() + "Final"));
1210  }
1211  else if (this->mesh().relaxField(this->name()))
1212  {
1213  relax(this->mesh().fieldRelaxationFactor(this->name()));
1214  }
1215 }
1216 
1217 
1218 template<class Type, template<class> class PatchField, class GeoMesh>
1221  bool final
1222 ) const
1223 {
1224  if (final)
1225  {
1226  return this->name() + "Final";
1227  }
1228  else
1229  {
1230  return this->name();
1231  }
1232 }
1233 
1234 
1235 template<class Type, template<class> class PatchField, class GeoMesh>
1238  Ostream& os
1239 ) const
1240 {
1241  os << "min/max(" << this->name() << ") = "
1242  << Foam::min(*this).value() << ", "
1243  << Foam::max(*this).value()
1244  << endl;
1245 }
1246 
1247 
1248 template<class Type, template<class> class PatchField, class GeoMesh>
1250 writeData(Ostream& os) const
1251 {
1252  os << *this;
1253  return os.good();
1254 }
1255 
1256 
1257 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1258 
1259 template<class Type, template<class> class PatchField, class GeoMesh>
1262 {
1264  (
1266  (
1267  this->name() + ".T()",
1268  this->mesh(),
1269  this->dimensions()
1270  )
1271  );
1272 
1273  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1274  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1275 
1276  return result;
1277 }
1278 
1279 
1280 template<class Type, template<class> class PatchField, class GeoMesh>
1281 Foam::tmp
1282 <
1284  <
1285  typename Foam::GeometricField<Type, PatchField, GeoMesh>::cmptType,
1286  PatchField,
1287  GeoMesh
1288  >
1289 >
1291 (
1292  const direction d
1293 ) const
1294 {
1296  (
1298  (
1299  this->name() + ".component(" + Foam::name(d) + ')',
1300  this->mesh(),
1301  this->dimensions()
1302  )
1303  );
1304 
1305  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1306  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1307 
1308  return Component;
1309 }
1310 
1311 
1312 template<class Type, template<class> class PatchField, class GeoMesh>
1314 (
1315  const direction d,
1316  const GeometricField
1317  <
1319  PatchField,
1320  GeoMesh
1321  >& gcf
1322 )
1323 {
1324  primitiveFieldRef().replace(d, gcf.primitiveField());
1325  boundaryFieldRef().replace(d, gcf.boundaryField());
1326 }
1327 
1328 
1329 template<class Type, template<class> class PatchField, class GeoMesh>
1332  const direction d,
1333  const dimensioned<cmptType>& ds
1334 )
1335 {
1336  primitiveFieldRef().replace(d, ds.value());
1337  boundaryFieldRef().replace(d, ds.value());
1338 }
1339 
1340 
1341 template<class Type, template<class> class PatchField, class GeoMesh>
1344  const dimensioned<Type>& dt
1345 )
1346 {
1347  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1348  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1349 }
1350 
1351 
1352 template<class Type, template<class> class PatchField, class GeoMesh>
1355  const dimensioned<Type>& dt
1356 )
1357 {
1358  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1359  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1360 }
1361 
1362 
1363 template<class Type, template<class> class PatchField, class GeoMesh>
1366  const dimensioned<Type>& minDt,
1367  const dimensioned<Type>& maxDt
1368 )
1369 {
1370  Foam::max(primitiveFieldRef(), primitiveField(), minDt.value());
1371  Foam::max(boundaryFieldRef(), boundaryField(), minDt.value());
1372  Foam::min(primitiveFieldRef(), primitiveField(), maxDt.value());
1373  Foam::min(boundaryFieldRef(), boundaryField(), maxDt.value());
1374 }
1375 
1376 
1377 template<class Type, template<class> class PatchField, class GeoMesh>
1379 {
1380  primitiveFieldRef().negate();
1381  boundaryFieldRef().negate();
1382 }
1383 
1384 
1385 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1386 
1387 template<class Type, template<class> class PatchField, class GeoMesh>
1388 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1391 )
1392 {
1393  if (this == &gf)
1394  {
1396  << "attempted assignment to self"
1397  << abort(FatalError);
1398  }
1399 
1400  checkField(*this, gf, "=");
1401 
1402  // Only assign field contents not ID
1403 
1404  ref() = gf();
1405  boundaryFieldRef() = gf.boundaryField();
1406 }
1407 
1408 
1409 template<class Type, template<class> class PatchField, class GeoMesh>
1410 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1413 )
1414 {
1415  if (this == &gf)
1416  {
1418  << "attempted assignment to self"
1419  << abort(FatalError);
1420  }
1421 
1422  checkField(*this, gf, "=");
1423 
1424  // Only assign field contents not ID
1425 
1426  ref() = move(gf());
1427  boundaryFieldRef() = move(gf.boundaryField());
1428 }
1429 
1430 
1431 template<class Type, template<class> class PatchField, class GeoMesh>
1432 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1435 )
1436 {
1437  if (this == &(tgf()))
1438  {
1440  << "attempted assignment to self"
1441  << abort(FatalError);
1442  }
1443 
1445 
1446  checkField(*this, gf, "=");
1447 
1448  // Only assign field contents not ID
1449 
1450  this->dimensions() = gf.dimensions();
1451 
1452  if (tgf.isTmp())
1453  {
1454  // Transfer the storage from the tmp
1455  primitiveFieldRef().transfer
1456  (
1457  const_cast<Field<Type>&>(gf.primitiveField())
1458  );
1459  }
1460  else
1461  {
1463  }
1464 
1465  boundaryFieldRef() = gf.boundaryField();
1466 
1467  tgf.clear();
1468 }
1469 
1470 
1471 template<class Type, template<class> class PatchField, class GeoMesh>
1472 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1474  const dimensioned<Type>& dt
1475 )
1476 {
1477  ref() = dt;
1478  boundaryFieldRef() = dt.value();
1479 }
1480 
1481 
1482 template<class Type, template<class> class PatchField, class GeoMesh>
1483 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1485  const zero&
1486 )
1487 {
1488  ref() = Zero;
1489  boundaryFieldRef() = Zero;
1490 }
1491 
1492 
1493 template<class Type, template<class> class PatchField, class GeoMesh>
1494 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1497 )
1498 {
1500 
1501  checkField(*this, gf, "==");
1502 
1503  // Only assign field contents not ID
1504 
1505  ref() = gf();
1506  boundaryFieldRef() == gf.boundaryField();
1507 
1508  tgf.clear();
1509 }
1510 
1511 
1512 template<class Type, template<class> class PatchField, class GeoMesh>
1513 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1515  const dimensioned<Type>& dt
1516 )
1517 {
1518  ref() = dt;
1519  boundaryFieldRef() == dt.value();
1520 }
1521 
1522 
1523 template<class Type, template<class> class PatchField, class GeoMesh>
1524 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1526  const zero&
1527 )
1528 {
1529  ref() = Zero;
1530  boundaryFieldRef() == Zero;
1531 }
1532 
1533 
1534 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1535  \
1536 template<class Type, template<class> class PatchField, class GeoMesh> \
1537 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1538 ( \
1539  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1540 ) \
1541 { \
1542  checkField(*this, gf, #op); \
1543  \
1544  ref() op gf(); \
1545  boundaryFieldRef() op gf.boundaryField(); \
1546 } \
1547  \
1548 template<class Type, template<class> class PatchField, class GeoMesh> \
1549 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1550 ( \
1551  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1552 ) \
1553 { \
1554  operator op(tgf()); \
1555  tgf.clear(); \
1556 } \
1557  \
1558 template<class Type, template<class> class PatchField, class GeoMesh> \
1559 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1560 ( \
1561  const dimensioned<TYPE>& dt \
1562 ) \
1563 { \
1564  ref() op dt; \
1565  boundaryFieldRef() op dt.value(); \
1566 }
1567 
1572 
1573 #undef COMPUTED_ASSIGNMENT
1574 
1575 
1576 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1577 
1578 template<class Type, template<class> class PatchField, class GeoMesh>
1579 Foam::Ostream& Foam::operator<<
1580 (
1581  Ostream& os,
1583 )
1584 {
1585  gf().writeData(os, "internalField");
1586  os << nl;
1587  gf.boundaryField().writeEntry("boundaryField", os);
1588 
1589  // Check state of IOstream
1590  os.check
1591  (
1592  "Ostream& operator<<(Ostream&, "
1593  "const GeometricField<Type, PatchField, GeoMesh>&)"
1594  );
1595 
1596  return (os);
1597 }
1598 
1599 
1600 template<class Type, template<class> class PatchField, class GeoMesh>
1601 Foam::Ostream& Foam::operator<<
1602 (
1603  Ostream& os,
1605 )
1606 {
1607  os << tgf();
1608  tgf.clear();
1609  return os;
1610 }
1611 
1612 
1613 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1614 
1615 #undef checkField
1616 
1617 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1618 
1619 #include "GeometricBoundaryField.C"
1620 #include "GeometricFieldFunctions.C"
1621 
1622 // ************************************************************************* //
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
UEqn relax()
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: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
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.
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()+
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 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
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
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.