Loading...
Searching...
No Matches
GeometricBoundaryField.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,2022 OpenFOAM Foundation
9 Copyright (C) 2016-2025 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
30#include "globalMeshData.H"
31#include "cyclicPolyPatch.H"
32#include "emptyPolyPatch.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36template<class Type, template<class> class PatchField, class GeoMesh>
37template<class CheckPatchFieldType>
38bool Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::checkConsistency
39(
40 const scalar tol,
41 const bool exitIfBad
42) const
43{
44 auto& bfld = this->constCast();
45
46 if (!bfld.size())
47 {
48 return true;
49 }
50
51 if (debug&2)
52 {
54 << " Checking boundary consistency for field "
55 << bfld[0].internalField().name() << endl;
56 }
57
58 // Store old values and states
59 List<Field<Type>> oldFields(bfld.size());
60 boolList oldUpdated(bfld.size());
61 boolList oldManipulated(bfld.size());
62
63 label nEvaluated(0);
64
65 for (auto& pfld : bfld)
66 {
67 if (isA<CheckPatchFieldType>(pfld))
68 {
69 const label patchi = pfld.patch().index();
70 oldFields[patchi] = pfld;
71 oldUpdated[patchi] = pfld.updated();
72 oldManipulated[patchi] = pfld.manipulatedMatrix();
73 ++nEvaluated;
74 }
75 }
76
77 if (!nEvaluated) return true; // Early termination
78
79 // Re-evaluate
80 {
81 const label startOfRequests = UPstream::nRequests();
82
83 nEvaluated = 0;
84
85 for (auto& pfld : bfld)
86 {
87 if (isA<CheckPatchFieldType>(pfld))
88 {
89 pfld.initEvaluate(UPstream::commsTypes::nonBlocking);
90 ++nEvaluated;
91 }
92 }
93
94 // Wait for outstanding requests (non-blocking)
95 UPstream::waitRequests(startOfRequests);
96
97 for (auto& pfld : bfld)
98 {
99 if (isA<CheckPatchFieldType>(pfld))
100 {
101 pfld.evaluate(UPstream::commsTypes::nonBlocking);
102 if (--nEvaluated == 0) break; // Early termination
103 }
104 }
105 }
106
107
108 // Check
109 bool allOk(true);
110
111 for (auto& pfld : bfld)
112 {
113 if (isA<CheckPatchFieldType>(pfld))
114 {
115 const label patchi = pfld.patch().index();
116 auto& oldPfld = oldFields[patchi];
117
118 bool localOk(true);
119
120 if (allOk)
121 {
122 // Only check for first failed patch
123 forAll(pfld, facei)
124 {
125 if (tol < Foam::mag(pfld[facei]-oldPfld[facei]))
126 {
127 allOk = false;
128 localOk = false;
129 break;
130 }
131 }
132 }
133
134 if (!localOk)
135 {
136 // Raise warning or error
137 OSstream& err =
138 (
139 exitIfBad
142 );
143
144 err << "Field "
145 << pfld.internalField().name()
146 << " is not evaluated?"
147 << " On patch " << pfld.patch().name()
148 << " type " << pfld.type()
149 << " : average of field = "
150 << average(oldPfld)
151 << ". Average of evaluated field = "
152 << average(pfld)
153 << ". Difference:" << average(pfld-oldPfld)
154 << ". Tolerance:" << tol << endl;
155
156 if (exitIfBad)
157 {
158 FatalError<< exit(FatalError);
159 }
160 }
161
162 // Restore patch field values and states
163 static_cast<Field<Type>&>(pfld) = std::move(oldPfld);
164 pfld.setUpdated(oldUpdated[patchi]);
165 pfld.setManipulated(oldManipulated[patchi]);
166 }
167 }
168
169 if (debug&2)
170 {
172 << " Result of checking for field "
173 << bfld[0].internalField().name() << " : " << allOk << endl;
175
176 return allOk;
177}
178
179
180template<class Type, template<class> class PatchField, class GeoMesh>
182(
183 const Internal& iField,
184 const dictionary& dict
185)
186{
187 // DebugInFunction << nl;
188
189 // Clear the boundary field if already initialised
190 this->clear();
191
192 this->resize(bmesh_.size());
193
194 label nUnset = this->size();
195
196 // 1. Handle explicit patch names. Note that there can be only one explicit
197 // patch name since is key of dictionary.
198
199 for (const entry& dEntry : dict)
200 {
201 const auto* subdict = dEntry.dictPtr();
202
203 if (subdict && dEntry.keyword().isLiteral())
204 {
205 const label patchi = bmesh_.findPatchID(dEntry.keyword());
206
207 if (patchi != -1)
208 {
209 this->set
210 (
211 patchi,
212 PatchField<Type>::New
213 (
214 bmesh_[patchi],
215 iField,
216 *subdict
217 )
218 );
219 --nUnset;
220 }
221 }
222 }
223
224 if (nUnset == 0)
226 return;
227 }
228
229
230 // 2. Patch-groups. (using non-wild card entries of dictionaries)
231 // (patchnames already matched above)
232 // Note: in reverse order of entries in the dictionary (last
233 // patchGroups wins). This is so it is consistent with dictionary wildcard
234 // behaviour
235 for (auto iter = dict.crbegin(); iter != dict.crend(); ++iter)
236 {
237 const entry& dEntry = *iter;
238 const auto* subdict = dEntry.dictPtr();
239
240 if (subdict && dEntry.keyword().isLiteral())
241 {
242 const labelList patchIds =
243 bmesh_.indices(dEntry.keyword(), true); // use patchGroups
244
245 for (const label patchi : patchIds)
246 {
247 if (!this->set(patchi))
248 {
249 this->set
250 (
251 patchi,
252 PatchField<Type>::New
253 (
254 bmesh_[patchi],
255 iField,
256 *subdict
257 )
258 );
259 }
261 }
262 }
263
264
265 // 3. Wildcard patch overrides
266 forAll(bmesh_, patchi)
267 {
268 if (!this->set(patchi))
270 if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
271 {
272 this->set
273 (
274 patchi,
275 PatchField<Type>::New
276 (
277 emptyPolyPatch::typeName,
278 bmesh_[patchi],
279 iField
280 )
281 );
282 }
283 else
284 {
285 const auto* subdict = dict.findDict(bmesh_[patchi].name());
286
287 if (subdict)
288 {
289 this->set
290 (
291 patchi,
292 PatchField<Type>::New
293 (
294 bmesh_[patchi],
295 iField,
296 *subdict
297 )
298 );
299 }
300 }
301 }
302 }
303
305 // Check for any unset patches
306 forAll(bmesh_, patchi)
307 {
308 if (!this->set(patchi))
310 if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
311 {
313 << "Cannot find patchField entry for cyclic "
314 << bmesh_[patchi].name() << endl
315 << "Is your field uptodate with split cyclics?" << endl
316 << "Run foamUpgradeCyclics to convert mesh and fields"
317 << " to split cyclics." << exit(FatalIOError);
318 }
319 else
320 {
322 << "Cannot find patchField entry for "
323 << bmesh_[patchi].name() << exit(FatalIOError);
324 }
325 }
327}
328
329
330// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
331
332template<class Type, template<class> class PatchField, class GeoMesh>
334(
335 const BoundaryMesh& bmesh
336)
338 FieldField<PatchField, Type>(bmesh.size()),
339 bmesh_(bmesh)
340{}
341
342
343template<class Type, template<class> class PatchField, class GeoMesh>
345(
346 const BoundaryMesh& bmesh,
347 const Internal& iField,
348 const word& patchFieldType
349)
350:
351 FieldField<PatchField, Type>(bmesh.size()),
352 bmesh_(bmesh)
353{
354 // DebugInFunction << nl;
355
356 forAll(bmesh_, patchi)
357 {
358 this->set
359 (
360 patchi,
361 PatchField<Type>::New
362 (
363 patchFieldType,
364 bmesh_[patchi],
365 iField
367 );
368 }
369}
370
371
372template<class Type, template<class> class PatchField, class GeoMesh>
374(
375 const BoundaryMesh& bmesh,
376 const Internal& iField,
377 const wordList& patchFieldTypes,
378 const wordList& constraintTypes
379)
380:
381 FieldField<PatchField, Type>(bmesh.size()),
382 bmesh_(bmesh)
383{
384 // DebugInFunction << nl;
385
386 if
387 (
388 patchFieldTypes.size() != this->size()
389 || (constraintTypes.size() && (constraintTypes.size() != this->size()))
390 )
391 {
392 FatalErrorInFunction
393 << "Incorrect number of patch type specifications given" << nl
394 << " Number of patches in mesh = " << bmesh.size()
395 << " number of patch type specifications = "
396 << patchFieldTypes.size()
397 << abort(FatalError);
398 }
399
400 if (constraintTypes.size())
401 {
402 forAll(bmesh_, patchi)
403 {
404 this->set
405 (
406 patchi,
407 PatchField<Type>::New
408 (
409 patchFieldTypes[patchi],
410 constraintTypes[patchi],
411 bmesh_[patchi],
412 iField
413 )
414 );
415 }
416 }
417 else
418 {
419 forAll(bmesh_, patchi)
420 {
421 this->set
422 (
423 patchi,
424 PatchField<Type>::New
425 (
426 patchFieldTypes[patchi],
427 bmesh_[patchi],
428 iField
429 )
430 );
431 }
432 }
433}
434
435
436template<class Type, template<class> class PatchField, class GeoMesh>
438(
439 const BoundaryMesh& bmesh,
440 const Internal& iField,
441 const PtrList<PatchField<Type>>& ptfl
442)
443:
444 FieldField<PatchField, Type>(bmesh.size()),
445 bmesh_(bmesh)
446{
447 // DebugInFunction << nl;
448
449 forAll(bmesh_, patchi)
451 this->set(patchi, ptfl[patchi].clone(iField));
452 }
453}
454
455
456template<class Type, template<class> class PatchField, class GeoMesh>
458(
459 const Internal& iField,
460 const GeometricBoundaryField<Type, PatchField, GeoMesh>& btf
461)
462:
463 FieldField<PatchField, Type>(btf.size()),
464 bmesh_(btf.bmesh_)
465{
466 // DebugInFunction << nl;
467
468 forAll(bmesh_, patchi)
470 this->set(patchi, btf[patchi].clone(iField));
471 }
472}
473
474
475template<class Type, template<class> class PatchField, class GeoMesh>
477(
478 const Internal& iField,
479 const GeometricBoundaryField<Type, PatchField, GeoMesh>& btf,
480 const labelList& patchIDs,
481 const word& patchFieldType
482)
483:
484 FieldField<PatchField, Type>(btf.size()),
485 bmesh_(btf.bmesh_)
486{
487 // DebugInFunction << nl;
488
489 for (const label patchi : patchIDs)
490 {
491 this->set
492 (
493 patchi,
494 PatchField<Type>::New
495 (
496 patchFieldType,
497 bmesh_[patchi],
498 iField
499 )
500 );
501 }
502
503 forAll(bmesh_, patchi)
504 {
505 if (!this->set(patchi))
506 {
507 this->set(patchi, btf[patchi].clone(iField));
508 }
509 }
510}
511
512
513template<class Type, template<class> class PatchField, class GeoMesh>
515(
517)
519 FieldField<PatchField, Type>(btf),
520 bmesh_(btf.bmesh_)
521{}
522
523
524template<class Type, template<class> class PatchField, class GeoMesh>
526(
527 const BoundaryMesh& bmesh,
528 const Internal& iField,
529 const dictionary& dict
530)
531:
532 FieldField<PatchField, Type>(bmesh.size()),
533 bmesh_(bmesh)
534{
535 readField(iField, dict);
536}
537
538
539// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
540
541template<class Type, template<class> class PatchField, class GeoMesh>
543{
544 // DebugInFunction << nl;
545
546 for (auto& pfld : *this)
548 pfld.updateCoeffs();
549 }
550}
551
552
553template<class Type, template<class> class PatchField, class GeoMesh>
555(
556 const UPstream::commsTypes commsType
557)
558{
559 if
560 (
563 )
564 {
565 const label startOfRequests = UPstream::nRequests();
566
567 for (auto& pfld : *this)
568 {
569 pfld.initEvaluate(commsType);
570 }
571
572 // Wait for outstanding requests (non-blocking)
573 UPstream::waitRequests(startOfRequests);
574
575 for (auto& pfld : *this)
576 {
577 pfld.evaluate(commsType);
578 }
579 }
580 else if (commsType == UPstream::commsTypes::scheduled)
581 {
582 const lduSchedule& patchSchedule =
583 bmesh_.mesh().globalData().patchSchedule();
584
585 for (const auto& schedEval : patchSchedule)
586 {
587 const label patchi = schedEval.patch;
588 auto& pfld = (*this)[patchi];
589
590 if (schedEval.init)
591 {
592 pfld.initEvaluate(commsType);
593 }
594 else
595 {
596 pfld.evaluate(commsType);
597 }
598 }
599 }
600 else
601 {
603 << "Unsupported communications type " << int(commsType) << nl
605 }
606}
607
608
609template<class Type, template<class> class PatchField, class GeoMesh>
610template<class UnaryPredicate>
612(
613 const UnaryPredicate& pred,
614 const UPstream::commsTypes commsType
615)
616{
617 if
618 (
621 )
622 {
623 label nEvaluated(0);
624 const label startOfRequests = UPstream::nRequests();
625
626 for (auto& pfld : *this)
627 {
628 if (pred(pfld))
629 {
630 pfld.initEvaluate(commsType);
631 ++nEvaluated;
632 }
633 }
634
635 // Wait for outstanding requests (non-blocking)
636 UPstream::waitRequests(startOfRequests);
637
638 if (!nEvaluated) return; // Early termination
639
640 for (auto& pfld : *this)
641 {
642 if (pred(pfld))
643 {
644 pfld.evaluate(commsType);
645 if (--nEvaluated == 0) break; // Early termination
646 }
647 }
648 }
649 else if (commsType == UPstream::commsTypes::scheduled)
650 {
651 const lduSchedule& patchSchedule =
652 bmesh_.mesh().globalData().patchSchedule();
653
654 for (const auto& schedEval : patchSchedule)
655 {
656 const label patchi = schedEval.patch;
657 auto& pfld = (*this)[patchi];
658
659 if (pred(pfld))
660 {
661 if (schedEval.init)
662 {
663 pfld.initEvaluate(commsType);
664 }
665 else
666 {
667 pfld.evaluate(commsType);
668 }
669 }
670 }
671 }
672 else
673 {
675 << "Unsupported communications type " << int(commsType) << nl
676 << exit(FatalError);
677 }
678}
679
680
681template<class Type, template<class> class PatchField, class GeoMesh>
683(
684 const labelUList& patchIDs
685)
686{
687 const label startOfRequests = UPstream::nRequests();
688
689 for (const label patchi : patchIDs)
690 {
691 auto& pfld = (*this)[patchi];
692
693 DebugInfo<< "Updating " << pfld.patch().name() << endl;
694
695 pfld.initEvaluate(UPstream::commsTypes::nonBlocking);
696 }
697
698 // Wait for outstanding requests (non-blocking)
699 UPstream::waitRequests(startOfRequests);
700
701 for (const label patchi : patchIDs)
702 {
703 auto& pfld = (*this)[patchi];
706 }
707}
708
709
710template<class Type, template<class> class PatchField, class GeoMesh>
712(
713 const UPstream::commsTypes commsType
714)
715{
716 // DebugInFunction << nl;
717
719 {
720 return;
721 }
722
723 if
724 (
727 )
728 {
729 const label startOfRequests = UPstream::nRequests();
730
731 for (auto& pfld : *this)
732 {
733 pfld.initEvaluateLocal(commsType);
734 }
735
736 // Wait for outstanding requests (non-blocking)
737 UPstream::waitRequests(startOfRequests);
738
739 for (auto& pfld : *this)
740 {
741 pfld.evaluateLocal(commsType);
742 }
743 }
744 else if (commsType == UPstream::commsTypes::scheduled)
745 {
746 const lduSchedule& patchSchedule =
747 bmesh_.mesh().globalData().patchSchedule();
748
749 for (const auto& schedEval : patchSchedule)
750 {
751 const label patchi = schedEval.patch;
752 auto& pfld = (*this)[patchi];
753
754 if (schedEval.init)
755 {
756 pfld.initEvaluateLocal(commsType);
757 }
758 else
759 {
760 pfld.evaluateLocal(commsType);
761 }
762 }
763 }
764 else
765 {
767 << "Unsupported communications type " << int(commsType) << nl
769 }
770}
771
772
773template<class Type, template<class> class PatchField, class GeoMesh>
774template<class CoupledPatchType>
776(
777 const UPstream::commsTypes commsType
778)
779{
780 if constexpr (std::is_void_v<CoupledPatchType>)
781 {
782 this->evaluate_if
783 (
784 [](const auto& pfld) { return pfld.coupled(); },
785 commsType
786 );
787 }
788 else
789 {
790 this->evaluate_if
791 (
792 [](const auto& pfld) -> bool
793 {
794 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
795 return (cpp && cpp->coupled());
796 },
797 commsType
798 );
799 }
800}
801
802
803template<class Type, template<class> class PatchField, class GeoMesh>
806{
807 const FieldField<PatchField, Type>& pff = *this;
808
809 wordList list(pff.size());
810
811 forAll(pff, patchi)
812 {
813 list[patchi] = pff[patchi].type();
814 }
816 return list;
817}
818
819
820template<class Type, template<class> class PatchField, class GeoMesh>
824{
826 (
828 *this
829 );
830
831 auto& result = tresult.ref();
832
833 forAll(result, patchi)
834 {
835 result[patchi] == this->operator[](patchi).patchInternalField();
836 }
838 return tresult;
839}
840
841
842template<class Type, template<class> class PatchField, class GeoMesh>
845{
846 LduInterfaceFieldPtrsList<Type> list(this->size());
847
848 forAll(list, patchi)
849 {
850 list.set
851 (
852 patchi,
853 isA<LduInterfaceField<Type>>(this->operator[](patchi))
854 );
855 }
857 return list;
858}
859
860
861template<class Type, template<class> class PatchField, class GeoMesh>
864scalarInterfaces() const
865{
866 lduInterfaceFieldPtrsList list(this->size());
867
868 forAll(list, patchi)
869 {
870 list.set
871 (
872 patchi,
873 isA<lduInterfaceField>(this->operator[](patchi))
874 );
876
877 return list;
878}
879
880
881template<class Type, template<class> class PatchField, class GeoMesh>
883(
884 const word& keyword,
885 Ostream& os
886) const
887{
888 os.beginBlock(keyword);
889 this->writeEntries(os);
890 os.endBlock();
891
892 os.check(FUNCTION_NAME);
893}
894
895
896template<class Type, template<class> class PatchField, class GeoMesh>
898(
899 Ostream& os
900) const
901{
902 for (const auto& pfld : *this)
903 {
904 os.beginBlock(pfld.patch().name());
905 os << pfld;
906 os.endBlock();
907 }
908}
909
910
911template<class Type, template<class> class PatchField, class GeoMesh>
913(
914) const
915{
916 // Dummy op - template specialisations provide logic (usually call
917 // to checkConsistency)
918 return true;
919}
920
921
922// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
923
924template<class Type, template<class> class PatchField, class GeoMesh>
926(
929{
931}
932
933
934template<class Type, template<class> class PatchField, class GeoMesh>
936(
939{
941}
942
943
944template<class Type, template<class> class PatchField, class GeoMesh>
946(
947 const Type& val
949{
951}
952
953
954template<class Type, template<class> class PatchField, class GeoMesh>
956(
958)
959{
961}
962
963
964template<class Type, template<class> class PatchField, class GeoMesh>
966(
968)
969{
970 forAll(*this, patchi)
972 this->operator[](patchi) == bf[patchi];
973 }
974}
975
976
977template<class Type, template<class> class PatchField, class GeoMesh>
979(
980 const FieldField<PatchField, Type>& bf
981)
982{
983 forAll(*this, patchi)
985 this->operator[](patchi) == bf[patchi];
986 }
987}
988
989
990template<class Type, template<class> class PatchField, class GeoMesh>
992(
993 const Type& val
994)
995{
996 for (auto& pfld : *this)
997 {
998 pfld == val;
999 }
1000}
1001
1002
1003// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
1004
1005template<class Type, template<class> class PatchField, class GeoMesh>
1006Foam::Ostream& Foam::operator<<
1007(
1008 Ostream& os,
1010)
1011{
1013 return os;
1014}
1015
1016
1017// ************************************************************************* //
if(patchID !=-1)
labelList patchIDs
static const this_type & null() noexcept
Return a null DimensionedField (reference to a nullObject).
static int localBoundaryConsistency() noexcept
Get flag for local boundary consistency checks.
Definition Field.H:133
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
tmp< FieldField< PatchField, Type > > clone() const
Definition FieldField.C:182
const Type & operator[](const labelPair &index) const
Definition FieldField.C:381
void operator=(const FieldField< Field, Type > &)
Copy assignment.
Definition FieldField.C:395
friend Ostream & operator(Ostream &, const FieldField< PatchField, Type > &)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
Generic GeometricBoundaryField class.
void evaluate_if(const UnaryPredicate &pred, const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions for patch fields matching the given predicate. Uses specified or default...
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set.
DimensionedField< Type, GeoMesh > Internal
The internal field type associated with the boundary fields.
void evaluateSelected(const labelUList &patchIDs)
Evaluate boundary conditions for selected patches. Uses non-blocking comms.
wordList types() const
Return a list of the patch types.
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set.
void evaluateCoupled(const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions on coupled patches of the given type, using specified or default comms.
void writeEntries(Ostream &os) const
Write dictionary entries of the individual boundary fields.
bool check() const
Helper: check if field has been evaluated. See instantiations.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
void updateCoeffs()
Update the boundary condition coefficients.
GeoMesh::BoundaryMesh BoundaryMesh
The boundary mesh type for the boundary fields.
GeometricBoundaryField(const BoundaryMesh &bmesh)
Construct from a BoundaryMesh, setting patches later.
void readField(const Internal &iField, const dictionary &dict)
Read the boundary field.
void evaluateLocal(const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions after change in local values. Uses specified or default comms.
tmp< GeometricBoundaryField > boundaryInternalField() const
Return boundary field of values neighbouring the boundary.
void evaluate(const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions for each patch field. Uses specified or default comms.
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const PatchField< Type > * set(const label i) const
Definition PtrList.H:171
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
commsTypes
Communications types.
Definition UPstream.H:81
@ scheduled
"scheduled" (MPI standard) : (MPI_Send, MPI_Recv)
Definition UPstream.H:83
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
label size() const noexcept
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
const keyType & keyword() const noexcept
Return keyword.
Definition entry.H:231
virtual const dictionary * dictPtr() const noexcept
Return pointer to dictionary, if entry is a dictionary, otherwise return nullptr.
Definition entry.H:290
bool isLiteral() const noexcept
The keyType is treated as literal, not as pattern.
Definition keyTypeI.H:91
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
labelList patchIds
patchWriters resize(patchIds.size())
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
surface1 clear()
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
#define PoutInFunction
Report using Foam::Pout with FUNCTION_NAME prefix.
#define FUNCTION_NAME
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition facAverage.C:40
type
Types of root.
Definition Roots.H:53
List< word > wordList
List of word.
Definition fileName.H:60
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition lduSchedule.H:46
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
UPtrList< const LduInterfaceField< Type > > LduInterfaceFieldPtrsList
Store lists of LduInterfaceField as a UPtrList.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299