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-2018 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>
37 (
39  const dictionary& dict
40 )
41 {
43 
44  // Clear the boundary field if already initialised
45  this->clear();
46 
47  this->setSize(bmesh_.size());
48 
49  label nUnset = this->size();
50 
51  // 1. Handle explicit patch names. Note that there can be only one explicit
52  // patch name since is key of dictionary.
53 
54  for (const entry& dEntry : dict)
55  {
56  if (dEntry.isDict() && dEntry.keyword().isLiteral())
57  {
58  const label patchi = bmesh_.findPatchID(dEntry.keyword());
59 
60  if (patchi != -1)
61  {
62  this->set
63  (
64  patchi,
66  (
67  bmesh_[patchi],
68  field,
69  dEntry.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 it is consistent with dictionary wildcard
87  // behaviour
88  for (auto iter = dict.crbegin(); iter != dict.crend(); ++iter)
89  {
90  const entry& dEntry = *iter;
91 
92  if (dEntry.isDict() && dEntry.keyword().isLiteral())
93  {
94  const labelList patchIds =
95  bmesh_.indices(dEntry.keyword(), true); // use patchGroups
96 
97  for (const label patchi : patchIds)
98  {
99  if (!this->set(patchi))
100  {
101  this->set
102  (
103  patchi,
105  (
106  bmesh_[patchi],
107  field,
108  dEntry.dict()
109  )
110  );
111  }
112  }
113  }
114  }
115 
116 
117  // 3. Wildcard patch overrides
118  forAll(bmesh_, patchi)
119  {
120  if (!this->set(patchi))
121  {
122  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
123  {
124  this->set
125  (
126  patchi,
128  (
129  emptyPolyPatch::typeName,
130  bmesh_[patchi],
131  field
132  )
133  );
134  }
135  else
136  {
137  bool found = dict.found(bmesh_[patchi].name());
138 
139  if (found)
140  {
141  this->set
142  (
143  patchi,
145  (
146  bmesh_[patchi],
147  field,
148  dict.subDict(bmesh_[patchi].name())
149  )
150  );
151  }
152  }
153  }
154  }
155 
156 
157  // Check for any unset patches
158  forAll(bmesh_, patchi)
159  {
160  if (!this->set(patchi))
161  {
162  if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
163  {
165  << "Cannot find patchField entry for cyclic "
166  << bmesh_[patchi].name() << endl
167  << "Is your field uptodate with split cyclics?" << endl
168  << "Run foamUpgradeCyclics to convert mesh and fields"
169  << " to split cyclics." << exit(FatalIOError);
170  }
171  else
172  {
174  << "Cannot find patchField entry for "
175  << bmesh_[patchi].name() << exit(FatalIOError);
176  }
177  }
178  }
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
183 
184 template<class Type, template<class> class PatchField, class GeoMesh>
186 Boundary
187 (
188  const BoundaryMesh& bmesh
189 )
190 :
191  FieldField<PatchField, Type>(bmesh.size()),
192  bmesh_(bmesh)
193 {}
194 
195 
196 template<class Type, template<class> class PatchField, class GeoMesh>
198 Boundary
199 (
200  const BoundaryMesh& bmesh,
202  const word& patchFieldType
203 )
204 :
205  FieldField<PatchField, Type>(bmesh.size()),
206  bmesh_(bmesh)
207 {
208  DebugInFunction << nl;
209 
210  forAll(bmesh_, patchi)
211  {
212  this->set
213  (
214  patchi,
216  (
217  patchFieldType,
218  bmesh_[patchi],
219  field
220  )
221  );
222  }
223 }
224 
225 
226 template<class Type, template<class> class PatchField, class GeoMesh>
228 Boundary
229 (
230  const BoundaryMesh& bmesh,
231  const DimensionedField<Type, GeoMesh>& field,
232  const wordList& patchFieldTypes,
233  const wordList& constraintTypes
234 )
235 :
236  FieldField<PatchField, Type>(bmesh.size()),
237  bmesh_(bmesh)
238 {
239  DebugInFunction << nl;
240 
241  if
242  (
243  patchFieldTypes.size() != this->size()
244  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
245  )
246  {
248  << "Incorrect number of patch type specifications given" << nl
249  << " Number of patches in mesh = " << bmesh.size()
250  << " number of patch type specifications = "
251  << patchFieldTypes.size()
252  << abort(FatalError);
253  }
254 
255  if (constraintTypes.size())
256  {
257  forAll(bmesh_, patchi)
258  {
259  this->set
260  (
261  patchi,
263  (
264  patchFieldTypes[patchi],
265  constraintTypes[patchi],
266  bmesh_[patchi],
267  field
268  )
269  );
270  }
271  }
272  else
273  {
274  forAll(bmesh_, patchi)
275  {
276  this->set
277  (
278  patchi,
280  (
281  patchFieldTypes[patchi],
282  bmesh_[patchi],
283  field
284  )
285  );
286  }
287  }
288 }
289 
290 
291 template<class Type, template<class> class PatchField, class GeoMesh>
293 Boundary
294 (
295  const BoundaryMesh& bmesh,
296  const DimensionedField<Type, GeoMesh>& field,
297  const PtrList<PatchField<Type>>& ptfl
298 )
299 :
300  FieldField<PatchField, Type>(bmesh.size()),
301  bmesh_(bmesh)
302 {
303  DebugInFunction << nl;
304 
305  forAll(bmesh_, patchi)
306  {
307  this->set(patchi, ptfl[patchi].clone(field));
308  }
309 }
310 
311 
312 template<class Type, template<class> class PatchField, class GeoMesh>
314 Boundary
315 (
316  const DimensionedField<Type, GeoMesh>& field,
317  const typename GeometricField<Type, PatchField, GeoMesh>::
318  Boundary& btf
319 )
320 :
321  FieldField<PatchField, Type>(btf.size()),
322  bmesh_(btf.bmesh_)
323 {
324  DebugInFunction << nl;
325 
326  forAll(bmesh_, patchi)
327  {
328  this->set(patchi, btf[patchi].clone(field));
329  }
330 }
331 
332 
333 template<class Type, template<class> class PatchField, class GeoMesh>
335 Boundary
336 (
337  const typename GeometricField<Type, PatchField, GeoMesh>::
338  Boundary& btf
339 )
340 :
341  FieldField<PatchField, Type>(btf),
342  bmesh_(btf.bmesh_)
343 {
344  DebugInFunction << nl;
345 }
346 
347 
348 template<class Type, template<class> class PatchField, class GeoMesh>
350 Boundary
351 (
352  const BoundaryMesh& bmesh,
353  const DimensionedField<Type, GeoMesh>& field,
354  const dictionary& dict
355 )
356 :
357  FieldField<PatchField, Type>(bmesh.size()),
358  bmesh_(bmesh)
359 {
360  readField(field, dict);
361 }
362 
363 
364 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
365 
366 template<class Type, template<class> class PatchField, class GeoMesh>
369 {
370  DebugInFunction << nl;
371 
372  forAll(*this, patchi)
373  {
374  this->operator[](patchi).updateCoeffs();
375  }
376 }
377 
378 
379 template<class Type, template<class> class PatchField, class GeoMesh>
382 {
383  DebugInFunction << nl;
384 
385  if
386  (
389  )
390  {
391  label nReq = Pstream::nRequests();
392 
393  forAll(*this, patchi)
394  {
395  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
396  }
397 
398  // Block for any outstanding requests
399  if
400  (
403  )
404  {
405  Pstream::waitRequests(nReq);
406  }
407 
408  forAll(*this, patchi)
409  {
410  this->operator[](patchi).evaluate(Pstream::defaultCommsType);
411  }
412  }
414  {
415  const lduSchedule& patchSchedule =
416  bmesh_.mesh().globalData().patchSchedule();
417 
418  forAll(patchSchedule, patchEvali)
419  {
420  if (patchSchedule[patchEvali].init)
421  {
422  this->operator[](patchSchedule[patchEvali].patch)
423  .initEvaluate(Pstream::commsTypes::scheduled);
424  }
425  else
426  {
427  this->operator[](patchSchedule[patchEvali].patch)
429  }
430  }
431  }
432  else
433  {
435  << "Unsuported communications type "
437  << exit(FatalError);
438  }
439 }
440 
441 
442 template<class Type, template<class> class PatchField, class GeoMesh>
445 types() const
446 {
447  const FieldField<PatchField, Type>& pff = *this;
448 
449  wordList list(pff.size());
450 
451  forAll(pff, patchi)
452  {
453  list[patchi] = pff[patchi].type();
454  }
455 
456  return list;
457 }
458 
459 
460 template<class Type, template<class> class PatchField, class GeoMesh>
464 {
466  BoundaryInternalField(*this);
467 
468  forAll(BoundaryInternalField, patchi)
469  {
470  BoundaryInternalField[patchi] ==
471  this->operator[](patchi).patchInternalField();
472  }
473 
474  return BoundaryInternalField;
475 }
476 
477 
478 template<class Type, template<class> class PatchField, class GeoMesh>
481 interfaces() const
482 {
483  LduInterfaceFieldPtrsList<Type> interfaces(this->size());
484 
485  forAll(interfaces, patchi)
486  {
487  if (isA<LduInterfaceField<Type>>(this->operator[](patchi)))
488  {
489  interfaces.set
490  (
491  patchi,
493  (
494  this->operator[](patchi)
495  )
496  );
497  }
498  }
499 
500  return interfaces;
501 }
502 
503 
504 template<class Type, template<class> class PatchField, class GeoMesh>
508 {
509  lduInterfaceFieldPtrsList interfaces(this->size());
510 
511  forAll(interfaces, patchi)
512  {
513  if (isA<lduInterfaceField>(this->operator[](patchi)))
514  {
515  interfaces.set
516  (
517  patchi,
518  &refCast<const lduInterfaceField>
519  (
520  this->operator[](patchi)
521  )
522  );
523  }
524  }
525 
526  return interfaces;
527 }
528 
529 
530 template<class Type, template<class> class PatchField, class GeoMesh>
532 writeEntry(const word& keyword, Ostream& os) const
533 {
534  os.beginBlock(keyword);
535  this->writeEntries(os);
536  os.endBlock();
537 
538  os.check(FUNCTION_NAME);
539 }
540 
541 
542 template<class Type, template<class> class PatchField, class GeoMesh>
545 {
546  forAll(*this, patchi)
547  {
548  os.beginBlock(this->operator[](patchi).patch().name());
549  os << this->operator[](patchi);
550  os.endBlock();
551  }
552 }
553 
554 
555 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
556 
557 template<class Type, template<class> class PatchField, class GeoMesh>
559 operator=
560 (
562  Boundary& bf
563 )
564 {
566 }
567 
568 
569 template<class Type, template<class> class PatchField, class GeoMesh>
571 operator=
572 (
573  const FieldField<PatchField, Type>& ptff
574 )
575 {
577 }
578 
579 
580 template<class Type, template<class> class PatchField, class GeoMesh>
582 operator=
583 (
584  const Type& t
585 )
586 {
588 }
589 
590 
591 template<class Type, template<class> class PatchField, class GeoMesh>
593 operator==
594 (
596  Boundary& bf
597 )
598 {
599  forAll((*this), patchi)
600  {
601  this->operator[](patchi) == bf[patchi];
602  }
603 }
604 
605 
606 template<class Type, template<class> class PatchField, class GeoMesh>
608 operator==
609 (
610  const FieldField<PatchField, Type>& ptff
611 )
612 {
613  forAll((*this), patchi)
614  {
615  this->operator[](patchi) == ptff[patchi];
616  }
617 }
618 
619 
620 template<class Type, template<class> class PatchField, class GeoMesh>
622 operator==
623 (
624  const Type& t
625 )
626 {
627  forAll((*this), patchi)
628  {
629  this->operator[](patchi) == t;
630  }
631 }
632 
633 
634 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
635 
636 template<class Type, template<class> class PatchField, class GeoMesh>
637 Foam::Ostream& Foam::operator<<
638 (
639  Ostream& os,
641  Boundary& bf
642 )
643 {
644  os << static_cast<const FieldField<PatchField, Type>&>(bf);
645  return os;
646 }
647 
648 
649 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::UPstream::commsTypes::blocking
setSize
points setSize(newPointi)
Foam::UPstream::commsTypes::nonBlocking
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::FieldField< PatchField, Type >
Foam::entry::keyword
const keyType & keyword() const
Return keyword.
Definition: entry.H:187
cyclicPolyPatch.H
globalMeshData.H
Foam::GeometricField::Boundary::Boundary
Boundary(const BoundaryMesh &bmesh)
Construct from a BoundaryMesh.
Definition: GeometricBoundaryField.C:187
Foam::GeometricField::Boundary::writeEntries
void writeEntries(Ostream &os) const
Write dictionary entries of the individual boundary fields.
Definition: GeometricBoundaryField.C:544
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::UPstream::waitRequests
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:162
Foam::UPtrList< const LduInterfaceField< Type > >::set
const const LduInterfaceField< Type > * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: UPtrListI.H:176
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::GeometricField::BoundaryMesh
GeoMesh::BoundaryMesh BoundaryMesh
Type of boundary mesh on which this.
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:273
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
Foam::GeometricField::Boundary::evaluate
void evaluate()
Evaluate boundary conditions.
Definition: GeometricBoundaryField.C:381
Foam::GeometricField::Boundary::updateCoeffs
void updateCoeffs()
Update the boundary condition coefficients.
Definition: GeometricBoundaryField.C:368
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
Foam::GeometricField::Boundary::interfaces
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
Definition: GeometricBoundaryField.C:481
commSchedule.H
patchIds
labelList patchIds
Definition: convertProcessorPatches.H:67
Foam::label
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:62
Foam::GeometricField::Boundary::readField
void readField(const Internal &field, const dictionary &dict)
Read the boundary field.
Definition: GeometricBoundaryField.C:37
Foam::entry::isDict
virtual bool isDict() const
Return true if this entry is a dictionary.
Definition: entry.H:222
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:356
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::UPtrList< const lduInterfaceField >
Foam::GeometricField::Boundary::types
wordList types() const
Return a list of the patch types.
Definition: GeometricBoundaryField.C:445
field
rDeltaTY field()
Foam::LduInterfaceField
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
Definition: LduInterfaceField.H:59
Foam::UPstream::commsTypeNames
static const Enum< commsTypes > commsTypeNames
Names of the communication types.
Definition: UPstream.H:74
Foam::UPstream::commsTypes::scheduled
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::FieldField::operator=
void operator=(const FieldField< Field, Type > &)
Copy assignment.
Definition: FieldField.C:291
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
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:121
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:137
Foam::entry::dict
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
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
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::UPstream::nRequests
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:152
Foam::GeometricField::Boundary
Definition: GeometricField.H:114
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:261
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
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::operator
friend Ostream & operator(Ostream &, const GeometricField< Type, PatchField, GeoMesh > &)
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::keyType::isLiteral
bool isLiteral() const
The keyType is treated as literal, not as pattern.
Definition: keyTypeI.H:122
Foam::GeometricField::Boundary::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
Definition: GeometricBoundaryField.C:532
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::GeometricField::Boundary::scalarInterfaces
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those.
Definition: GeometricBoundaryField.C:507
Foam::GeometricField::Boundary::boundaryInternalField
Boundary boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
Definition: GeometricBoundaryField.C:463