Fix tutorials: coupled/conjugateHeatFoam/conjugateCavity: fix Allrun file
[OpenFOAM-1.6-ext.git] / src / finiteVolume / fields / fvPatchFields / basic / sliced / slicedFvPatchField.C
blob5de12b3aaeb4ffcf140cd38ce4eeae86f78f31b2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "slicedFvPatchField.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
34 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
36 template<class Type>
37 slicedFvPatchField<Type>::slicedFvPatchField
39     const fvPatch& p,
40     const DimensionedField<Type, volMesh>& iF,
41     const Field<Type>& completeField
44     fvPatchField<Type>(p, iF, Field<Type>())
46     // Set the fvPatchField to a slice of the given complete field
47     UList<Type>::operator=(p.patchSlice(completeField));
51 template<class Type>
52 slicedFvPatchField<Type>::slicedFvPatchField
54     const fvPatch& p,
55     const DimensionedField<Type, volMesh>& iF
58     fvPatchField<Type>(p, iF, Field<Type>())
62 template<class Type>
63 slicedFvPatchField<Type>::slicedFvPatchField
65     const slicedFvPatchField<Type>& ptf,
66     const fvPatch& p,
67     const DimensionedField<Type, volMesh>& iF,
68     const fvPatchFieldMapper& mapper
71     fvPatchField<Type>(ptf, p, iF, mapper)
73     notImplemented
74     (
75         "slicedFvPatchField<Type>::"
76         "slicedFvPatchField(const slicedFvPatchField<Type>&, "
77         "const fvPatch&, const Field<Type>&, const fvPatchFieldMapper&)"
78     );
82 template<class Type>
83 slicedFvPatchField<Type>::slicedFvPatchField
85     const fvPatch& p,
86     const DimensionedField<Type, volMesh>& iF,
87     const dictionary& dict
90     fvPatchField<Type>(p, iF, dict)
92     notImplemented
93     (
94         "slicedFvPatchField<Type>::"
95         "slicedFvPatchField(const Field<Type>&, const dictionary&)"
96     );
100 template<class Type>
101 slicedFvPatchField<Type>::slicedFvPatchField
103     const slicedFvPatchField<Type>& ptf,
104     const DimensionedField<Type, volMesh>& iF
107     fvPatchField<Type>(ptf.patch(), iF, Field<Type>())
109     // Transfer the slice from the argument
110     UList<Type>::operator=(ptf);
113 template<class Type>
114 tmp<fvPatchField<Type> > slicedFvPatchField<Type>::clone() const
116     return tmp<fvPatchField<Type> >
117     (
118         new slicedFvPatchField<Type>(*this)
119     );
123 template<class Type>
124 slicedFvPatchField<Type>::slicedFvPatchField
126     const slicedFvPatchField<Type>& ptf
129     fvPatchField<Type>
130     (
131         ptf.patch(),
132         ptf.dimensionedInternalField(),
133         Field<Type>()
134     )
136     // Transfer the slice from the argument
137     UList<Type>::operator=(ptf);
141 template<class Type>
142 tmp<fvPatchField<Type> > slicedFvPatchField<Type>::clone
144     const DimensionedField<Type, volMesh>& iF
145 ) const
147     return tmp<fvPatchField<Type> >
148     (
149         new slicedFvPatchField<Type>(*this, iF)
150     );
154 template<class Type>
155 slicedFvPatchField<Type>::~slicedFvPatchField<Type>()
157     // Set the fvPatchField storage pointer to NULL before its destruction
158     // to protect the field it a slice of.
159     UList<Type>::operator=(UList<Type>(NULL, 0));
163 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
165 template<class Type>
166 tmp<Field<Type> > slicedFvPatchField<Type>::snGrad() const
168     notImplemented
169     (
170         "slicedFvPatchField<Type>::"
171         "snGrad()"
172     );
174     return Field<Type>::null();
178 template<class Type>
179 void slicedFvPatchField<Type>::updateCoeffs()
181     notImplemented
182     (
183         "slicedFvPatchField<Type>::"
184         "updateCoeffs()"
185     );
189 template<class Type>
190 tmp<Field<Type> > slicedFvPatchField<Type>::patchInternalField() const
192     notImplemented
193     (
194         "slicedFvPatchField<Type>::"
195         "patchInternalField()"
196     );
198     return Field<Type>::null();
202 template<class Type>
203 tmp<Field<Type> > slicedFvPatchField<Type>::patchNeighbourField
205     const Field<Type>& iField
206 ) const
208     notImplemented
209     (
210         "slicedFvPatchField<Type>::"
211         "patchNeighbourField(const DimensionedField<Type, volMesh>& iField)"
212     );
214     return Field<Type>::null();
218 template<class Type>
219 tmp<Field<Type> > slicedFvPatchField<Type>::patchNeighbourField() const
221     notImplemented
222     (
223         "slicedFvPatchField<Type>::"
224         "patchNeighbourField()"
225     );
227     return Field<Type>::null();
231 template<class Type>
232 tmp<Field<Type> > slicedFvPatchField<Type>::valueInternalCoeffs
234     const tmp<scalarField>&
235 ) const
237     notImplemented
238     (
239         "slicedFvPatchField<Type>::"
240         "valueInternalCoeffs(const tmp<scalarField>&)"
241     );
243     return Field<Type>::null();
247 template<class Type>
248 tmp<Field<Type> > slicedFvPatchField<Type>::valueBoundaryCoeffs
250     const tmp<scalarField>&
251 ) const
253     notImplemented
254     (
255         "slicedFvPatchField<Type>::"
256         "valueBoundaryCoeffs(const tmp<scalarField>&)"
257     );
259     return Field<Type>::null();
263 template<class Type>
264 tmp<Field<Type> > slicedFvPatchField<Type>::gradientInternalCoeffs() const
266     notImplemented
267     (
268         "slicedFvPatchField<Type>::"
269         "gradientInternalCoeffs()"
270     );
272     return Field<Type>::null();
276 template<class Type>
277 tmp<Field<Type> > slicedFvPatchField<Type>::gradientBoundaryCoeffs() const
279     notImplemented
280     (
281         "slicedFvPatchField<Type>::"
282         "gradientBoundaryCoeffs()"
283     );
285     return Field<Type>::null();
289 template<class Type>
290 void slicedFvPatchField<Type>::write(Ostream& os) const
292     fvPatchField<Type>::write(os);
293     this->writeEntry("value", os);
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 } // End namespace Foam
301 // ************************************************************************* //