initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / fields / fvPatchFields / basic / sliced / slicedFvPatchField.C
blob5cf94bb064804fcc840c0eb900ce9a7c42e7d3db
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 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)
60     notImplemented
61     (
62         "slicedFvPatchField<Type>::"
63         "slicedFvPatchField(const fvPatch&, const Field<Type>&)"
64     );
68 template<class Type>
69 slicedFvPatchField<Type>::slicedFvPatchField
71     const slicedFvPatchField<Type>& ptf,
72     const fvPatch& p,
73     const DimensionedField<Type, volMesh>& iF,
74     const fvPatchFieldMapper& mapper
77     fvPatchField<Type>(ptf, p, iF, mapper)
79     notImplemented
80     (
81         "slicedFvPatchField<Type>::"
82         "slicedFvPatchField(const slicedFvPatchField<Type>&, "
83         "const fvPatch&, const Field<Type>&, const fvPatchFieldMapper&)"
84     );
88 template<class Type>
89 slicedFvPatchField<Type>::slicedFvPatchField
91     const fvPatch& p,
92     const DimensionedField<Type, volMesh>& iF,
93     const dictionary& dict
96     fvPatchField<Type>(p, iF, dict)
98     notImplemented
99     (
100         "slicedFvPatchField<Type>::"
101         "slicedFvPatchField(const Field<Type>&, const dictionary&)"
102     );
106 template<class Type>
107 slicedFvPatchField<Type>::slicedFvPatchField
109     const slicedFvPatchField<Type>& ptf,
110     const DimensionedField<Type, volMesh>& iF
113     fvPatchField<Type>(ptf.patch(), iF, Field<Type>())
115     // Transfer the slice from the argument
116     UList<Type>::operator=(ptf);
119 template<class Type>
120 tmp<fvPatchField<Type> > slicedFvPatchField<Type>::clone() const
122     return tmp<fvPatchField<Type> >
123     (
124         new slicedFvPatchField<Type>(*this)
125     );
129 template<class Type>
130 slicedFvPatchField<Type>::slicedFvPatchField
132     const slicedFvPatchField<Type>& ptf
135     fvPatchField<Type>
136     (
137         ptf.patch(),
138         ptf.dimensionedInternalField(),
139         Field<Type>()
140     )
142     // Transfer the slice from the argument
143     UList<Type>::operator=(ptf);
147 template<class Type>
148 tmp<fvPatchField<Type> > slicedFvPatchField<Type>::clone
150     const DimensionedField<Type, volMesh>& iF
151 ) const
153     return tmp<fvPatchField<Type> >
154     (
155         new slicedFvPatchField<Type>(*this, iF)
156     );
160 template<class Type>
161 slicedFvPatchField<Type>::~slicedFvPatchField<Type>()
163     // Set the fvPatchField storage pointer to NULL before its destruction
164     // to protect the field it a slice of.
165     UList<Type>::operator=(UList<Type>(NULL, 0));
169 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
171 template<class Type>
172 tmp<Field<Type> > slicedFvPatchField<Type>::snGrad() const
174     notImplemented
175     (
176         "slicedFvPatchField<Type>::"
177         "snGrad()"
178     );
180     return Field<Type>::null();
184 template<class Type>
185 void slicedFvPatchField<Type>::updateCoeffs()
187     notImplemented
188     (
189         "slicedFvPatchField<Type>::"
190         "updateCoeffs()"
191     );
195 template<class Type>
196 tmp<Field<Type> > slicedFvPatchField<Type>::patchInternalField() const
198     notImplemented
199     (
200         "slicedFvPatchField<Type>::"
201         "patchInternalField()"
202     );
204     return Field<Type>::null();
208 template<class Type>
209 tmp<Field<Type> > slicedFvPatchField<Type>::patchNeighbourField
211     const Field<Type>& iField
212 ) const
214     notImplemented
215     (
216         "slicedFvPatchField<Type>::"
217         "patchNeighbourField(const DimensionedField<Type, volMesh>& iField)"
218     );
220     return Field<Type>::null();
224 template<class Type>
225 tmp<Field<Type> > slicedFvPatchField<Type>::patchNeighbourField() const
227     notImplemented
228     (
229         "slicedFvPatchField<Type>::"
230         "patchNeighbourField()"
231     );
233     return Field<Type>::null();
237 template<class Type>
238 tmp<Field<Type> > slicedFvPatchField<Type>::valueInternalCoeffs
240     const tmp<scalarField>&
241 ) const
243     notImplemented
244     (
245         "slicedFvPatchField<Type>::"
246         "valueInternalCoeffs(const tmp<scalarField>&)"
247     );
249     return Field<Type>::null();
253 template<class Type>
254 tmp<Field<Type> > slicedFvPatchField<Type>::valueBoundaryCoeffs
256     const tmp<scalarField>&
257 ) const
259     notImplemented
260     (
261         "slicedFvPatchField<Type>::"
262         "valueBoundaryCoeffs(const tmp<scalarField>&)"
263     );
265     return Field<Type>::null();
269 template<class Type>
270 tmp<Field<Type> > slicedFvPatchField<Type>::gradientInternalCoeffs() const
272     notImplemented
273     (
274         "slicedFvPatchField<Type>::"
275         "gradientInternalCoeffs()"
276     );
278     return Field<Type>::null();
282 template<class Type>
283 tmp<Field<Type> > slicedFvPatchField<Type>::gradientBoundaryCoeffs() const
285     notImplemented
286     (
287         "slicedFvPatchField<Type>::"
288         "gradientBoundaryCoeffs()"
289     );
291     return Field<Type>::null();
295 template<class Type>
296 void slicedFvPatchField<Type>::write(Ostream& os) const
298     fvPatchField<Type>::write(os);
299     this->writeEntry("value", os);
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 } // End namespace Foam
307 // ************************************************************************* //