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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "GeometricField.H"
30#include "Time.H"
31#include "demandDrivenData.H"
32#include "dictionary.H"
33#include "localIOdictionary.H"
34#include "data.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38#define checkField(gf1, gf2, op) \
39if ((gf1).mesh() != (gf2).mesh()) \
40{ \
41 FatalErrorInFunction \
42 << "different mesh for fields " \
43 << (gf1).name() << " and " << (gf2).name() \
44 << " during operation " << op \
45 << abort(FatalError); \
46}
47
48
49// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50
51template<class Type, template<class> class PatchField, class GeoMesh>
53(
54 const dictionary& dict
55)
56{
57 Internal::readField(dict, "internalField");
58
59 boundaryField_.readField(*this, dict.subDict("boundaryField"));
60
61 Type refLevel;
62
63 if (dict.readIfPresent("referenceLevel", refLevel))
64 {
65 Field<Type>::operator+=(refLevel);
66
67 forAll(boundaryField_, patchi)
68 {
69 boundaryField_[patchi] == boundaryField_[patchi] + refLevel;
70 }
71 }
72}
73
74
75template<class Type, template<class> class PatchField, class GeoMesh>
77{
78 const localIOdictionary dict
79 (
80 IOobject
81 (
82 this->name(),
83 this->instance(),
84 this->local(),
85 this->db(),
86 IOobject::MUST_READ,
87 IOobject::NO_WRITE,
88 false
89 ),
90 typeName
91 );
92
93 this->close();
94
96}
97
98
99template<class Type, template<class> class PatchField, class GeoMesh>
101{
102 if
103 (
104 this->readOpt() == IOobject::MUST_READ
105 || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
106 )
107 {
109 << "read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED"
110 << " suggests that a read constructor for field " << this->name()
111 << " would be more appropriate." << endl;
112 }
113 else if
114 (
115 this->readOpt() == IOobject::READ_IF_PRESENT
116 && this->template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
117 (
118 true
119 )
120 )
121 {
122 readFields();
123
124 // Check compatibility between field and mesh
125 if (this->size() != GeoMesh::size(this->mesh()))
126 {
127 FatalIOErrorInFunction(this->readStream(typeName))
128 << " number of field elements = " << this->size()
129 << " number of mesh elements = "
130 << GeoMesh::size(this->mesh())
131 << exit(FatalIOError);
132 }
133
134 readOldTimeIfPresent();
135
136 return true;
137 }
138
139 return false;
140}
141
142
143template<class Type, template<class> class PatchField, class GeoMesh>
145{
146 // Read the old time field if present
147 IOobject field0
148 (
149 this->name() + "_0",
150 this->time().timeName(),
151 this->db(),
152 IOobject::READ_IF_PRESENT,
153 IOobject::AUTO_WRITE,
154 this->registerObject()
155 );
156
157 if
158 (
159 field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
160 (
161 true
162 )
163 )
164 {
166 << "Reading old time level for field" << nl << this->info() << endl;
167
169 (
170 field0,
171 this->mesh()
172 );
173
174 // Ensure the old time field oriented flag is set to the parent's state
175 // Note: required for backwards compatibility in case of restarting from
176 // an old run where the oriented state may not have been set
177 field0Ptr_->oriented() = this->oriented();
178
179 field0Ptr_->timeIndex_ = timeIndex_ - 1;
180
181 if (!field0Ptr_->readOldTimeIfPresent())
182 {
183 field0Ptr_->oldTime();
184 }
185
186 return true;
187 }
189 return false;
190}
191
192
193// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
194
195template<class Type, template<class> class PatchField, class GeoMesh>
197(
198 const IOobject& io,
199 const Mesh& mesh,
200 const dimensionSet& ds,
201 const word& patchFieldType
202)
203:
204 Internal(io, mesh, ds, false),
205 timeIndex_(this->time().timeIndex()),
206 field0Ptr_(nullptr),
207 fieldPrevIterPtr_(nullptr),
208 boundaryField_(mesh.boundary(), *this, patchFieldType)
209{
211 << "Creating temporary" << nl << this->info() << endl;
212
215
216
217template<class Type, template<class> class PatchField, class GeoMesh>
219(
220 const IOobject& io,
221 const Mesh& mesh,
222 const dimensionSet& ds,
223 const wordList& patchFieldTypes,
224 const wordList& actualPatchTypes
225)
226:
227 Internal(io, mesh, ds, false),
228 timeIndex_(this->time().timeIndex()),
229 field0Ptr_(nullptr),
230 fieldPrevIterPtr_(nullptr),
231 boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
232{
234 << "Creating temporary" << nl << this->info() << endl;
235
237}
238
239
240template<class Type, template<class> class PatchField, class GeoMesh>
242(
243 const IOobject& io,
244 const Mesh& mesh,
245 const dimensioned<Type>& dt,
246 const word& patchFieldType
247)
248:
249 Internal(io, mesh, dt, false),
250 timeIndex_(this->time().timeIndex()),
251 field0Ptr_(nullptr),
252 fieldPrevIterPtr_(nullptr),
253 boundaryField_(mesh.boundary(), *this, patchFieldType)
254{
256 << "Creating temporary" << nl << this->info() << endl;
257
258 boundaryField_ == dt.value();
259
261}
262
263
264template<class Type, template<class> class PatchField, class GeoMesh>
267 const IOobject& io,
268 const Mesh& mesh,
269 const dimensioned<Type>& dt,
270 const wordList& patchFieldTypes,
271 const wordList& actualPatchTypes
272)
273:
274 Internal(io, mesh, dt, false),
275 timeIndex_(this->time().timeIndex()),
276 field0Ptr_(nullptr),
277 fieldPrevIterPtr_(nullptr),
278 boundaryField_(mesh.boundary(), *this, patchFieldTypes, actualPatchTypes)
279{
281 << "Creating temporary" << nl << this->info() << endl;
282
283 boundaryField_ == dt.value();
284
287
288
289template<class Type, template<class> class PatchField, class GeoMesh>
291(
292 const IOobject& io,
293 const Internal& diField,
294 const PtrList<PatchField<Type>>& ptfl
295)
296:
297 Internal(io, diField),
298 timeIndex_(this->time().timeIndex()),
299 field0Ptr_(nullptr),
300 fieldPrevIterPtr_(nullptr),
301 boundaryField_(this->mesh().boundary(), *this, ptfl)
304 << "Copy construct from components" << nl << this->info() << endl;
305
307}
309
310template<class Type, template<class> class PatchField, class GeoMesh>
312(
313 const IOobject& io,
314 Internal&& diField,
315 const PtrList<PatchField<Type>>& ptfl
316)
317:
318 Internal(io, std::move(diField)),
319 timeIndex_(this->time().timeIndex()),
320 field0Ptr_(nullptr),
321 fieldPrevIterPtr_(nullptr),
322 boundaryField_(this->mesh().boundary(), *this, ptfl)
323{
325 << "Move construct from components" << nl << this->info() << endl;
326
329
330
331template<class Type, template<class> class PatchField, class GeoMesh>
333(
334 const IOobject& io,
335 const tmp<Internal>& tfield,
336 const PtrList<PatchField<Type>>& ptfl
337)
338:
339 Internal(io, tfield),
340 timeIndex_(this->time().timeIndex()),
341 field0Ptr_(nullptr),
342 fieldPrevIterPtr_(nullptr),
343 boundaryField_(this->mesh().boundary(), *this, ptfl)
344{
346 << "Construct from tmp internalField" << nl << this->info() << endl;
347
349}
350
352template<class Type, template<class> class PatchField, class GeoMesh>
354(
355 const Internal& diField,
356 const PtrList<PatchField<Type>>& ptfl
357)
358:
359 Internal(diField),
360 timeIndex_(this->time().timeIndex()),
361 field0Ptr_(nullptr),
362 fieldPrevIterPtr_(nullptr),
363 boundaryField_(this->mesh().boundary(), *this, ptfl)
364{
366 << "Copy construct from components" << nl << this->info() << endl;
367
370
371
372template<class Type, template<class> class PatchField, class GeoMesh>
374(
375 Internal&& diField,
376 const PtrList<PatchField<Type>>& ptfl
377)
379 Internal(std::move(diField)),
380 timeIndex_(this->time().timeIndex()),
381 field0Ptr_(nullptr),
382 fieldPrevIterPtr_(nullptr),
383 boundaryField_(this->mesh().boundary(), *this, ptfl)
384{
386 << "Move construct from components" << nl << this->info() << endl;
387
389}
390
391
392template<class Type, template<class> class PatchField, class GeoMesh>
394(
395 const IOobject& io,
396 const Mesh& mesh,
397 const dimensionSet& ds,
398 const Field<Type>& iField,
399 const word& patchFieldType
400)
401:
402 Internal(io, mesh, ds, iField),
403 timeIndex_(this->time().timeIndex()),
404 field0Ptr_(nullptr),
405 fieldPrevIterPtr_(nullptr),
406 boundaryField_(mesh.boundary(), *this, patchFieldType)
407{
409 << "Copy construct from internal field" << nl << this->info() << endl;
410
411 readIfPresent();
412}
413
414
415template<class Type, template<class> class PatchField, class GeoMesh>
417(
418 const IOobject& io,
419 const Mesh& mesh,
420 const dimensionSet& ds,
421 Field<Type>&& iField,
422 const word& patchFieldType
423)
424:
425 Internal(io, mesh, ds, std::move(iField)),
426 timeIndex_(this->time().timeIndex()),
427 field0Ptr_(nullptr),
428 fieldPrevIterPtr_(nullptr),
429 boundaryField_(mesh.boundary(), *this, patchFieldType)
430{
432 << "Move construct from internal field" << nl << this->info() << endl;
433
434 readIfPresent();
435}
436
437
438template<class Type, template<class> class PatchField, class GeoMesh>
440(
441 const IOobject& io,
442 const Mesh& mesh,
443 const dimensionSet& ds,
444 const Field<Type>& iField,
445 const PtrList<PatchField<Type>>& ptfl
446)
447:
448 Internal(io, mesh, ds, iField),
449 timeIndex_(this->time().timeIndex()),
450 field0Ptr_(nullptr),
451 fieldPrevIterPtr_(nullptr),
452 boundaryField_(mesh.boundary(), *this, ptfl)
453{
455 << "Copy construct from components" << nl << this->info() << endl;
456
457 readIfPresent();
458}
459
460
461template<class Type, template<class> class PatchField, class GeoMesh>
463(
464 const IOobject& io,
465 const Mesh& mesh,
466 const dimensionSet& ds,
467 Field<Type>&& iField,
468 const PtrList<PatchField<Type>>& ptfl
469)
470:
471 Internal(io, mesh, ds, std::move(iField)),
472 timeIndex_(this->time().timeIndex()),
473 field0Ptr_(nullptr),
474 fieldPrevIterPtr_(nullptr),
475 boundaryField_(mesh.boundary(), *this, ptfl)
476{
478 << "Move construct from components" << nl << this->info() << endl;
479
481}
482
483
484template<class Type, template<class> class PatchField, class GeoMesh>
486(
487 const IOobject& io,
488 const Mesh& mesh,
489 const dimensionSet& ds,
490 const tmp<Field<Type>>& tfield,
491 const PtrList<PatchField<Type>>& ptfl
492)
493:
494 Internal(io, mesh, ds, tfield),
495 timeIndex_(this->time().timeIndex()),
496 field0Ptr_(nullptr),
497 fieldPrevIterPtr_(nullptr),
498 boundaryField_(mesh.boundary(), *this, ptfl)
499{
501 << "Construct from tmp internalField" << nl << this->info() << endl;
502
504}
505
506
507template<class Type, template<class> class PatchField, class GeoMesh>
509(
510 const IOobject& io,
511 const Mesh& mesh,
512 const bool readOldTime
513)
514:
515 Internal(io, mesh, dimless, false),
516 timeIndex_(this->time().timeIndex()),
517 field0Ptr_(nullptr),
518 fieldPrevIterPtr_(nullptr),
519 boundaryField_(mesh.boundary())
520{
521 readFields();
523 // Check compatibility between field and mesh
524
525 if (this->size() != GeoMesh::size(this->mesh()))
526 {
527 FatalIOErrorInFunction(this->readStream(typeName))
528 << " number of field elements = " << this->size()
529 << " number of mesh elements = " << GeoMesh::size(this->mesh())
530 << exit(FatalIOError);
531 }
533 if (readOldTime)
534 {
535 readOldTimeIfPresent();
536 }
537
539 << "Finishing read-construction" << nl << this->info() << endl;
540}
542
543template<class Type, template<class> class PatchField, class GeoMesh>
545(
546 const IOobject& io,
547 const Mesh& mesh,
548 const dictionary& dict
549)
550:
551 Internal(io, mesh, dimless, false),
552 timeIndex_(this->time().timeIndex()),
553 field0Ptr_(nullptr),
554 fieldPrevIterPtr_(nullptr),
555 boundaryField_(mesh.boundary())
558
559 // Check compatibility between field and mesh
560
561 if (this->size() != GeoMesh::size(this->mesh()))
564 << " number of field elements = " << this->size()
565 << " number of mesh elements = " << GeoMesh::size(this->mesh())
566 << exit(FatalIOError);
567 }
568
570 << "Finishing dictionary-construct" << nl << this->info() << endl;
572
573
574template<class Type, template<class> class PatchField, class GeoMesh>
576(
578)
579:
581 timeIndex_(gf.timeIndex()),
582 field0Ptr_(nullptr),
583 fieldPrevIterPtr_(nullptr),
584 boundaryField_(*this, gf.boundaryField_)
585{
587 << "Copy construct" << nl << this->info() << endl;
588
589 if (gf.field0Ptr_)
590 {
592 (
593 *gf.field0Ptr_
594 );
595 }
596
597 this->writeOpt(IOobject::NO_WRITE);
598}
599
600
601template<class Type, template<class> class PatchField, class GeoMesh>
603(
606:
607 Internal(tgf.constCast(), tgf.movable()),
608 timeIndex_(tgf().timeIndex()),
609 field0Ptr_(nullptr),
610 fieldPrevIterPtr_(nullptr),
611 boundaryField_(*this, tgf().boundaryField_)
612{
614 << "Constructing from tmp" << nl << this->info() << endl;
615
616 this->writeOpt(IOobject::NO_WRITE);
617
618 tgf.clear();
619}
621
622template<class Type, template<class> class PatchField, class GeoMesh>
624(
625 const IOobject& io,
627)
628:
629 Internal(io, gf),
630 timeIndex_(gf.timeIndex()),
631 field0Ptr_(nullptr),
632 fieldPrevIterPtr_(nullptr),
633 boundaryField_(*this, gf.boundaryField_)
636 << "Copy construct, resetting IO params" << nl
637 << this->info() << endl;
639 if (!readIfPresent() && gf.field0Ptr_)
642 (
643 io.name() + "_0",
644 *gf.field0Ptr_
645 );
648
650template<class Type, template<class> class PatchField, class GeoMesh>
653 const IOobject& io,
657 Internal(io, tgf.constCast(), tgf.movable()),
658 timeIndex_(tgf().timeIndex()),
659 field0Ptr_(nullptr),
660 fieldPrevIterPtr_(nullptr),
661 boundaryField_(*this, tgf().boundaryField_)
662{
664 << "Constructing from tmp resetting IO params" << nl
665 << this->info() << endl;
666
667 tgf.clear();
668
670}
671
672
673template<class Type, template<class> class PatchField, class GeoMesh>
675(
676 const word& newName,
678)
679:
680 Internal(newName, gf),
681 timeIndex_(gf.timeIndex()),
682 field0Ptr_(nullptr),
683 fieldPrevIterPtr_(nullptr),
684 boundaryField_(*this, gf.boundaryField_)
685{
687 << "Copy construct, resetting name" << nl
688 << this->info() << endl;
689
690 if (!readIfPresent() && gf.field0Ptr_)
691 {
693 (
694 newName + "_0",
695 *gf.field0Ptr_
696 );
697 }
698}
699
700
701template<class Type, template<class> class PatchField, class GeoMesh>
703(
704 const word& newName,
706)
707:
708 Internal(newName, tgf.constCast(), tgf.movable()),
709 timeIndex_(tgf().timeIndex()),
710 field0Ptr_(nullptr),
711 fieldPrevIterPtr_(nullptr),
712 boundaryField_(*this, tgf().boundaryField_)
713{
715 << "Constructing from tmp resetting name" << nl
716 << this->info() << endl;
717
718 tgf.clear();
719}
720
721
722template<class Type, template<class> class PatchField, class GeoMesh>
724(
725 const IOobject& io,
727 const word& patchFieldType
728)
729:
730 Internal(io, gf),
731 timeIndex_(gf.timeIndex()),
732 field0Ptr_(nullptr),
733 fieldPrevIterPtr_(nullptr),
734 boundaryField_(this->mesh().boundary(), *this, patchFieldType)
735{
737 << "Copy construct, resetting IO params" << nl
738 << this->info() << endl;
739
740 boundaryField_ == gf.boundaryField_;
741
742 if (!readIfPresent() && gf.field0Ptr_)
743 {
745 (
746 io.name() + "_0",
747 *gf.field0Ptr_
748 );
749 }
750}
751
752
753template<class Type, template<class> class PatchField, class GeoMesh>
755(
756 const IOobject& io,
758 const wordList& patchFieldTypes,
759 const wordList& actualPatchTypes
760)
761:
762 Internal(io, gf),
763 timeIndex_(gf.timeIndex()),
764 field0Ptr_(nullptr),
765 fieldPrevIterPtr_(nullptr),
766 boundaryField_
767 (
768 this->mesh().boundary(),
769 *this,
770 patchFieldTypes,
771 actualPatchTypes
772 )
773{
775 << "Copy construct, resetting IO params and patch types" << nl
776 << this->info() << endl;
777
778 boundaryField_ == gf.boundaryField_;
779
780 if (!readIfPresent() && gf.field0Ptr_)
781 {
783 (
784 io.name() + "_0",
785 *gf.field0Ptr_
786 );
787 }
788}
789
790
791template<class Type, template<class> class PatchField, class GeoMesh>
793(
794 const IOobject& io,
796 const labelList& patchIDs,
797 const word& patchFieldType
798)
799:
800 Internal(io, gf),
801 timeIndex_(gf.timeIndex()),
802 field0Ptr_(nullptr),
803 fieldPrevIterPtr_(nullptr),
804 boundaryField_(*this, gf.boundaryField_, patchIDs, patchFieldType)
805{
807 << "Copy construct, resetting IO params and setting patchFieldType "
808 << "for patchIDs" << nl
809 << this->info() << endl;
810
811 if (!readIfPresent() && gf.field0Ptr_)
812 {
814 (
815 io.name() + "_0",
816 *gf.field0Ptr_
817 );
818 }
819}
820
821
822template<class Type, template<class> class PatchField, class GeoMesh>
824(
825 const IOobject& io,
827 const wordList& patchFieldTypes,
828 const wordList& actualPatchTypes
829)
830:
831 Internal(io, tgf.constCast(), tgf.movable()),
832 timeIndex_(tgf().timeIndex()),
833 field0Ptr_(nullptr),
834 fieldPrevIterPtr_(nullptr),
835 boundaryField_
836 (
837 this->mesh().boundary(),
838 *this,
839 patchFieldTypes,
840 actualPatchTypes
841 )
842{
844 << "Constructing from tmp resetting IO params and patch types" << nl
845 << this->info() << endl;
846
847 boundaryField_ == tgf().boundaryField_;
848
849 tgf.clear();
850}
851
852
853template<class Type, template<class> class PatchField, class GeoMesh>
856{
858}
859
860
861// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
862
863template<class Type, template<class> class PatchField, class GeoMesh>
865{
866 deleteDemandDrivenData(field0Ptr_);
867 deleteDemandDrivenData(fieldPrevIterPtr_);
868}
869
870
871// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
872
873template<class Type, template<class> class PatchField, class GeoMesh>
874typename
877(
878 const bool updateAccessTime
879)
880{
881 if (updateAccessTime)
882 {
883 this->setUpToDate();
884 storeOldTimes();
885 }
886 return *this;
887}
888
889
890template<class Type, template<class> class PatchField, class GeoMesh>
891typename
894(
895 const bool updateAccessTime
896)
897{
898 if (updateAccessTime)
899 {
900 this->setUpToDate();
901 storeOldTimes();
902 }
903 return *this;
904}
905
906
907template<class Type, template<class> class PatchField, class GeoMesh>
908typename
911(
912 const bool updateAccessTime
913)
914{
915 if (updateAccessTime)
916 {
917 this->setUpToDate();
918 storeOldTimes();
919 }
920 return boundaryField_;
921}
922
923
924template<class Type, template<class> class PatchField, class GeoMesh>
926{
927 if
928 (
929 field0Ptr_
930 && timeIndex_ != this->time().timeIndex()
931 && !this->name().ends_with("_0")
932 )
933 {
934 storeOldTime();
935 timeIndex_ = this->time().timeIndex();
936 }
937
938 // Correct time index
939 //timeIndex_ = this->time().timeIndex();
940}
941
942
943template<class Type, template<class> class PatchField, class GeoMesh>
945{
946 if (field0Ptr_)
947 {
948 field0Ptr_->storeOldTime();
949
951 << "Storing old time field for field" << nl << this->info() << endl;
952
953 *field0Ptr_ == *this;
954 field0Ptr_->timeIndex_ = timeIndex_;
955
956 if (field0Ptr_->field0Ptr_)
957 {
958 field0Ptr_->writeOpt(this->writeOpt());
959 }
960 }
961}
962
963
964template<class Type, template<class> class PatchField, class GeoMesh>
966{
967 if (field0Ptr_)
968 {
969 return field0Ptr_->nOldTimes() + 1;
970 }
971
972 return 0;
973}
974
975
976template<class Type, template<class> class PatchField, class GeoMesh>
979{
980 if (!field0Ptr_)
981 {
983 (
985 (
986 this->name() + "_0",
987 this->time().timeName(),
988 this->db(),
991 this->registerObject()
992 ),
993 *this
994 );
995
996 if (debug)
997 {
999 << "created old time field " << field0Ptr_->info() << endl;
1000
1001 if (debug&2)
1002 {
1004 }
1005 }
1006 }
1007 else
1008 {
1009 storeOldTimes();
1010 }
1011
1012 return *field0Ptr_;
1013}
1014
1015
1016template<class Type, template<class> class PatchField, class GeoMesh>
1019{
1020 static_cast<const GeometricField<Type, PatchField, GeoMesh>&>(*this)
1021 .oldTime();
1022
1023 return *field0Ptr_;
1024}
1025
1026
1027template<class Type, template<class> class PatchField, class GeoMesh>
1029{
1030 if (!fieldPrevIterPtr_)
1031 {
1033 << "Allocating previous iteration field" << nl
1034 << this->info() << endl;
1035
1036 fieldPrevIterPtr_ = new GeometricField<Type, PatchField, GeoMesh>
1037 (
1038 this->name() + "PrevIter",
1039 *this
1040 );
1041 }
1042 else
1043 {
1044 *fieldPrevIterPtr_ == *this;
1045 }
1046}
1047
1048
1049template<class Type, template<class> class PatchField, class GeoMesh>
1052{
1053 if (!fieldPrevIterPtr_)
1054 {
1056 << "previous iteration field" << endl << this->info() << endl
1057 << " not stored."
1058 << " Use field.storePrevIter() at start of iteration."
1059 << abort(FatalError);
1060 }
1061
1062 return *fieldPrevIterPtr_;
1063}
1064
1065
1066template<class Type, template<class> class PatchField, class GeoMesh>
1069{
1070 this->setUpToDate();
1071 storeOldTimes();
1072 boundaryField_.evaluate();
1073}
1074
1075
1076template<class Type, template<class> class PatchField, class GeoMesh>
1078{
1079 // Search all boundary conditions, if any are
1080 // fixed-value or mixed (Robin) do not set reference level for solution.
1081
1082 bool needRef = true;
1083
1084 forAll(boundaryField_, patchi)
1085 {
1086 if (boundaryField_[patchi].fixesValue())
1087 {
1088 needRef = false;
1089 break;
1090 }
1091 }
1092
1093 reduce(needRef, andOp<bool>());
1094
1095 return needRef;
1096}
1097
1098
1099template<class Type, template<class> class PatchField, class GeoMesh>
1101{
1103 << "Relaxing" << nl << this->info() << " by " << alpha << endl;
1104
1105 operator==(prevIter() + alpha*(*this - prevIter()));
1106}
1107
1108
1109template<class Type, template<class> class PatchField, class GeoMesh>
1111{
1112 word name = this->name();
1113
1114 if
1115 (
1116 this->mesh().data::template getOrDefault<bool>
1117 (
1118 "finalIteration",
1119 false
1120 )
1121 )
1122 {
1123 name += "Final";
1124 }
1125
1126 if (this->mesh().relaxField(name))
1127 {
1128 relax(this->mesh().fieldRelaxationFactor(name));
1129 }
1130}
1131
1132
1133template<class Type, template<class> class PatchField, class GeoMesh>
1135(
1136 bool final
1137) const
1138{
1139 if (final)
1140 {
1141 return this->name() + "Final";
1142 }
1143
1144 return this->name();
1145}
1146
1147
1148template<class Type, template<class> class PatchField, class GeoMesh>
1150(
1151 Ostream& os
1152) const
1153{
1154 MinMax<Type> range = Foam::minMax(*this).value();
1155
1156 os << "min/max(" << this->name() << ") = "
1157 << range.min() << ", " << range.max() << endl;
1158}
1159
1160
1161template<class Type, template<class> class PatchField, class GeoMesh>
1163writeData(Ostream& os) const
1164{
1165 os << *this;
1166 return os.good();
1167}
1168
1169
1170// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1171
1172template<class Type, template<class> class PatchField, class GeoMesh>
1175{
1177 (
1178 IOobject
1179 (
1180 this->name() + ".T()",
1181 this->instance(),
1182 this->db()
1183 ),
1184 this->mesh(),
1185 this->dimensions()
1186 );
1187
1188 Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1189 Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1190
1191 return tresult;
1192}
1193
1194
1195template<class Type, template<class> class PatchField, class GeoMesh>
1197<
1199 <
1201 PatchField,
1202 GeoMesh
1203 >
1204>
1206(
1207 const direction d
1208) const
1209{
1211 (
1212 IOobject
1213 (
1214 this->name() + ".component(" + Foam::name(d) + ')',
1215 this->instance(),
1216 this->db()
1217 ),
1218 this->mesh(),
1219 this->dimensions()
1220 );
1221
1222 Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1223 Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1224
1225 return tresult;
1226}
1227
1228
1229template<class Type, template<class> class PatchField, class GeoMesh>
1231(
1232 const direction d,
1233 const GeometricField
1234 <
1236 PatchField,
1237 GeoMesh
1238 >& gcf
1239)
1240{
1241 primitiveFieldRef().replace(d, gcf.primitiveField());
1242 boundaryFieldRef().replace(d, gcf.boundaryField());
1243}
1244
1245
1246template<class Type, template<class> class PatchField, class GeoMesh>
1248(
1249 const direction d,
1250 const dimensioned<cmptType>& ds
1251)
1252{
1253 primitiveFieldRef().replace(d, ds.value());
1254 boundaryFieldRef().replace(d, ds.value());
1255}
1256
1257
1258template<class Type, template<class> class PatchField, class GeoMesh>
1260(
1261 const dimensioned<Type>& dt
1262)
1263{
1264 Foam::min(primitiveFieldRef(), primitiveField(), dt.value());
1265 Foam::min(boundaryFieldRef(), boundaryField(), dt.value());
1266}
1267
1268
1269template<class Type, template<class> class PatchField, class GeoMesh>
1271(
1272 const dimensioned<Type>& dt
1273)
1274{
1275 Foam::max(primitiveFieldRef(), primitiveField(), dt.value());
1276 Foam::max(boundaryFieldRef(), boundaryField(), dt.value());
1277}
1278
1279
1280template<class Type, template<class> class PatchField, class GeoMesh>
1282(
1284)
1285{
1286 Foam::clip(primitiveFieldRef(), primitiveField(), range.value());
1287 Foam::clip(boundaryFieldRef(), boundaryField(), range.value());
1288}
1289
1290
1291template<class Type, template<class> class PatchField, class GeoMesh>
1293(
1294 const dimensioned<Type>& minVal,
1295 const dimensioned<Type>& maxVal
1296)
1297{
1298 MinMax<Type> range(minVal.value(), maxVal.value());
1299
1300 Foam::clip(primitiveFieldRef(), primitiveField(), range);
1301 Foam::clip(boundaryFieldRef(), boundaryField(), range);
1302}
1303
1304
1305template<class Type, template<class> class PatchField, class GeoMesh>
1307(
1308 const dimensioned<Type>& minVal,
1309 const dimensioned<Type>& maxVal
1310)
1311{
1312 this->clip(minVal, maxVal);
1313}
1314
1315
1316template<class Type, template<class> class PatchField, class GeoMesh>
1318{
1319 primitiveFieldRef().negate();
1320 boundaryFieldRef().negate();
1321}
1322
1323
1324template<class Type, template<class> class PatchField, class GeoMesh>
1326{
1327 primitiveFieldRef().normalise();
1328 boundaryFieldRef().normalise();
1329}
1330
1331
1332// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1333
1334template<class Type, template<class> class PatchField, class GeoMesh>
1336(
1338)
1339{
1340 if (this == &gf)
1341 {
1342 return; // Self-assignment is a no-op
1343 }
1344
1345 checkField(*this, gf, "=");
1346
1347 // Only assign field contents not ID
1348
1349 ref() = gf();
1350 boundaryFieldRef() = gf.boundaryField();
1351}
1352
1353
1354template<class Type, template<class> class PatchField, class GeoMesh>
1356(
1358)
1359{
1360 const auto& gf = tgf();
1361
1362 if (this == &gf)
1363 {
1364 return; // Self-assignment is a no-op
1365 }
1366
1367 checkField(*this, gf, "=");
1368
1369 // Only assign field contents not ID
1370
1371 this->dimensions() = gf.dimensions();
1372 this->oriented() = gf.oriented();
1373
1374 if (tgf.movable())
1375 {
1376 // Transfer storage from the tmp
1377 primitiveFieldRef().transfer(tgf.constCast().primitiveFieldRef());
1378 }
1379 else
1380 {
1382 }
1383
1384 boundaryFieldRef() = gf.boundaryField();
1385
1386 tgf.clear();
1387}
1388
1389
1390template<class Type, template<class> class PatchField, class GeoMesh>
1392(
1393 const dimensioned<Type>& dt
1394)
1395{
1396 ref() = dt;
1397 boundaryFieldRef() = dt.value();
1398}
1399
1400
1401template<class Type, template<class> class PatchField, class GeoMesh>
1403(
1405)
1406{
1407 const auto& gf = tgf();
1408
1409 checkField(*this, gf, "==");
1410
1411 // Only assign field contents not ID
1412
1413 ref() = gf();
1414 boundaryFieldRef() == gf.boundaryField();
1415
1416 tgf.clear();
1417}
1418
1419
1420template<class Type, template<class> class PatchField, class GeoMesh>
1422(
1423 const dimensioned<Type>& dt
1424)
1425{
1426 ref() = dt;
1427 boundaryFieldRef() == dt.value();
1428}
1429
1430
1431#define COMPUTED_ASSIGNMENT(TYPE, op) \
1432 \
1433template<class Type, template<class> class PatchField, class GeoMesh> \
1434void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1435( \
1436 const GeometricField<TYPE, PatchField, GeoMesh>& gf \
1437) \
1438{ \
1439 checkField(*this, gf, #op); \
1440 \
1441 ref() op gf(); \
1442 boundaryFieldRef() op gf.boundaryField(); \
1443} \
1444 \
1445template<class Type, template<class> class PatchField, class GeoMesh> \
1446void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1447( \
1448 const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \
1449) \
1450{ \
1451 operator op(tgf()); \
1452 tgf.clear(); \
1453} \
1454 \
1455template<class Type, template<class> class PatchField, class GeoMesh> \
1456void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \
1457( \
1458 const dimensioned<TYPE>& dt \
1459) \
1460{ \
1461 ref() op dt; \
1462 boundaryFieldRef() op dt.value(); \
1463}
1464
1469
1470#undef COMPUTED_ASSIGNMENT
1471
1472
1473// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1474
1475template<class Type, template<class> class PatchField, class GeoMesh>
1476Foam::Ostream& Foam::operator<<
1477(
1478 Ostream& os,
1480)
1481{
1482 gf().writeData(os, "internalField");
1483 os << nl;
1484 gf.boundaryField().writeEntry("boundaryField", os);
1485
1487 return os;
1488}
1489
1490
1491template<class Type, template<class> class PatchField, class GeoMesh>
1492Foam::Ostream& Foam::operator<<
1493(
1494 Ostream& os,
1496)
1497{
1498 os << tgf();
1499 tgf.clear();
1500 return os;
1501}
1502
1503
1504// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1505
1506#undef checkField
1507
1508// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1509
1510#include "GeometricFieldNew.C"
1512
1513// ************************************************************************* //
#define checkField(df1, df2, op)
#define COMPUTED_ASSIGNMENT(TYPE, op)
scalar range
const dimensionSet & dimensions() const
Return dimensions.
const orientedType & oriented() const noexcept
Return oriented type.
Generic templated field type.
Definition: Field.H:82
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Generic GeometricBoundaryField class.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
Generic GeometricField class.
void storeOldTime() const
Store the old-time field.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
void maxMin(const dimensioned< Type > &minVal, const dimensioned< Type > &maxVal)
Deprecated(2019-01) identical to clip()
void relax()
Relax field (for steady-state solution).
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
GeoMesh::Mesh Mesh
The mesh type for the GeometricField.
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void clip(const dimensioned< MinMax< Type > > &range)
Clip the field to be bounded within the specified range.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
bool needReference() const
Does the field need a reference level for solution.
virtual ~GeometricField()
Destructor.
void negate()
Negate the field inplace. See notes in Field.
void storePrevIter() const
Store the field as the previous iteration value.
label nOldTimes() const
Return the number of old time fields stored.
word select(bool final) const
Select the final iteration parameters if `final' is true.
void storeOldTimes() const
Store the old-time fields.
void correctBoundaryConditions()
Correct boundary field.
void normalise()
Normalise the field inplace. See notes in Field.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
InfoProxy< IOobject > info() const
Return info proxy, for printing information to a stream.
Definition: IOobject.H:690
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition: MinMax.H:128
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:97
reference ref() const
A reference to the entry (Error if not found)
Definition: dictionary.H:225
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Generic dimensioned Type class.
const Type & value() const
Return const reference to value.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
static void printStack(Ostream &os)
Helper function to print a stack.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
static const complex min
complex (-VGREAT,-VGREAT)
Definition: complex.H:292
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
UEqn relax()
rDeltaT ref()
faceListList boundary
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:278
const dimensionSet dimless
Dimensionless.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void deleteDemandDrivenData(DataPtr &dataPtr)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
propsDict readIfPresent("fields", acceptFields)
conserve primitiveFieldRef()+