GeometricField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 "data.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 #define checkField(gf1, gf2, op) \
35 if ((gf1).mesh() != (gf2).mesh()) \
36 { \
37  FatalErrorInFunction \
38  << "different mesh for fields " \
39  << (gf1).name() << " and " << (gf2).name() \
40  << " during operatrion " << op \
41  << abort(FatalError); \
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 template<class Type, template<class> class PatchField, class GeoMesh>
49 (
50  const dictionary& dict
51 )
52 {
53  Internal::readField(dict, "internalField");
54 
55  boundaryField_.readField(*this, dict.subDict("boundaryField"));
56 
57  if (dict.found("referenceLevel"))
58  {
59  Type fieldAverage(pTraits<Type>(dict.lookup("referenceLevel")));
60 
61  Field<Type>::operator+=(fieldAverage);
62 
63  forAll(boundaryField_, patchi)
64  {
65  boundaryField_[patchi] == boundaryField_[patchi] + fieldAverage;
66  }
67  }
68 }
69 
70 
71 template<class Type, template<class> class PatchField, class GeoMesh>
73 {
74  const IOdictionary dict
75  (
76  IOobject
77  (
78  this->name(),
79  this->time().timeName(),
80  this->db(),
81  IOobject::NO_READ,
82  IOobject::NO_WRITE,
83  false
84  ),
85  this->readStream(typeName)
86  );
87 
88  this->close();
89 
90  readFields(dict);
91 }
92 
93 
94 template<class Type, template<class> class PatchField, class GeoMesh>
96 {
97  if
98  (
99  this->readOpt() == IOobject::MUST_READ
100  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
101  )
102  {
104  << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
105  << " suggests that a read constructor for field " << this->name()
106  << " would be more appropriate." << endl;
107  }
108  else if (this->readOpt() == IOobject::READ_IF_PRESENT && this->headerOk())
109  {
110  readFields();
111 
112  // Check compatibility between field and mesh
113  if (this->size() != GeoMesh::size(this->mesh()))
114  {
115  FatalIOErrorInFunction(this->readStream(typeName))
116  << " number of field elements = " << this->size()
117  << " number of mesh elements = "
118  << GeoMesh::size(this->mesh())
119  << exit(FatalIOError);
120  }
121 
122  readOldTimeIfPresent();
123 
124  return true;
125  }
126 
127  return false;
128 }
129 
130 
131 template<class Type, template<class> class PatchField, class GeoMesh>
133 {
134  // Read the old time field if present
135  IOobject field0
136  (
137  this->name() + "_0",
138  this->time().timeName(),
139  this->db(),
140  IOobject::READ_IF_PRESENT,
141  IOobject::AUTO_WRITE,
142  this->registerObject()
143  );
144 
145  if (field0.headerOk())
146  {
147  if (debug)
148  {
149  InfoInFunction << "Reading old time level for field"
150  << endl << this->info() << endl;
151  }
152 
153  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
154  (
155  field0,
156  this->mesh()
157  );
158 
159  field0Ptr_->timeIndex_ = timeIndex_ - 1;
160 
161  if (!field0Ptr_->readOldTimeIfPresent())
162  {
163  field0Ptr_->oldTime();
164  }
165 
166  return true;
167  }
168 
169  return false;
170 }
171 
172 
173 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
174 
175 template<class Type, template<class> class PatchField, class GeoMesh>
177 (
178  const IOobject& io,
179  const Mesh& mesh,
180  const dimensionSet& ds,
181  const word& patchFieldType
182 )
183 :
184  Internal(io, mesh, ds, false),
185  timeIndex_(this->time().timeIndex()),
186  field0Ptr_(NULL),
187  fieldPrevIterPtr_(NULL),
188  boundaryField_(mesh.boundary(), *this, patchFieldType)
189 {
190  if (debug)
191  {
192  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
193  }
194 
195  readIfPresent();
196 }
197 
198 
199 template<class Type, template<class> class PatchField, class GeoMesh>
201 (
202  const IOobject& io,
203  const Mesh& mesh,
204  const dimensionSet& ds,
205  const wordList& patchFieldTypes,
206  const wordList& actualPatchTypes
207 )
208 :
209  Internal(io, mesh, ds, false),
210  timeIndex_(this->time().timeIndex()),
211  field0Ptr_(NULL),
212  fieldPrevIterPtr_(NULL),
213  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
214 {
215  if (debug)
216  {
217  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
218  }
219 
220  readIfPresent();
221 }
222 
223 
224 template<class Type, template<class> class PatchField, class GeoMesh>
226 (
227  const IOobject& io,
228  const Mesh& mesh,
229  const dimensioned<Type>& dt,
230  const word& patchFieldType
231 )
232 :
233  Internal(io, mesh, dt, false),
234  timeIndex_(this->time().timeIndex()),
235  field0Ptr_(NULL),
236  fieldPrevIterPtr_(NULL),
237  boundaryField_(mesh.boundary(), *this, patchFieldType)
238 {
239  if (debug)
240  {
241  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
242  }
243 
244  boundaryField_ == dt.value();
245 
246  readIfPresent();
247 }
248 
249 
250 template<class Type, template<class> class PatchField, class GeoMesh>
252 (
253  const IOobject& io,
254  const Mesh& mesh,
255  const dimensioned<Type>& dt,
256  const wordList& patchFieldTypes,
257  const wordList& actualPatchTypes
258 )
259 :
260  Internal(io, mesh, dt, false),
261  timeIndex_(this->time().timeIndex()),
262  field0Ptr_(NULL),
263  fieldPrevIterPtr_(NULL),
264  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
265 {
266  if (debug)
267  {
268  InfoInFunction << "Creating temporary" << endl << this->info() << endl;
269  }
270 
271  boundaryField_ == dt.value();
272 
273  readIfPresent();
274 }
275 
276 
277 template<class Type, template<class> class PatchField, class GeoMesh>
279 (
280  const IOobject& io,
281  const Internal& diField,
282  const PtrList<PatchField<Type>>& ptfl
283 )
284 :
285  Internal(io, diField),
286  timeIndex_(this->time().timeIndex()),
287  field0Ptr_(NULL),
288  fieldPrevIterPtr_(NULL),
289  boundaryField_(this->mesh().boundary(), *this, ptfl)
290 {
291  if (debug)
292  {
294  << "Constructing from components" << endl << this->info() << endl;
295  }
296 
297  readIfPresent();
298 }
299 
300 
301 template<class Type, template<class> class PatchField, class GeoMesh>
303 (
304  const IOobject& io,
305  const Mesh& mesh,
306  const dimensionSet& ds,
307  const Field<Type>& iField,
308  const PtrList<PatchField<Type>>& ptfl
309 )
310 :
311  Internal(io, mesh, ds, iField),
312  timeIndex_(this->time().timeIndex()),
313  field0Ptr_(NULL),
314  fieldPrevIterPtr_(NULL),
315  boundaryField_(mesh.boundary(), *this, ptfl)
316 {
317  if (debug)
318  {
320  << "Constructing from components" << endl << this->info() << endl;
321  }
322 
323  readIfPresent();
324 }
325 
326 
327 template<class Type, template<class> class PatchField, class GeoMesh>
329 (
330  const IOobject& io,
331  const Mesh& mesh
332 )
333 :
334  Internal(io, mesh, dimless, false),
335  timeIndex_(this->time().timeIndex()),
336  field0Ptr_(NULL),
337  fieldPrevIterPtr_(NULL),
338  boundaryField_(mesh.boundary())
339 {
340  readFields();
341 
342  // Check compatibility between field and mesh
343 
344  if (this->size() != GeoMesh::size(this->mesh()))
345  {
346  FatalIOErrorInFunction(this->readStream(typeName))
347  << " number of field elements = " << this->size()
348  << " number of mesh elements = " << GeoMesh::size(this->mesh())
349  << exit(FatalIOError);
350  }
351 
352  readOldTimeIfPresent();
353 
354  if (debug)
355  {
357  << "Finishing read-construction of" << endl << this->info() << endl;
358  }
359 }
360 
361 
362 template<class Type, template<class> class PatchField, class GeoMesh>
364 (
365  const IOobject& io,
366  const Mesh& mesh,
367  const dictionary& dict
368 )
369 :
370  Internal(io, mesh, dimless, false),
371  timeIndex_(this->time().timeIndex()),
372  field0Ptr_(NULL),
373  fieldPrevIterPtr_(NULL),
374  boundaryField_(mesh.boundary())
375 {
376  readFields(dict);
377 
378  // Check compatibility between field and mesh
379 
380  if (this->size() != GeoMesh::size(this->mesh()))
381  {
383  << " number of field elements = " << this->size()
384  << " number of mesh elements = " << GeoMesh::size(this->mesh())
385  << exit(FatalIOError);
386  }
387 
388  readOldTimeIfPresent();
389 
390  if (debug)
391  {
393  << "Finishing dictionary-construct of "
394  << endl << this->info() << endl;
395  }
396 }
397 
398 
399 template<class Type, template<class> class PatchField, class GeoMesh>
401 (
403 )
404 :
405  Internal(gf),
406  timeIndex_(gf.timeIndex()),
407  field0Ptr_(NULL),
408  fieldPrevIterPtr_(NULL),
409  boundaryField_(*this, gf.boundaryField_)
410 {
411  if (debug)
412  {
414  << "Constructing as copy" << endl << this->info() << endl;
415  }
416 
417  if (gf.field0Ptr_)
418  {
420  (
421  *gf.field0Ptr_
422  );
423  }
424 
425  this->writeOpt() = IOobject::NO_WRITE;
426 }
427 
428 
429 #ifndef NoConstructFromTmp
430 template<class Type, template<class> class PatchField, class GeoMesh>
432 (
434 )
435 :
436  Internal
437  (
438  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
439  tgf.isTmp()
440  ),
441  timeIndex_(tgf().timeIndex()),
442  field0Ptr_(NULL),
443  fieldPrevIterPtr_(NULL),
444  boundaryField_(*this, tgf().boundaryField_)
445 {
446  if (debug)
447  {
449  << "Constructing from tmp" << endl << this->info() << endl;
450  }
451 
452  this->writeOpt() = IOobject::NO_WRITE;
453 
454  tgf.clear();
455 }
456 #endif
457 
458 
459 template<class Type, template<class> class PatchField, class GeoMesh>
461 (
462  const IOobject& io,
464 )
465 :
466  Internal(io, gf),
467  timeIndex_(gf.timeIndex()),
468  field0Ptr_(NULL),
469  fieldPrevIterPtr_(NULL),
470  boundaryField_(*this, gf.boundaryField_)
471 {
472  if (debug)
473  {
475  << "Constructing as copy resetting IO params"
476  << endl << this->info() << endl;
477  }
478 
479  if (!readIfPresent() && gf.field0Ptr_)
480  {
482  (
483  io.name() + "_0",
484  *gf.field0Ptr_
485  );
486  }
487 }
488 
489 
490 #ifndef NoConstructFromTmp
491 template<class Type, template<class> class PatchField, class GeoMesh>
493 (
494  const IOobject& io,
496 )
497 :
498  Internal
499  (
500  io,
501  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
502  tgf.isTmp()
503  ),
504  timeIndex_(tgf().timeIndex()),
505  field0Ptr_(NULL),
506  fieldPrevIterPtr_(NULL),
507  boundaryField_(*this, tgf().boundaryField_)
508 {
509  if (debug)
510  {
512  << "Constructing from tmp resetting IO params"
513  << endl << this->info() << endl;
514  }
515 
516  tgf.clear();
517 
518  readIfPresent();
519 }
520 #endif
521 
522 
523 template<class Type, template<class> class PatchField, class GeoMesh>
525 (
526  const word& newName,
528 )
529 :
530  Internal(newName, gf),
531  timeIndex_(gf.timeIndex()),
532  field0Ptr_(NULL),
533  fieldPrevIterPtr_(NULL),
534  boundaryField_(*this, gf.boundaryField_)
535 {
536  if (debug)
537  {
539  << "Constructing as copy resetting name"
540  << endl << this->info() << endl;
541  }
542 
543  if (!readIfPresent() && gf.field0Ptr_)
544  {
546  (
547  newName + "_0",
548  *gf.field0Ptr_
549  );
550  }
551 }
552 
553 
554 #ifndef NoConstructFromTmp
555 template<class Type, template<class> class PatchField, class GeoMesh>
557 (
558  const word& newName,
560 )
561 :
562  Internal
563  (
564  newName,
565  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
566  tgf.isTmp()
567  ),
568  timeIndex_(tgf().timeIndex()),
569  field0Ptr_(NULL),
570  fieldPrevIterPtr_(NULL),
571  boundaryField_(*this, tgf().boundaryField_)
572 {
573  if (debug)
574  {
576  << "Constructing from tmp resetting name"
577  << endl << this->info() << endl;
578  }
579 
580  tgf.clear();
581 }
582 #endif
583 
584 
585 template<class Type, template<class> class PatchField, class GeoMesh>
587 (
588  const IOobject& io,
590  const word& patchFieldType
591 )
592 :
593  Internal(io, gf),
594  timeIndex_(gf.timeIndex()),
595  field0Ptr_(NULL),
596  fieldPrevIterPtr_(NULL),
597  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
598 {
599  if (debug)
600  {
602  << "Constructing as copy resetting IO params"
603  << endl << this->info() << endl;
604  }
605 
606  boundaryField_ == gf.boundaryField_;
607 
608  if (!readIfPresent() && gf.field0Ptr_)
609  {
611  (
612  io.name() + "_0",
613  *gf.field0Ptr_
614  );
615  }
616 }
617 
618 
619 template<class Type, template<class> class PatchField, class GeoMesh>
621 (
622  const IOobject& io,
624  const wordList& patchFieldTypes,
625  const wordList& actualPatchTypes
626 
627 )
628 :
629  Internal(io, gf),
630  timeIndex_(gf.timeIndex()),
631  field0Ptr_(NULL),
632  fieldPrevIterPtr_(NULL),
633  boundaryField_
634  (
635  this->mesh().boundary(),
636  *this,
637  patchFieldTypes,
638  actualPatchTypes
639  )
640 {
641  if (debug)
642  {
644  << "Constructing as copy resetting IO params and patch types"
645  << endl << this->info() << endl;
646  }
647 
648  boundaryField_ == gf.boundaryField_;
649 
650  if (!readIfPresent() && gf.field0Ptr_)
651  {
653  (
654  io.name() + "_0",
655  *gf.field0Ptr_
656  );
657  }
658 }
659 
660 
661 #ifndef NoConstructFromTmp
662 template<class Type, template<class> class PatchField, class GeoMesh>
664 (
665  const IOobject& io,
667  const wordList& patchFieldTypes,
668  const wordList& actualPatchTypes
669 )
670 :
671  Internal
672  (
673  io,
674  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
675  tgf.isTmp()
676  ),
677  timeIndex_(tgf().timeIndex()),
678  field0Ptr_(NULL),
679  fieldPrevIterPtr_(NULL),
680  boundaryField_
681  (
682  this->mesh().boundary(),
683  *this,
684  patchFieldTypes,
685  actualPatchTypes
686  )
687 {
688  if (debug)
689  {
691  << "Constructing from tmp resetting IO params and patch types"
692  << endl << this->info() << endl;
693  }
694 
695  boundaryField_ == tgf().boundaryField_;
696 
697  tgf.clear();
698 }
699 #endif
700 
701 
702 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
703 
704 template<class Type, template<class> class PatchField, class GeoMesh>
706 {
707  deleteDemandDrivenData(field0Ptr_);
708  deleteDemandDrivenData(fieldPrevIterPtr_);
709 }
710 
711 
712 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
713 
714 template<class Type, template<class> class PatchField, class GeoMesh>
715 typename
718 {
719  this->setUpToDate();
720  storeOldTimes();
721  return *this;
722 }
723 
724 
725 template<class Type, template<class> class PatchField, class GeoMesh>
726 typename
729 {
730  this->setUpToDate();
731  storeOldTimes();
732  return *this;
733 }
734 
735 
736 template<class Type, template<class> class PatchField, class GeoMesh>
737 typename
740 {
741  this->setUpToDate();
742  storeOldTimes();
743  return boundaryField_;
744 }
745 
746 
747 template<class Type, template<class> class PatchField, class GeoMesh>
749 {
750  if
751  (
752  field0Ptr_
753  && timeIndex_ != this->time().timeIndex()
754  && !(
755  this->name().size() > 2
756  && this->name()(this->name().size()-2, 2) == "_0"
757  )
758  )
759  {
760  storeOldTime();
761  }
762 
763  // Correct time index
764  timeIndex_ = this->time().timeIndex();
765 }
766 
767 
768 template<class Type, template<class> class PatchField, class GeoMesh>
770 {
771  if (field0Ptr_)
772  {
773  field0Ptr_->storeOldTime();
774 
775  if (debug)
776  {
778  << "Storing old time field for field" << endl
779  << this->info() << endl;
780  }
781 
782  *field0Ptr_ == *this;
783  field0Ptr_->timeIndex_ = timeIndex_;
784 
785  if (field0Ptr_->field0Ptr_)
786  {
787  field0Ptr_->writeOpt() = this->writeOpt();
788  }
789  }
790 }
791 
792 
793 template<class Type, template<class> class PatchField, class GeoMesh>
795 {
796  if (field0Ptr_)
797  {
798  return field0Ptr_->nOldTimes() + 1;
799  }
800  else
801  {
802  return 0;
803  }
804 }
805 
806 
807 template<class Type, template<class> class PatchField, class GeoMesh>
810 {
811  if (!field0Ptr_)
812  {
814  (
815  IOobject
816  (
817  this->name() + "_0",
818  this->time().timeName(),
819  this->db(),
820  IOobject::NO_READ,
821  IOobject::NO_WRITE,
822  this->registerObject()
823  ),
824  *this
825  );
826  }
827  else
828  {
829  storeOldTimes();
830  }
831 
832  return *field0Ptr_;
833 }
834 
835 
836 template<class Type, template<class> class PatchField, class GeoMesh>
839 {
840  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
841  .oldTime();
842 
843  return *field0Ptr_;
844 }
845 
846 
847 template<class Type, template<class> class PatchField, class GeoMesh>
849 {
850  if (!fieldPrevIterPtr_)
851  {
852  if (debug)
853  {
855  << "Allocating previous iteration field" << endl
856  << this->info() << endl;
857  }
858 
859  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
860  (
861  this->name() + "PrevIter",
862  *this
863  );
864  }
865  else
866  {
867  *fieldPrevIterPtr_ == *this;
868  }
869 }
870 
871 
872 template<class Type, template<class> class PatchField, class GeoMesh>
875 {
876  if (!fieldPrevIterPtr_)
877  {
879  << "previous iteration field" << endl << this->info() << endl
880  << " not stored."
881  << " Use field.storePrevIter() at start of iteration."
882  << abort(FatalError);
883  }
884 
885  return *fieldPrevIterPtr_;
886 }
887 
888 
889 template<class Type, template<class> class PatchField, class GeoMesh>
892 {
893  this->setUpToDate();
894  storeOldTimes();
895  boundaryField_.evaluate();
896 }
897 
898 
899 template<class Type, template<class> class PatchField, class GeoMesh>
901 {
902  // Search all boundary conditions, if any are
903  // fixed-value or mixed (Robin) do not set reference level for solution.
904 
905  bool needRef = true;
906 
907  forAll(boundaryField_, patchi)
908  {
909  if (boundaryField_[patchi].fixesValue())
910  {
911  needRef = false;
912  break;
913  }
914  }
915 
916  reduce(needRef, andOp<bool>());
917 
918  return needRef;
919 }
920 
921 
922 template<class Type, template<class> class PatchField, class GeoMesh>
924 {
925  if (debug)
926  {
928  << "Relaxing" << endl << this->info() << " by " << alpha << endl;
929  }
930 
931  operator==(prevIter() + alpha*(*this - prevIter()));
932 }
933 
934 
935 template<class Type, template<class> class PatchField, class GeoMesh>
937 {
938  word name = this->name();
939 
940  if
941  (
942  this->mesh().data::template lookupOrDefault<bool>
943  (
944  "finalIteration",
945  false
946  )
947  )
948  {
949  name += "Final";
950  }
951 
952  if (this->mesh().relaxField(name))
953  {
954  relax(this->mesh().fieldRelaxationFactor(name));
955  }
956 }
957 
958 
959 template<class Type, template<class> class PatchField, class GeoMesh>
961 (
962  bool final
963 ) const
964 {
965  if (final)
966  {
967  return this->name() + "Final";
968  }
969  else
970  {
971  return this->name();
972  }
973 }
974 
975 
976 template<class Type, template<class> class PatchField, class GeoMesh>
978 (
979  Ostream& os
980 ) const
981 {
982  os << "min/max(" << this->name() << ") = "
983  << Foam::min(*this).value() << ", "
984  << Foam::max(*this).value()
985  << endl;
986 }
987 
988 
989 template<class Type, template<class> class PatchField, class GeoMesh>
991 writeData(Ostream& os) const
992 {
993  os << *this;
994  return os.good();
995 }
996 
997 
998 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
999 
1000 template<class Type, template<class> class PatchField, class GeoMesh>
1003 {
1005  (
1007  (
1008  IOobject
1009  (
1010  this->name() + ".T()",
1011  this->instance(),
1012  this->db()
1013  ),
1014  this->mesh(),
1015  this->dimensions()
1016  )
1017  );
1018 
1019  Foam::T(result.ref().primitiveFieldRef(), primitiveField());
1020  Foam::T(result.ref().boundaryFieldRef(), boundaryField());
1021 
1022  return result;
1023 }
1024 
1025 
1026 template<class Type, template<class> class PatchField, class GeoMesh>
1027 Foam::tmp
1028 <
1030  <
1031  typename Foam::GeometricField<Type, PatchField, GeoMesh>::cmptType,
1032  PatchField,
1033  GeoMesh
1034  >
1035 >
1037 (
1038  const direction d
1039 ) const
1040 {
1042  (
1044  (
1045  IOobject
1046  (
1047  this->name() + ".component(" + Foam::name(d) + ')',
1048  this->instance(),
1049  this->db()
1050  ),
1051  this->mesh(),
1052  this->dimensions()
1053  )
1054  );
1055 
1056  Foam::component(Component.ref().primitiveFieldRef(), primitiveField(), d);
1057  Foam::component(Component.ref().boundaryFieldRef(), boundaryField(), d);
1058 
1059  return Component;
1060 }
1061 
1062 
1063 template<class Type, template<class> class PatchField, class GeoMesh>
1065 (
1066  const direction d,
1067  const GeometricField
1068  <
1070  PatchField,
1071  GeoMesh
1072  >& gcf
1073 )
1074 {
1075  primitiveFieldRef().replace(d, gcf.primitiveField());
1076  boundaryFieldRef().replace(d, gcf.boundaryField());
1077 }
1078 
1079 
1080 template<class Type, template<class> class PatchField, class GeoMesh>
1083  const direction d,
1084  const dimensioned<cmptType>& ds
1085 )
1086 {
1087  primitiveFieldRef().replace(d, ds.value());
1088  boundaryFieldRef().replace(d, ds.value());
1089 }
1090 
1091 
1092 template<class Type, template<class> class PatchField, class GeoMesh>
1095  const dimensioned<Type>& dt
1096 )
1097 {
1098  Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1099  Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1100 }
1101 
1102 
1103 template<class Type, template<class> class PatchField, class GeoMesh>
1106  const dimensioned<Type>& dt
1107 )
1108 {
1109  Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1110  Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1111 }
1112 
1113 
1114 template<class Type, template<class> class PatchField, class GeoMesh>
1116 {
1117  primitiveFieldRef().negate();
1118  boundaryFieldRef().negate();
1119 }
1120 
1121 
1122 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1123 
1124 template<class Type, template<class> class PatchField, class GeoMesh>
1125 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1128 )
1129 {
1130  if (this == &gf)
1131  {
1133  << "attempted assignment to self"
1134  << abort(FatalError);
1135  }
1136 
1137  checkField(*this, gf, "=");
1138 
1139  // Only assign field contents not ID
1140 
1141  ref() = gf();
1142  boundaryFieldRef() = gf.boundaryField();
1143 }
1144 
1145 
1146 template<class Type, template<class> class PatchField, class GeoMesh>
1147 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1150 )
1151 {
1152  if (this == &(tgf()))
1153  {
1155  << "attempted assignment to self"
1156  << abort(FatalError);
1157  }
1158 
1160 
1161  checkField(*this, gf, "=");
1162 
1163  // Only assign field contents not ID
1164 
1165  this->dimensions() = gf.dimensions();
1166 
1167  // Transfer the storage from the tmp
1168  primitiveFieldRef().transfer
1169  (
1170  const_cast<Field<Type>&>(gf.primitiveField())
1171  );
1172 
1173  boundaryFieldRef() = gf.boundaryField();
1174 
1175  tgf.clear();
1176 }
1177 
1178 
1179 template<class Type, template<class> class PatchField, class GeoMesh>
1180 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator=
1182  const dimensioned<Type>& dt
1183 )
1184 {
1185  ref() = dt;
1186  boundaryFieldRef() = dt.value();
1187 }
1188 
1189 
1190 template<class Type, template<class> class PatchField, class GeoMesh>
1191 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1194 )
1195 {
1197 
1198  checkField(*this, gf, "==");
1199 
1200  // Only assign field contents not ID
1201 
1202  ref() = gf();
1203  boundaryFieldRef() == gf.boundaryField();
1204 
1205  tgf.clear();
1206 }
1207 
1208 
1209 template<class Type, template<class> class PatchField, class GeoMesh>
1210 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator==
1212  const dimensioned<Type>& dt
1213 )
1214 {
1215  ref() = dt;
1216  boundaryFieldRef() == dt.value();
1217 }
1218 
1219 
1220 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1221  \
1222 template<class Type, template<class> class PatchField, class GeoMesh> \
1223 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1224 ( \
1225  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1226 ) \
1227 { \
1228  checkField(*this, gf, #op); \
1229  \
1230  ref() op gf(); \
1231  boundaryFieldRef() op gf.boundaryField(); \
1232 } \
1233  \
1234 template<class Type, template<class> class PatchField, class GeoMesh> \
1235 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1236 ( \
1237  const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1238 ) \
1239 { \
1240  operator op(tgf()); \
1241  tgf.clear(); \
1242 } \
1243  \
1244 template<class Type, template<class> class PatchField, class GeoMesh> \
1245 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1246 ( \
1247  const dimensioned<TYPE>& dt \
1248 ) \
1249 { \
1250  ref() op dt; \
1251  boundaryFieldRef() op dt.value(); \
1252 }
1253 
1258 
1259 #undef COMPUTED_ASSIGNMENT
1260 
1261 
1262 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1263 
1264 template<class Type, template<class> class PatchField, class GeoMesh>
1265 Foam::Ostream& Foam::operator<<
1266 (
1267  Ostream& os,
1269 )
1270 {
1271  gf().writeData(os, "internalField");
1272  os << nl;
1273  gf.boundaryField().writeEntry("boundaryField", os);
1274 
1275  // Check state of IOstream
1276  os.check
1277  (
1278  "Ostream& operator<<(Ostream&, "
1279  "const GeometricField<Type, PatchField, GeoMesh>&)"
1280  );
1281 
1282  return (os);
1283 }
1284 
1285 
1286 template<class Type, template<class> class PatchField, class GeoMesh>
1287 Foam::Ostream& Foam::operator<<
1288 (
1289  Ostream& os,
1291 )
1292 {
1293  os << tgf();
1294  tgf.clear();
1295  return os;
1296 }
1297 
1298 
1299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1300 
1301 #undef checkField
1302 
1303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1304 
1305 #include "GeometricBoundaryField.C"
1306 #include "GeometricFieldFunctions.C"
1307 
1308 // ************************************************************************* //
bool needReference() const
Does the field need a reference level for solution.
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
#define COMPUTED_ASSIGNMENT(TYPE, op)
dictionary dict
void storeOldTime() const
Store the old-time field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
#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
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
UEqn relax()
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Field< Type >::cmptType cmptType
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Type & value() const
Return const reference to value.
Generic GeometricField class.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Generic dimensioned Type class.
conserve primitiveFieldRef()+
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
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
Dimension set for the base types.
Definition: dimensionSet.H:118
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Pre-declare SubField and related Field type.
Definition: Field.H:57
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
label timeIndex() const
Return the time index of the field.
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
word timeName
Definition: getTimeIndex.H:3
void min(const dimensioned< Type > &)
faceListList boundary(nPatches)
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet & dimensions() const
Return dimensions.
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
label nOldTimes() const
Return the number of old time fields stored.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:262
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
Internal & ref()
Return a reference to the dimensioned internal field.
Template functions to aid in the implementation of demand driven data.
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:62
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...
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
void correctBoundaryConditions()
Correct boundary field.
rDeltaT ref()
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void deleteDemandDrivenData(DataPtr &dataPtr)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
void storePrevIter() const
Store the field as the previous iteration value.
void relax()
Relax field (for steady-state solution).
const word & name() const
Return name.
Definition: IOobject.H:260
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.
void storeOldTimes() const
Store the old-time fields.