GeometricBoundaryField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "GeometricBoundaryField.H"
27 #include "emptyPolyPatch.H"
28 #include "commSchedule.H"
29 #include "globalMeshData.H"
30 #include "cyclicPolyPatch.H"
31 #include "processorPolyPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class Type, template<class> class PatchField, class GeoMesh>
37 (
39  const dictionary& dict
40 )
41 {
42  // Clear the boundary field if already initialised
43  this->clear();
44 
45  this->setSize(bmesh_.size());
46 
48  {
50  }
51 
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  forAllConstIter(dictionary, dict, iter)
58  {
59  if (iter().isDict() && !iter().keyword().isPattern())
60  {
61  label patchi = bmesh_.findPatchID(iter().keyword());
62 
63  if (patchi != -1)
64  {
65  this->set
66  (
67  patchi,
69  (
70  bmesh_[patchi],
71  field,
72  iter().dict()
73  )
74  );
75  nUnset--;
76  }
77  }
78  }
79 
80  if (nUnset == 0)
81  {
82  return;
83  }
84 
85 
86  // 2. Patch-groups. (using non-wild card entries of dictionaries)
87  // (patchnames already matched above)
88  // Note: in reverse order of entries in the dictionary (last
89  // patchGroups wins). This is so is consistent with dictionary wildcard
90  // behaviour
91  if (dict.size())
92  {
93  for
94  (
96  iter != dict.rend();
97  ++iter
98  )
99  {
100  const entry& e = iter();
101 
102  if (e.isDict() && !e.keyword().isPattern())
103  {
104  const labelList patchIDs = bmesh_.findIndices
105  (
106  wordRe(e.keyword()),
107  true // use patchGroups
108  );
109 
110  forAll(patchIDs, i)
111  {
112  label patchi = patchIDs[i];
113 
114  if (!this->set(patchi))
115  {
116  this->set
117  (
118  patchi,
120  (
121  bmesh_[patchi],
122  field,
123  e.dict()
124  )
125  );
126  }
127  }
128  }
129  }
130  }
131 
132 
133  // 3. Wildcard patch overrides
134  forAll(bmesh_, patchi)
135  {
136  if (!this->set(patchi))
137  {
138  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
139  {
140  this->set
141  (
142  patchi,
144  (
145  emptyPolyPatch::typeName,
146  bmesh_[patchi],
147  field
148  )
149  );
150  }
151  else
152  {
153  bool found = dict.found(bmesh_[patchi].name());
154 
155  if (found)
156  {
157  this->set
158  (
159  patchi,
161  (
162  bmesh_[patchi],
163  field,
164  dict.subDict(bmesh_[patchi].name())
165  )
166  );
167  }
168  }
169  }
170  }
171 
172  // Check for any unset patches
173  forAll(bmesh_, patchi)
174  {
175  if (!this->set(patchi))
176  {
177  if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
178  {
180  (
181  dict
182  ) << "Cannot find patchField entry for cyclic "
183  << bmesh_[patchi].name() << endl << exit(FatalIOError);
184  }
185  else
186  {
188  (
189  dict
190  ) << "Cannot find patchField entry for "
191  << bmesh_[patchi].name() << exit(FatalIOError);
192  }
193  }
194  }
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
199 
200 template<class Type, template<class> class PatchField, class GeoMesh>
202 (
203  const BoundaryMesh& bmesh
204 )
205 :
206  FieldField<PatchField, Type>(bmesh.size()),
207  bmesh_(bmesh)
208 {}
209 
210 
211 template<class Type, template<class> class PatchField, class GeoMesh>
213 (
214  const BoundaryMesh& bmesh,
215  const DimensionedField<Type, GeoMesh>& field,
216  const word& patchFieldType
217 )
218 :
219  FieldField<PatchField, Type>(bmesh.size()),
220  bmesh_(bmesh)
221 {
223  {
224  InfoInFunction << endl;
225  }
226 
227  forAll(bmesh_, patchi)
228  {
229  this->set
230  (
231  patchi,
233  (
234  patchFieldType,
235  bmesh_[patchi],
236  field
237  )
238  );
239  }
240 }
241 
242 
243 template<class Type, template<class> class PatchField, class GeoMesh>
245 (
246  const BoundaryMesh& bmesh,
247  const DimensionedField<Type, GeoMesh>& field,
248  const wordList& patchFieldTypes,
249  const wordList& constraintTypes
250 )
251 :
252  FieldField<PatchField, Type>(bmesh.size()),
253  bmesh_(bmesh)
254 {
256  {
257  InfoInFunction << endl;
258  }
259 
260  if
261  (
262  patchFieldTypes.size() != this->size()
263  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
264  )
265  {
267  << "Incorrect number of patch type specifications given" << nl
268  << " Number of patches in mesh = " << bmesh.size()
269  << " number of patch type specifications = "
270  << patchFieldTypes.size()
271  << abort(FatalError);
272  }
273 
274  if (constraintTypes.size())
275  {
276  forAll(bmesh_, patchi)
277  {
278  this->set
279  (
280  patchi,
282  (
283  patchFieldTypes[patchi],
284  constraintTypes[patchi],
285  bmesh_[patchi],
286  field
287  )
288  );
289  }
290  }
291  else
292  {
293  forAll(bmesh_, patchi)
294  {
295  this->set
296  (
297  patchi,
299  (
300  patchFieldTypes[patchi],
301  bmesh_[patchi],
302  field
303  )
304  );
305  }
306  }
307 }
308 
309 
310 template<class Type, template<class> class PatchField, class GeoMesh>
312 (
313  const BoundaryMesh& bmesh,
314  const DimensionedField<Type, GeoMesh>& field,
315  const PtrList<PatchField<Type>>& ptfl
316 )
317 :
318  FieldField<PatchField, Type>(bmesh.size()),
319  bmesh_(bmesh)
320 {
322  {
323  InfoInFunction << endl;
324  }
325 
326  forAll(bmesh_, patchi)
327  {
328  this->set(patchi, ptfl[patchi].clone(field));
329  }
330 }
331 
332 
333 template<class Type, template<class> class PatchField, class GeoMesh>
335 (
336  const DimensionedField<Type, GeoMesh>& field,
338 )
339 :
341  bmesh_(btf.bmesh_)
342 {
344  {
345  InfoInFunction << endl;
346  }
347 
348  forAll(bmesh_, patchi)
349  {
350  this->set(patchi, btf[patchi].clone(field));
351  }
352 }
353 
354 
355 template<class Type, template<class> class PatchField, class GeoMesh>
357 (
358  const BoundaryMesh& bmesh,
359  const DimensionedField<Type, GeoMesh>& field,
360  const dictionary& dict
361 )
362 :
363  FieldField<PatchField, Type>(bmesh.size()),
364  bmesh_(bmesh)
365 {
366  readField(field, dict);
367 }
368 
369 
370 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
371 
372 template<class Type, template<class> class PatchField, class GeoMesh>
374 {
376  {
377  InfoInFunction << endl;
378  }
379 
380  forAll(*this, patchi)
381  {
382  this->operator[](patchi).updateCoeffs();
383  }
384 }
385 
386 
387 template<class Type, template<class> class PatchField, class GeoMesh>
389 {
391  {
392  InfoInFunction << endl;
393  }
394 
395  if
396  (
397  Pstream::defaultCommsType == Pstream::commsTypes::blocking
398  || Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
399  )
400  {
401  label nReq = Pstream::nRequests();
402 
403  forAll(*this, patchi)
404  {
405  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
406  }
407 
408  // Block for any outstanding requests
409  if
410  (
411  Pstream::parRun()
412  && Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
413  )
414  {
415  Pstream::waitRequests(nReq);
416  }
417 
418  forAll(*this, patchi)
419  {
420  this->operator[](patchi).evaluate(Pstream::defaultCommsType);
421  }
422  }
423  else if (Pstream::defaultCommsType == Pstream::commsTypes::scheduled)
424  {
425  const lduSchedule& patchSchedule =
426  bmesh_.mesh().globalData().patchSchedule();
427 
428  forAll(patchSchedule, patchEvali)
429  {
430  if (patchSchedule[patchEvali].init)
431  {
432  this->operator[](patchSchedule[patchEvali].patch)
433  .initEvaluate(Pstream::commsTypes::scheduled);
434  }
435  else
436  {
437  this->operator[](patchSchedule[patchEvali].patch)
438  .evaluate(Pstream::commsTypes::scheduled);
439  }
440  }
441  }
442  else
443  {
445  << "Unsupported communications type "
446  << Pstream::commsTypeNames[Pstream::defaultCommsType]
447  << exit(FatalError);
448  }
449 }
450 
451 
452 template<class Type, template<class> class PatchField, class GeoMesh>
455 {
456  const FieldField<PatchField, Type>& pff = *this;
457 
458  wordList Types(pff.size());
459 
460  forAll(pff, patchi)
461  {
462  Types[patchi] = pff[patchi].type();
463  }
464 
465  return Types;
466 }
467 
468 
469 template<class Type, template<class> class PatchField, class GeoMesh>
473 {
475  (
477  (
479  *this
480  )
481  );
482 
484 
485  forAll(*this, patchi)
486  {
487  result[patchi] == this->operator[](patchi).patchInternalField();
488  }
489 
490  return tresult;
491 }
492 
493 
494 template<class Type, template<class> class PatchField, class GeoMesh>
498 {
500  (
502  (
504  *this
505  )
506  );
507 
509 
510  if
511  (
512  Pstream::defaultCommsType == Pstream::commsTypes::blocking
513  || Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
514  )
515  {
516  const label nReq = Pstream::nRequests();
517 
518  forAll(*this, patchi)
519  {
520  if (this->operator[](patchi).coupled())
521  {
522  this->operator[](patchi)
523  .initPatchNeighbourField(Pstream::defaultCommsType);
524  }
525  }
526 
527  // Block for any outstanding requests
528  if
529  (
530  Pstream::parRun()
531  && Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
532  )
533  {
534  Pstream::waitRequests(nReq);
535  }
536 
537  forAll(*this, patchi)
538  {
539  if (this->operator[](patchi).coupled())
540  {
541  result[patchi] =
542  this->operator[](patchi)
543  .patchNeighbourField(Pstream::defaultCommsType);
544  }
545  }
546  }
547  else if (Pstream::defaultCommsType == Pstream::commsTypes::scheduled)
548  {
549  const lduSchedule& patchSchedule =
550  bmesh_.mesh().globalData().patchSchedule();
551 
552  forAll(patchSchedule, patchEvali)
553  {
554  if (this->operator[](patchSchedule[patchEvali].patch).coupled())
555  {
556  if (patchSchedule[patchEvali].init)
557  {
558  this->operator[](patchSchedule[patchEvali].patch)
559  .initPatchNeighbourField(Pstream::defaultCommsType);
560  }
561  else
562  {
563  result[patchSchedule[patchEvali].patch] =
564  this->operator[](patchSchedule[patchEvali].patch)
565  .patchNeighbourField(Pstream::defaultCommsType);
566  }
567  }
568  }
569  }
570  else
571  {
573  << "Unsupported communications type "
574  << Pstream::commsTypeNames[Pstream::defaultCommsType]
575  << exit(FatalError);
576  }
577 
578  return tresult;
579 }
580 
581 
582 template<class Type, template<class> class PatchField, class GeoMesh>
585 {
586  LduInterfaceFieldPtrsList<Type> interfaces(this->size());
587 
588  forAll(interfaces, patchi)
589  {
590  if (isA<LduInterfaceField<Type>>(this->operator[](patchi)))
591  {
592  interfaces.set
593  (
594  patchi,
596  (
597  this->operator[](patchi)
598  )
599  );
600  }
601  }
602 
603  return interfaces;
604 }
605 
606 
607 template<class Type, template<class> class PatchField, class GeoMesh>
611 {
612  lduInterfaceFieldPtrsList interfaces(this->size());
613 
614  forAll(interfaces, patchi)
615  {
616  if (isA<lduInterfaceField>(this->operator[](patchi)))
617  {
618  interfaces.set
619  (
620  patchi,
621  &refCast<const lduInterfaceField>
622  (
623  this->operator[](patchi)
624  )
625  );
626  }
627  }
628 
629  return interfaces;
630 }
631 
632 
633 template<class Type, template<class> class PatchField, class GeoMesh>
635 (
636  const DimensionedField<Type, GeoMesh>& field,
638 )
639 {
640  // Reset the number of patches in case the decomposition changed
641  this->setSize(btf.size());
642 
643  const polyBoundaryMesh& pbm = field.mesh()().boundaryMesh();
644 
645  forAll(*this, patchi)
646  {
647  // Construct new processor patch fields in case the decomposition
648  // changed
649  if (isA<processorPolyPatch>(pbm[patchi]))
650  {
651  this->set(patchi, btf[patchi].clone(bmesh_[patchi], field));
652  }
653  else
654  {
655  this->operator[](patchi).reset(btf[patchi]);
656  }
657  }
658 }
659 
660 
661 
662 template<class Type, template<class> class PatchField, class GeoMesh>
664 (
665  const word& keyword,
666  Ostream& os
667 ) const
668 {
669  os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
670 
671  forAll(*this, patchi)
672  {
673  os << indent << this->operator[](patchi).patch().name() << nl
674  << indent << token::BEGIN_BLOCK << nl
675  << incrIndent << this->operator[](patchi) << decrIndent
676  << indent << token::END_BLOCK << endl;
677  }
678 
679  os << decrIndent << token::END_BLOCK << endl;
680 
681  // Check state of IOstream
682  os.check
683  (
684  "GeometricBoundaryField<Type, PatchField, GeoMesh>::"
685  "writeEntry(const word& keyword, Ostream& os) const"
686  );
687 }
688 
689 
690 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
691 
692 template<class Type, template<class> class PatchField, class GeoMesh>
693 void Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::operator=
694 (
696 )
697 {
699 }
700 
701 
702 template<class Type, template<class> class PatchField, class GeoMesh>
703 void Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::operator=
704 (
706 )
707 {
709 }
710 
711 
712 template<class Type, template<class> class PatchField, class GeoMesh>
713 void Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::operator=
714 (
715  const FieldField<PatchField, Type>& ptff
716 )
717 {
719 }
720 
721 
722 template<class Type, template<class> class PatchField, class GeoMesh>
723 template<template<class> class OtherPatchField>
724 void Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::operator=
725 (
727 )
728 {
730 }
731 
732 
733 template<class Type, template<class> class PatchField, class GeoMesh>
735 operator=
736 (
737  const Type& t
738 )
739 {
741 }
742 
743 
744 template<class Type, template<class> class PatchField, class GeoMesh>
745 void Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::operator==
746 (
748 )
749 {
750  forAll(*this, patchi)
751  {
752  this->operator[](patchi) == bf[patchi];
753  }
754 }
755 
756 
757 template<class Type, template<class> class PatchField, class GeoMesh>
759 operator==
760 (
761  const FieldField<PatchField, Type>& ptff
762 )
763 {
764  forAll(*this, patchi)
765  {
766  this->operator[](patchi) == ptff[patchi];
767  }
768 }
769 
770 
771 template<class Type, template<class> class PatchField, class GeoMesh>
772 template<template<class> class OtherPatchField>
774 operator==
775 (
777 )
778 {
779  forAll(*this, patchi)
780  {
781  this->operator[](patchi) == ptff[patchi];
782  }
783 }
784 
785 
786 template<class Type, template<class> class PatchField, class GeoMesh>
787 void Foam::GeometricBoundaryField<Type, PatchField, GeoMesh>::operator==
788 (
789  const Type& t
790 )
791 {
792  forAll(*this, patchi)
793  {
794  this->operator[](patchi) == t;
795  }
796 }
797 
798 
799 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
800 
801 template<class Type, template<class> class PatchField, class GeoMesh>
802 Foam::Ostream& Foam::operator<<
803 (
804  Ostream& os,
806 )
807 {
808  os << static_cast<const FieldField<PatchField, Type>&>(bf);
809  return os;
810 }
811 
812 
813 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
const_reverse_iterator rbegin() const
Definition: UILList.H:349
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tUEqn clear()
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
const keyType & keyword() const
Return keyword.
Definition: entry.H:123
tmp< GeometricBoundaryField > boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void readField(const Internal &field, const dictionary &dict)
Read the boundary field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those.
wordList types() const
Return a list of the patch field types.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
Generic GeometricField class.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
points setSize(newPointi)
void evaluate()
Evaluate boundary conditions.
pointMesh ::BoundaryMesh BoundaryMesh
Type of boundary mesh on which this boundary is instantiated.
A class for handling words, derived from string.
Definition: word.H:59
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
void reset(const Internal &, const GeometricBoundaryField< Type, PatchField, GeoMesh > &)
Reset the boundary field contents to the given field.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:74
void updateCoeffs()
Update the boundary condition coefficients.
virtual bool isDict() const
Return true if this entry is a dictionary.
Definition: entry.H:156
Foam::polyBoundaryMesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const const_reverse_iterator & rend() const
Definition: UILList.H:354
static const char nl
Definition: Ostream.H:260
bool set(const label) const
Is element set.
Definition: UPtrListI.H:78
const Mesh & mesh() const
Return mesh.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
GeometricBoundaryField(const BoundaryMesh &)
Construct from a BoundaryMesh.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
An STL-conforming const_reverse_iterator.
Definition: UILList.H:299
Generic GeometricBoundaryField class.
Definition: fvMesh.H:80
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
bool found
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
tmp< FieldField< PatchField, Type > > clone() const
Clone.
Definition: FieldField.C:188
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:65
bool isPattern() const
Should be treated as a match rather than a literal string.
Definition: keyTypeI.H:97
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.