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-2019 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  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 
170  // Check for any unset patches
171  forAll(bmesh_, patchi)
172  {
173  if (!this->set(patchi))
174  {
175  if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
176  {
178  (
179  dict
180  ) << "Cannot find patchField entry for cyclic "
181  << bmesh_[patchi].name() << endl
182  << "Is your field uptodate with split cyclics?" << endl
183  << "Run foamUpgradeCyclics to convert mesh and fields"
184  << " to split cyclics." << exit(FatalIOError);
185  }
186  else
187  {
189  (
190  dict
191  ) << "Cannot find patchField entry for "
192  << bmesh_[patchi].name() << exit(FatalIOError);
193  }
194  }
195  }
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
200 
201 template<class Type, template<class> class PatchField, class GeoMesh>
203 Boundary
204 (
205  const BoundaryMesh& bmesh
206 )
207 :
208  FieldField<PatchField, Type>(bmesh.size()),
209  bmesh_(bmesh)
210 {}
211 
212 
213 template<class Type, template<class> class PatchField, class GeoMesh>
215 Boundary
216 (
217  const BoundaryMesh& bmesh,
218  const DimensionedField<Type, GeoMesh>& field,
219  const word& patchFieldType
220 )
221 :
222  FieldField<PatchField, Type>(bmesh.size()),
223  bmesh_(bmesh)
224 {
225  if (debug)
226  {
227  InfoInFunction << endl;
228  }
229 
230  forAll(bmesh_, patchi)
231  {
232  this->set
233  (
234  patchi,
236  (
237  patchFieldType,
238  bmesh_[patchi],
239  field
240  )
241  );
242  }
243 }
244 
245 
246 template<class Type, template<class> class PatchField, class GeoMesh>
248 Boundary
249 (
250  const BoundaryMesh& bmesh,
251  const DimensionedField<Type, GeoMesh>& field,
252  const wordList& patchFieldTypes,
253  const wordList& constraintTypes
254 )
255 :
256  FieldField<PatchField, Type>(bmesh.size()),
257  bmesh_(bmesh)
258 {
259  if (debug)
260  {
261  InfoInFunction << endl;
262  }
263 
264  if
265  (
266  patchFieldTypes.size() != this->size()
267  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
268  )
269  {
271  << "Incorrect number of patch type specifications given" << nl
272  << " Number of patches in mesh = " << bmesh.size()
273  << " number of patch type specifications = "
274  << patchFieldTypes.size()
275  << abort(FatalError);
276  }
277 
278  if (constraintTypes.size())
279  {
280  forAll(bmesh_, patchi)
281  {
282  this->set
283  (
284  patchi,
286  (
287  patchFieldTypes[patchi],
288  constraintTypes[patchi],
289  bmesh_[patchi],
290  field
291  )
292  );
293  }
294  }
295  else
296  {
297  forAll(bmesh_, patchi)
298  {
299  this->set
300  (
301  patchi,
303  (
304  patchFieldTypes[patchi],
305  bmesh_[patchi],
306  field
307  )
308  );
309  }
310  }
311 }
312 
313 
314 template<class Type, template<class> class PatchField, class GeoMesh>
316 Boundary
317 (
318  const BoundaryMesh& bmesh,
319  const DimensionedField<Type, GeoMesh>& field,
320  const PtrList<PatchField<Type>>& ptfl
321 )
322 :
323  FieldField<PatchField, Type>(bmesh.size()),
324  bmesh_(bmesh)
325 {
326  if (debug)
327  {
328  InfoInFunction << endl;
329  }
330 
331  forAll(bmesh_, patchi)
332  {
333  this->set(patchi, ptfl[patchi].clone(field));
334  }
335 }
336 
337 
338 template<class Type, template<class> class PatchField, class GeoMesh>
340 Boundary
341 (
342  const DimensionedField<Type, GeoMesh>& field,
344  Boundary& btf
345 )
346 :
348  bmesh_(btf.bmesh_)
349 {
350  if (debug)
351  {
352  InfoInFunction << endl;
353  }
354 
355  forAll(bmesh_, patchi)
356  {
357  this->set(patchi, btf[patchi].clone(field));
358  }
359 }
360 
361 
362 template<class Type, template<class> class PatchField, class GeoMesh>
364 Boundary
365 (
367  Boundary& btf
368 )
369 :
371  bmesh_(btf.bmesh_)
372 {
373  if (debug)
374  {
375  InfoInFunction << endl;
376  }
377 }
378 
379 
380 template<class Type, template<class> class PatchField, class GeoMesh>
382 Boundary
383 (
385  Boundary&& btf
386 )
387 :
388  FieldField<PatchField, Type>(move(btf)),
389  bmesh_(btf.bmesh_)
390 {
391  if (debug)
392  {
393  InfoInFunction << endl;
394  }
395 }
396 
397 
398 template<class Type, template<class> class PatchField, class GeoMesh>
400 Boundary
401 (
402  const BoundaryMesh& bmesh,
403  const DimensionedField<Type, GeoMesh>& field,
404  const dictionary& dict
405 )
406 :
407  FieldField<PatchField, Type>(bmesh.size()),
408  bmesh_(bmesh)
409 {
410  readField(field, dict);
411 }
412 
413 
414 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
415 
416 template<class Type, template<class> class PatchField, class GeoMesh>
419 {
420  if (debug)
421  {
422  InfoInFunction << endl;
423  }
424 
425  forAll(*this, patchi)
426  {
427  this->operator[](patchi).updateCoeffs();
428  }
429 }
430 
431 
432 template<class Type, template<class> class PatchField, class GeoMesh>
435 {
436  if (debug)
437  {
438  InfoInFunction << endl;
439  }
440 
441  if
442  (
445  )
446  {
447  label nReq = Pstream::nRequests();
448 
449  forAll(*this, patchi)
450  {
451  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
452  }
453 
454  // Block for any outstanding requests
455  if
456  (
459  )
460  {
461  Pstream::waitRequests(nReq);
462  }
463 
464  forAll(*this, patchi)
465  {
467  }
468  }
470  {
471  const lduSchedule& patchSchedule =
472  bmesh_.mesh().globalData().patchSchedule();
473 
474  forAll(patchSchedule, patchEvali)
475  {
476  if (patchSchedule[patchEvali].init)
477  {
478  this->operator[](patchSchedule[patchEvali].patch)
479  .initEvaluate(Pstream::commsTypes::scheduled);
480  }
481  else
482  {
483  this->operator[](patchSchedule[patchEvali].patch)
485  }
486  }
487  }
488  else
489  {
491  << "Unsupported communications type "
493  << exit(FatalError);
494  }
495 }
496 
497 
498 template<class Type, template<class> class PatchField, class GeoMesh>
501 types() const
502 {
503  const FieldField<PatchField, Type>& pff = *this;
504 
505  wordList Types(pff.size());
506 
507  forAll(pff, patchi)
508  {
509  Types[patchi] = pff[patchi].type();
510  }
511 
512  return Types;
513 }
514 
515 
516 template<class Type, template<class> class PatchField, class GeoMesh>
520 {
522  BoundaryInternalField(*this);
523 
524  forAll(BoundaryInternalField, patchi)
525  {
526  BoundaryInternalField[patchi] ==
527  this->operator[](patchi).patchInternalField();
528  }
529 
530  return BoundaryInternalField;
531 }
532 
533 
534 template<class Type, template<class> class PatchField, class GeoMesh>
537 interfaces() const
538 {
540 
541  forAll(interfaces, patchi)
542  {
543  if (isA<LduInterfaceField<Type>>(this->operator[](patchi)))
544  {
545  interfaces.set
546  (
547  patchi,
549  (
550  this->operator[](patchi)
551  )
552  );
553  }
554  }
555 
556  return interfaces;
557 }
558 
559 
560 template<class Type, template<class> class PatchField, class GeoMesh>
564 {
566 
567  forAll(interfaces, patchi)
568  {
569  if (isA<lduInterfaceField>(this->operator[](patchi)))
570  {
571  interfaces.set
572  (
573  patchi,
574  &refCast<const lduInterfaceField>
575  (
576  this->operator[](patchi)
577  )
578  );
579  }
580  }
581 
582  return interfaces;
583 }
584 
585 
586 template<class Type, template<class> class PatchField, class GeoMesh>
588 writeEntry(const word& keyword, Ostream& os) const
589 {
590  os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
591 
592  forAll(*this, patchi)
593  {
594  os << indent << this->operator[](patchi).patch().name() << nl
595  << indent << token::BEGIN_BLOCK << nl
596  << incrIndent << this->operator[](patchi) << decrIndent
597  << indent << token::END_BLOCK << endl;
598  }
599 
600  os << decrIndent << token::END_BLOCK << endl;
601 
602  // Check state of IOstream
603  os.check
604  (
605  "GeometricField<Type, PatchField, GeoMesh>::Boundary::"
606  "writeEntry(const word& keyword, Ostream& os) const"
607  );
608 }
609 
610 
611 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
612 
613 template<class Type, template<class> class PatchField, class GeoMesh>
615 operator=
616 (
618  Boundary& bf
619 )
620 {
622 }
623 
624 
625 template<class Type, template<class> class PatchField, class GeoMesh>
627 operator=
628 (
630  Boundary&& bf
631 )
632 {
634 }
635 
636 
637 template<class Type, template<class> class PatchField, class GeoMesh>
639 operator=
640 (
641  const FieldField<PatchField, Type>& ptff
642 )
643 {
645 }
646 
647 
648 template<class Type, template<class> class PatchField, class GeoMesh>
650 operator=
651 (
652  const Type& t
653 )
654 {
656 }
657 
658 
659 template<class Type, template<class> class PatchField, class GeoMesh>
661 operator==
662 (
664  Boundary& bf
665 )
666 {
667  forAll((*this), patchi)
668  {
669  this->operator[](patchi) == bf[patchi];
670  }
671 }
672 
673 
674 template<class Type, template<class> class PatchField, class GeoMesh>
676 operator==
677 (
678  const FieldField<PatchField, Type>& ptff
679 )
680 {
681  forAll((*this), patchi)
682  {
683  this->operator[](patchi) == ptff[patchi];
684  }
685 }
686 
687 
688 template<class Type, template<class> class PatchField, class GeoMesh>
690 operator==
691 (
692  const Type& t
693 )
694 {
695  forAll((*this), patchi)
696  {
697  this->operator[](patchi) == t;
698  }
699 }
700 
701 
702 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
703 
704 template<class Type, template<class> class PatchField, class GeoMesh>
705 Foam::Ostream& Foam::operator<<
706 (
707  Ostream& os,
709  Boundary& bf
710 )
711 {
712  os << static_cast<const FieldField<PatchField, Type>&>(bf);
713  return os;
714 }
715 
716 
717 // ************************************************************************* //
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:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Boundary(const BoundaryMesh &)
Construct from a BoundaryMesh.
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:226
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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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.
tUEqn clear()
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:163
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:256
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 types.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
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)
volScalarField & e
Elementary charge.
Definition: createFields.H:11
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
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:53
const const_reverse_iterator & rend() const
Definition: UILList.H:354
static const char nl
Definition: Ostream.H:265
bool set(const label) const
Is element set.
Definition: UPtrListI.H:78
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
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:331
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:233
Boundary boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
bool found
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:76
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.