BUGFIX: Seg-fault in multiphaseInterFoam. Author: Henrik Rusche. Merge: Hrvoje Jasak
[foam-extend-3.2.git] / src / finiteVolume / fields / fvPatchFields / basic / generic / genericFvPatchField.C
blob85cc6b6494e2cfe2aeac2d4ad0d74731ddc38ea7
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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                "    or link the boundary-condition into libfoamUtil.so"
82             << exit(FatalIOError);
83     }
85     for
86     (
87         dictionary::const_iterator iter = dict_.begin();
88         iter != dict_.end();
89         ++iter
90     )
91     {
92         if (iter().keyword() != "type" && iter().keyword() != "value")
93         {
94             if
95             (
96                 iter().isStream()
97              && iter().stream().size()
98             )
99             {
100                 ITstream& is = iter().stream();
102                 // Read first token
103                 token firstToken(is);
105                 if
106                 (
107                     firstToken.isWord()
108                  && firstToken.wordToken() == "nonuniform"
109                 )
110                 {
111                     token fieldToken(is);
113                     if (!fieldToken.isCompound())
114                     {
115                         if
116                         (
117                             fieldToken.isLabel()
118                          && fieldToken.labelToken() == 0
119                         )
120                         {
121                             scalarFields_.insert
122                             (
123                                 iter().keyword(),
124                                 new scalarField(0)
125                             );
126                         }
127                         else
128                         {
129                             FatalIOErrorIn
130                             (
131                                 "genericFvPatchField<Type>::genericFvPatchField"
132                                 "(const fvPatch&, const Field<Type>&, "
133                                 "const dictionary&)",
134                                 dict
135                             )   << "\n    token following 'nonuniform' "
136                                   "is not a compound"
137                                 << "\n    on patch " << this->patch().name()
138                                 << " of field "
139                                 << this->dimensionedInternalField().name()
140                                 << " in file "
141                                 << this->dimensionedInternalField().objectPath()
142                             << exit(FatalIOError);
143                         }
144                     }
145                     else if
146                     (
147                         fieldToken.compoundToken().type()
148                      == token::Compound<List<scalar> >::typeName
149                     )
150                     {
151                         scalarField* fPtr = new scalarField;
152                         fPtr->transfer
153                         (
154                             dynamicCast<token::Compound<List<scalar> > >
155                             (
156                                 fieldToken.transferCompoundToken()
157                             )
158                         );
160                         if (fPtr->size() != this->size())
161                         {
162                             FatalIOErrorIn
163                             (
164                                 "genericFvPatchField<Type>::genericFvPatchField"
165                                 "(const fvPatch&, const Field<Type>&, "
166                                 "const dictionary&)",
167                                 dict
168                             )   << "\n    size of field " << iter().keyword()
169                                 << " (" << fPtr->size() << ')'
170                                 << " is not the same size as the patch ("
171                                 << this->size() << ')'
172                                 << "\n    on patch " << this->patch().name()
173                                 << " of field "
174                                 << this->dimensionedInternalField().name()
175                                 << " in file "
176                                 << this->dimensionedInternalField().objectPath()
177                                 << exit(FatalIOError);
178                         }
180                         scalarFields_.insert(iter().keyword(), fPtr);
181                     }
182                     else if
183                     (
184                         fieldToken.compoundToken().type()
185                      == token::Compound<List<vector> >::typeName
186                     )
187                     {
188                         vectorField* fPtr = new vectorField;
189                         fPtr->transfer
190                         (
191                             dynamicCast<token::Compound<List<vector> > >
192                             (
193                                 fieldToken.transferCompoundToken()
194                             )
195                         );
197                         if (fPtr->size() != this->size())
198                         {
199                             FatalIOErrorIn
200                             (
201                                 "genericFvPatchField<Type>::genericFvPatchField"
202                                 "(const fvPatch&, const Field<Type>&, "
203                                 "const dictionary&)",
204                                 dict
205                             )   << "\n    size of field " << iter().keyword()
206                                 << " (" << fPtr->size() << ')'
207                                 << " is not the same size as the patch ("
208                                 << this->size() << ')'
209                                 << "\n    on patch " << this->patch().name()
210                                 << " of field "
211                                 << this->dimensionedInternalField().name()
212                                 << " in file "
213                                 << this->dimensionedInternalField().objectPath()
214                                 << exit(FatalIOError);
215                         }
217                         vectorFields_.insert(iter().keyword(), fPtr);
218                     }
219                     else if
220                     (
221                         fieldToken.compoundToken().type()
222                      == token::Compound<List<sphericalTensor> >::typeName
223                     )
224                     {
225                         sphericalTensorField* fPtr = new sphericalTensorField;
226                         fPtr->transfer
227                         (
228                             dynamicCast
229                             <
230                                 token::Compound<List<sphericalTensor> >
231                             >
232                             (
233                                 fieldToken.transferCompoundToken()
234                             )
235                         );
237                         if (fPtr->size() != this->size())
238                         {
239                             FatalIOErrorIn
240                             (
241                                 "genericFvPatchField<Type>::genericFvPatchField"
242                                 "(const fvPatch&, const Field<Type>&, "
243                                 "const dictionary&)",
244                                 dict
245                             )   << "\n    size of field " << iter().keyword()
246                                 << " (" << fPtr->size() << ')'
247                                 << " is not the same size as the patch ("
248                                 << this->size() << ')'
249                                 << "\n    on patch " << this->patch().name()
250                                 << " of field "
251                                 << this->dimensionedInternalField().name()
252                                 << " in file "
253                                 << this->dimensionedInternalField().objectPath()
254                                 << exit(FatalIOError);
255                         }
257                         sphericalTensorFields_.insert(iter().keyword(), fPtr);
258                     }
259                     else if
260                     (
261                         fieldToken.compoundToken().type()
262                      == token::Compound<List<symmTensor> >::typeName
263                     )
264                     {
265                         symmTensorField* fPtr = new symmTensorField;
266                         fPtr->transfer
267                         (
268                             dynamicCast
269                             <
270                                 token::Compound<List<symmTensor> >
271                             >
272                             (
273                                 fieldToken.transferCompoundToken()
274                             )
275                         );
277                         if (fPtr->size() != this->size())
278                         {
279                             FatalIOErrorIn
280                             (
281                                 "genericFvPatchField<Type>::genericFvPatchField"
282                                 "(const fvPatch&, const Field<Type>&, "
283                                 "const dictionary&)",
284                                 dict
285                             )   << "\n    size of field " << iter().keyword()
286                                 << " (" << fPtr->size() << ')'
287                                 << " is not the same size as the patch ("
288                                 << this->size() << ')'
289                                 << "\n    on patch " << this->patch().name()
290                                 << " of field "
291                                 << this->dimensionedInternalField().name()
292                                 << " in file "
293                                 << this->dimensionedInternalField().objectPath()
294                                 << exit(FatalIOError);
295                         }
297                         symmTensorFields_.insert(iter().keyword(), fPtr);
298                     }
299                     else if
300                     (
301                         fieldToken.compoundToken().type()
302                      == token::Compound<List<tensor> >::typeName
303                     )
304                     {
305                         tensorField* fPtr = new tensorField;
306                         fPtr->transfer
307                         (
308                             dynamicCast<token::Compound<List<tensor> > >
309                             (
310                                 fieldToken.transferCompoundToken()
311                             )
312                         );
314                         if (fPtr->size() != this->size())
315                         {
316                             FatalIOErrorIn
317                             (
318                                 "genericFvPatchField<Type>::genericFvPatchField"
319                                 "(const fvPatch&, const Field<Type>&, "
320                                 "const dictionary&)",
321                                 dict
322                             )   << "\n    size of field " << iter().keyword()
323                                 << " (" << fPtr->size() << ')'
324                                 << " is not the same size as the patch ("
325                                 << this->size() << ')'
326                                 << "\n    on patch " << this->patch().name()
327                                 << " of field "
328                                 << this->dimensionedInternalField().name()
329                                 << " in file "
330                                 << this->dimensionedInternalField().objectPath()
331                                 << exit(FatalIOError);
332                         }
334                         tensorFields_.insert(iter().keyword(), fPtr);
335                     }
336                     else
337                     {
338                         FatalIOErrorIn
339                         (
340                             "genericFvPatchField<Type>::genericFvPatchField"
341                             "(const fvPatch&, const Field<Type>&, "
342                             "const dictionary&)",
343                             dict
344                         )   << "\n    compound " << fieldToken.compoundToken()
345                             << " not supported"
346                             << "\n    on patch " << this->patch().name()
347                             << " of field "
348                             << this->dimensionedInternalField().name()
349                             << " in file "
350                             << this->dimensionedInternalField().objectPath()
351                             << exit(FatalIOError);
352                     }
353                 }
354                 else if
355                 (
356                     firstToken.isWord()
357                  && firstToken.wordToken() == "uniform"
358                 )
359                 {
360                     token fieldToken(is);
362                     if (!fieldToken.isPunctuation())
363                     {
364                         scalarFields_.insert
365                         (
366                             iter().keyword(),
367                             new scalarField
368                             (
369                                 this->size(),
370                                 fieldToken.number()
371                             )
372                         );
373                     }
374                     else
375                     {
376                         // Read as scalarList.
377                         is.putBack(fieldToken);
379                         scalarList l(is);
381                         if (l.size() == vector::nComponents)
382                         {
383                             vector vs(l[0], l[1], l[2]);
385                             vectorFields_.insert
386                             (
387                                 iter().keyword(),
388                                 new vectorField(this->size(), vs)
389                             );
390                         }
391                         else if (l.size() == sphericalTensor::nComponents)
392                         {
393                             sphericalTensor vs(l[0]);
395                             sphericalTensorFields_.insert
396                             (
397                                 iter().keyword(),
398                                 new sphericalTensorField(this->size(), vs)
399                             );
400                         }
401                         else if (l.size() == symmTensor::nComponents)
402                         {
403                             symmTensor vs(l[0], l[1], l[2], l[3], l[4], l[5]);
405                             symmTensorFields_.insert
406                             (
407                                 iter().keyword(),
408                                 new symmTensorField(this->size(), vs)
409                             );
410                         }
411                         else if (l.size() == tensor::nComponents)
412                         {
413                             tensor vs
414                             (
415                                 l[0], l[1], l[2],
416                                 l[3], l[4], l[5],
417                                 l[6], l[7], l[8]
418                             );
420                             tensorFields_.insert
421                             (
422                                 iter().keyword(),
423                                 new tensorField(this->size(), vs)
424                             );
425                         }
426                         else
427                         {
428                             FatalIOErrorIn
429                             (
430                                 "genericFvPatchField<Type>::genericFvPatchField"
431                                 "(const fvPatch&, const Field<Type>&, "
432                                 "const dictionary&)",
433                                 dict
434                             )   << "\n    unrecognised native type " << l
435                                 << "\n    on patch " << this->patch().name()
436                                 << " of field "
437                                 << this->dimensionedInternalField().name()
438                                 << " in file "
439                                 << this->dimensionedInternalField().objectPath()
440                                 << exit(FatalIOError);
441                         }
442                     }
443                 }
444             }
445         }
446     }
450 template<class Type>
451 Foam::genericFvPatchField<Type>::genericFvPatchField
453     const genericFvPatchField<Type>& ptf,
454     const fvPatch& p,
455     const DimensionedField<Type, volMesh>& iF,
456     const fvPatchFieldMapper& mapper
459     calculatedFvPatchField<Type>(ptf, p, iF, mapper),
460     actualTypeName_(ptf.actualTypeName_),
461     dict_(ptf.dict_)
463     for
464     (
465         HashPtrTable<scalarField>::const_iterator iter =
466             ptf.scalarFields_.begin();
467         iter != ptf.scalarFields_.end();
468         ++iter
469     )
470     {
471         scalarFields_.insert(iter.key(), new scalarField(*iter(), mapper));
472     }
474     for
475     (
476         HashPtrTable<vectorField>::const_iterator iter =
477             ptf.vectorFields_.begin();
478         iter != ptf.vectorFields_.end();
479         ++iter
480     )
481     {
482         vectorFields_.insert(iter.key(), new vectorField(*iter(), mapper));
483     }
485     for
486     (
487         HashPtrTable<sphericalTensorField>::const_iterator iter =
488             ptf.sphericalTensorFields_.begin();
489         iter != ptf.sphericalTensorFields_.end();
490         ++iter
491     )
492     {
493         sphericalTensorFields_.insert
494         (
495             iter.key(),
496             new sphericalTensorField(*iter(), mapper)
497         );
498     }
500     for
501     (
502         HashPtrTable<symmTensorField>::const_iterator iter =
503             ptf.symmTensorFields_.begin();
504         iter != ptf.symmTensorFields_.end();
505         ++iter
506     )
507     {
508         symmTensorFields_.insert
509         (
510             iter.key(),
511             new symmTensorField(*iter(), mapper)
512         );
513     }
515     for
516     (
517         HashPtrTable<symmTensor4thOrderField>::const_iterator iter =
518             ptf.symmTensor4thOrderFields_.begin();
519         iter != ptf.symmTensor4thOrderFields_.end();
520         ++iter
521     )
522     {
523         symmTensor4thOrderFields_.insert
524         (
525             iter.key(),
526             new symmTensor4thOrderField(*iter(), mapper)
527         );
528     }
530     for
531     (
532         HashPtrTable<diagTensorField>::const_iterator iter =
533             ptf.diagTensorFields_.begin();
534         iter != ptf.diagTensorFields_.end();
535         ++iter
536     )
537     {
538         diagTensorFields_.insert
539         (
540             iter.key(),
541             new diagTensorField(*iter(), mapper)
542         );
543     }
545     for
546     (
547         HashPtrTable<tensorField>::const_iterator iter =
548             ptf.tensorFields_.begin();
549         iter != ptf.tensorFields_.end();
550         ++iter
551     )
552     {
553         tensorFields_.insert(iter.key(), new tensorField(*iter(), mapper));
554     }
558 template<class Type>
559 Foam::genericFvPatchField<Type>::genericFvPatchField
561     const genericFvPatchField<Type>& ptf
564     calculatedFvPatchField<Type>(ptf),
565     actualTypeName_(ptf.actualTypeName_),
566     dict_(ptf.dict_),
567     scalarFields_(ptf.scalarFields_),
568     vectorFields_(ptf.vectorFields_),
569     sphericalTensorFields_(ptf.sphericalTensorFields_),
570     symmTensorFields_(ptf.symmTensorFields_),
571     tensorFields_(ptf.tensorFields_)
575 template<class Type>
576 Foam::genericFvPatchField<Type>::genericFvPatchField
578     const genericFvPatchField<Type>& ptf,
579     const DimensionedField<Type, volMesh>& iF
582     calculatedFvPatchField<Type>(ptf, iF),
583     actualTypeName_(ptf.actualTypeName_),
584     dict_(ptf.dict_),
585     scalarFields_(ptf.scalarFields_),
586     vectorFields_(ptf.vectorFields_),
587     sphericalTensorFields_(ptf.sphericalTensorFields_),
588     symmTensorFields_(ptf.symmTensorFields_),
589     tensorFields_(ptf.tensorFields_)
593 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
595 template<class Type>
596 void Foam::genericFvPatchField<Type>::autoMap
598     const fvPatchFieldMapper& m
601     calculatedFvPatchField<Type>::autoMap(m);
603     for
604     (
605         HashPtrTable<scalarField>::iterator iter = scalarFields_.begin();
606         iter != scalarFields_.end();
607         ++iter
608     )
609     {
610         iter()->autoMap(m);
611     }
613     for
614     (
615         HashPtrTable<vectorField>::iterator iter = vectorFields_.begin();
616         iter != vectorFields_.end();
617         ++iter
618     )
619     {
620         iter()->autoMap(m);
621     }
623     for
624     (
625         HashPtrTable<sphericalTensorField>::iterator iter =
626             sphericalTensorFields_.begin();
627         iter != sphericalTensorFields_.end();
628         ++iter
629     )
630     {
631         iter()->autoMap(m);
632     }
634     for
635     (
636         HashPtrTable<symmTensorField>::iterator iter =
637             symmTensorFields_.begin();
638         iter != symmTensorFields_.end();
639         ++iter
640     )
641     {
642         iter()->autoMap(m);
643     }
645     for
646     (
647         HashPtrTable<tensorField>::iterator iter = tensorFields_.begin();
648         iter != tensorFields_.end();
649         ++iter
650     )
651     {
652         iter()->autoMap(m);
653     }
657 template<class Type>
658 void Foam::genericFvPatchField<Type>::rmap
660     const fvPatchField<Type>& ptf,
661     const labelList& addr
664     calculatedFvPatchField<Type>::rmap(ptf, addr);
666     const genericFvPatchField<Type>& dptf =
667         refCast<const genericFvPatchField<Type> >(ptf);
669     for
670     (
671         HashPtrTable<scalarField>::iterator iter = scalarFields_.begin();
672         iter != scalarFields_.end();
673         ++iter
674     )
675     {
676         HashPtrTable<scalarField>::const_iterator dptfIter =
677             dptf.scalarFields_.find(iter.key());
679         if (dptfIter != dptf.scalarFields_.end())
680         {
681             iter()->rmap(*dptfIter(), addr);
682         }
683     }
685     for
686     (
687         HashPtrTable<vectorField>::iterator iter = vectorFields_.begin();
688         iter != vectorFields_.end();
689         ++iter
690     )
691     {
692         HashPtrTable<vectorField>::const_iterator dptfIter =
693             dptf.vectorFields_.find(iter.key());
695         if (dptfIter != dptf.vectorFields_.end())
696         {
697             iter()->rmap(*dptfIter(), addr);
698         }
699     }
701     for
702     (
703         HashPtrTable<sphericalTensorField>::iterator iter =
704             sphericalTensorFields_.begin();
705         iter != sphericalTensorFields_.end();
706         ++iter
707     )
708     {
709         HashPtrTable<sphericalTensorField>::const_iterator dptfIter =
710             dptf.sphericalTensorFields_.find(iter.key());
712         if (dptfIter != dptf.sphericalTensorFields_.end())
713         {
714             iter()->rmap(*dptfIter(), addr);
715         }
716     }
718     for
719     (
720         HashPtrTable<symmTensorField>::iterator iter =
721             symmTensorFields_.begin();
722         iter != symmTensorFields_.end();
723         ++iter
724     )
725     {
726         HashPtrTable<symmTensorField>::const_iterator dptfIter =
727             dptf.symmTensorFields_.find(iter.key());
729         if (dptfIter != dptf.symmTensorFields_.end())
730         {
731             iter()->rmap(*dptfIter(), addr);
732         }
733     }
735     for
736     (
737         HashPtrTable<tensorField>::iterator iter = tensorFields_.begin();
738         iter != tensorFields_.end();
739         ++iter
740     )
741     {
742         HashPtrTable<tensorField>::const_iterator dptfIter =
743             dptf.tensorFields_.find(iter.key());
745         if (dptfIter != dptf.tensorFields_.end())
746         {
747             iter()->rmap(*dptfIter(), addr);
748         }
749     }
753 template<class Type>
754 Foam::tmp<Foam::Field<Type> >
755 Foam::genericFvPatchField<Type>::valueInternalCoeffs
757     const tmp<scalarField>&
758 ) const
760     FatalErrorIn
761     (
762         "genericFvPatchField<Type>::"
763         "valueInternalCoeffs(const tmp<scalarField>&) const"
764     )   << "\n    "
765            "valueInternalCoeffs cannot be called for a genericFvPatchField"
766            " (actual type " << actualTypeName_ << ")"
767         << "\n    on patch " << this->patch().name()
768         << " of field " << this->dimensionedInternalField().name()
769         << " in file " << this->dimensionedInternalField().objectPath()
770         << "\n    You are probably trying to solve for a field with a "
771            "generic boundary condition."
772         << exit(FatalError);
774     return *this;
778 template<class Type>
779 Foam::tmp<Foam::Field<Type> >
780 Foam::genericFvPatchField<Type>::valueBoundaryCoeffs
782     const tmp<scalarField>&
783 ) const
785     FatalErrorIn
786     (
787         "genericFvPatchField<Type>::"
788         "valueBoundaryCoeffs(const tmp<scalarField>&) const"
789     )   << "\n    "
790            "valueBoundaryCoeffs cannot be called for a genericFvPatchField"
791            " (actual type " << actualTypeName_ << ")"
792         << "\n    on patch " << this->patch().name()
793         << " of field " << this->dimensionedInternalField().name()
794         << " in file " << this->dimensionedInternalField().objectPath()
795         << "\n    You are probably trying to solve for a field with a "
796            "generic boundary condition."
797         << exit(FatalError);
799     return *this;
803 template<class Type>
804 Foam::tmp<Foam::Field<Type> >
805 Foam::genericFvPatchField<Type>::gradientInternalCoeffs() const
807     FatalErrorIn
808     (
809         "genericFvPatchField<Type>::"
810         "gradientInternalCoeffs() const"
811     )   << "\n    "
812            "gradientInternalCoeffs cannot be called for a genericFvPatchField"
813            " (actual type " << actualTypeName_ << ")"
814         << "\n    on patch " << this->patch().name()
815         << " of field " << this->dimensionedInternalField().name()
816         << " in file " << this->dimensionedInternalField().objectPath()
817         << "\n    You are probably trying to solve for a field with a "
818            "generic boundary condition."
819         << exit(FatalError);
821     return *this;
824 template<class Type>
825 Foam::tmp<Foam::Field<Type> >
826 Foam::genericFvPatchField<Type>::gradientBoundaryCoeffs() const
828     FatalErrorIn
829     (
830         "genericFvPatchField<Type>::"
831         "gradientBoundaryCoeffs() const"
832     )   << "\n    "
833            "gradientBoundaryCoeffs cannot be called for a genericFvPatchField"
834            " (actual type " << actualTypeName_ << ")"
835         << "\n    on patch " << this->patch().name()
836         << " of field " << this->dimensionedInternalField().name()
837         << " in file " << this->dimensionedInternalField().objectPath()
838         << "\n    You are probably trying to solve for a field with a "
839            "generic boundary condition."
840         << exit(FatalError);
842     return *this;
846 template<class Type>
847 void Foam::genericFvPatchField<Type>::write(Ostream& os) const
849     os.writeKeyword("type") << actualTypeName_ << token::END_STATEMENT << nl;
851     for
852     (
853         dictionary::const_iterator iter = dict_.begin();
854         iter != dict_.end();
855         ++iter
856     )
857     {
858         if (iter().keyword() != "type" && iter().keyword() != "value")
859         {
860             if
861             (
862                 iter().isStream()
863              && iter().stream().size()
864              && iter().stream()[0].isWord()
865              && iter().stream()[0].wordToken() == "nonuniform"
866             )
867             {
868                 if (scalarFields_.found(iter().keyword()))
869                 {
870                     scalarFields_.find(iter().keyword())()
871                         ->writeEntry(iter().keyword(), os);
872                 }
873                 else if (vectorFields_.found(iter().keyword()))
874                 {
875                     vectorFields_.find(iter().keyword())()
876                         ->writeEntry(iter().keyword(), os);
877                 }
878                 else if (sphericalTensorFields_.found(iter().keyword()))
879                 {
880                     sphericalTensorFields_.find(iter().keyword())()
881                         ->writeEntry(iter().keyword(), os);
882                 }
883                 else if (symmTensorFields_.found(iter().keyword()))
884                 {
885                     symmTensorFields_.find(iter().keyword())()
886                         ->writeEntry(iter().keyword(), os);
887                 }
888                 else if (tensorFields_.found(iter().keyword()))
889                 {
890                     tensorFields_.find(iter().keyword())()
891                         ->writeEntry(iter().keyword(), os);
892                 }
893             }
894             else
895             {
896                iter().write(os);
897             }
898         }
899     }
901     this->writeEntry("value", os);
905 // ************************************************************************* //