Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / genericPatchFields / genericFvPatchField / genericFvPatchField.C
blob7587a7be485fb564867f5bd67bf1084464da578b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2009-2011 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "genericFvPatchField.H"
27 #include "fvPatchFieldMapper.H"
29 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
31 template<class Type>
32 Foam::genericFvPatchField<Type>::genericFvPatchField
34     const fvPatch& p,
35     const DimensionedField<Type, volMesh>& iF
38     calculatedFvPatchField<Type>(p, iF)
40     FatalErrorIn
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);
52 template<class Type>
53 Foam::genericFvPatchField<Type>::genericFvPatchField
55     const fvPatch& p,
56     const DimensionedField<Type, volMesh>& iF,
57     const dictionary& dict
60     calculatedFvPatchField<Type>(p, iF, dict, false),
61     actualTypeName_(dict.lookup("type")),
62     dict_(dict)
64     if (!dict.found("value"))
65     {
66         FatalIOErrorIn
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     }
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();
96                 // Read first token
97                 token firstToken(is);
99                 if
100                 (
101                     firstToken.isWord()
102                  && firstToken.wordToken() == "nonuniform"
103                 )
104                 {
105                     token fieldToken(is);
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                         {
123                             FatalIOErrorIn
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                         (
148                             dynamicCast<token::Compound<List<scalar> > >
149                             (
150                                 fieldToken.transferCompoundToken()
151                             )
152                         );
154                         if (fPtr->size() != this->size())
155                         {
156                             FatalIOErrorIn
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                         }
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                         (
185                             dynamicCast<token::Compound<List<vector> > >
186                             (
187                                 fieldToken.transferCompoundToken()
188                             )
189                         );
191                         if (fPtr->size() != this->size())
192                         {
193                             FatalIOErrorIn
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                         }
211                         vectorFields_.insert(iter().keyword(), fPtr);
212                     }
213                     else if
214                     (
215                         fieldToken.compoundToken().type()
216                      == token::Compound<List<sphericalTensor> >::typeName
217                     )
218                     {
219                         sphericalTensorField* fPtr = new sphericalTensorField;
220                         fPtr->transfer
221                         (
222                             dynamicCast
223                             <
224                                 token::Compound<List<sphericalTensor> >
225                             >
226                             (
227                                 fieldToken.transferCompoundToken()
228                             )
229                         );
231                         if (fPtr->size() != this->size())
232                         {
233                             FatalIOErrorIn
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                         }
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                         (
262                             dynamicCast
263                             <
264                                 token::Compound<List<symmTensor> >
265                             >
266                             (
267                                 fieldToken.transferCompoundToken()
268                             )
269                         );
271                         if (fPtr->size() != this->size())
272                         {
273                             FatalIOErrorIn
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                         }
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                         (
302                             dynamicCast<token::Compound<List<tensor> > >
303                             (
304                                 fieldToken.transferCompoundToken()
305                             )
306                         );
308                         if (fPtr->size() != this->size())
309                         {
310                             FatalIOErrorIn
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                         }
328                         tensorFields_.insert(iter().keyword(), fPtr);
329                     }
330                     else
331                     {
332                         FatalIOErrorIn
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);
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);
373                         scalarList l(is);
375                         if (l.size() == vector::nComponents)
376                         {
377                             vector vs(l[0], l[1], l[2]);
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]);
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]);
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                             );
414                             tensorFields_.insert
415                             (
416                                 iter().keyword(),
417                                 new tensorField(this->size(), vs)
418                             );
419                         }
420                         else
421                         {
422                             FatalIOErrorIn
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     }
444 template<class Type>
445 Foam::genericFvPatchField<Type>::genericFvPatchField
447     const genericFvPatchField<Type>& ptf,
448     const fvPatch& p,
449     const DimensionedField<Type, volMesh>& iF,
450     const fvPatchFieldMapper& mapper
453     calculatedFvPatchField<Type>(ptf, p, iF, mapper),
454     actualTypeName_(ptf.actualTypeName_),
455     dict_(ptf.dict_)
457     forAllConstIter
458     (
459         HashPtrTable<scalarField>,
460         ptf.scalarFields_,
461         iter
462     )
463     {
464         scalarFields_.insert
465         (
466             iter.key(),
467             new scalarField(*iter(), mapper)
468         );
469     }
471     forAllConstIter
472     (
473         HashPtrTable<vectorField>,
474         ptf.vectorFields_,
475         iter
476     )
477     {
478         vectorFields_.insert
479         (
480             iter.key(),
481             new vectorField(*iter(), mapper)
482         );
483     }
485     forAllConstIter
486     (
487         HashPtrTable<sphericalTensorField>,
488         ptf.sphericalTensorFields_,
489         iter
490     )
491     {
492         sphericalTensorFields_.insert
493         (
494             iter.key(),
495             new sphericalTensorField(*iter(), mapper)
496         );
497     }
499     forAllConstIter
500     (
501         HashPtrTable<symmTensorField>,
502         ptf.symmTensorFields_,
503         iter
504     )
505     {
506         symmTensorFields_.insert
507         (
508             iter.key(),
509             new symmTensorField(*iter(), mapper)
510         );
511     }
513     forAllConstIter
514     (
515         HashPtrTable<tensorField>,
516         ptf.tensorFields_,
517         iter
518     )
519     {
520         tensorFields_.insert
521         (
522             iter.key(),
523             new tensorField(*iter(), mapper)
524         );
525     }
529 template<class Type>
530 Foam::genericFvPatchField<Type>::genericFvPatchField
532     const genericFvPatchField<Type>& ptf
535     calculatedFvPatchField<Type>(ptf),
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_)
546 template<class Type>
547 Foam::genericFvPatchField<Type>::genericFvPatchField
549     const genericFvPatchField<Type>& ptf,
550     const DimensionedField<Type, volMesh>& iF
553     calculatedFvPatchField<Type>(ptf, iF),
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_)
564 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
566 template<class Type>
567 void Foam::genericFvPatchField<Type>::autoMap
569     const fvPatchFieldMapper& m
572     calculatedFvPatchField<Type>::autoMap(m);
574     forAllIter
575     (
576         HashPtrTable<scalarField>,
577         scalarFields_,
578         iter
579     )
580     {
581         iter()->autoMap(m);
582     }
584     forAllIter
585     (
586         HashPtrTable<vectorField>,
587         vectorFields_,
588         iter
589     )
590     {
591         iter()->autoMap(m);
592     }
594     forAllIter
595     (
596         HashPtrTable<sphericalTensorField>,
597         sphericalTensorFields_,
598         iter
599     )
600     {
601         iter()->autoMap(m);
602     }
604     forAllIter
605     (
606         HashPtrTable<symmTensorField>,
607         symmTensorFields_,
608         iter
609     )
610     {
611         iter()->autoMap(m);
612     }
614     forAllIter
615     (
616         HashPtrTable<tensorField>,
617         tensorFields_,
618         iter
619     )
620     {
621         iter()->autoMap(m);
622     }
626 template<class Type>
627 void Foam::genericFvPatchField<Type>::rmap
629     const fvPatchField<Type>& ptf,
630     const labelList& addr
633     calculatedFvPatchField<Type>::rmap(ptf, addr);
635     const genericFvPatchField<Type>& dptf =
636         refCast<const genericFvPatchField<Type> >(ptf);
638     forAllIter
639     (
640         HashPtrTable<scalarField>,
641         scalarFields_,
642         iter
643     )
644     {
645         HashPtrTable<scalarField>::const_iterator dptfIter =
646             dptf.scalarFields_.find(iter.key());
648         if (dptfIter != dptf.scalarFields_.end())
649         {
650             iter()->rmap(*dptfIter(), addr);
651         }
652     }
654     forAllIter
655     (
656         HashPtrTable<vectorField>,
657         vectorFields_,
658         iter
659     )
660     {
661         HashPtrTable<vectorField>::const_iterator dptfIter =
662             dptf.vectorFields_.find(iter.key());
664         if (dptfIter != dptf.vectorFields_.end())
665         {
666             iter()->rmap(*dptfIter(), addr);
667         }
668     }
670     forAllIter
671     (
672         HashPtrTable<sphericalTensorField>,
673         sphericalTensorFields_,
674         iter
675     )
676     {
677         HashPtrTable<sphericalTensorField>::const_iterator dptfIter =
678             dptf.sphericalTensorFields_.find(iter.key());
680         if (dptfIter != dptf.sphericalTensorFields_.end())
681         {
682             iter()->rmap(*dptfIter(), addr);
683         }
684     }
686     forAllIter
687     (
688         HashPtrTable<symmTensorField>,
689         symmTensorFields_,
690         iter
691     )
692     {
693         HashPtrTable<symmTensorField>::const_iterator dptfIter =
694             dptf.symmTensorFields_.find(iter.key());
696         if (dptfIter != dptf.symmTensorFields_.end())
697         {
698             iter()->rmap(*dptfIter(), addr);
699         }
700     }
702     forAllIter
703     (
704         HashPtrTable<tensorField>,
705         tensorFields_,
706         iter
707     )
708     {
709         HashPtrTable<tensorField>::const_iterator dptfIter =
710             dptf.tensorFields_.find(iter.key());
712         if (dptfIter != dptf.tensorFields_.end())
713         {
714             iter()->rmap(*dptfIter(), addr);
715         }
716     }
720 template<class Type>
721 Foam::tmp<Foam::Field<Type> >
722 Foam::genericFvPatchField<Type>::valueInternalCoeffs
724     const tmp<scalarField>&
725 ) const
727     FatalErrorIn
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);
741     return *this;
745 template<class Type>
746 Foam::tmp<Foam::Field<Type> >
747 Foam::genericFvPatchField<Type>::valueBoundaryCoeffs
749     const tmp<scalarField>&
750 ) const
752     FatalErrorIn
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);
766     return *this;
770 template<class Type>
771 Foam::tmp<Foam::Field<Type> >
772 Foam::genericFvPatchField<Type>::gradientInternalCoeffs() const
774     FatalErrorIn
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);
788     return *this;
791 template<class Type>
792 Foam::tmp<Foam::Field<Type> >
793 Foam::genericFvPatchField<Type>::gradientBoundaryCoeffs() const
795     FatalErrorIn
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);
809     return *this;
813 template<class Type>
814 void Foam::genericFvPatchField<Type>::write(Ostream& os) const
816     os.writeKeyword("type") << actualTypeName_ << token::END_STATEMENT << nl;
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     }
863     this->writeEntry("value", os);
867 // ************************************************************************* //