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 | Copyright (C) 2011-2015 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  {
46  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
47  "GeometricBoundaryField::readField"
48  "("
49  "const DimensionedField<Type, GeoMesh>&, "
50  "const dictionary&"
51  ")"
52  << endl;
53  }
54 
55 
56  label nUnset = this->size();
57 
58  // 1. Handle explicit patch names. Note that there can be only one explicit
59  // patch name since is key of dictionary.
61  {
62  if (iter().isDict() && !iter().keyword().isPattern())
63  {
64  label patchi = bmesh_.findPatchID(iter().keyword());
65 
66  if (patchi != -1)
67  {
68  this->set
69  (
70  patchi,
72  (
73  bmesh_[patchi],
74  field,
75  iter().dict()
76  )
77  );
78  nUnset--;
79  }
80  }
81  }
82 
83  if (nUnset == 0)
84  {
85  return;
86  }
87 
88 
89  // 2. Patch-groups. (using non-wild card entries of dictionaries)
90  // (patchnames already matched above)
91  // Note: in reverse order of entries in the dictionary (last
92  // patchGroups wins). This is so is consistent with dictionary wildcard
93  // behaviour
94  if (dict.size())
95  {
96  for
97  (
99  iter != dict.rend();
100  ++iter
101  )
102  {
103  const entry& e = iter();
104 
105  if (e.isDict() && !e.keyword().isPattern())
106  {
107  const labelList patchIDs = bmesh_.findIndices
108  (
109  e.keyword(),
110  true // use patchGroups
111  );
112 
113  forAll(patchIDs, i)
114  {
115  label patchi = patchIDs[i];
116 
117  if (!this->set(patchi))
118  {
119  this->set
120  (
121  patchi,
123  (
124  bmesh_[patchi],
125  field,
126  e.dict()
127  )
128  );
129  }
130  }
131  }
132  }
133  }
134 
135 
136  // 3. Wildcard patch overrides
137  forAll(bmesh_, patchi)
138  {
139  if (!this->set(patchi))
140  {
141  if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
142  {
143  this->set
144  (
145  patchi,
147  (
148  emptyPolyPatch::typeName,
149  bmesh_[patchi],
150  field
151  )
152  );
153  }
154  else
155  {
156  bool found = dict.found(bmesh_[patchi].name());
157 
158  if (found)
159  {
160  this->set
161  (
162  patchi,
164  (
165  bmesh_[patchi],
166  field,
167  dict.subDict(bmesh_[patchi].name())
168  )
169  );
170  }
171  }
172  }
173  }
174 
175 
176  // Check for any unset patches
177  forAll(bmesh_, patchi)
178  {
179  if (!this->set(patchi))
180  {
181  if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
182  {
184  (
185  dict
186  ) << "Cannot find patchField entry for cyclic "
187  << bmesh_[patchi].name() << endl
188  << "Is your field uptodate with split cyclics?" << endl
189  << "Run foamUpgradeCyclics to convert mesh and fields"
190  << " to split cyclics." << exit(FatalIOError);
191  }
192  else
193  {
195  (
196  dict
197  ) << "Cannot find patchField entry for "
198  << bmesh_[patchi].name() << exit(FatalIOError);
199  }
200  }
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
207 template<class Type, template<class> class PatchField, class GeoMesh>
210 (
211  const BoundaryMesh& bmesh
212 )
213 :
214  FieldField<PatchField, Type>(bmesh.size()),
215  bmesh_(bmesh)
216 {}
217 
218 
219 template<class Type, template<class> class PatchField, class GeoMesh>
222 (
223  const BoundaryMesh& bmesh,
224  const DimensionedField<Type, GeoMesh>& field,
225  const word& patchFieldType
226 )
227 :
228  FieldField<PatchField, Type>(bmesh.size()),
229  bmesh_(bmesh)
230 {
231  if (debug)
232  {
233  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
234  "GeometricBoundaryField::"
235  "GeometricBoundaryField(const BoundaryMesh&, "
236  "const DimensionedField<Type>&, const word&)"
237  << endl;
238  }
239 
241  {
242  this->set
243  (
244  patchi,
246  (
247  patchFieldType,
248  bmesh_[patchi],
249  field
250  )
251  );
252  }
253 }
254 
255 
256 template<class Type, template<class> class PatchField, class GeoMesh>
259 (
260  const BoundaryMesh& bmesh,
261  const DimensionedField<Type, GeoMesh>& field,
262  const wordList& patchFieldTypes,
263  const wordList& constraintTypes
264 )
265 :
266  FieldField<PatchField, Type>(bmesh.size()),
267  bmesh_(bmesh)
268 {
269  if (debug)
270  {
271  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
272  "GeometricBoundaryField::"
273  "GeometricBoundaryField"
274  "("
275  "const BoundaryMesh&, "
276  "const DimensionedField<Type>&, "
277  "const wordList&, "
278  "const wordList&"
279  ")"
280  << endl;
281  }
282 
283  if
284  (
285  patchFieldTypes.size() != this->size()
286  || (constraintTypes.size() && (constraintTypes.size() != this->size()))
287  )
288  {
290  << "Incorrect number of patch type specifications given" << nl
291  << " Number of patches in mesh = " << bmesh.size()
292  << " number of patch type specifications = "
293  << patchFieldTypes.size()
294  << abort(FatalError);
295  }
296 
297  if (constraintTypes.size())
298  {
300  {
301  this->set
302  (
303  patchi,
305  (
306  patchFieldTypes[patchi],
307  constraintTypes[patchi],
308  bmesh_[patchi],
309  field
310  )
311  );
312  }
313  }
314  else
315  {
317  {
318  this->set
319  (
320  patchi,
322  (
323  patchFieldTypes[patchi],
324  bmesh_[patchi],
325  field
326  )
327  );
328  }
329  }
330 }
331 
332 
333 template<class Type, template<class> class PatchField, class GeoMesh>
336 (
337  const BoundaryMesh& bmesh,
338  const DimensionedField<Type, GeoMesh>& field,
339  const PtrList<PatchField<Type> >& ptfl
340 )
341 :
342  FieldField<PatchField, Type>(bmesh.size()),
343  bmesh_(bmesh)
344 {
345  if (debug)
346  {
347  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
348  "GeometricBoundaryField::"
349  "GeometricBoundaryField"
350  "("
351  "const BoundaryMesh&, "
352  "const DimensionedField<Type, GeoMesh>&, "
353  "const PtrLIst<PatchField<Type> >&"
354  ")"
355  << endl;
356  }
357 
359  {
360  this->set(patchi, ptfl[patchi].clone(field));
361  }
362 }
363 
364 
365 template<class Type, template<class> class PatchField, class GeoMesh>
368 (
369  const DimensionedField<Type, GeoMesh>& field,
370  const typename GeometricField<Type, PatchField, GeoMesh>::
371  GeometricBoundaryField& btf
372 )
373 :
374  FieldField<PatchField, Type>(btf.size()),
375  bmesh_(btf.bmesh_)
376 {
377  if (debug)
378  {
379  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
380  "GeometricBoundaryField::"
381  "GeometricBoundaryField"
382  "("
383  "const DimensionedField<Type, GeoMesh>&, "
384  "const typename GeometricField<Type, PatchField, GeoMesh>::"
385  "GeometricBoundaryField&"
386  ")"
387  << endl;
388  }
389 
391  {
392  this->set(patchi, btf[patchi].clone(field));
393  }
394 }
395 
396 
397 // Construct as copy
398 // Dangerous because Field may be set to a field which gets deleted.
399 // Need new type of GeometricBoundaryField, one which IS part of a geometric
400 // field for which snGrad etc. may be called and a free standing
401 // GeometricBoundaryField for which such operations are unavailable.
402 template<class Type, template<class> class PatchField, class GeoMesh>
405 (
406  const typename GeometricField<Type, PatchField, GeoMesh>::
407  GeometricBoundaryField& btf
408 )
409 :
410  FieldField<PatchField, Type>(btf),
411  bmesh_(btf.bmesh_)
412 {
413  if (debug)
414  {
415  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
416  "GeometricBoundaryField::"
417  "GeometricBoundaryField"
418  "("
419  "const GeometricField<Type, PatchField, GeoMesh>::"
420  "GeometricBoundaryField&"
421  ")"
422  << endl;
423  }
424 }
425 
426 
427 template<class Type, template<class> class PatchField, class GeoMesh>
430 (
431  const BoundaryMesh& bmesh,
432  const DimensionedField<Type, GeoMesh>& field,
433  const dictionary& dict
434 )
435 :
436  FieldField<PatchField, Type>(bmesh.size()),
437  bmesh_(bmesh)
438 {
439  readField(field, dict);
440 }
441 
442 
443 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
444 
445 template<class Type, template<class> class PatchField, class GeoMesh>
448 {
449  if (debug)
450  {
451  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
452  "GeometricBoundaryField::"
453  "updateCoeffs()" << endl;
454  }
455 
456  forAll(*this, patchi)
457  {
458  this->operator[](patchi).updateCoeffs();
459  }
460 }
461 
462 
463 template<class Type, template<class> class PatchField, class GeoMesh>
466 {
467  if (debug)
468  {
469  Info<< "GeometricField<Type, PatchField, GeoMesh>::"
470  "GeometricBoundaryField::"
471  "evaluate()" << endl;
472  }
473 
474  if
475  (
478  )
479  {
480  label nReq = Pstream::nRequests();
481 
482  forAll(*this, patchi)
483  {
484  this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
485  }
486 
487  // Block for any outstanding requests
488  if
489  (
492  )
493  {
494  Pstream::waitRequests(nReq);
495  }
496 
497  forAll(*this, patchi)
498  {
499  this->operator[](patchi).evaluate(Pstream::defaultCommsType);
500  }
501  }
503  {
504  const lduSchedule& patchSchedule =
505  bmesh_.mesh().globalData().patchSchedule();
506 
507  forAll(patchSchedule, patchEvali)
508  {
509  if (patchSchedule[patchEvali].init)
510  {
511  this->operator[](patchSchedule[patchEvali].patch)
512  .initEvaluate(Pstream::scheduled);
513  }
514  else
515  {
516  this->operator[](patchSchedule[patchEvali].patch)
517  .evaluate(Pstream::scheduled);
518  }
519  }
520  }
521  else
522  {
524  << "Unsuported communications type "
526  << exit(FatalError);
527  }
528 }
529 
530 
531 template<class Type, template<class> class PatchField, class GeoMesh>
534 types() const
535 {
536  const FieldField<PatchField, Type>& pff = *this;
537 
538  wordList Types(pff.size());
539 
540  forAll(pff, patchi)
541  {
542  Types[patchi] = pff[patchi].type();
543  }
544 
545  return Types;
546 }
547 
548 
549 template<class Type, template<class> class PatchField, class GeoMesh>
553 {
555  BoundaryInternalField(*this);
556 
557  forAll(BoundaryInternalField, patchi)
558  {
559  BoundaryInternalField[patchi] ==
560  this->operator[](patchi).patchInternalField();
561  }
562 
563  return BoundaryInternalField;
564 }
565 
566 
567 template<class Type, template<class> class PatchField, class GeoMesh>
570 interfaces() const
571 {
572  LduInterfaceFieldPtrsList<Type> interfaces(this->size());
573 
574  forAll(interfaces, patchi)
575  {
576  if (isA<LduInterfaceField<Type> >(this->operator[](patchi)))
577  {
578  interfaces.set
579  (
580  patchi,
582  (
583  this->operator[](patchi)
584  )
585  );
586  }
587  }
588 
589  return interfaces;
590 }
591 
592 
593 template<class Type, template<class> class PatchField, class GeoMesh>
597 {
598  lduInterfaceFieldPtrsList interfaces(this->size());
599 
600  forAll(interfaces, patchi)
601  {
602  if (isA<lduInterfaceField>(this->operator[](patchi)))
603  {
604  interfaces.set
605  (
606  patchi,
607  &refCast<const lduInterfaceField>
608  (
609  this->operator[](patchi)
610  )
611  );
612  }
613  }
614 
615  return interfaces;
616 }
617 
618 
619 template<class Type, template<class> class PatchField, class GeoMesh>
621 writeEntry(const word& keyword, Ostream& os) const
622 {
623  os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
624 
625  forAll(*this, patchi)
626  {
627  os << indent << this->operator[](patchi).patch().name() << nl
628  << indent << token::BEGIN_BLOCK << nl
629  << incrIndent << this->operator[](patchi) << decrIndent
630  << indent << token::END_BLOCK << endl;
631  }
632 
633  os << decrIndent << token::END_BLOCK << endl;
634 
635  // Check state of IOstream
636  os.check
637  (
638  "GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::"
639  "writeEntry(const word& keyword, Ostream& os) const"
640  );
641 }
642 
643 
644 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
645 
646 template<class Type, template<class> class PatchField, class GeoMesh>
648 operator=
649 (
651  GeometricBoundaryField& bf
652 )
653 {
655 }
656 
657 
658 template<class Type, template<class> class PatchField, class GeoMesh>
660 operator=
661 (
662  const FieldField<PatchField, Type>& ptff
663 )
664 {
666 }
667 
668 
669 template<class Type, template<class> class PatchField, class GeoMesh>
671 operator=
672 (
673  const Type& t
674 )
675 {
677 }
678 
679 
680 // Forced assignments
681 template<class Type, template<class> class PatchField, class GeoMesh>
683 operator==
684 (
685  const typename GeometricField<Type, PatchField, GeoMesh>::
686  GeometricBoundaryField& bf
687 )
688 {
689  forAll((*this), patchI)
690  {
691  this->operator[](patchI) == bf[patchI];
692  }
693 }
694 
695 
696 template<class Type, template<class> class PatchField, class GeoMesh>
698 operator==
699 (
700  const FieldField<PatchField, Type>& ptff
701 )
702 {
703  forAll((*this), patchI)
704  {
705  this->operator[](patchI) == ptff[patchI];
706  }
707 }
708 
709 
710 template<class Type, template<class> class PatchField, class GeoMesh>
712 operator==
713 (
714  const Type& t
715 )
716 {
717  forAll((*this), patchI)
718  {
719  this->operator[](patchI) == t;
720  }
721 }
722 
723 
724 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
725 
726 template<class Type, template<class> class PatchField, class GeoMesh>
727 Foam::Ostream& Foam::operator<<
728 (
729  Ostream& os,
730  const typename GeometricField<Type, PatchField, GeoMesh>::
731  GeometricBoundaryField& bf
732 )
733 {
734  os << static_cast<const FieldField<PatchField, Type>&>(bf);
735  return os;
736 }
737 
738 
739 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
Foam::GeometricField::GeometricBoundaryField::interfaces
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
Definition: GeometricBoundaryField.C:570
Foam::GeometricField::GeometricBoundaryField::evaluate
void evaluate()
Evaluate boundary conditions.
Definition: GeometricBoundaryField.C:465
setSize
points setSize(newPointi)
Foam::GeometricField::GeometricBoundaryField::updateCoeffs
void updateCoeffs()
Update the boundary condition coefficients.
Definition: GeometricBoundaryField.C:447
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
clear
UEqn clear()
Foam::UPstream::scheduled
@ scheduled
Definition: UPstream.H:67
cyclicPolyPatch.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
globalMeshData.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::UPstream::commsTypeNames
static const NamedEnum< commsTypes, 3 > commsTypeNames
Definition: UPstream.H:71
Foam::UPstream::waitRequests
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:106
Foam::FatalIOError
IOerror FatalIOError
Foam::GeometricField::GeometricBoundaryField::readField
void readField(const DimensionedField< Type, GeoMesh > &field, const dictionary &dict)
Read the boundary field.
Definition: GeometricBoundaryField.C:34
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::GeometricField::BoundaryMesh
GeoMesh::BoundaryMesh BoundaryMesh
Definition: GeometricField.H:98
Foam::FieldField< PatchField, Type >::clone
tmp< FieldField< Field, Type > > clone() const
Clone.
Definition: FieldField.C:184
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:252
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
Foam::GeometricField::GeometricBoundaryField::bmesh_
const BoundaryMesh & bmesh_
Reference to BoundaryMesh for which this field is defined.
Definition: GeometricField.H:112
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::GeometricField::GeometricBoundaryField::scalarInterfaces
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those.
Definition: GeometricBoundaryField.C:596
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
commSchedule.H
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
Foam::UPstream::nonBlocking
@ nonBlocking
Definition: UPstream.H:68
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:59
Foam::IDLList
Intrusive doubly-linked list.
Definition: IDLList.H:47
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::GeometricField::GeometricBoundaryField::GeometricBoundaryField
GeometricBoundaryField(const BoundaryMesh &)
Construct from a BoundaryMesh.
Foam::UPtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:53
Foam::LduInterfaceField
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
Definition: LduInterfaceField.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FieldField::operator=
void operator=(const FieldField< Field, Type > &)
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::GeometricField::GeometricBoundaryField::types
wordList types() const
Return a list of the patch types.
Definition: GeometricBoundaryField.C:534
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::UPstream::nRequests
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:96
Foam::token::END_BLOCK
@ END_BLOCK
Definition: token.H:105
Foam::isA
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::UPtrList< const LduInterfaceField< Type > >::set
bool set(const label) const
Is element set.
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::token::BEGIN_BLOCK
@ BEGIN_BLOCK
Definition: token.H:104
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::LduInterfaceFieldPtrsList
Definition: LduInterfaceFieldPtrsList.H:48
Foam::GeometricField::GeometricBoundaryField::boundaryInternalField
GeometricBoundaryField boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
Definition: GeometricBoundaryField.C:552
Foam::GeometricField::operator
friend Ostream & operator(Ostream &, const GeometricField< Type, PatchField, GeoMesh > &)
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
Foam::GeometricField::GeometricBoundaryField::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
Definition: GeometricBoundaryField.C:621
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51