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-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
30#include "globalMeshData.H"
31#include "cyclicPolyPatch.H"
32#include "emptyPolyPatch.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36template<class Type, template<class> class PatchField, class GeoMesh>
38(
40 const dictionary& dict
41)
42{
47
48 // Clear the boundary field if already initialised
49 this->clear();
50
51 this->resize(bmesh_.size());
52
53 label nUnset = this->size();
54
55 // 1. Handle explicit patch names. Note that there can be only one explicit
56 // patch name since is key of dictionary.
57
58 for (const entry& dEntry : dict)
59 {
60 if (dEntry.isDict() && dEntry.keyword().isLiteral())
61 {
62 const label patchi = bmesh_.findPatchID(dEntry.keyword());
63
64 if (patchi != -1)
65 {
66 this->set
67 (
68 patchi,
69 PatchField<Type>::New
70 (
71 bmesh_[patchi],
72 field,
73 dEntry.dict()
74 )
75 );
76 nUnset--;
77 }
78 }
79 }
80
81 if (nUnset == 0)
82 {
83 return;
84 }
85
86
87 // 2. Patch-groups. (using non-wild card entries of dictionaries)
88 // (patchnames already matched above)
89 // Note: in reverse order of entries in the dictionary (last
90 // patchGroups wins). This is so it is consistent with dictionary wildcard
91 // behaviour
92 for (auto iter = dict.crbegin(); iter != dict.crend(); ++iter)
93 {
94 const entry& dEntry = *iter;
95
96 if (dEntry.isDict() && dEntry.keyword().isLiteral())
97 {
98 const labelList patchIds =
99 bmesh_.indices(dEntry.keyword(), true); // use patchGroups
100
101 for (const label patchi : patchIds)
102 {
103 if (!this->set(patchi))
104 {
105 this->set
106 (
107 patchi,
108 PatchField<Type>::New
109 (
110 bmesh_[patchi],
111 field,
112 dEntry.dict()
113 )
114 );
115 }
116 }
118 }
119
120
121 // 3. Wildcard patch overrides
122 forAll(bmesh_, patchi)
123 {
124 if (!this->set(patchi))
126 if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
127 {
128 this->set
129 (
130 patchi,
131 PatchField<Type>::New
132 (
133 emptyPolyPatch::typeName,
134 bmesh_[patchi],
135 field
136 )
137 );
138 }
139 else
140 {
141 bool found = dict.found(bmesh_[patchi].name());
142
143 if (found)
144 {
145 this->set
146 (
147 patchi,
148 PatchField<Type>::New
150 bmesh_[patchi],
151 field,
152 dict.subDict(bmesh_[patchi].name())
153 )
154 );
155 }
156 }
157 }
158 }
159
161 // Check for any unset patches
162 forAll(bmesh_, patchi)
163 {
164 if (!this->set(patchi))
165 {
166 if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
169 << "Cannot find patchField entry for cyclic "
170 << bmesh_[patchi].name() << endl
171 << "Is your field uptodate with split cyclics?" << endl
172 << "Run foamUpgradeCyclics to convert mesh and fields"
173 << " to split cyclics." << exit(FatalIOError);
175 else
176 {
178 << "Cannot find patchField entry for "
179 << bmesh_[patchi].name() << exit(FatalIOError);
181 }
182 }
183}
185
186// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
187
188template<class Type, template<class> class PatchField, class GeoMesh>
190(
191 const BoundaryMesh& bmesh
192)
193:
194 FieldField<PatchField, Type>(bmesh.size()),
195 bmesh_(bmesh)
196{}
197
198
199template<class Type, template<class> class PatchField, class GeoMesh>
201(
202 const BoundaryMesh& bmesh,
204 const word& patchFieldType
205)
207 FieldField<PatchField, Type>(bmesh.size()),
208 bmesh_(bmesh)
209{
214
215 forAll(bmesh_, patchi)
216 {
217 this->set
218 (
219 patchi,
220 PatchField<Type>::New
221 (
222 patchFieldType,
223 bmesh_[patchi],
224 field
225 )
226 );
227 }
228}
229
230
231template<class Type, template<class> class PatchField, class GeoMesh>
233(
234 const BoundaryMesh& bmesh,
236 const wordList& patchFieldTypes,
237 const wordList& constraintTypes
238)
239:
240 FieldField<PatchField, Type>(bmesh.size()),
241 bmesh_(bmesh)
242{
247
248 if
249 (
250 patchFieldTypes.size() != this->size()
251 || (constraintTypes.size() && (constraintTypes.size() != this->size()))
252 )
253 {
255 << "Incorrect number of patch type specifications given" << nl
256 << " Number of patches in mesh = " << bmesh.size()
257 << " number of patch type specifications = "
258 << patchFieldTypes.size()
259 << abort(FatalError);
260 }
261
262 if (constraintTypes.size())
263 {
264 forAll(bmesh_, patchi)
265 {
266 this->set
267 (
268 patchi,
269 PatchField<Type>::New
270 (
271 patchFieldTypes[patchi],
272 constraintTypes[patchi],
273 bmesh_[patchi],
274 field
275 )
276 );
277 }
278 }
279 else
280 {
281 forAll(bmesh_, patchi)
282 {
283 this->set
284 (
285 patchi,
286 PatchField<Type>::New
287 (
288 patchFieldTypes[patchi],
289 bmesh_[patchi],
290 field
291 )
292 );
293 }
294 }
295}
296
297
298template<class Type, template<class> class PatchField, class GeoMesh>
300(
301 const BoundaryMesh& bmesh,
303 const PtrList<PatchField<Type>>& ptfl
304)
305:
306 FieldField<PatchField, Type>(bmesh.size()),
307 bmesh_(bmesh)
308{
313
314 forAll(bmesh_, patchi)
315 {
316 this->set(patchi, ptfl[patchi].clone(field));
317 }
318}
319
320
321template<class Type, template<class> class PatchField, class GeoMesh>
323(
326)
327:
328 FieldField<PatchField, Type>(btf.size()),
329 bmesh_(btf.bmesh_)
330{
335
336 forAll(bmesh_, patchi)
337 {
338 this->set(patchi, btf[patchi].clone(field));
339 }
340}
341
342
343template<class Type, template<class> class PatchField, class GeoMesh>
345(
348 const labelList& patchIDs,
349 const word& patchFieldType
350)
351:
352 FieldField<PatchField, Type>(btf.size()),
353 bmesh_(btf.bmesh_)
354{
359
360 for (const label patchi : patchIDs)
361 {
362 this->set
363 (
364 patchi,
365 PatchField<Type>::New
366 (
367 patchFieldType,
368 bmesh_[patchi],
369 field
370 )
371 );
372 }
373
374 forAll(bmesh_, patchi)
375 {
376 if (!this->set(patchi))
377 {
378 this->set(patchi, btf[patchi].clone(field));
379 }
380 }
381}
382
383
384template<class Type, template<class> class PatchField, class GeoMesh>
386(
388)
389:
390 FieldField<PatchField, Type>(btf),
391 bmesh_(btf.bmesh_)
392{
397}
398
399
400template<class Type, template<class> class PatchField, class GeoMesh>
402(
403 const BoundaryMesh& bmesh,
405 const dictionary& dict
406)
407:
408 FieldField<PatchField, Type>(bmesh.size()),
409 bmesh_(bmesh)
410{
412}
413
414
415// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
416
417template<class Type, template<class> class PatchField, class GeoMesh>
419{
424
425 for (auto& pfld : *this)
426 {
427 pfld.updateCoeffs();
428 }
429}
430
431
432template<class Type, template<class> class PatchField, class GeoMesh>
434{
439
441
442 if
443 (
446 )
447 {
448 const label startOfRequests = UPstream::nRequests();
449
450 for (auto& pfld : *this)
451 {
452 pfld.initEvaluate(commsType);
453 }
454
455 // Wait for outstanding requests
456 if
457 (
460 )
461 {
462 UPstream::waitRequests(startOfRequests);
463 }
464
465 for (auto& pfld : *this)
466 {
467 pfld.evaluate(commsType);
468 }
469 }
470 else if (commsType == UPstream::commsTypes::scheduled)
471 {
472 const lduSchedule& patchSchedule =
473 bmesh_.mesh().globalData().patchSchedule();
474
475 for (const auto& schedEval : patchSchedule)
476 {
477 const label patchi = schedEval.patch;
478 auto& pfld = (*this)[patchi];
479
480 if (schedEval.init)
481 {
482 pfld.initEvaluate(commsType);
483 }
484 else
485 {
486 pfld.evaluate(commsType);
487 }
488 }
489 }
490 else
491 {
493 << "Unsupported communications type "
494 << UPstream::commsTypeNames[commsType]
495 << exit(FatalError);
496 }
497}
498
499
500template<class Type, template<class> class PatchField, class GeoMesh>
501template<class CoupledPatchType>
503{
508
510
511 if
512 (
515 )
516 {
517 const label startOfRequests = UPstream::nRequests();
518
519 for (auto& pfld : *this)
520 {
521 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
522
523 if (cpp && cpp->coupled())
524 {
525 pfld.initEvaluate(commsType);
526 }
527 }
528
529 // Wait for outstanding requests
530 if
531 (
534 )
535 {
536 UPstream::waitRequests(startOfRequests);
537 }
538
539 for (auto& pfld : *this)
540 {
541 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
542
543 if (cpp && cpp->coupled())
544 {
545 pfld.evaluate(commsType);
546 }
547 }
548 }
549 else if (commsType == UPstream::commsTypes::scheduled)
550 {
551 const lduSchedule& patchSchedule =
552 bmesh_.mesh().globalData().patchSchedule();
553
554 for (const auto& schedEval : patchSchedule)
555 {
556 const label patchi = schedEval.patch;
557 auto& pfld = (*this)[patchi];
558
559 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
560
561 if (cpp && cpp->coupled())
562 {
563 if (schedEval.init)
564 {
565 pfld.initEvaluate(commsType);
566 }
567 else
568 {
569 pfld.evaluate(commsType);
570 }
571 }
572 }
573 }
574 else
575 {
577 << "Unsupported communications type "
578 << UPstream::commsTypeNames[commsType]
579 << exit(FatalError);
580 }
581}
582
583
584template<class Type, template<class> class PatchField, class GeoMesh>
587{
588 const FieldField<PatchField, Type>& pff = *this;
589
590 wordList list(pff.size());
591
592 forAll(pff, patchi)
593 {
594 list[patchi] = pff[patchi].type();
595 }
596
597 return list;
598}
599
600
601template<class Type, template<class> class PatchField, class GeoMesh>
605{
607
608 forAll(result, patchi)
609 {
610 result[patchi] == this->operator[](patchi).patchInternalField();
611 }
612
613 return result;
614}
615
616
617template<class Type, template<class> class PatchField, class GeoMesh>
620{
621 LduInterfaceFieldPtrsList<Type> list(this->size());
622
623 forAll(list, patchi)
624 {
625 const auto* lduPtr =
626 isA<LduInterfaceField<Type>>(this->operator[](patchi));
627
628 if (lduPtr)
629 {
630 list.set(patchi, lduPtr);
631 }
632 }
633
634 return list;
635}
636
637
638template<class Type, template<class> class PatchField, class GeoMesh>
641scalarInterfaces() const
642{
643 lduInterfaceFieldPtrsList list(this->size());
644
645 forAll(list, patchi)
646 {
647 const auto* lduPtr =
648 isA<lduInterfaceField>(this->operator[](patchi));
649
650 if (lduPtr)
651 {
652 list.set(patchi, lduPtr);
653 }
654 }
655
656 return list;
657}
658
659
660template<class Type, template<class> class PatchField, class GeoMesh>
662(
663 const word& keyword,
664 Ostream& os
665) const
666{
667 os.beginBlock(keyword);
668 this->writeEntries(os);
669 os.endBlock();
670
672}
673
674
675template<class Type, template<class> class PatchField, class GeoMesh>
677(
678 Ostream& os
679) const
680{
681 for (const auto& pfld : *this)
682 {
683 os.beginBlock(pfld.patch().name());
684 os << pfld;
685 os.endBlock();
686 }
687}
688
689
690// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
691
692template<class Type, template<class> class PatchField, class GeoMesh>
694(
696)
697{
699}
700
701
702template<class Type, template<class> class PatchField, class GeoMesh>
704(
706)
707{
709}
710
711
712template<class Type, template<class> class PatchField, class GeoMesh>
714(
715 const Type& val
716)
717{
719}
720
721
722template<class Type, template<class> class PatchField, class GeoMesh>
724(
726)
727{
728 forAll(*this, patchi)
729 {
730 this->operator[](patchi) == bf[patchi];
731 }
732}
733
734
735template<class Type, template<class> class PatchField, class GeoMesh>
737(
739)
740{
741 forAll(*this, patchi)
742 {
743 this->operator[](patchi) == bf[patchi];
744 }
745}
746
747
748template<class Type, template<class> class PatchField, class GeoMesh>
750(
751 const Type& val
752)
753{
754 forAll(*this, patchi)
755 {
756 this->operator[](patchi) == val;
757 }
758}
759
760
761// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
762
763template<class Type, template<class> class PatchField, class GeoMesh>
764Foam::Ostream& Foam::operator<<
765(
766 Ostream& os,
768)
769{
770 os << static_cast<const FieldField<PatchField, Type>&>(bf);
771 return os;
772}
773
774
775// ************************************************************************* //
bool found
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
tmp< FieldField< PatchField, Type > > clone() const
Clone.
Definition: FieldField.C:189
void operator=(const FieldField< Field, Type > &)
Copy assignment.
Definition: FieldField.C:301
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Generic GeometricBoundaryField class.
lduInterfaceFieldPtrsList scalarInterfaces() const
void readField(const DimensionedField< Type, GeoMesh > &field, const dictionary &dict)
Read the boundary field.
wordList types() const
Return a list of the patch types.
void evaluateCoupled()
Evaluate boundary conditions on a subset of coupled patches.
void evaluate()
Evaluate boundary conditions.
LduInterfaceFieldPtrsList< Type > interfaces() const
void writeEntries(Ostream &os) const
Write dictionary entries of the individual boundary fields.
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 boundaryInternalField() const
Return boundary field of values neighbouring the boundary.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:105
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
@ nonBlocking
"nonBlocking"
static const Enum< commsTypes > commsTypeNames
Names of the communication types.
Definition: UPstream.H:74
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:90
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:100
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
const T * set(const label i) const
Definition: UPtrList.H:248
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
virtual bool isDict() const noexcept
Return true if this entry is a dictionary.
Definition: entry.H:233
const keyType & keyword() const noexcept
Return keyword.
Definition: entry.H:195
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
bool isLiteral() const noexcept
The keyType is treated as literal, not as pattern.
Definition: keyTypeI.H:98
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
labelList patchIds
patchWriters resize(patchIds.size())
patchWriters clear()
#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)
#define FUNCTION_NAME
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
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
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333