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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 operation " << 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 {
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 
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  Info<< "Reading old time level for field"
150  << endl << this->info() << endl;
151  }
152 
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  DimensionedField<Type, GeoMesh>(io, mesh, ds, false),
185  timeIndex_(this->time().timeIndex()),
186  field0Ptr_(NULL),
187  fieldPrevIterPtr_(NULL),
188  boundaryField_(mesh.boundary(), *this, patchFieldType)
189 {
190 // Info << "construct From " <<endl;
191 // Info << "foam::GeometricField<Type, PatchField, GeoMesh>::GeometricField"<< endl;
192 // Info << "("<< endl;
193 // Info << " const IOobject& io,"<< endl;
194 // Info << " const Mesh& mesh,"<< endl;
195 // Info << " const dimensionSet& ds,"<< endl;
196 // Info << " const word& patchFieldType"<< endl;
197 // Info << ")"<< endl;
198 
199  if (debug)
200  {
201  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
202  "creating temporary"
203  << endl << this->info() << endl;
204 }
205 
206  readIfPresent();
207 }
208 
209 
210 template<class Type, template<class> class PatchField, class GeoMesh>
212 (
213  const IOobject& io,
214  const Mesh& mesh,
215  const dimensionSet& ds,
216  const wordList& patchFieldTypes,
217  const wordList& actualPatchTypes
218 )
219 :
220  DimensionedField<Type, GeoMesh>(io, mesh, ds, false),
221  timeIndex_(this->time().timeIndex()),
222  field0Ptr_(NULL),
223  fieldPrevIterPtr_(NULL),
224  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
225 {
226 // Info << "construct From " <<endl;
227 // Info << "foam::GeometricField<Type, PatchField, GeoMesh>::GeometricField"<< endl;
228 // Info << "("<< endl;
229 // Info << " const IOobject& io,"<< endl;
230 // Info << " const Mesh& mesh,"<< endl;
231 // Info << " const dimensionSet& ds,"<< endl;
232 // Info << " const word& patchFieldType"<< endl;
233 // Info << " const wordList& actualPatchTypes"<< endl;
234 // Info << ")"<< endl;
235 //
236  if (debug)
237  {
238  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
239  "creating temporary"
240  << endl << this->info() << endl;
241  }
242 
243  readIfPresent();
244 }
245 
246 
247 template<class Type, template<class> class PatchField, class GeoMesh>
249 (
250  const IOobject& io,
251  const Mesh& mesh,
252  const dimensioned<Type>& dt,
253  const word& patchFieldType
254 )
255 :
256  DimensionedField<Type, GeoMesh>(io, mesh, dt, false),
257  timeIndex_(this->time().timeIndex()),
258  field0Ptr_(NULL),
259  fieldPrevIterPtr_(NULL),
260  boundaryField_(mesh.boundary(), *this, patchFieldType)
261 {
262 // Info << "construct From 3" <<endl;
263 // Info << "foam::GeometricField<Type, PatchField, GeoMesh>::GeometricField"<< endl;
264 // Info << "("<< endl;
265 // Info << " const IOobject& io,"<< endl;
266 // Info << " const Mesh& mesh,"<< endl;
267 // Info << " const dimensionSet& ds,"<< endl;
268 // Info << " const word& patchFieldType"<< endl;
269 // Info << ")"<< endl;
270 //
271  if (debug)
272  {
273  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
274  "creating temporary"
275  << endl << this->info() << endl;
276  }
277 
278  boundaryField_ == dt.value();
279 
280  readIfPresent();
281 }
282 
283 
284 template<class Type, template<class> class PatchField, class GeoMesh>
286 (
287  const IOobject& io,
288  const Mesh& mesh,
289  const dimensioned<Type>& dt,
290  const wordList& patchFieldTypes,
291  const wordList& actualPatchTypes
292 )
293 :
294  DimensionedField<Type, GeoMesh>(io, mesh, dt, false),
295  timeIndex_(this->time().timeIndex()),
296  field0Ptr_(NULL),
297  fieldPrevIterPtr_(NULL),
298  boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
299 {
300  // Info << "construct From 4" <<endl;
301  if (debug)
302  {
303  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
304  "creating temporary"
305  << endl << this->info() << endl;
306  }
307 
308  boundaryField_ == dt.value();
309 
310  readIfPresent();
311 }
312 
313 
314 template<class Type, template<class> class PatchField, class GeoMesh>
316 (
317  const IOobject& io,
318  const Mesh& mesh,
319  const dimensionSet& ds,
320  const Field<Type>& iField,
321  const PtrList<PatchField<Type> >& ptfl
322 )
323 :
324  DimensionedField<Type, GeoMesh>(io, mesh, ds, iField),
325  timeIndex_(this->time().timeIndex()),
326  field0Ptr_(NULL),
327  fieldPrevIterPtr_(NULL),
328  boundaryField_(mesh.boundary(), *this, ptfl)
329 {
330  // Info << "construct From 5" <<endl;
331  // if (debug)
332  // {
333  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
334  "constructing from components"
335  << 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  const bool readOldTime
348 )
349 :
350  DimensionedField<Type, GeoMesh>(io, mesh, dimless, false),
351  timeIndex_(this->time().timeIndex()),
352  field0Ptr_(NULL),
353  fieldPrevIterPtr_(NULL),
354  boundaryField_(mesh.boundary())
355 {
356 // Info << "construct From 5.5" <<endl;
357 // Info << "foam::GeometricField<Type, PatchField, GeoMesh>::GeometricField"<< endl;
358 // Info << "("<< endl;
359 // Info << " const IOobject& io,"<< endl;
360 // Info << " const Mesh& mesh,"<< endl;
361 // Info << " const boold readOldTime"<< endl;
362 // Info << ")"<< endl;
363 
364  readFields();
365 
366  // Check compatibility between field and mesh
367 
368  if (this->size() != GeoMesh::size(this->mesh()))
369  {
370  FatalIOErrorInFunction(this->readStream(typeName))
371  << " number of field elements = " << this->size()
372  << " number of mesh elements = " << GeoMesh::size(this->mesh())
373  << exit(FatalIOError);
374  }
375 
376  if (readOldTime)
377  {
378  readOldTimeIfPresent();
379  }
380 
381  if (debug)
382  {
383  Info<< "Finishing read-construct of "
384  "GeometricField<Type, PatchField, GeoMesh>"
385  << endl << this->info() << endl;
386  }
387 }
388 
389 
390 template<class Type, template<class> class PatchField, class GeoMesh>
392 (
393  const IOobject& io,
394  const Mesh& mesh,
395  const dictionary& dict
396 )
397 :
398  DimensionedField<Type, GeoMesh>(io, mesh, dimless, false),
399  timeIndex_(this->time().timeIndex()),
400  field0Ptr_(NULL),
401  fieldPrevIterPtr_(NULL),
402  boundaryField_(mesh.boundary())
403 {
404 // Info << "construct From 6" <<endl;
405  readFields(dict);
406 
407  // Check compatibility between field and mesh
408 
409  if (this->size() != GeoMesh::size(this->mesh()))
410  {
412  << " number of field elements = " << this->size()
413  << " number of mesh elements = " << GeoMesh::size(this->mesh())
414  << exit(FatalIOError);
415  }
416 
417  if (debug)
418  {
419  Info<< "Finishing dictionary-construct of "
420  "GeometricField<Type, PatchField, GeoMesh>"
421  << endl << this->info() << endl;
422  }
423 }
424 
425 
426 template<class Type, template<class> class PatchField, class GeoMesh>
428 (
429  const GeometricField<Type, PatchField, GeoMesh>& gf
430 )
431 :
432  DimensionedField<Type, GeoMesh>(gf),
433  timeIndex_(gf.timeIndex()),
434  field0Ptr_(NULL),
435  fieldPrevIterPtr_(NULL),
436  boundaryField_(*this, gf.boundaryField_)
437 {
438 // Info << "construct From 7" <<endl;
439  if (debug)
440  {
441  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
442  "constructing as copy"
443  << endl << this->info() << endl;
444  }
445 
446  if (gf.field0Ptr_)
447  {
448  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
449  (
450  *gf.field0Ptr_
451  );
452  }
453 
454  this->writeOpt() = IOobject::NO_WRITE;
455 }
456 
457 
458 #ifndef NoConstructFromTmp
459 template<class Type, template<class> class PatchField, class GeoMesh>
461 (
462  const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
463 )
464 :
465  DimensionedField<Type, GeoMesh>
466  (
467  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
468  tgf.isTmp()
469  ),
470  timeIndex_(tgf().timeIndex()),
471  field0Ptr_(NULL),
472  fieldPrevIterPtr_(NULL),
473  boundaryField_(*this, tgf().boundaryField_)
474 {
475 // Info << "construct From 8" <<endl;
476  if (debug)
477  {
478  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
479  "constructing as copy"
480  << endl << this->info() << endl;
481  }
482 
483  this->writeOpt() = IOobject::NO_WRITE;
484 
485  tgf.clear();
486 }
487 #endif
488 
489 
490 template<class Type, template<class> class PatchField, class GeoMesh>
492 (
493  const IOobject& io,
494  const GeometricField<Type, PatchField, GeoMesh>& gf
495 )
496 :
497  DimensionedField<Type, GeoMesh>(io, gf),
498  timeIndex_(gf.timeIndex()),
499  field0Ptr_(NULL),
500  fieldPrevIterPtr_(NULL),
501  boundaryField_(*this, gf.boundaryField_)
502 {
503 // Info << "construct From 9" <<endl;
504  if (debug)
505  {
506  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
507  "constructing as copy resetting IO params"
508  << endl << this->info() << endl;
509  }
510 
511  if (!readIfPresent() && gf.field0Ptr_)
512  {
513  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
514  (
515  io.name() + "_0",
516  *gf.field0Ptr_
517  );
518  }
519 }
520 
521 
522 #ifndef NoConstructFromTmp
523 template<class Type, template<class> class PatchField, class GeoMesh>
525 (
526  const IOobject& io,
527  const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
528 )
529 :
530  DimensionedField<Type, GeoMesh>
531  (
532  io,
533  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
534  tgf.isTmp()
535  ),
536  timeIndex_(tgf().timeIndex()),
537  field0Ptr_(NULL),
538  fieldPrevIterPtr_(NULL),
539  boundaryField_(*this, tgf().boundaryField_)
540 {
541 // Info << "construct From 10" <<endl;
542  if (debug)
543  {
544  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
545  "constructing from tmp resetting IO params"
546  << endl << this->info() << endl;
547  }
548 
549  tgf.clear();
550 
551  readIfPresent();
552 }
553 #endif
554 
555 
556 template<class Type, template<class> class PatchField, class GeoMesh>
558 (
559  const word& newName,
560  const GeometricField<Type, PatchField, GeoMesh>& gf
561 )
562 :
563  DimensionedField<Type, GeoMesh>(newName, gf),
564  timeIndex_(gf.timeIndex()),
565  field0Ptr_(NULL),
566  fieldPrevIterPtr_(NULL),
567  boundaryField_(*this, gf.boundaryField_)
568 {
569 // Info << "construct From 11" <<endl;
570  if (debug)
571  {
572  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
573  "constructing as copy resetting name"
574  << endl << this->info() << endl;
575  }
576 
577  if (!readIfPresent() && gf.field0Ptr_)
578  {
579  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
580  (
581  newName + "_0",
582  *gf.field0Ptr_
583  );
584  }
585 }
586 
587 
588 #ifndef NoConstructFromTmp
589 template<class Type, template<class> class PatchField, class GeoMesh>
591 (
592  const word& newName,
593  const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
594 )
595 :
596  DimensionedField<Type, GeoMesh>
597  (
598  newName,
599  const_cast<GeometricField<Type, PatchField, GeoMesh>&>(tgf()),
600  tgf.isTmp()
601  ),
602  timeIndex_(tgf().timeIndex()),
603  field0Ptr_(NULL),
604  fieldPrevIterPtr_(NULL),
605  boundaryField_(*this, tgf().boundaryField_)
606 {
607 // Info << "construct From 12" <<endl;
608  if (debug)
609  {
610  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
611  "constructing from tmp resetting name"
612  << endl << this->info() << endl;
613  }
614 
615  tgf.clear();
616 }
617 #endif
618 
619 
620 template<class Type, template<class> class PatchField, class GeoMesh>
622 (
623  const IOobject& io,
624  const GeometricField<Type, PatchField, GeoMesh>& gf,
625  const word& patchFieldType
626 )
627 :
628  DimensionedField<Type, GeoMesh>(io, gf),
629  timeIndex_(gf.timeIndex()),
630  field0Ptr_(NULL),
631  fieldPrevIterPtr_(NULL),
632  boundaryField_(this->mesh().boundary(), *this, patchFieldType)
633 {
634 // Info << "construct From 13" <<endl;
635  if (debug)
636  {
637  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
638  "constructing as copy resetting IO params and patch type"
639  << endl << this->info() << endl;
640  }
641 
642  boundaryField_ == gf.boundaryField_;
643 
644  if (!readIfPresent() && gf.field0Ptr_)
645  {
646  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
647  (
648  io.name() + "_0",
649  *gf.field0Ptr_
650  );
651  }
652 }
653 
654 
655 template<class Type, template<class> class PatchField, class GeoMesh>
657 (
658  const IOobject& io,
659  const GeometricField<Type, PatchField, GeoMesh>& gf,
660  const wordList& patchFieldTypes,
661  const wordList& actualPatchTypes
662 
663 )
664 :
665  DimensionedField<Type, GeoMesh>(io, gf),
666  timeIndex_(gf.timeIndex()),
667  field0Ptr_(NULL),
668  fieldPrevIterPtr_(NULL),
669  boundaryField_
670  (
671  this->mesh().boundary(),
672  *this,
673  patchFieldTypes,
674  actualPatchTypes
675  )
676 {
677 // Info << "construct From 14" <<endl;
678  if (debug)
679  {
680  Info<< "GeometricField<Type, PatchField, GeoMesh>::GeometricField : "
681  "constructing as copy resetting IO params and patch types"
682  << endl << this->info() << endl;
683  }
684 
685  boundaryField_ == gf.boundaryField_;
686 
687  if (!readIfPresent() && gf.field0Ptr_)
688  {
689  field0Ptr_ = new GeometricField<Type, PatchField, GeoMesh>
690  (
691  io.name() + "_0",
692  *gf.field0Ptr_
693  );
694  }
695 }
696 
697 
698 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
699 
700 template<class Type, template<class> class PatchField, class GeoMesh>
702 {
703  deleteDemandDrivenData(field0Ptr_);
704  deleteDemandDrivenData(fieldPrevIterPtr_);
705 }
706 
707 
708 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
709 
710 template<class Type, template<class> class PatchField, class GeoMesh>
711 typename
714 {
715  this->setUpToDate();
716  storeOldTimes();
717  return *this;
718 }
719 
720 
721 template<class Type, template<class> class PatchField, class GeoMesh>
722 typename
725 {
726  this->setUpToDate();
727  storeOldTimes();
728  return *this;
729 }
730 
731 
732 template<class Type, template<class> class PatchField, class GeoMesh>
733 typename
736 {
737  this->setUpToDate();
738  storeOldTimes();
739  return boundaryField_;
740 }
741 
742 
743 template<class Type, template<class> class PatchField, class GeoMesh>
745 {
746  if
747  (
748  field0Ptr_
749  && timeIndex_ != this->time().timeIndex()
750  && !(
751  this->name().size() > 2
752  && this->name()(this->name().size()-2, 2) == "_0"
753  )
754  )
755  {
756  storeOldTime();
757  }
758 
759  // Correct time index
760  timeIndex_ = this->time().timeIndex();
761 }
762 
763 
764 template<class Type, template<class> class PatchField, class GeoMesh>
766 {
767  if (field0Ptr_)
768  {
769  field0Ptr_->storeOldTime();
770 
771  if (debug)
772  {
773  Info<< "Storing old time field for field" << endl
774  << this->info() << endl;
775  }
776 
777  *field0Ptr_ == *this;
778  field0Ptr_->timeIndex_ = timeIndex_;
779 
780  if (field0Ptr_->field0Ptr_)
781  {
782  field0Ptr_->writeOpt() = this->writeOpt();
783  }
784  }
785 }
786 
787 
788 template<class Type, template<class> class PatchField, class GeoMesh>
790 {
791  if (field0Ptr_)
792  {
793  return field0Ptr_->nOldTimes() + 1;
794  }
795  else
796  {
797  return 0;
798  }
799 }
800 
801 
802 template<class Type, template<class> class PatchField, class GeoMesh>
805 {
806  if (!field0Ptr_)
807  {
809  (
810  IOobject
811  (
812  this->name() + "_0",
813  this->time().timeName(),
814  this->db(),
815  IOobject::NO_READ,
816  IOobject::NO_WRITE,
817  this->registerObject()
818  ),
819  *this
820  );
821  }
822  else
823  {
824  storeOldTimes();
825  }
826 
827  return *field0Ptr_;
828 }
829 
830 
831 template<class Type, template<class> class PatchField, class GeoMesh>
834 {
835  static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
836  .oldTime();
837 
838  return *field0Ptr_;
839 }
840 
841 
842 template<class Type, template<class> class PatchField, class GeoMesh>
844 {
845  if (!fieldPrevIterPtr_)
846  {
847  if (debug)
848  {
849  Info<< "Allocating previous iteration field" << endl
850  << this->info() << endl;
851  }
852 
853  fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
854  (
855  this->name() + "PrevIter",
856  *this
857  );
858  }
859  else
860  {
861  *fieldPrevIterPtr_ == *this;
862  }
863 }
864 
865 
866 template<class Type, template<class> class PatchField, class GeoMesh>
869 {
870  if (!fieldPrevIterPtr_)
871  {
873  << "previous iteration field" << endl << this->info() << endl
874  << " not stored."
875  << " Use field.storePrevIter() at start of iteration."
876  << abort(FatalError);
877  }
878 
879  return *fieldPrevIterPtr_;
880 }
881 
882 
883 template<class Type, template<class> class PatchField, class GeoMesh>
886 {
887  this->setUpToDate();
888  storeOldTimes();
889  boundaryField_.evaluate();
890 }
891 
892 
893 template<class Type, template<class> class PatchField, class GeoMesh>
895 {
896  // Search all boundary conditions, if any are
897  // fixed-value or mixed (Robin) do not set reference level for solution.
898 
899  bool needRef = true;
900 
901  forAll(boundaryField_, patchi)
902  {
903  if (boundaryField_[patchi].fixesValue())
904  {
905  needRef = false;
906  break;
907  }
908  }
909 
910  reduce(needRef, andOp<bool>());
911 
912  return needRef;
913 }
914 
915 
916 template<class Type, template<class> class PatchField, class GeoMesh>
918 {
919  if (debug)
920  {
922  << "Relaxing" << endl << this->info() << " by " << alpha << endl;
923  }
924 
925  operator==(prevIter() + alpha*(*this - prevIter()));
926 }
927 
928 
929 template<class Type, template<class> class PatchField, class GeoMesh>
931 {
932  word name = this->name();
933 
934  if
935  (
936  this->mesh().data::template lookupOrDefault<bool>
937  (
938  "finalIteration",
939  false
940  )
941  )
942  {
943  name += "Final";
944  }
945 
946  if (this->mesh().relaxField(name))
947  {
948  relax(this->mesh().fieldRelaxationFactor(name));
949  }
950 }
951 
952 
953 template<class Type, template<class> class PatchField, class GeoMesh>
955 (
956  bool final
957 ) const
958 {
959  if (final)
960  {
961  return this->name() + "Final";
962  }
963  else
964  {
965  return this->name();
966  }
967 }
968 
969 
970 template<class Type, template<class> class PatchField, class GeoMesh>
972 (
973  Ostream& os
974 ) const
975 {
976  os << "min/max(" << this->name() << ") = "
977  << Foam::min(*this).value() << ", "
978  << Foam::max(*this).value()
979  << endl;
980 }
981 
982 
983 template<class Type, template<class> class PatchField, class GeoMesh>
985 writeData(Ostream& os) const
986 {
987  os << *this;
988  return os.good();
989 }
990 
991 
992 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
993 
994 template<class Type, template<class> class PatchField, class GeoMesh>
997 {
999  (
1001  (
1002  IOobject
1003  (
1004  this->name() + ".T()",
1005  this->instance(),
1006  this->db()
1007  ),
1008  this->mesh(),
1009  this->dimensions()
1010  )
1011  );
1012 
1013  Foam::T(result().internalField(), internalField());
1014  Foam::T(result().boundaryField(), boundaryField());
1015 
1016  return result;
1017 }
1018 
1019 
1020 template<class Type, template<class> class PatchField, class GeoMesh>
1021 Foam::tmp
1022 <
1024  <
1026  PatchField,
1027  GeoMesh
1028  >
1029 >
1031 (
1032  const direction d
1033 ) const
1034 {
1036  (
1038  (
1039  IOobject
1040  (
1041  this->name() + ".component(" + Foam::name(d) + ')',
1042  this->instance(),
1043  this->db()
1044  ),
1045  this->mesh(),
1046  this->dimensions()
1047  )
1048  );
1049 
1050  Foam::component(Component().internalField(), internalField(), d);
1051  Foam::component(Component().boundaryField(), boundaryField(), d);
1052 
1053  return Component;
1054 }
1055 
1056 
1057 template<class Type, template<class> class PatchField, class GeoMesh>
1059 (
1060  const direction d,
1061  const GeometricField
1062  <
1063  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1064  PatchField,
1065  GeoMesh
1066  >& gcf
1067 )
1068 {
1069  internalField().replace(d, gcf.internalField());
1070  boundaryField().replace(d, gcf.boundaryField());
1071 }
1072 
1073 
1074 template<class Type, template<class> class PatchField, class GeoMesh>
1076 (
1077  const direction d,
1078  const dimensioned<cmptType>& ds
1079 )
1080 {
1081  internalField().replace(d, ds.value());
1082  boundaryField().replace(d, ds.value());
1083 }
1084 
1085 
1086 template<class Type, template<class> class PatchField, class GeoMesh>
1088 (
1089  const dimensioned<Type>& dt
1090 )
1091 {
1092  Foam::max(internalField(), internalField(), dt.value());
1093  Foam::max(boundaryField(), boundaryField(), dt.value());
1094 }
1095 
1096 
1097 template<class Type, template<class> class PatchField, class GeoMesh>
1099 (
1100  const dimensioned<Type>& dt
1101 )
1102 {
1103  Foam::min(internalField(), internalField(), dt.value());
1104  Foam::min(boundaryField(), boundaryField(), dt.value());
1105 }
1106 
1107 
1108 template<class Type, template<class> class PatchField, class GeoMesh>
1110 {
1111  internalField().negate();
1112  boundaryField().negate();
1113 }
1114 
1115 
1116 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1117 
1118 template<class Type, template<class> class PatchField, class GeoMesh>
1120 (
1122 )
1123 {
1124  if (this == &gf)
1125  {
1127  << "attempted assignment to self"
1128  << abort(FatalError);
1129  }
1130 
1131  checkField(*this, gf, "=");
1132 
1133  // only equate field contents not ID
1134 
1135  dimensionedInternalField() = gf.dimensionedInternalField();
1136  boundaryField() = gf.boundaryField();
1137 }
1138 
1139 
1140 template<class Type, template<class> class PatchField, class GeoMesh>
1142 (
1143  const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
1144 )
1145 {
1146  if (this == &(tgf()))
1147  {
1149  << "attempted assignment to self"
1150  << abort(FatalError);
1151  }
1152 
1153  const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
1154 
1155  checkField(*this, gf, "=");
1156 
1157  // only equate field contents not ID
1158 
1159  this->dimensions() = gf.dimensions();
1160 
1161  // This is dodgy stuff, don't try it at home.
1162  internalField().transfer
1163  (
1164  const_cast<Field<Type>&>(gf.internalField())
1165  );
1166 
1167  boundaryField() = gf.boundaryField();
1168 
1169  tgf.clear();
1170 }
1171 
1172 
1173 template<class Type, template<class> class PatchField, class GeoMesh>
1175 (
1176  const dimensioned<Type>& dt
1177 )
1178 {
1179  dimensionedInternalField() = dt;
1180  boundaryField() = dt.value();
1181 }
1182 
1183 
1184 template<class Type, template<class> class PatchField, class GeoMesh>
1186 (
1187  const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
1188 )
1189 {
1190  const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
1191 
1192  checkField(*this, gf, "==");
1193 
1194  // only equate field contents not ID
1195 
1196  dimensionedInternalField() = gf.dimensionedInternalField();
1197  boundaryField() == gf.boundaryField();
1198 
1199  tgf.clear();
1200 }
1201 
1202 
1203 template<class Type, template<class> class PatchField, class GeoMesh>
1205 (
1206  const dimensioned<Type>& dt
1207 )
1208 {
1209  dimensionedInternalField() = dt;
1210  boundaryField() == dt.value();
1211 }
1212 
1213 
1214 #define COMPUTED_ASSIGNMENT(TYPE, op) \
1215  \
1216 template<class Type, template<class> class PatchField, class GeoMesh> \
1217 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1218 ( \
1219  const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1220 ) \
1221 { \
1222  checkField(*this, gf, #op); \
1223  \
1224  dimensionedInternalField() op gf.dimensionedInternalField(); \
1225  boundaryField() op gf.boundaryField(); \
1226 } \
1227  \
1228 template<class Type, template<class> class PatchField, class GeoMesh> \
1229 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1230 ( \
1231  const tmp<GeometricField<TYPE, PatchField, GeoMesh> >& tgf \
1232 ) \
1233 { \
1234  operator op(tgf()); \
1235  tgf.clear(); \
1236 } \
1237  \
1238 template<class Type, template<class> class PatchField, class GeoMesh> \
1239 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1240 ( \
1241  const dimensioned<TYPE>& dt \
1242 ) \
1243 { \
1244  dimensionedInternalField() op dt; \
1245  boundaryField() op dt.value(); \
1246 }
1247 
1248 COMPUTED_ASSIGNMENT(Type, +=)
1249 COMPUTED_ASSIGNMENT(Type, -=)
1250 COMPUTED_ASSIGNMENT(scalar, *=)
1251 COMPUTED_ASSIGNMENT(scalar, /=)
1252 
1253 #undef COMPUTED_ASSIGNMENT
1254 
1255 
1256 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1257 
1258 template<class Type, template<class> class PatchField, class GeoMesh>
1259 Foam::Ostream& Foam::operator<<
1260 (
1261  Ostream& os,
1262  const GeometricField<Type, PatchField, GeoMesh>& gf
1263 )
1264 {
1265  gf.dimensionedInternalField().writeData(os, "internalField");
1266  os << nl;
1267  gf.boundaryField().writeEntry("boundaryField", os);
1268 
1269  // Check state of IOstream
1270  os.check
1271  (
1272  "Ostream& operator<<(Ostream&, "
1273  "const GeometricField<Type, PatchField, GeoMesh>&)"
1274  );
1275 
1276  return (os);
1277 }
1278 
1279 
1280 template<class Type, template<class> class PatchField, class GeoMesh>
1281 Foam::Ostream& Foam::operator<<
1282 (
1283  Ostream& os,
1284  const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
1285 )
1286 {
1287  os << tgf();
1288  tgf.clear();
1289  return os;
1290 }
1291 
1292 
1293 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1294 
1295 #undef checkField
1296 
1297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1298 
1299 #include "GeometricBoundaryField.C"
1300 #include "GeometricFieldFunctions.C"
1301 
1302 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
InfoInFunction
#define InfoInFunction
Report a information message using Foam::Info.
Definition: messageStream.H:281
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
relax
p relax()
Foam::GeometricField::dimensionedInternalField
DimensionedInternalField & dimensionedInternalField()
Return dimensioned internal field.
Definition: GeometricField.C:713
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::GeometricField::writeData
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Definition: GeometricField.C:985
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
readField
void readField(const IOobject &io, const fvMesh &mesh, const label i, PtrList< GeoField > &fields)
Definition: redistributePar.C:557
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
dimensionedInternalField
rDeltaT dimensionedInternalField()
COMPUTED_ASSIGNMENT
#define COMPUTED_ASSIGNMENT(TYPE, op)
Definition: GeometricField.C:1214
boundary
faceListList boundary(nPatches)
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:804
Foam::GeometricField::readFields
void readFields()
Read the field - create the field dictionary on-the-fly.
Definition: GeometricField.C:72
checkField
#define checkField(gf1, gf2, op)
Definition: GeometricField.C:34
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::Field< Type >::operator+=
void operator+=(const UList< Type > &)
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::GeometricField::select
word select(bool final) const
Select the final iteration parameters if `final' is true.
Definition: GeometricField.C:955
Foam::GeometricField::cmptType
Field< Type >::cmptType cmptType
Definition: GeometricField.H:273
Foam::GeometricField::needReference
bool needReference() const
Does the field need a reference level for solution.
Definition: GeometricField.C:894
Foam::GeometricField::nOldTimes
label nOldTimes() const
Return the number of old time fields stored.
Definition: GeometricField.C:789
Foam::GeometricField::readIfPresent
bool readIfPresent()
Read from file if it is present.
Definition: GeometricField.C:95
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::GeometricField::min
void min(const dimensioned< Type > &)
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::GeometricField::storeOldTimes
void storeOldTimes() const
Store the old-time fields.
Definition: GeometricField.C:744
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::GeometricField::negate
void negate()
Definition: GeometricField.C:1109
Foam::label
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
Foam::Field< Type >
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::GeometricField::writeMinMax
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
Definition: GeometricField.C:972
Foam::GeometricField::readOldTimeIfPresent
bool readOldTimeIfPresent()
Read old time field from file if it is present.
Definition: GeometricField.C:132
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::GeometricField::prevIter
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
Definition: GeometricField.C:868
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::GeometricField::storeOldTime
void storeOldTime() const
Store the old-time field.
Definition: GeometricField.C:765
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::GeometricField::T
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Definition: GeometricField.C:996
internalField
conserve internalField()+
GeometricField.H
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::GeometricField::max
void max(const dimensioned< Type > &)
Foam::GeometricField::storePrevIter
void storePrevIter() const
Store the field as the previous iteration value.
Definition: GeometricField.C:843
Foam::andOp
Definition: ops.H:177
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::GeometricField::GeometricField
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch type.
Foam::GeometricField::relax
void relax()
Relax field (for steady-state solution).
Definition: GeometricField.C:930
readFields
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Definition: redistributePar.C:589
patchi
label patchi
Definition: getPatchFieldScalar.H:1
dictionary.H
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::GeometricField::~GeometricField
virtual ~GeometricField()
Destructor.
Definition: GeometricField.C:701
Foam::GeometricField::replace
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
data.H
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
GeometricBoundaryField.C
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
GeometricFieldFunctions.C
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51