genericFvPatchField.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 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 "genericFvPatchField.H"
27 #include "fvPatchFieldMapper.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const fvPatch& p,
36 )
37 :
39 {
41  (
42  "genericFvPatchField<Type>::genericFvPatchField"
43  "(const fvPatch& p, const DimensionedField<Type, volMesh>& iF)"
44  ) << "Not Implemented\n "
45  << "Trying to construct an genericFvPatchField on patch "
46  << this->patch().name()
47  << " of field " << this->dimensionedInternalField().name()
48  << abort(FatalError);
49 }
50 
51 
52 template<class Type>
54 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
60  calculatedFvPatchField<Type>(p, iF, dict, false),
61  actualTypeName_(dict.lookup("type")),
62  dict_(dict)
63 {
64  if (!dict.found("value"))
65  {
67  (
68  "genericFvPatchField<Type>::genericFvPatchField"
69  "(const fvPatch&, const Field<Type>&, const dictionary&)",
70  dict
71  ) << "\n Cannot find 'value' entry"
72  << " on patch " << this->patch().name()
73  << " of field " << this->dimensionedInternalField().name()
74  << " in file " << this->dimensionedInternalField().objectPath()
75  << nl
76  << " which is required to set the"
77  " values of the generic patch field." << nl
78  << " (Actual type " << actualTypeName_ << ")" << nl
79  << "\n Please add the 'value' entry to the write function "
80  "of the user-defined boundary-condition\n"
81  << exit(FatalIOError);
82  }
83 
84  forAllConstIter(dictionary, dict_, iter)
85  {
86  if (iter().keyword() != "type" && iter().keyword() != "value")
87  {
88  if
89  (
90  iter().isStream()
91  && iter().stream().size()
92  )
93  {
94  ITstream& is = iter().stream();
95 
96  // Read first token
97  token firstToken(is);
98 
99  if
100  (
101  firstToken.isWord()
102  && firstToken.wordToken() == "nonuniform"
103  )
104  {
105  token fieldToken(is);
106 
107  if (!fieldToken.isCompound())
108  {
109  if
110  (
111  fieldToken.isLabel()
112  && fieldToken.labelToken() == 0
113  )
114  {
115  scalarFields_.insert
116  (
117  iter().keyword(),
118  new scalarField(0)
119  );
120  }
121  else
122  {
124  (
125  "genericFvPatchField<Type>::genericFvPatchField"
126  "(const fvPatch&, const Field<Type>&, "
127  "const dictionary&)",
128  dict
129  ) << "\n token following 'nonuniform' "
130  "is not a compound"
131  << "\n on patch " << this->patch().name()
132  << " of field "
133  << this->dimensionedInternalField().name()
134  << " in file "
135  << this->dimensionedInternalField().objectPath()
136  << exit(FatalIOError);
137  }
138  }
139  else if
140  (
141  fieldToken.compoundToken().type()
142  == token::Compound<List<scalar> >::typeName
143  )
144  {
145  scalarField* fPtr = new scalarField;
146  fPtr->transfer
147  (
149  (
150  fieldToken.transferCompoundToken(is)
151  )
152  );
153 
154  if (fPtr->size() != this->size())
155  {
157  (
158  "genericFvPatchField<Type>::genericFvPatchField"
159  "(const fvPatch&, const Field<Type>&, "
160  "const dictionary&)",
161  dict
162  ) << "\n size of field " << iter().keyword()
163  << " (" << fPtr->size() << ')'
164  << " is not the same size as the patch ("
165  << this->size() << ')'
166  << "\n on patch " << this->patch().name()
167  << " of field "
168  << this->dimensionedInternalField().name()
169  << " in file "
170  << this->dimensionedInternalField().objectPath()
171  << exit(FatalIOError);
172  }
173 
174  scalarFields_.insert(iter().keyword(), fPtr);
175  }
176  else if
177  (
178  fieldToken.compoundToken().type()
179  == token::Compound<List<vector> >::typeName
180  )
181  {
182  vectorField* fPtr = new vectorField;
183  fPtr->transfer
184  (
186  (
187  fieldToken.transferCompoundToken(is)
188  )
189  );
190 
191  if (fPtr->size() != this->size())
192  {
194  (
195  "genericFvPatchField<Type>::genericFvPatchField"
196  "(const fvPatch&, const Field<Type>&, "
197  "const dictionary&)",
198  dict
199  ) << "\n size of field " << iter().keyword()
200  << " (" << fPtr->size() << ')'
201  << " is not the same size as the patch ("
202  << this->size() << ')'
203  << "\n on patch " << this->patch().name()
204  << " of field "
205  << this->dimensionedInternalField().name()
206  << " in file "
207  << this->dimensionedInternalField().objectPath()
208  << exit(FatalIOError);
209  }
210 
211  vectorFields_.insert(iter().keyword(), fPtr);
212  }
213  else if
214  (
215  fieldToken.compoundToken().type()
217  )
218  {
220  fPtr->transfer
221  (
223  <
225  >
226  (
227  fieldToken.transferCompoundToken(is)
228  )
229  );
230 
231  if (fPtr->size() != this->size())
232  {
234  (
235  "genericFvPatchField<Type>::genericFvPatchField"
236  "(const fvPatch&, const Field<Type>&, "
237  "const dictionary&)",
238  dict
239  ) << "\n size of field " << iter().keyword()
240  << " (" << fPtr->size() << ')'
241  << " is not the same size as the patch ("
242  << this->size() << ')'
243  << "\n on patch " << this->patch().name()
244  << " of field "
245  << this->dimensionedInternalField().name()
246  << " in file "
247  << this->dimensionedInternalField().objectPath()
248  << exit(FatalIOError);
249  }
250 
251  sphericalTensorFields_.insert(iter().keyword(), fPtr);
252  }
253  else if
254  (
255  fieldToken.compoundToken().type()
256  == token::Compound<List<symmTensor> >::typeName
257  )
258  {
259  symmTensorField* fPtr = new symmTensorField;
260  fPtr->transfer
261  (
263  <
265  >
266  (
267  fieldToken.transferCompoundToken(is)
268  )
269  );
270 
271  if (fPtr->size() != this->size())
272  {
274  (
275  "genericFvPatchField<Type>::genericFvPatchField"
276  "(const fvPatch&, const Field<Type>&, "
277  "const dictionary&)",
278  dict
279  ) << "\n size of field " << iter().keyword()
280  << " (" << fPtr->size() << ')'
281  << " is not the same size as the patch ("
282  << this->size() << ')'
283  << "\n on patch " << this->patch().name()
284  << " of field "
285  << this->dimensionedInternalField().name()
286  << " in file "
287  << this->dimensionedInternalField().objectPath()
288  << exit(FatalIOError);
289  }
290 
291  symmTensorFields_.insert(iter().keyword(), fPtr);
292  }
293  else if
294  (
295  fieldToken.compoundToken().type()
296  == token::Compound<List<tensor> >::typeName
297  )
298  {
299  tensorField* fPtr = new tensorField;
300  fPtr->transfer
301  (
303  (
304  fieldToken.transferCompoundToken(is)
305  )
306  );
307 
308  if (fPtr->size() != this->size())
309  {
311  (
312  "genericFvPatchField<Type>::genericFvPatchField"
313  "(const fvPatch&, const Field<Type>&, "
314  "const dictionary&)",
315  dict
316  ) << "\n size of field " << iter().keyword()
317  << " (" << fPtr->size() << ')'
318  << " is not the same size as the patch ("
319  << this->size() << ')'
320  << "\n on patch " << this->patch().name()
321  << " of field "
322  << this->dimensionedInternalField().name()
323  << " in file "
324  << this->dimensionedInternalField().objectPath()
325  << exit(FatalIOError);
326  }
327 
328  tensorFields_.insert(iter().keyword(), fPtr);
329  }
330  else
331  {
333  (
334  "genericFvPatchField<Type>::genericFvPatchField"
335  "(const fvPatch&, const Field<Type>&, "
336  "const dictionary&)",
337  dict
338  ) << "\n compound " << fieldToken.compoundToken()
339  << " not supported"
340  << "\n on patch " << this->patch().name()
341  << " of field "
342  << this->dimensionedInternalField().name()
343  << " in file "
344  << this->dimensionedInternalField().objectPath()
345  << exit(FatalIOError);
346  }
347  }
348  else if
349  (
350  firstToken.isWord()
351  && firstToken.wordToken() == "uniform"
352  )
353  {
354  token fieldToken(is);
355 
356  if (!fieldToken.isPunctuation())
357  {
358  scalarFields_.insert
359  (
360  iter().keyword(),
361  new scalarField
362  (
363  this->size(),
364  fieldToken.number()
365  )
366  );
367  }
368  else
369  {
370  // Read as scalarList.
371  is.putBack(fieldToken);
372 
373  scalarList l(is);
374 
375  if (l.size() == vector::nComponents)
376  {
377  vector vs(l[0], l[1], l[2]);
378 
379  vectorFields_.insert
380  (
381  iter().keyword(),
382  new vectorField(this->size(), vs)
383  );
384  }
385  else if (l.size() == sphericalTensor::nComponents)
386  {
387  sphericalTensor vs(l[0]);
388 
389  sphericalTensorFields_.insert
390  (
391  iter().keyword(),
392  new sphericalTensorField(this->size(), vs)
393  );
394  }
395  else if (l.size() == symmTensor::nComponents)
396  {
397  symmTensor vs(l[0], l[1], l[2], l[3], l[4], l[5]);
398 
399  symmTensorFields_.insert
400  (
401  iter().keyword(),
402  new symmTensorField(this->size(), vs)
403  );
404  }
405  else if (l.size() == tensor::nComponents)
406  {
407  tensor vs
408  (
409  l[0], l[1], l[2],
410  l[3], l[4], l[5],
411  l[6], l[7], l[8]
412  );
413 
414  tensorFields_.insert
415  (
416  iter().keyword(),
417  new tensorField(this->size(), vs)
418  );
419  }
420  else
421  {
423  (
424  "genericFvPatchField<Type>::genericFvPatchField"
425  "(const fvPatch&, const Field<Type>&, "
426  "const dictionary&)",
427  dict
428  ) << "\n unrecognised native type " << l
429  << "\n on patch " << this->patch().name()
430  << " of field "
431  << this->dimensionedInternalField().name()
432  << " in file "
433  << this->dimensionedInternalField().objectPath()
434  << exit(FatalIOError);
435  }
436  }
437  }
438  }
439  }
440  }
441 }
442 
443 
444 template<class Type>
446 (
447  const genericFvPatchField<Type>& ptf,
448  const fvPatch& p,
450  const fvPatchFieldMapper& mapper
451 )
452 :
453  calculatedFvPatchField<Type>(ptf, p, iF, mapper),
454  actualTypeName_(ptf.actualTypeName_),
455  dict_(ptf.dict_)
456 {
458  (
460  ptf.scalarFields_,
461  iter
462  )
463  {
464  scalarFields_.insert
465  (
466  iter.key(),
467  new scalarField(*iter(), mapper)
468  );
469  }
470 
472  (
474  ptf.vectorFields_,
475  iter
476  )
477  {
478  vectorFields_.insert
479  (
480  iter.key(),
481  new vectorField(*iter(), mapper)
482  );
483  }
484 
486  (
488  ptf.sphericalTensorFields_,
489  iter
490  )
491  {
492  sphericalTensorFields_.insert
493  (
494  iter.key(),
495  new sphericalTensorField(*iter(), mapper)
496  );
497  }
498 
500  (
502  ptf.symmTensorFields_,
503  iter
504  )
505  {
506  symmTensorFields_.insert
507  (
508  iter.key(),
509  new symmTensorField(*iter(), mapper)
510  );
511  }
512 
514  (
516  ptf.tensorFields_,
517  iter
518  )
519  {
520  tensorFields_.insert
521  (
522  iter.key(),
523  new tensorField(*iter(), mapper)
524  );
525  }
526 }
527 
528 
529 template<class Type>
531 (
532  const genericFvPatchField<Type>& ptf
533 )
534 :
536  actualTypeName_(ptf.actualTypeName_),
537  dict_(ptf.dict_),
538  scalarFields_(ptf.scalarFields_),
539  vectorFields_(ptf.vectorFields_),
540  sphericalTensorFields_(ptf.sphericalTensorFields_),
541  symmTensorFields_(ptf.symmTensorFields_),
542  tensorFields_(ptf.tensorFields_)
543 {}
544 
545 
546 template<class Type>
548 (
549  const genericFvPatchField<Type>& ptf,
551 )
552 :
554  actualTypeName_(ptf.actualTypeName_),
555  dict_(ptf.dict_),
556  scalarFields_(ptf.scalarFields_),
557  vectorFields_(ptf.vectorFields_),
558  sphericalTensorFields_(ptf.sphericalTensorFields_),
559  symmTensorFields_(ptf.symmTensorFields_),
560  tensorFields_(ptf.tensorFields_)
561 {}
562 
563 
564 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
565 
566 template<class Type>
568 (
569  const fvPatchFieldMapper& m
570 )
571 {
573 
574  forAllIter
575  (
577  scalarFields_,
578  iter
579  )
580  {
581  iter()->autoMap(m);
582  }
583 
584  forAllIter
585  (
587  vectorFields_,
588  iter
589  )
590  {
591  iter()->autoMap(m);
592  }
593 
594  forAllIter
595  (
597  sphericalTensorFields_,
598  iter
599  )
600  {
601  iter()->autoMap(m);
602  }
603 
604  forAllIter
605  (
607  symmTensorFields_,
608  iter
609  )
610  {
611  iter()->autoMap(m);
612  }
613 
614  forAllIter
615  (
617  tensorFields_,
618  iter
619  )
620  {
621  iter()->autoMap(m);
622  }
623 }
624 
625 
626 template<class Type>
628 (
629  const fvPatchField<Type>& ptf,
630  const labelList& addr
631 )
632 {
634 
635  const genericFvPatchField<Type>& dptf =
636  refCast<const genericFvPatchField<Type> >(ptf);
637 
638  forAllIter
639  (
641  scalarFields_,
642  iter
643  )
644  {
646  dptf.scalarFields_.find(iter.key());
647 
648  if (dptfIter != dptf.scalarFields_.end())
649  {
650  iter()->rmap(*dptfIter(), addr);
651  }
652  }
653 
654  forAllIter
655  (
657  vectorFields_,
658  iter
659  )
660  {
662  dptf.vectorFields_.find(iter.key());
663 
664  if (dptfIter != dptf.vectorFields_.end())
665  {
666  iter()->rmap(*dptfIter(), addr);
667  }
668  }
669 
670  forAllIter
671  (
673  sphericalTensorFields_,
674  iter
675  )
676  {
678  dptf.sphericalTensorFields_.find(iter.key());
679 
680  if (dptfIter != dptf.sphericalTensorFields_.end())
681  {
682  iter()->rmap(*dptfIter(), addr);
683  }
684  }
685 
686  forAllIter
687  (
689  symmTensorFields_,
690  iter
691  )
692  {
694  dptf.symmTensorFields_.find(iter.key());
695 
696  if (dptfIter != dptf.symmTensorFields_.end())
697  {
698  iter()->rmap(*dptfIter(), addr);
699  }
700  }
701 
702  forAllIter
703  (
705  tensorFields_,
706  iter
707  )
708  {
710  dptf.tensorFields_.find(iter.key());
711 
712  if (dptfIter != dptf.tensorFields_.end())
713  {
714  iter()->rmap(*dptfIter(), addr);
715  }
716  }
717 }
718 
719 
720 template<class Type>
723 (
724  const tmp<scalarField>&
725 ) const
726 {
728  (
729  "genericFvPatchField<Type>::"
730  "valueInternalCoeffs(const tmp<scalarField>&) const"
731  ) << "\n "
732  "valueInternalCoeffs cannot be called for a genericFvPatchField"
733  " (actual type " << actualTypeName_ << ")"
734  << "\n on patch " << this->patch().name()
735  << " of field " << this->dimensionedInternalField().name()
736  << " in file " << this->dimensionedInternalField().objectPath()
737  << "\n You are probably trying to solve for a field with a "
738  "generic boundary condition."
739  << exit(FatalError);
740 
741  return *this;
742 }
743 
744 
745 template<class Type>
748 (
749  const tmp<scalarField>&
750 ) const
751 {
753  (
754  "genericFvPatchField<Type>::"
755  "valueBoundaryCoeffs(const tmp<scalarField>&) const"
756  ) << "\n "
757  "valueBoundaryCoeffs cannot be called for a genericFvPatchField"
758  " (actual type " << actualTypeName_ << ")"
759  << "\n on patch " << this->patch().name()
760  << " of field " << this->dimensionedInternalField().name()
761  << " in file " << this->dimensionedInternalField().objectPath()
762  << "\n You are probably trying to solve for a field with a "
763  "generic boundary condition."
764  << exit(FatalError);
765 
766  return *this;
767 }
768 
769 
770 template<class Type>
773 {
775  (
776  "genericFvPatchField<Type>::"
777  "gradientInternalCoeffs() const"
778  ) << "\n "
779  "gradientInternalCoeffs cannot be called for a genericFvPatchField"
780  " (actual type " << actualTypeName_ << ")"
781  << "\n on patch " << this->patch().name()
782  << " of field " << this->dimensionedInternalField().name()
783  << " in file " << this->dimensionedInternalField().objectPath()
784  << "\n You are probably trying to solve for a field with a "
785  "generic boundary condition."
786  << exit(FatalError);
787 
788  return *this;
789 }
790 
791 template<class Type>
794 {
796  (
797  "genericFvPatchField<Type>::"
798  "gradientBoundaryCoeffs() const"
799  ) << "\n "
800  "gradientBoundaryCoeffs cannot be called for a genericFvPatchField"
801  " (actual type " << actualTypeName_ << ")"
802  << "\n on patch " << this->patch().name()
803  << " of field " << this->dimensionedInternalField().name()
804  << " in file " << this->dimensionedInternalField().objectPath()
805  << "\n You are probably trying to solve for a field with a "
806  "generic boundary condition."
807  << exit(FatalError);
808 
809  return *this;
810 }
811 
812 
813 template<class Type>
815 {
816  os.writeKeyword("type") << actualTypeName_ << token::END_STATEMENT << nl;
817 
818  forAllConstIter(dictionary, dict_, iter)
819  {
820  if (iter().keyword() != "type" && iter().keyword() != "value")
821  {
822  if
823  (
824  iter().isStream()
825  && iter().stream().size()
826  && iter().stream()[0].isWord()
827  && iter().stream()[0].wordToken() == "nonuniform"
828  )
829  {
830  if (scalarFields_.found(iter().keyword()))
831  {
832  scalarFields_.find(iter().keyword())()
833  ->writeEntry(iter().keyword(), os);
834  }
835  else if (vectorFields_.found(iter().keyword()))
836  {
837  vectorFields_.find(iter().keyword())()
838  ->writeEntry(iter().keyword(), os);
839  }
840  else if (sphericalTensorFields_.found(iter().keyword()))
841  {
842  sphericalTensorFields_.find(iter().keyword())()
843  ->writeEntry(iter().keyword(), os);
844  }
845  else if (symmTensorFields_.found(iter().keyword()))
846  {
847  symmTensorFields_.find(iter().keyword())()
848  ->writeEntry(iter().keyword(), os);
849  }
850  else if (tensorFields_.found(iter().keyword()))
851  {
852  tensorFields_.find(iter().keyword())()
853  ->writeEntry(iter().keyword(), os);
854  }
855  }
856  else
857  {
858  iter().write(os);
859  }
860  }
861  }
862 
863  this->writeEntry("value", os);
864 }
865 
866 
867 // ************************************************************************* //
scalar number() const
Definition: tokenI.H:345
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
bool isWord() const
Definition: tokenI.H:221
bool isLabel() const
Definition: tokenI.H:262
Input token stream.
Definition: ITstream.H:49
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
An STL-conforming const_iterator.
Definition: HashTable.H:470
genericFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Templated 3D symmetric tensor derived from VectorSpace adding construction from 6 components...
Definition: SymmTensor.H:53
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
#define forAllIter(Container, container, iter)
Definition: UList.H:440
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool isPunctuation() const
Definition: tokenI.H:203
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component, element access using th ii() member function and the inner-product (dot-product) and outer-product operators.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const word & wordToken() const
Definition: tokenI.H:226
const compound & compoundToken() const
Definition: tokenI.H:367
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
dictionary dict
tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
static const char nl
Definition: Ostream.H:260
IOerror FatalIOError
void putBack(const token &)
Put back token.
Definition: Istream.C:30
volVectorField vectorField(fieldObject, mesh)
volScalarField & p
Definition: createFields.H:51
This boundary condition provides a generic version of the calculated condition, useful as a fallback ...
A templated class for holding compound tokens.
Definition: token.H:213
A token holds items read from Istream.
Definition: token.H:67
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
compound & transferCompoundToken(const Istream &is)
Definition: token.C:93
label labelToken() const
Definition: tokenI.H:267
bool isCompound() const
Definition: tokenI.H:362
error FatalError
tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
To & dynamicCast(From &r)
Reference type cast template function,.
Definition: typeInfo.H:85
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual void write(Ostream &) const
Write.
volScalarField scalarField(fieldObject, mesh)
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
A class for managing temporary objects.
Definition: PtrList.H:118
Field< sphericalTensor > sphericalTensorField
Specialisation of Field<T> for sphericalTensor.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50