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 "GeometricFieldFwd.H"
28 #include "emptyPolyPatch.H"
29 #include "processorPolyPatch.H"
30 #include "commSchedule.H"
31 #include "globalMeshData.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class Type, class GeoMesh, template<class> class PrimitiveField>
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  // Construct a list of entry pointers for each patch
53  UPtrList<const entry> patchEntries(this->size());
54 
55  // 1. Explicit patch names
57  {
58  if (iter().isDict() && !iter().keyword().isPattern())
59  {
60  const label patchi = bmesh_.findIndex(iter().keyword());
61 
62  if (patchi != -1)
63  {
64  patchEntries.set(patchi, &(iter()));
65  }
66  }
67  }
68 
69  // 2. Patch-groups
70  // Note: This is done in reverse order of the entries in the dictionary,
71  // so that it is consistent with dictionary wildcard behaviour.
72  if (dict.size())
73  {
74  for
75  (
77  iter != dict.rend();
78  ++iter
79  )
80  {
81  if (iter().isDict() && !iter().keyword().isPattern())
82  {
83  const labelList patchIDs =
84  bmesh_.findIndices(wordRe(iter().keyword()), true);
85 
86  forAll(patchIDs, i)
87  {
88  const label patchi = patchIDs[i];
89 
90  if (!patchEntries.set(patchi))
91  {
92  patchEntries.set(patchi, &(iter()));
93  }
94  }
95  }
96  }
97  }
98 
99  // 3. Empty patches
100  // These take precedence over wildcards
101  // (... apparently. Why not wedges and/or other constraints too?)
102  forAll(bmesh_, patchi)
103  {
104  if (!patchEntries.set(patchi))
105  {
106  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
107  {
108  patchEntries.set(patchi, NullObjectPtr<entry>());
109  }
110  }
111  }
112 
113  // 4. Wildcards
114  forAll(bmesh_, patchi)
115  {
116  if (!patchEntries.set(patchi))
117  {
118  const entry* ePtr =
119  dict.lookupEntryPtr(bmesh_[patchi].name(), false, true);
120 
121  if (ePtr)
122  {
123  patchEntries.set(patchi, ePtr);
124  }
125  }
126  }
127 
128  // Construct all the patches in order
129  forAll(bmesh_, patchi)
130  {
131  if (patchEntries.set(patchi) && !isNull(patchEntries(patchi)))
132  {
133  this->set
134  (
135  patchi,
136  Patch::New
137  (
138  bmesh_[patchi],
139  field,
140  patchEntries[patchi].dict()
141  )
142  );
143  }
144  else if (patchEntries.set(patchi) && isNull(patchEntries[patchi]))
145  {
146  this->set
147  (
148  patchi,
149  Patch::New
150  (
151  emptyPolyPatch::typeName,
152  bmesh_[patchi],
153  field
154  )
155  );
156  }
157  else
158  {
160  << "Cannot find patchField entry for "
161  << bmesh_[patchi].name() << exit(FatalIOError);
162  }
163  }
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168 
169 template<class Type, class GeoMesh, template<class> class PrimitiveField>
172 (
173  const BoundaryMesh& bmesh
174 )
175 :
177  bmesh_(bmesh)
178 {}
179 
180 
181 template<class Type, class GeoMesh, template<class> class PrimitiveField>
184 (
185  const BoundaryMesh& bmesh,
188 )
189 :
191  bmesh_(bmesh)
192 {
194  {
195  InfoInFunction << endl;
196  }
197 
198  forAll(bmesh_, patchi)
199  {
200  this->set
201  (
202  patchi,
204  (
206  bmesh_[patchi],
207  field
208  )
209  );
210  }
211 }
212 
213 
214 template<class Type, class GeoMesh, template<class> class PrimitiveField>
217 (
218  const BoundaryMesh& bmesh,
220  const wordList& patchFieldTypes,
221  const wordList& constraintTypes
222 )
223 :
225  bmesh_(bmesh)
226 {
228  {
229  InfoInFunction << endl;
230  }
231 
232  if
233  (
234  patchFieldTypes.size() != this->size()
235  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
236  )
237  {
239  << "Incorrect number of patch type specifications given" << nl
240  << " Number of patches in mesh = " << bmesh.size()
241  << " number of patch type specifications = "
242  << patchFieldTypes.size()
243  << abort(FatalError);
244  }
245 
246  if (constraintTypes.size())
247  {
248  forAll(bmesh_, patchi)
249  {
250  this->set
251  (
252  patchi,
253  Patch::New
254  (
255  patchFieldTypes[patchi],
256  constraintTypes[patchi],
257  bmesh_[patchi],
258  field
259  )
260  );
261  }
262  }
263  else
264  {
265  forAll(bmesh_, patchi)
266  {
267  this->set
268  (
269  patchi,
270  Patch::New
271  (
272  patchFieldTypes[patchi],
273  bmesh_[patchi],
274  field
275  )
276  );
277  }
278  }
279 }
280 
281 
282 template<class Type, class GeoMesh, template<class> class PrimitiveField>
285 (
286  const BoundaryMesh& bmesh,
288  const PtrList<Patch>& ptfl
289 )
290 :
291  FieldField<GeoMesh::template PatchField, Type>(bmesh.size()),
292  bmesh_(bmesh)
293 {
295  {
296  InfoInFunction << endl;
297  }
298 
299  forAll(bmesh_, patchi)
300  {
301  this->set(patchi, ptfl[patchi].clone(field));
302  }
303 }
304 
305 
306 template<class Type, class GeoMesh, template<class> class PrimitiveField>
309 (
312 )
313 :
314  FieldField<GeoMesh::template PatchField, Type>(btf.size()),
315  bmesh_(btf.bmesh_)
316 {
318  {
319  InfoInFunction << endl;
320  }
321 
322  forAll(bmesh_, patchi)
323  {
324  this->set(patchi, btf[patchi].clone(field));
325  }
326 }
327 
328 
329 template<class Type, class GeoMesh, template<class> class PrimitiveField>
332 (
333  const BoundaryMesh& bmesh,
335  const dictionary& dict
336 )
337 :
338  FieldField<GeoMesh::template PatchField, Type>(bmesh.size()),
339  bmesh_(bmesh)
340 {
341  readField(field, dict);
342 }
343 
344 
345 template<class Type, class GeoMesh, template<class> class PrimitiveField>
346 template<template<class> class PrimitiveField2>
349 (
352 )
353 :
354  FieldField<GeoMesh::template PatchField, Type>(btf.size()),
355  bmesh_(btf.bmesh_)
356 {
358  {
359  InfoInFunction << endl;
360  }
361 
362  forAll(bmesh_, patchi)
363  {
364  this->set(patchi, btf[patchi].clone(field));
365  }
366 }
367 
368 
369 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
370 
371 template<class Type, class GeoMesh, template<class> class PrimitiveField>
373 {
375  {
376  InfoInFunction << endl;
377  }
378 
379  forAll(*this, patchi)
380  {
381  this->operator[](patchi).updateCoeffs();
382  }
383 }
384 
385 
386 template<class Type, class GeoMesh, template<class> class PrimitiveField>
388 {
390  {
391  InfoInFunction << endl;
392  }
393 
394  if
395  (
398  )
399  {
400  label nReq = Pstream::nRequests();
401 
402  forAll(*this, patchi)
403  {
404  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
405  }
406 
407  // Block for any outstanding requests
408  if
409  (
412  )
413  {
414  Pstream::waitRequests(nReq);
415  }
416 
417  forAll(*this, patchi)
418  {
419  this->operator[](patchi).evaluate(Pstream::defaultCommsType);
420  }
421  }
423  {
424  const lduSchedule& patchSchedule =
425  bmesh_.mesh().globalData().patchSchedule();
426 
427  forAll(patchSchedule, patchEvali)
428  {
429  if (patchSchedule[patchEvali].init)
430  {
431  this->operator[](patchSchedule[patchEvali].patch)
432  .initEvaluate(Pstream::commsTypes::scheduled);
433  }
434  else
435  {
436  this->operator[](patchSchedule[patchEvali].patch)
438  }
439  }
440  }
441  else
442  {
444  << "Unsupported communications type "
446  << exit(FatalError);
447  }
448 }
449 
450 
451 template<class Type, class GeoMesh, template<class> class PrimitiveField>
454 {
456 
457  wordList Types(pff.size());
458 
459  forAll(pff, patchi)
460  {
461  Types[patchi] = pff[patchi].type();
462  }
463 
464  return Types;
465 }
466 
467 
468 template<class Type, class GeoMesh, template<class> class PrimitiveField>
471 boundaryInternalField() const
472 {
473  typedef
475  resultType;
476 
477  tmp<resultType> tresult(new resultType(Internal::null(), *this));
478  resultType& result = tresult.ref();
479 
480  forAll(*this, patchi)
481  {
482  result[patchi] == this->operator[](patchi).patchInternalField();
483  }
484 
485  return tresult;
486 }
487 
488 
489 template<class Type, class GeoMesh, template<class> class PrimitiveField>
493 {
494  typedef
496  resultType;
497 
498  tmp<resultType> tresult(new resultType(Internal::null(), *this));
499  resultType& result = tresult.ref();
500 
501  if
502  (
505  )
506  {
507  const label nReq = Pstream::nRequests();
508 
509  forAll(*this, patchi)
510  {
511  if (this->operator[](patchi).coupled())
512  {
513  this->operator[](patchi)
514  .initPatchNeighbourField(Pstream::defaultCommsType);
515  }
516  }
517 
518  // Block for any outstanding requests
519  if
520  (
523  )
524  {
525  Pstream::waitRequests(nReq);
526  }
527 
528  forAll(*this, patchi)
529  {
530  if (this->operator[](patchi).coupled())
531  {
532  result[patchi] =
533  this->operator[](patchi)
534  .patchNeighbourField(Pstream::defaultCommsType);
535  }
536  }
537  }
539  {
540  const lduSchedule& patchSchedule =
541  bmesh_.mesh().globalData().patchSchedule();
542 
543  forAll(patchSchedule, patchEvali)
544  {
545  if (this->operator[](patchSchedule[patchEvali].patch).coupled())
546  {
547  if (patchSchedule[patchEvali].init)
548  {
549  this->operator[](patchSchedule[patchEvali].patch)
550  .initPatchNeighbourField(Pstream::defaultCommsType);
551  }
552  else
553  {
554  result[patchSchedule[patchEvali].patch] =
555  this->operator[](patchSchedule[patchEvali].patch)
556  .patchNeighbourField(Pstream::defaultCommsType);
557  }
558  }
559  }
560  }
561  else
562  {
564  << "Unsupported communications type "
566  << exit(FatalError);
567  }
568 
569  return tresult;
570 }
571 
572 
573 template<class Type, class GeoMesh, template<class> class PrimitiveField>
576 {
577  LduInterfaceFieldPtrsList<Type> interfaces(this->size());
578 
579  forAll(interfaces, patchi)
580  {
581  if (isA<LduInterfaceField<Type>>(this->operator[](patchi)))
582  {
583  interfaces.set
584  (
585  patchi,
587  (
588  this->operator[](patchi)
589  )
590  );
591  }
592  }
593 
594  return interfaces;
595 }
596 
597 
598 template<class Type, class GeoMesh, template<class> class PrimitiveField>
601 scalarInterfaces() const
602 {
603  lduInterfaceFieldPtrsList interfaces(this->size());
604 
605  forAll(interfaces, patchi)
606  {
607  if (isA<lduInterfaceField>(this->operator[](patchi)))
608  {
609  interfaces.set
610  (
611  patchi,
612  &refCast<const lduInterfaceField>
613  (
614  this->operator[](patchi)
615  )
616  );
617  }
618  }
619 
620  return interfaces;
621 }
622 
623 
624 template<class Type, class GeoMesh, template<class> class PrimitiveField>
626 (
628 )
629 {
630  // Reset the number of patches in case the decomposition changed
631  this->setSize(btf.size());
632 
633  const polyBoundaryMesh& pbm = bmesh_.mesh().mesh().boundaryMesh();
634 
635  forAll(*this, patchi)
636  {
637  // Construct new processor patch fields in case the decomposition
638  // changed
639  if (isA<processorPolyPatch>(pbm[patchi]))
640  {
641  this->set
642  (
643  patchi,
644  btf[patchi].clone
645  (
646  bmesh_[patchi],
647  this->operator[](0).internalField()
648  )
649  );
650  }
651  else
652  {
653  this->operator[](patchi).reset(btf[patchi]);
654  }
655  }
656 }
657 
658 
659 
660 template<class Type, class GeoMesh, template<class> class PrimitiveField>
662 (
663  const word& keyword,
664  Ostream& os
665 ) const
666 {
667  os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
668 
669  forAll(*this, patchi)
670  {
671  os << indent << this->operator[](patchi).patch().name() << nl
672  << indent << token::BEGIN_BLOCK << nl
673  << incrIndent << this->operator[](patchi) << decrIndent
674  << indent << token::END_BLOCK << endl;
675  }
676 
677  os << decrIndent << token::END_BLOCK << endl;
678 
679  // Check state of IOstream
680  os.check
681  (
682  "GeometricBoundaryField<Type, GeoMesh, PrimitiveField>::"
683  "writeEntry(const word& keyword, Ostream& os) const"
684  );
685 }
686 
687 
688 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
689 
690 template<class Type, class GeoMesh, template<class> class PrimitiveField>
693 {
695 }
696 
697 
698 template<class Type, class GeoMesh, template<class> class PrimitiveField>
700 operator=(GeometricBoundaryField<Type, GeoMesh, PrimitiveField>&& bf)
701 {
703 }
704 
705 
706 template<class Type, class GeoMesh, template<class> class PrimitiveField>
709 {
711 }
712 
713 
714 template<class Type, class GeoMesh, template<class> class PrimitiveField>
715 template<template<class> class OtherPatchField>
718 {
720 }
721 
722 
723 template<class Type, class GeoMesh, template<class> class PrimitiveField>
725 operator=(const Type& t)
726 {
728 }
729 
730 
731 template<class Type, class GeoMesh, template<class> class PrimitiveField>
734 {
735  forAll(*this, patchi)
736  {
737  this->operator[](patchi) == bf[patchi];
738  }
739 }
740 
741 
742 template<class Type, class GeoMesh, template<class> class PrimitiveField>
745 {
746  forAll(*this, patchi)
747  {
748  this->operator[](patchi) == ptff[patchi];
749  }
750 }
751 
752 
753 template<class Type, class GeoMesh, template<class> class PrimitiveField>
754 template<template<class> class OtherPatchField>
757 {
758  forAll(*this, patchi)
759  {
760  this->operator[](patchi) == ptff[patchi];
761  }
762 }
763 
764 
765 template<class Type, class GeoMesh, template<class> class PrimitiveField>
767 operator==(const Type& t)
768 {
769  forAll(*this, patchi)
770  {
771  this->operator[](patchi) == t;
772  }
773 }
774 
775 
776 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
777 
778 template<class Type, class GeoMesh, template<class> class PrimitiveField>
779 Foam::Ostream& Foam::operator<<
780 (
781  Ostream& os,
783 )
784 {
785  os <<
786  static_cast<const FieldField<GeoMesh::template PatchField, Type>&>(bf);
787  return os;
788 }
789 
790 
791 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
tmp< FieldField< GeoMesh::template PatchField, Type > > clone() const
Clone.
Definition: FieldField.C:188
void operator=(const FieldField< Field, Type > &)
Definition: FieldField.C:290
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricBoundaryField class.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those.
friend class GeometricBoundaryField
Declare friendship with other geometric boundary fields.
void readField(const Internal &field, const dictionary &dict)
Read the boundary field.
wordList types() const
Return a list of the patch field types.
void evaluate()
Evaluate boundary conditions.
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
void operator==(const GeometricBoundaryField &)
Forced assignment to.
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 &)
Reset the boundary field contents to the given field.
void operator=(const GeometricBoundaryField &)
Assignment operator.
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:297
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
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
An STL-conforming const_reverse_iterator.
Definition: UILList.H:302
static const NamedEnum< commsTypes, 3 > commsTypeNames
Definition: UPstream.H:71
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Foam::polyBoundaryMesh.
const polyMesh & mesh() const
Return the mesh reference.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
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:197
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
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:242
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:134
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:171
word patchFieldType(const PatchField &pf)
T clone(const T &t)
Definition: List.H:55
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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:228
static const char nl
Definition: Ostream.H:267
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