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-2024 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 "commSchedule.H"
28 #include "globalMeshData.H"
29 #include "emptyPolyPatch.H"
30 #include "processorPolyPatch.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template<class Type, template<class> class PatchField, class GeoMesh>
36 (
38  const dictionary& dict
39 )
40 {
41  // Clear the boundary field if already initialised
42  this->clear();
43 
44  this->setSize(bmesh_.size());
45 
47  {
49  }
50 
51  // Construct a list of entry pointers for each patch
52  UPtrList<const entry> patchEntries(this->size());
53 
54  // 1. Explicit patch names
56  {
57  if (iter().isDict() && !iter().keyword().isPattern())
58  {
59  const label patchi = bmesh_.findIndex(iter().keyword());
60 
61  if (patchi != -1)
62  {
63  patchEntries.set(patchi, &(iter()));
64  }
65  }
66  }
67 
68  // 2. Patch-groups
69  // Note: This is done in reverse order of the entries in the dictionary,
70  // so that it is consistent with dictionary wildcard behaviour.
71  if (dict.size())
72  {
73  for
74  (
76  iter != dict.rend();
77  ++iter
78  )
79  {
80  if (iter().isDict() && !iter().keyword().isPattern())
81  {
82  const labelList patchIDs =
83  bmesh_.findIndices(wordRe(iter().keyword()), true);
84 
85  forAll(patchIDs, i)
86  {
87  const label patchi = patchIDs[i];
88 
89  if (!patchEntries.set(patchi))
90  {
91  patchEntries.set(patchi, &(iter()));
92  }
93  }
94  }
95  }
96  }
97 
98  // 3. Empty patches
99  // These take precedence over wildcards
100  // (... apparently. Why not wedges and/or other constraints too?)
101  forAll(bmesh_, patchi)
102  {
103  if (!patchEntries.set(patchi))
104  {
105  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
106  {
107  patchEntries.set(patchi, NullObjectPtr<entry>());
108  }
109  }
110  }
111 
112  // 4. Wildcards
113  forAll(bmesh_, patchi)
114  {
115  if (!patchEntries.set(patchi))
116  {
117  const entry* ePtr =
118  dict.lookupEntryPtr(bmesh_[patchi].name(), false, true);
119 
120  if (ePtr)
121  {
122  patchEntries.set(patchi, ePtr);
123  }
124  }
125  }
126 
127  // Construct all the patches in order
128  forAll(bmesh_, patchi)
129  {
130  if (patchEntries.set(patchi) && !isNull(patchEntries(patchi)))
131  {
132  this->set
133  (
134  patchi,
136  (
137  bmesh_[patchi],
138  field,
139  patchEntries[patchi].dict()
140  )
141  );
142  }
143  else if (patchEntries.set(patchi) && isNull(patchEntries[patchi]))
144  {
145  this->set
146  (
147  patchi,
149  (
150  emptyPolyPatch::typeName,
151  bmesh_[patchi],
152  field
153  )
154  );
155  }
156  else
157  {
159  << "Cannot find patchField entry for "
160  << bmesh_[patchi].name() << exit(FatalIOError);
161  }
162  }
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167 
168 template<class Type, template<class> class PatchField, class GeoMesh>
170 (
171  const BoundaryMesh& bmesh
172 )
173 :
174  FieldField<PatchField, Type>(bmesh.size()),
175  bmesh_(bmesh)
176 {}
177 
178 
179 template<class Type, template<class> class PatchField, class GeoMesh>
181 (
182  const BoundaryMesh& bmesh,
183  const DimensionedField<Type, GeoMesh>& field,
184  const word& patchFieldType
185 )
186 :
187  FieldField<PatchField, Type>(bmesh.size()),
188  bmesh_(bmesh)
189 {
191  {
192  InfoInFunction << endl;
193  }
194 
195  forAll(bmesh_, patchi)
196  {
197  this->set
198  (
199  patchi,
201  (
203  bmesh_[patchi],
204  field
205  )
206  );
207  }
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 wordList& patchFieldTypes,
217  const wordList& constraintTypes
218 )
219 :
220  FieldField<PatchField, Type>(bmesh.size()),
221  bmesh_(bmesh)
222 {
224  {
225  InfoInFunction << endl;
226  }
227 
228  if
229  (
230  patchFieldTypes.size() != this->size()
231  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
232  )
233  {
235  << "Incorrect number of patch type specifications given" << nl
236  << " Number of patches in mesh = " << bmesh.size()
237  << " number of patch type specifications = "
238  << patchFieldTypes.size()
239  << abort(FatalError);
240  }
241 
242  if (constraintTypes.size())
243  {
244  forAll(bmesh_, patchi)
245  {
246  this->set
247  (
248  patchi,
250  (
251  patchFieldTypes[patchi],
252  constraintTypes[patchi],
253  bmesh_[patchi],
254  field
255  )
256  );
257  }
258  }
259  else
260  {
261  forAll(bmesh_, patchi)
262  {
263  this->set
264  (
265  patchi,
267  (
268  patchFieldTypes[patchi],
269  bmesh_[patchi],
270  field
271  )
272  );
273  }
274  }
275 }
276 
277 
278 template<class Type, template<class> class PatchField, class GeoMesh>
280 (
281  const BoundaryMesh& bmesh,
282  const DimensionedField<Type, GeoMesh>& field,
283  const PtrList<PatchField<Type>>& ptfl
284 )
285 :
286  FieldField<PatchField, Type>(bmesh.size()),
287  bmesh_(bmesh)
288 {
289  if (GeometricField<Type, PatchField, GeoMesh>::debug)
290  {
291  InfoInFunction << endl;
292  }
293 
294  forAll(bmesh_, patchi)
295  {
296  this->set(patchi, ptfl[patchi].clone(field));
297  }
298 }
299 
300 
301 template<class Type, template<class> class PatchField, class GeoMesh>
303 (
304  const DimensionedField<Type, GeoMesh>& field,
305  const GeometricBoundaryField<Type, PatchField, GeoMesh>& btf
306 )
307 :
308  FieldField<PatchField, Type>(btf.size()),
309  bmesh_(btf.bmesh_)
310 {
311  if (GeometricField<Type, PatchField, GeoMesh>::debug)
312  {
313  InfoInFunction << endl;
314  }
315 
316  forAll(bmesh_, patchi)
317  {
318  this->set(patchi, btf[patchi].clone(field));
319  }
320 }
321 
322 
323 template<class Type, template<class> class PatchField, class GeoMesh>
325 (
326  const BoundaryMesh& bmesh,
327  const DimensionedField<Type, GeoMesh>& field,
328  const dictionary& dict
329 )
330 :
331  FieldField<PatchField, Type>(bmesh.size()),
332  bmesh_(bmesh)
333 {
334  readField(field, dict);
335 }
336 
337 
338 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
339 
340 template<class Type, template<class> class PatchField, class GeoMesh>
342 {
344  {
345  InfoInFunction << endl;
346  }
347 
348  forAll(*this, patchi)
349  {
350  this->operator[](patchi).updateCoeffs();
351  }
352 }
353 
354 
355 template<class Type, template<class> class PatchField, class GeoMesh>
357 {
359  {
360  InfoInFunction << endl;
361  }
362 
363  if
364  (
365  Pstream::defaultCommsType == Pstream::commsTypes::blocking
366  || Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
367  )
368  {
369  label nReq = Pstream::nRequests();
370 
371  forAll(*this, patchi)
372  {
373  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
374  }
375 
376  // Block for any outstanding requests
377  if
378  (
379  Pstream::parRun()
380  && Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
381  )
382  {
383  Pstream::waitRequests(nReq);
384  }
385 
386  forAll(*this, patchi)
387  {
388  this->operator[](patchi).evaluate(Pstream::defaultCommsType);
389  }
390  }
391  else if (Pstream::defaultCommsType == Pstream::commsTypes::scheduled)
392  {
393  const lduSchedule& patchSchedule =
394  bmesh_.mesh().globalData().patchSchedule();
395 
396  forAll(patchSchedule, patchEvali)
397  {
398  if (patchSchedule[patchEvali].init)
399  {
400  this->operator[](patchSchedule[patchEvali].patch)
401  .initEvaluate(Pstream::commsTypes::scheduled);
402  }
403  else
404  {
405  this->operator[](patchSchedule[patchEvali].patch)
406  .evaluate(Pstream::commsTypes::scheduled);
407  }
408  }
409  }
410  else
411  {
413  << "Unsupported communications type "
414  << Pstream::commsTypeNames[Pstream::defaultCommsType]
415  << exit(FatalError);
416  }
417 }
418 
419 
420 template<class Type, template<class> class PatchField, class GeoMesh>
423 {
424  const FieldField<PatchField, Type>& pff = *this;
425 
426  wordList Types(pff.size());
427 
428  forAll(pff, patchi)
429  {
430  Types[patchi] = pff[patchi].type();
431  }
432 
433  return Types;
434 }
435 
436 
437 template<class Type, template<class> class PatchField, class GeoMesh>
440 boundaryInternalField() const
441 {
443  (
445  (
447  *this
448  )
449  );
450 
452 
453  forAll(*this, patchi)
454  {
455  result[patchi] == this->operator[](patchi).patchInternalField();
456  }
457 
458  return tresult;
459 }
460 
461 
462 template<class Type, template<class> class PatchField, class GeoMesh>
466 {
468  (
470  (
472  *this
473  )
474  );
475 
477 
478  if
479  (
480  Pstream::defaultCommsType == Pstream::commsTypes::blocking
481  || Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
482  )
483  {
484  const label nReq = Pstream::nRequests();
485 
486  forAll(*this, patchi)
487  {
488  if (this->operator[](patchi).coupled())
489  {
490  this->operator[](patchi)
491  .initPatchNeighbourField(Pstream::defaultCommsType);
492  }
493  }
494 
495  // Block for any outstanding requests
496  if
497  (
498  Pstream::parRun()
499  && Pstream::defaultCommsType == Pstream::commsTypes::nonBlocking
500  )
501  {
502  Pstream::waitRequests(nReq);
503  }
504 
505  forAll(*this, patchi)
506  {
507  if (this->operator[](patchi).coupled())
508  {
509  result[patchi] =
510  this->operator[](patchi)
511  .patchNeighbourField(Pstream::defaultCommsType);
512  }
513  }
514  }
515  else if (Pstream::defaultCommsType == Pstream::commsTypes::scheduled)
516  {
517  const lduSchedule& patchSchedule =
518  bmesh_.mesh().globalData().patchSchedule();
519 
520  forAll(patchSchedule, patchEvali)
521  {
522  if (this->operator[](patchSchedule[patchEvali].patch).coupled())
523  {
524  if (patchSchedule[patchEvali].init)
525  {
526  this->operator[](patchSchedule[patchEvali].patch)
527  .initPatchNeighbourField(Pstream::defaultCommsType);
528  }
529  else
530  {
531  result[patchSchedule[patchEvali].patch] =
532  this->operator[](patchSchedule[patchEvali].patch)
533  .patchNeighbourField(Pstream::defaultCommsType);
534  }
535  }
536  }
537  }
538  else
539  {
541  << "Unsupported communications type "
542  << Pstream::commsTypeNames[Pstream::defaultCommsType]
543  << exit(FatalError);
544  }
545 
546  return tresult;
547 }
548 
549 
550 template<class Type, template<class> class PatchField, class GeoMesh>
553 {
554  LduInterfaceFieldPtrsList<Type> interfaces(this->size());
555 
556  forAll(interfaces, patchi)
557  {
558  if (isA<LduInterfaceField<Type>>(this->operator[](patchi)))
559  {
560  interfaces.set
561  (
562  patchi,
564  (
565  this->operator[](patchi)
566  )
567  );
568  }
569  }
570 
571  return interfaces;
572 }
573 
574 
575 template<class Type, template<class> class PatchField, class GeoMesh>
578 scalarInterfaces() const
579 {
580  lduInterfaceFieldPtrsList interfaces(this->size());
581 
582  forAll(interfaces, patchi)
583  {
584  if (isA<lduInterfaceField>(this->operator[](patchi)))
585  {
586  interfaces.set
587  (
588  patchi,
589  &refCast<const lduInterfaceField>
590  (
591  this->operator[](patchi)
592  )
593  );
594  }
595  }
596 
597  return interfaces;
598 }
599 
600 
601 template<class Type, template<class> class PatchField, class GeoMesh>
603 (
605 )
606 {
607  // Reset the number of patches in case the decomposition changed
608  this->setSize(btf.size());
609 
610  const polyBoundaryMesh& pbm = bmesh_.mesh().mesh().boundaryMesh();
611 
612  forAll(*this, patchi)
613  {
614  // Construct new processor patch fields in case the decomposition
615  // changed
616  if (isA<processorPolyPatch>(pbm[patchi]))
617  {
618  this->set
619  (
620  patchi,
621  btf[patchi].clone
622  (
623  bmesh_[patchi],
624  this->operator[](0).internalField()
625  )
626  );
627  }
628  else
629  {
630  this->operator[](patchi).reset(btf[patchi]);
631  }
632  }
633 }
634 
635 
636 
637 template<class Type, template<class> class PatchField, class GeoMesh>
639 (
640  const word& keyword,
641  Ostream& os
642 ) const
643 {
644  os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
645 
646  forAll(*this, patchi)
647  {
648  os << indent << this->operator[](patchi).patch().name() << nl
649  << indent << token::BEGIN_BLOCK << nl
650  << incrIndent << this->operator[](patchi) << decrIndent
651  << indent << token::END_BLOCK << endl;
652  }
653 
654  os << decrIndent << token::END_BLOCK << endl;
655 
656  // Check state of IOstream
657  os.check
658  (
659  "GeometricBoundaryField<Type, PatchField, GeoMesh>::"
660  "writeEntry(const word& keyword, Ostream& os) const"
661  );
662 }
663 
664 
665 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
666 
667 template<class Type, template<class> class PatchField, class GeoMesh>
669 (
671 )
672 {
674 }
675 
676 
677 template<class Type, template<class> class PatchField, class GeoMesh>
679 (
680  GeometricBoundaryField<Type, PatchField, GeoMesh>&& bf
681 )
682 {
684 }
685 
686 
687 template<class Type, template<class> class PatchField, class GeoMesh>
689 (
690  const FieldField<PatchField, Type>& ptff
691 )
692 {
694 }
695 
696 
697 template<class Type, template<class> class PatchField, class GeoMesh>
698 template<template<class> class OtherPatchField>
700 (
702 )
703 {
705 }
706 
707 
708 template<class Type, template<class> class PatchField, class GeoMesh>
710 operator=
711 (
712  const Type& t
713 )
714 {
716 }
717 
718 
719 template<class Type, template<class> class PatchField, class GeoMesh>
721 (
723 )
724 {
725  forAll(*this, patchi)
726  {
727  this->operator[](patchi) == bf[patchi];
728  }
729 }
730 
731 
732 template<class Type, template<class> class PatchField, class GeoMesh>
734 operator==
735 (
736  const FieldField<PatchField, Type>& ptff
737 )
738 {
739  forAll(*this, patchi)
740  {
741  this->operator[](patchi) == ptff[patchi];
742  }
743 }
744 
745 
746 template<class Type, template<class> class PatchField, class GeoMesh>
747 template<template<class> class OtherPatchField>
749 operator==
750 (
752 )
753 {
754  forAll(*this, patchi)
755  {
756  this->operator[](patchi) == ptff[patchi];
757  }
758 }
759 
760 
761 template<class Type, template<class> class PatchField, class GeoMesh>
763 (
764  const Type& t
765 )
766 {
767  forAll(*this, patchi)
768  {
769  this->operator[](patchi) == t;
770  }
771 }
772 
773 
774 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
775 
776 template<class Type, template<class> class PatchField, class GeoMesh>
777 Foam::Ostream& Foam::operator<<
778 (
779  Ostream& os,
781 )
782 {
783  os << static_cast<const FieldField<PatchField, Type>&>(bf);
784  return os;
785 }
786 
787 
788 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic field type.
Definition: FieldField.H:77
void operator=(const FieldField< PatchField, Type > &)
Definition: FieldField.C:290
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricBoundaryField class.
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those.
void readField(const Internal &field, const dictionary &dict)
Read the boundary field.
wordList types() const
Return a list of the patch field types.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
void evaluate()
Evaluate boundary conditions.
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
void updateCoeffs()
Update the boundary condition coefficients.
GeoMesh::BoundaryMesh BoundaryMesh
Type of boundary mesh on which this boundary is instantiated.
tmp< GeometricBoundaryField > boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
void reset(const GeometricBoundaryField< Type, PatchField, GeoMesh > &)
Reset the boundary field contents to the given field.
GeometricBoundaryField(const BoundaryMesh &)
Construct from a BoundaryMesh.
Generic GeometricField class.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An STL-conforming const_reverse_iterator.
Definition: UILList.H:302
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Foam::polyBoundaryMesh.
const polyMesh & mesh() const
Return the mesh reference.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:77
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
tUEqn clear()
#define InfoInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:166
word patchFieldType(const PatchField &pf)
T clone(const T &t)
Definition: List.H:55
IOerror FatalIOError
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
points setSize(newPointi)
dictionary dict