GeometricBoundaryField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "emptyPolyPatch.H"
30 #include "commSchedule.H"
31 #include "globalMeshData.H"
32 #include "cyclicPolyPatch.H"
33 
34 template<class Type, template<class> class PatchField, class GeoMesh>
36 (
38  const dictionary& dict
39 )
40 {
42 
43  // Clear the boundary field if already initialised
44  this->clear();
45 
46  this->setSize(bmesh_.size());
47 
48  label nUnset = this->size();
49 
50  // 1. Handle explicit patch names. Note that there can be only one explicit
51  // patch name since is key of dictionary.
52 
53  for (const entry& dEntry : dict)
54  {
55  if (dEntry.isDict() && dEntry.keyword().isLiteral())
56  {
57  const label patchi = bmesh_.findPatchID(dEntry.keyword());
58 
59  if (patchi != -1)
60  {
61  this->set
62  (
63  patchi,
65  (
66  bmesh_[patchi],
67  field,
68  dEntry.dict()
69  )
70  );
71  nUnset--;
72  }
73  }
74  }
75 
76  if (nUnset == 0)
77  {
78  return;
79  }
80 
81 
82  // 2. Patch-groups. (using non-wild card entries of dictionaries)
83  // (patchnames already matched above)
84  // Note: in reverse order of entries in the dictionary (last
85  // patchGroups wins). This is so it is consistent with dictionary wildcard
86  // behaviour
87  for (auto iter = dict.crbegin(); iter != dict.crend(); ++iter)
88  {
89  const entry& dEntry = *iter;
90 
91  if (dEntry.isDict() && dEntry.keyword().isLiteral())
92  {
93  const labelList patchIds =
94  bmesh_.indices(dEntry.keyword(), true); // use patchGroups
95 
96  for (const label patchi : patchIds)
97  {
98  if (!this->set(patchi))
99  {
100  this->set
101  (
102  patchi,
104  (
105  bmesh_[patchi],
106  field,
107  dEntry.dict()
108  )
109  );
110  }
111  }
112  }
113  }
114 
115 
116  // 3. Wildcard patch overrides
117  forAll(bmesh_, patchi)
118  {
119  if (!this->set(patchi))
120  {
121  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
122  {
123  this->set
124  (
125  patchi,
127  (
128  emptyPolyPatch::typeName,
129  bmesh_[patchi],
130  field
131  )
132  );
133  }
134  else
135  {
136  bool found = dict.found(bmesh_[patchi].name());
137 
138  if (found)
139  {
140  this->set
141  (
142  patchi,
144  (
145  bmesh_[patchi],
146  field,
147  dict.subDict(bmesh_[patchi].name())
148  )
149  );
150  }
151  }
152  }
153  }
154 
155 
156  // Check for any unset patches
157  forAll(bmesh_, patchi)
158  {
159  if (!this->set(patchi))
160  {
161  if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
162  {
164  << "Cannot find patchField entry for cyclic "
165  << bmesh_[patchi].name() << endl
166  << "Is your field uptodate with split cyclics?" << endl
167  << "Run foamUpgradeCyclics to convert mesh and fields"
168  << " to split cyclics." << exit(FatalIOError);
169  }
170  else
171  {
173  << "Cannot find patchField entry for "
174  << bmesh_[patchi].name() << exit(FatalIOError);
175  }
176  }
177  }
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
182 
183 template<class Type, template<class> class PatchField, class GeoMesh>
185 (
186  const BoundaryMesh& bmesh
187 )
188 :
189  FieldField<PatchField, Type>(bmesh.size()),
190  bmesh_(bmesh)
191 {}
192 
193 
194 template<class Type, template<class> class PatchField, class GeoMesh>
196 (
197  const BoundaryMesh& bmesh,
199  const word& patchFieldType
200 )
201 :
202  FieldField<PatchField, Type>(bmesh.size()),
203  bmesh_(bmesh)
204 {
205  DebugInFunction << nl;
206 
207  forAll(bmesh_, patchi)
208  {
209  this->set
210  (
211  patchi,
213  (
214  patchFieldType,
215  bmesh_[patchi],
216  field
217  )
218  );
219  }
220 }
221 
222 
223 template<class Type, template<class> class PatchField, class GeoMesh>
225 (
226  const BoundaryMesh& bmesh,
227  const DimensionedField<Type, GeoMesh>& field,
228  const wordList& patchFieldTypes,
229  const wordList& constraintTypes
230 )
231 :
232  FieldField<PatchField, Type>(bmesh.size()),
233  bmesh_(bmesh)
234 {
235  DebugInFunction << nl;
236 
237  if
238  (
239  patchFieldTypes.size() != this->size()
240  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
241  )
242  {
244  << "Incorrect number of patch type specifications given" << nl
245  << " Number of patches in mesh = " << bmesh.size()
246  << " number of patch type specifications = "
247  << patchFieldTypes.size()
248  << abort(FatalError);
249  }
250 
251  if (constraintTypes.size())
252  {
253  forAll(bmesh_, patchi)
254  {
255  this->set
256  (
257  patchi,
259  (
260  patchFieldTypes[patchi],
261  constraintTypes[patchi],
262  bmesh_[patchi],
263  field
264  )
265  );
266  }
267  }
268  else
269  {
270  forAll(bmesh_, patchi)
271  {
272  this->set
273  (
274  patchi,
276  (
277  patchFieldTypes[patchi],
278  bmesh_[patchi],
279  field
280  )
281  );
282  }
283  }
284 }
285 
286 
287 template<class Type, template<class> class PatchField, class GeoMesh>
289 (
290  const BoundaryMesh& bmesh,
291  const DimensionedField<Type, GeoMesh>& field,
292  const PtrList<PatchField<Type>>& ptfl
293 )
294 :
295  FieldField<PatchField, Type>(bmesh.size()),
296  bmesh_(bmesh)
297 {
298  DebugInFunction << nl;
299 
300  forAll(bmesh_, patchi)
301  {
302  this->set(patchi, ptfl[patchi].clone(field));
303  }
304 }
305 
306 
307 template<class Type, template<class> class PatchField, class GeoMesh>
309 (
310  const DimensionedField<Type, GeoMesh>& field,
311  const typename GeometricField<Type, PatchField, GeoMesh>::Boundary& btf
312 )
313 :
314  FieldField<PatchField, Type>(btf.size()),
315  bmesh_(btf.bmesh_)
316 {
317  DebugInFunction << nl;
318 
319  forAll(bmesh_, patchi)
320  {
321  this->set(patchi, btf[patchi].clone(field));
322  }
323 }
324 
325 
326 template<class Type, template<class> class PatchField, class GeoMesh>
328 (
329  const DimensionedField<Type, GeoMesh>& field,
330  const typename GeometricField<Type, PatchField, GeoMesh>::Boundary& btf,
331  const labelList& patchIDs,
332  const word& patchFieldType
333 )
334 :
335  FieldField<PatchField, Type>(btf.size()),
336  bmesh_(btf.bmesh_)
337 {
338  DebugInFunction << nl;
339 
340  for (const label patchi : patchIDs)
341  {
342  this->set
343  (
344  patchi,
346  (
347  patchFieldType,
348  bmesh_[patchi],
349  field
350  )
351  );
352  }
353 
354  forAll(bmesh_, patchi)
355  {
356  if (!this->set(patchi))
357  {
358  this->set(patchi, btf[patchi].clone(field));
359  }
360  }
361 }
362 
363 
364 template<class Type, template<class> class PatchField, class GeoMesh>
366 (
367  const typename GeometricField<Type, PatchField, GeoMesh>::Boundary& btf
368 )
369 :
370  FieldField<PatchField, Type>(btf),
371  bmesh_(btf.bmesh_)
372 {
373  DebugInFunction << nl;
374 }
375 
376 
377 template<class Type, template<class> class PatchField, class GeoMesh>
379 (
380  const BoundaryMesh& bmesh,
381  const DimensionedField<Type, GeoMesh>& field,
382  const dictionary& dict
383 )
384 :
385  FieldField<PatchField, Type>(bmesh.size()),
386  bmesh_(bmesh)
387 {
388  readField(field, dict);
389 }
390 
391 
392 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
393 
394 template<class Type, template<class> class PatchField, class GeoMesh>
396 {
397  DebugInFunction << nl;
398 
399  forAll(*this, patchi)
400  {
401  this->operator[](patchi).updateCoeffs();
402  }
403 }
404 
405 
406 template<class Type, template<class> class PatchField, class GeoMesh>
408 {
409  DebugInFunction << nl;
410 
411  if
412  (
415  )
416  {
417  const label nReq = Pstream::nRequests();
418 
419  forAll(*this, patchi)
420  {
421  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
422  }
423 
424  // Block for any outstanding requests
425  if
426  (
429  )
430  {
431  Pstream::waitRequests(nReq);
432  }
433 
434  forAll(*this, patchi)
435  {
436  this->operator[](patchi).evaluate(Pstream::defaultCommsType);
437  }
438  }
440  {
441  const lduSchedule& patchSchedule =
442  bmesh_.mesh().globalData().patchSchedule();
443 
444  forAll(patchSchedule, patchEvali)
445  {
446  if (patchSchedule[patchEvali].init)
447  {
448  this->operator[](patchSchedule[patchEvali].patch)
449  .initEvaluate(Pstream::commsTypes::scheduled);
450  }
451  else
452  {
453  this->operator[](patchSchedule[patchEvali].patch)
455  }
456  }
457  }
458  else
459  {
461  << "Unsupported communications type "
463  << exit(FatalError);
464  }
465 }
466 
467 
468 template<class Type, template<class> class PatchField, class GeoMesh>
471 {
472  const FieldField<PatchField, Type>& pff = *this;
473 
474  wordList list(pff.size());
475 
476  forAll(pff, patchi)
477  {
478  list[patchi] = pff[patchi].type();
479  }
480 
481  return list;
482 }
483 
484 
485 template<class Type, template<class> class PatchField, class GeoMesh>
489 {
491  BoundaryInternalField(*this);
492 
493  forAll(BoundaryInternalField, patchi)
494  {
495  BoundaryInternalField[patchi] ==
496  this->operator[](patchi).patchInternalField();
497  }
498 
499  return BoundaryInternalField;
500 }
501 
502 
503 template<class Type, template<class> class PatchField, class GeoMesh>
506 {
507  LduInterfaceFieldPtrsList<Type> list(this->size());
508 
509  forAll(list, patchi)
510  {
511  const LduInterfaceField<Type>* lduPtr =
512  isA<LduInterfaceField<Type>>(this->operator[](patchi));
513 
514  if (lduPtr)
515  {
516  list.set(patchi, lduPtr);
517  }
518  }
519 
520  return list;
521 }
522 
523 
524 template<class Type, template<class> class PatchField, class GeoMesh>
528 {
529  lduInterfaceFieldPtrsList list(this->size());
530 
531  forAll(list, patchi)
532  {
533  const lduInterfaceField* lduPtr =
534  isA<lduInterfaceField>(this->operator[](patchi));
535 
536  if (lduPtr)
537  {
538  list.set(patchi, lduPtr);
539  }
540  }
541 
542  return list;
543 }
544 
545 
546 template<class Type, template<class> class PatchField, class GeoMesh>
548 (
549  const word& keyword,
550  Ostream& os
551 ) const
552 {
553  os.beginBlock(keyword);
554  this->writeEntries(os);
555  os.endBlock();
556 
557  os.check(FUNCTION_NAME);
558 }
559 
560 
561 template<class Type, template<class> class PatchField, class GeoMesh>
563 (
564  Ostream& os
565 ) const
566 {
567  forAll(*this, patchi)
568  {
569  os.beginBlock(this->operator[](patchi).patch().name());
570  os << this->operator[](patchi);
571  os.endBlock();
572  }
573 }
574 
575 
576 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
577 
578 template<class Type, template<class> class PatchField, class GeoMesh>
580 (
582  Boundary& bf
583 )
584 {
586 }
587 
588 
589 template<class Type, template<class> class PatchField, class GeoMesh>
591 (
592  const FieldField<PatchField, Type>& ptff
593 )
594 {
596 }
597 
598 
599 template<class Type, template<class> class PatchField, class GeoMesh>
601 (
602  const Type& t
603 )
604 {
606 }
607 
608 
609 template<class Type, template<class> class PatchField, class GeoMesh>
611 (
613  Boundary& bf
614 )
615 {
616  forAll(*this, patchi)
617  {
618  this->operator[](patchi) == bf[patchi];
619  }
620 }
621 
622 
623 template<class Type, template<class> class PatchField, class GeoMesh>
625 (
626  const FieldField<PatchField, Type>& ptff
627 )
628 {
629  forAll(*this, patchi)
630  {
631  this->operator[](patchi) == ptff[patchi];
632  }
633 }
634 
635 
636 template<class Type, template<class> class PatchField, class GeoMesh>
638 (
639  const Type& t
640 )
641 {
642  forAll(*this, patchi)
643  {
644  this->operator[](patchi) == t;
645  }
646 }
647 
648 
649 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
650 
651 template<class Type, template<class> class PatchField, class GeoMesh>
652 Foam::Ostream& Foam::operator<<
653 (
654  Ostream& os,
656  Boundary& bf
657 )
658 {
659  os << static_cast<const FieldField<PatchField, Type>&>(bf);
660  return os;
661 }
662 
663 
664 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
setSize
points setSize(newPointi)
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::FieldField< PatchField, Type >
cyclicPolyPatch.H
globalMeshData.H
Foam::GeometricField::Boundary::Boundary
Boundary(const BoundaryMesh &bmesh)
Construct from a BoundaryMesh, setting patches later.
Definition: GeometricBoundaryField.C:185
Foam::entry::isDict
virtual bool isDict() const noexcept
Return true if this entry is a dictionary.
Definition: entry.H:233
Foam::GeometricField::Boundary::writeEntries
void writeEntries(Ostream &os) const
Write dictionary entries of the individual boundary fields.
Definition: GeometricBoundaryField.C:563
Foam::UPstream::waitRequests
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:262
Foam::entry::keyword
const keyType & keyword() const noexcept
Return keyword.
Definition: entry.H:195
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::GeometricField::BoundaryMesh
GeoMesh::BoundaryMesh BoundaryMesh
Definition: GeometricField.H:104
Foam::FieldField< PatchField, Type >::clone
tmp< FieldField< PatchField, Type > > clone() const
Clone.
Definition: FieldField.C:189
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
Foam::GeometricField::Boundary::evaluate
void evaluate()
Evaluate boundary conditions.
Definition: GeometricBoundaryField.C:407
Foam::GeometricField::Boundary::updateCoeffs
void updateCoeffs()
Update the boundary condition coefficients.
Definition: GeometricBoundaryField.C:395
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::GeometricField::Boundary::interfaces
LduInterfaceFieldPtrsList< Type > interfaces() const
Definition: GeometricBoundaryField.C:505
commSchedule.H
patchIds
labelList patchIds
Definition: convertProcessorPatches.H:67
Foam::GeometricField::Boundary::readField
void readField(const Internal &field, const dictionary &dict)
Read the boundary field.
Definition: GeometricBoundaryField.C:36
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
init
mesh init(true)
Foam::UPtrList< const lduInterfaceField >
Foam::GeometricField::Boundary::types
wordList types() const
Return a list of the patch types.
Definition: GeometricBoundaryField.C:470
field
rDeltaTY field()
Foam::LduInterfaceField
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
Definition: LduInterfaceField.H:58
Foam::UPstream::commsTypeNames
static const Enum< commsTypes > commsTypeNames
Names of the communication types.
Definition: UPstream.H:77
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FieldField::operator=
void operator=(const FieldField< Field, Type > &)
Copy assignment.
Definition: FieldField.C:291
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:48
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::entry::dict
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
Foam::lduInterfaceField
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
Definition: lduInterfaceField.H:58
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
found
bool found
Definition: TABSMDCalcMethod2.H:32
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPtrList< const LduInterfaceField< Type > >::set
const const LduInterfaceField< Type > * set(const label i) const
Definition: UPtrList.H:176
Foam::UPstream::commsTypes::nonBlocking
Foam::UPstream::nRequests
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:252
Foam::GeometricField::Boundary
The boundary fields.
Definition: GeometricField.H:115
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::UPstream::commsTypes::scheduled
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::LduInterfaceFieldPtrsList
Definition: LduInterfaceFieldPtrsList.H:50
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::keyType::isLiteral
bool isLiteral() const noexcept
The keyType is treated as literal, not as pattern.
Definition: keyTypeI.H:98
Foam::GeometricField::Boundary::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
Definition: GeometricBoundaryField.C:548
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::UPstream::commsTypes::blocking
Foam::GeometricField::Boundary::scalarInterfaces
lduInterfaceFieldPtrsList scalarInterfaces() const
Definition: GeometricBoundaryField.C:527
Foam::GeometricField::Boundary::boundaryInternalField
Boundary boundaryInternalField() const
Definition: GeometricBoundaryField.C:488