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-2020 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 "emptyPolyPatch.H"
27 #include "commSchedule.H"
28 #include "globalMeshData.H"
29 #include "cyclicPolyPatch.H"
30 
31 template<class Type, template<class> class PatchField, class GeoMesh>
34 (
36  const dictionary& dict
37 )
38 {
39  // Clear the boundary field if already initialised
40  this->clear();
41 
42  this->setSize(bmesh_.size());
43 
44  if (debug)
45  {
47  }
48 
49 
50  label nUnset = this->size();
51 
52  // 1. Handle explicit patch names. Note that there can be only one explicit
53  // patch name since is key of dictionary.
54  forAllConstIter(dictionary, dict, iter)
55  {
56  if (iter().isDict() && !iter().keyword().isPattern())
57  {
58  label patchi = bmesh_.findPatchID(iter().keyword());
59 
60  if (patchi != -1)
61  {
62  this->set
63  (
64  patchi,
66  (
67  bmesh_[patchi],
68  field,
69  iter().dict()
70  )
71  );
72  nUnset--;
73  }
74  }
75  }
76 
77  if (nUnset == 0)
78  {
79  return;
80  }
81 
82 
83  // 2. Patch-groups. (using non-wild card entries of dictionaries)
84  // (patchnames already matched above)
85  // Note: in reverse order of entries in the dictionary (last
86  // patchGroups wins). This is so is consistent with dictionary wildcard
87  // behaviour
88  if (dict.size())
89  {
90  for
91  (
93  iter != dict.rend();
94  ++iter
95  )
96  {
97  const entry& e = iter();
98 
99  if (e.isDict() && !e.keyword().isPattern())
100  {
101  const labelList patchIDs = bmesh_.findIndices
102  (
103  wordRe(e.keyword()),
104  true // use patchGroups
105  );
106 
107  forAll(patchIDs, i)
108  {
109  label patchi = patchIDs[i];
110 
111  if (!this->set(patchi))
112  {
113  this->set
114  (
115  patchi,
117  (
118  bmesh_[patchi],
119  field,
120  e.dict()
121  )
122  );
123  }
124  }
125  }
126  }
127  }
128 
129 
130  // 3. Wildcard patch overrides
131  forAll(bmesh_, patchi)
132  {
133  if (!this->set(patchi))
134  {
135  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
136  {
137  this->set
138  (
139  patchi,
141  (
142  emptyPolyPatch::typeName,
143  bmesh_[patchi],
144  field
145  )
146  );
147  }
148  else
149  {
150  bool found = dict.found(bmesh_[patchi].name());
151 
152  if (found)
153  {
154  this->set
155  (
156  patchi,
158  (
159  bmesh_[patchi],
160  field,
161  dict.subDict(bmesh_[patchi].name())
162  )
163  );
164  }
165  }
166  }
167  }
168 
169  // Check for any unset patches
170  forAll(bmesh_, patchi)
171  {
172  if (!this->set(patchi))
173  {
174  if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
175  {
177  (
178  dict
179  ) << "Cannot find patchField entry for cyclic "
180  << bmesh_[patchi].name() << endl << exit(FatalIOError);
181  }
182  else
183  {
185  (
186  dict
187  ) << "Cannot find patchField entry for "
188  << bmesh_[patchi].name() << exit(FatalIOError);
189  }
190  }
191  }
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
196 
197 template<class Type, template<class> class PatchField, class GeoMesh>
199 Boundary
200 (
201  const BoundaryMesh& bmesh
202 )
203 :
204  FieldField<PatchField, Type>(bmesh.size()),
205  bmesh_(bmesh)
206 {}
207 
208 
209 template<class Type, template<class> class PatchField, class GeoMesh>
211 Boundary
212 (
213  const BoundaryMesh& bmesh,
214  const DimensionedField<Type, GeoMesh>& field,
215  const word& patchFieldType
216 )
217 :
218  FieldField<PatchField, Type>(bmesh.size()),
219  bmesh_(bmesh)
220 {
221  if (debug)
222  {
223  InfoInFunction << endl;
224  }
225 
226  forAll(bmesh_, patchi)
227  {
228  this->set
229  (
230  patchi,
232  (
233  patchFieldType,
234  bmesh_[patchi],
235  field
236  )
237  );
238  }
239 }
240 
241 
242 template<class Type, template<class> class PatchField, class GeoMesh>
244 Boundary
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 {
255  if (debug)
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 Boundary
313 (
314  const BoundaryMesh& bmesh,
315  const DimensionedField<Type, GeoMesh>& field,
316  const PtrList<PatchField<Type>>& ptfl
317 )
318 :
319  FieldField<PatchField, Type>(bmesh.size()),
320  bmesh_(bmesh)
321 {
322  if (debug)
323  {
324  InfoInFunction << endl;
325  }
326 
327  forAll(bmesh_, patchi)
328  {
329  this->set(patchi, ptfl[patchi].clone(field));
330  }
331 }
332 
333 
334 template<class Type, template<class> class PatchField, class GeoMesh>
336 Boundary
337 (
338  const DimensionedField<Type, GeoMesh>& field,
340  Boundary& btf
341 )
342 :
344  bmesh_(btf.bmesh_)
345 {
346  if (debug)
347  {
348  InfoInFunction << endl;
349  }
350 
351  forAll(bmesh_, patchi)
352  {
353  this->set(patchi, btf[patchi].clone(field));
354  }
355 }
356 
357 
358 template<class Type, template<class> class PatchField, class GeoMesh>
360 Boundary
361 (
362  const BoundaryMesh& bmesh,
363  const DimensionedField<Type, GeoMesh>& field,
364  const dictionary& dict
365 )
366 :
367  FieldField<PatchField, Type>(bmesh.size()),
368  bmesh_(bmesh)
369 {
370  readField(field, dict);
371 }
372 
373 
374 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
375 
376 template<class Type, template<class> class PatchField, class GeoMesh>
379 {
380  if (debug)
381  {
382  InfoInFunction << endl;
383  }
384 
385  forAll(*this, patchi)
386  {
387  this->operator[](patchi).updateCoeffs();
388  }
389 }
390 
391 
392 template<class Type, template<class> class PatchField, class GeoMesh>
395 {
396  if (debug)
397  {
398  InfoInFunction << endl;
399  }
400 
401  if
402  (
405  )
406  {
407  label nReq = Pstream::nRequests();
408 
409  forAll(*this, patchi)
410  {
411  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
412  }
413 
414  // Block for any outstanding requests
415  if
416  (
419  )
420  {
421  Pstream::waitRequests(nReq);
422  }
423 
424  forAll(*this, patchi)
425  {
427  }
428  }
430  {
431  const lduSchedule& patchSchedule =
432  bmesh_.mesh().globalData().patchSchedule();
433 
434  forAll(patchSchedule, patchEvali)
435  {
436  if (patchSchedule[patchEvali].init)
437  {
438  this->operator[](patchSchedule[patchEvali].patch)
439  .initEvaluate(Pstream::commsTypes::scheduled);
440  }
441  else
442  {
443  this->operator[](patchSchedule[patchEvali].patch)
445  }
446  }
447  }
448  else
449  {
451  << "Unsupported communications type "
453  << exit(FatalError);
454  }
455 }
456 
457 
458 template<class Type, template<class> class PatchField, class GeoMesh>
461 types() const
462 {
463  const FieldField<PatchField, Type>& pff = *this;
464 
465  wordList Types(pff.size());
466 
467  forAll(pff, patchi)
468  {
469  Types[patchi] = pff[patchi].type();
470  }
471 
472  return Types;
473 }
474 
475 
476 template<class Type, template<class> class PatchField, class GeoMesh>
480 {
482  BoundaryInternalField(*this);
483 
484  forAll(BoundaryInternalField, patchi)
485  {
486  BoundaryInternalField[patchi] ==
487  this->operator[](patchi).patchInternalField();
488  }
489 
490  return BoundaryInternalField;
491 }
492 
493 
494 template<class Type, template<class> class PatchField, class GeoMesh>
497 interfaces() const
498 {
500 
501  forAll(interfaces, patchi)
502  {
503  if (isA<LduInterfaceField<Type>>(this->operator[](patchi)))
504  {
505  interfaces.set
506  (
507  patchi,
509  (
510  this->operator[](patchi)
511  )
512  );
513  }
514  }
515 
516  return interfaces;
517 }
518 
519 
520 template<class Type, template<class> class PatchField, class GeoMesh>
524 {
526 
527  forAll(interfaces, patchi)
528  {
529  if (isA<lduInterfaceField>(this->operator[](patchi)))
530  {
531  interfaces.set
532  (
533  patchi,
534  &refCast<const lduInterfaceField>
535  (
536  this->operator[](patchi)
537  )
538  );
539  }
540  }
541 
542  return interfaces;
543 }
544 
545 
546 template<class Type, template<class> class PatchField, class GeoMesh>
548 writeEntry(const word& keyword, Ostream& os) const
549 {
550  os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
551 
552  forAll(*this, patchi)
553  {
554  os << indent << this->operator[](patchi).patch().name() << nl
555  << indent << token::BEGIN_BLOCK << nl
556  << incrIndent << this->operator[](patchi) << decrIndent
557  << indent << token::END_BLOCK << endl;
558  }
559 
560  os << decrIndent << token::END_BLOCK << endl;
561 
562  // Check state of IOstream
563  os.check
564  (
565  "GeometricField<Type, PatchField, GeoMesh>::Boundary::"
566  "writeEntry(const word& keyword, Ostream& os) const"
567  );
568 }
569 
570 
571 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
572 
573 template<class Type, template<class> class PatchField, class GeoMesh>
575 operator=
576 (
578  Boundary& bf
579 )
580 {
582 }
583 
584 
585 template<class Type, template<class> class PatchField, class GeoMesh>
587 operator=
588 (
590  Boundary&& bf
591 )
592 {
594 }
595 
596 
597 template<class Type, template<class> class PatchField, class GeoMesh>
599 operator=
600 (
601  const FieldField<PatchField, Type>& ptff
602 )
603 {
605 }
606 
607 
608 template<class Type, template<class> class PatchField, class GeoMesh>
610 operator=
611 (
612  const Type& t
613 )
614 {
616 }
617 
618 
619 template<class Type, template<class> class PatchField, class GeoMesh>
621 operator==
622 (
624  Boundary& bf
625 )
626 {
627  forAll((*this), patchi)
628  {
629  this->operator[](patchi) == bf[patchi];
630  }
631 }
632 
633 
634 template<class Type, template<class> class PatchField, class GeoMesh>
636 operator==
637 (
638  const FieldField<PatchField, Type>& ptff
639 )
640 {
641  forAll((*this), patchi)
642  {
643  this->operator[](patchi) == ptff[patchi];
644  }
645 }
646 
647 
648 template<class Type, template<class> class PatchField, class GeoMesh>
650 operator==
651 (
652  const Type& t
653 )
654 {
655  forAll((*this), patchi)
656  {
657  this->operator[](patchi) == t;
658  }
659 }
660 
661 
662 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
663 
664 template<class Type, template<class> class PatchField, class GeoMesh>
665 Foam::Ostream& Foam::operator<<
666 (
667  Ostream& os,
669  Boundary& bf
670 )
671 {
672  os << static_cast<const FieldField<PatchField, Type>&>(bf);
673  return os;
674 }
675 
676 
677 // ************************************************************************* //
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:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Boundary(const BoundaryMesh &)
Construct from a BoundaryMesh.
tUEqn clear()
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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
const keyType & keyword() const
Return keyword.
Definition: entry.H:123
GeoMesh::BoundaryMesh BoundaryMesh
Type of boundary mesh on which this.
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
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
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void evaluate()
Evaluate boundary conditions.
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.
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
Generic GeometricField class.
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
wordList types() const
Return a list of the patch field types.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
void readField(const Internal &field, const dictionary &dict)
Read the boundary field.
void operator=(const FieldField< Field, Type > &)
Definition: FieldField.C:290
void updateCoeffs()
Update the boundary condition coefficients.
points setSize(newPointi)
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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
virtual bool isDict() const
Return true if this entry is a dictionary.
Definition: entry.H:156
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
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
static const NamedEnum< commsTypes, 3 > commsTypeNames
Definition: UPstream.H:71
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
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
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
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
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 scalarInterfaces() const
Return a list of pointers for each patch field with only those.
friend Ostream & operator(Ostream &, const FieldField< PatchField, Type > &)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
Boundary boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
bool found
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
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.