slicedFvPatchField.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-2013 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 "slicedFvPatchField.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
40  const Field<Type>& completeField
41 )
42 :
44 {
45  // Set the fvPatchField to a slice of the given complete field
46  UList<Type>::operator=(p.patchSlice(completeField));
47 }
48 
49 
50 template<class Type>
52 (
53  const fvPatch& p,
55 )
56 :
58 {}
59 
60 
61 template<class Type>
63 (
64  const slicedFvPatchField<Type>& ptf,
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  fvPatchField<Type>(ptf, p, iF, mapper)
71 {
73  (
74  "slicedFvPatchField<Type>::"
75  "slicedFvPatchField(const slicedFvPatchField<Type>&, "
76  "const fvPatch&, const Field<Type>&, const fvPatchFieldMapper&)"
77  );
78 }
79 
80 
81 template<class Type>
83 (
84  const fvPatch& p,
86  const dictionary& dict
87 )
88 :
90 {
92  (
93  "slicedFvPatchField<Type>::"
94  "slicedFvPatchField(const Field<Type>&, const dictionary&)"
95  );
96 }
97 
98 
99 template<class Type>
101 (
102  const slicedFvPatchField<Type>& ptf,
104 )
105 :
106  fvPatchField<Type>(ptf.patch(), iF, Field<Type>())
107 {
108  // Transfer the slice from the argument
110 }
111 
112 template<class Type>
114 {
115  return tmp<fvPatchField<Type> >
116  (
117  new slicedFvPatchField<Type>(*this)
118  );
119 }
120 
121 
122 template<class Type>
124 (
125  const slicedFvPatchField<Type>& ptf
126 )
127 :
129  (
130  ptf.patch(),
132  Field<Type>()
133  )
134 {
135  // Transfer the slice from the argument
137 }
138 
139 
140 template<class Type>
142 (
144 ) const
145 {
146  return tmp<fvPatchField<Type> >
147  (
148  new slicedFvPatchField<Type>(*this, iF)
149  );
150 }
151 
152 
153 template<class Type>
154 slicedFvPatchField<Type>::~slicedFvPatchField<Type>()
155 {
156  // Set the fvPatchField storage pointer to NULL before its destruction
157  // to protect the field it a slice of.
159 }
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
164 template<class Type>
166 {
168  (
169  "slicedFvPatchField<Type>::"
170  "snGrad()"
171  );
172 
173  return Field<Type>::null();
174 }
175 
176 
177 template<class Type>
179 {
181  (
182  "slicedFvPatchField<Type>::"
183  "updateCoeffs()"
184  );
185 }
186 
187 
188 template<class Type>
190 {
192  (
193  "slicedFvPatchField<Type>::"
194  "patchInternalField()"
195  );
196 
197  return Field<Type>::null();
198 }
199 
200 
201 template<class Type>
203 {
205  (
206  "slicedFvPatchField<Type>::"
207  "patchInternalField(Field<Type>&)"
208  );
209 }
210 
211 
212 template<class Type>
214 (
215  const Field<Type>& iField
216 ) const
217 {
219  (
220  "slicedFvPatchField<Type>::"
221  "patchNeighbourField(const DimensionedField<Type, volMesh>& iField)"
222  );
223 
224  return Field<Type>::null();
225 }
226 
227 
228 template<class Type>
230 {
232  (
233  "slicedFvPatchField<Type>::"
234  "patchNeighbourField()"
235  );
236 
237  return Field<Type>::null();
238 }
239 
240 
241 template<class Type>
243 (
244  const tmp<scalarField>&
245 ) const
246 {
248  (
249  "slicedFvPatchField<Type>::"
250  "valueInternalCoeffs(const tmp<scalarField>&)"
251  );
252 
253  return Field<Type>::null();
254 }
255 
256 
257 template<class Type>
259 (
260  const tmp<scalarField>&
261 ) const
262 {
264  (
265  "slicedFvPatchField<Type>::"
266  "valueBoundaryCoeffs(const tmp<scalarField>&)"
267  );
268 
269  return Field<Type>::null();
270 }
271 
272 
273 template<class Type>
275 {
277  (
278  "slicedFvPatchField<Type>::"
279  "gradientInternalCoeffs()"
280  );
281 
282  return Field<Type>::null();
283 }
284 
285 
286 template<class Type>
288 {
290  (
291  "slicedFvPatchField<Type>::"
292  "gradientBoundaryCoeffs()"
293  );
294 
295  return Field<Type>::null();
296 }
297 
298 
299 template<class Type>
301 {
303  this->writeEntry("value", os);
304 }
305 
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 } // End namespace Foam
310 
311 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
Specialization of fvPatchField which creates the underlying fvPatchField as a slice of the given comp...
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Foam::fvPatchFieldMapper.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:354
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
static const Field< Type > & null()
Return a null field.
Definition: Field.H:100
dictionary dict
volScalarField & p
Definition: createFields.H:51
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
void operator=(const T &)
Assignment of all entries to the given value.
Definition: UList.C:70
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:307
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField of the values on the patch or on the.
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
slicedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &)
Construct from patch, internal field and field to slice.
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
A class for managing temporary objects.
Definition: PtrList.H:118