Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
finiteVolume
fields
fvPatchFields
basic
sliced
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) 1991-2010 OpenCFD Ltd.
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 <
finiteVolume/slicedFvPatchField.H
>
27
28
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29
30
namespace
Foam
31
{
32
33
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34
35
template
<
class
Type>
36
slicedFvPatchField<Type>::slicedFvPatchField
37
(
38
const
fvPatch
&
p
,
39
const
DimensionedField<Type, volMesh>
& iF,
40
const
Field<Type>
& completeField
41
)
42
:
43
fvPatchField<Type>
(
p
, iF,
Field<Type>
())
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>
51
slicedFvPatchField<Type>::slicedFvPatchField
52
(
53
const
fvPatch
&
p
,
54
const
DimensionedField<Type, volMesh>
& iF
55
)
56
:
57
fvPatchField<Type>
(
p
, iF,
Field<Type>
())
58
{}
59
60
61
template
<
class
Type>
62
slicedFvPatchField<Type>::slicedFvPatchField
63
(
64
const
slicedFvPatchField<Type>
& ptf,
65
const
fvPatch
&
p
,
66
const
DimensionedField<Type, volMesh>
& iF,
67
const
fvPatchFieldMapper
& mapper
68
)
69
:
70
fvPatchField<Type>
(ptf,
p
, iF, mapper)
71
{
72
notImplemented
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>
82
slicedFvPatchField<Type>::slicedFvPatchField
83
(
84
const
fvPatch
& p,
85
const
DimensionedField<Type, volMesh>
& iF,
86
const
dictionary
& dict
87
)
88
:
89
fvPatchField<Type>
(
p
, iF, dict)
90
{
91
notImplemented
92
(
93
"slicedFvPatchField<Type>::"
94
"slicedFvPatchField(const Field<Type>&, const dictionary&)"
95
);
96
}
97
98
99
template
<
class
Type>
100
slicedFvPatchField<Type>::slicedFvPatchField
101
(
102
const
slicedFvPatchField<Type>
& ptf,
103
const
DimensionedField<Type, volMesh>
& iF
104
)
105
:
106
fvPatchField<Type>
(ptf.
patch
(), iF,
Field<Type>
())
107
{
108
// Transfer the slice from the argument
109
UList<Type>::operator=
(ptf);
110
}
111
112
template
<
class
Type>
113
tmp<fvPatchField<Type>
>
slicedFvPatchField<Type>::clone
()
const
114
{
115
return
tmp<fvPatchField<Type>
>
116
(
117
new
slicedFvPatchField<Type>
(*this)
118
);
119
}
120
121
122
template
<
class
Type>
123
slicedFvPatchField<Type>::slicedFvPatchField
124
(
125
const
slicedFvPatchField<Type>
& ptf
126
)
127
:
128
fvPatchField<Type>
129
(
130
ptf.
patch
(),
131
ptf.
dimensionedInternalField
(),
132
Field<Type>
()
133
)
134
{
135
// Transfer the slice from the argument
136
UList<Type>::operator=
(ptf);
137
}
138
139
140
template
<
class
Type>
141
tmp<fvPatchField<Type>
>
slicedFvPatchField<Type>::clone
142
(
143
const
DimensionedField<Type, volMesh>
& iF
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.
158
UList<Type>::operator=
(
UList<Type>
(NULL, 0));
159
}
160
161
162
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163
164
template
<
class
Type>
165
tmp<Field<Type>
>
slicedFvPatchField<Type>::snGrad
()
const
166
{
167
notImplemented
168
(
169
"slicedFvPatchField<Type>::"
170
"snGrad()"
171
);
172
173
return
Field<Type>::null
();
174
}
175
176
177
template
<
class
Type>
178
void
slicedFvPatchField<Type>::updateCoeffs
()
179
{
180
notImplemented
181
(
182
"slicedFvPatchField<Type>::"
183
"updateCoeffs()"
184
);
185
}
186
187
188
template
<
class
Type>
189
tmp<Field<Type>
>
slicedFvPatchField<Type>::patchInternalField
()
const
190
{
191
notImplemented
192
(
193
"slicedFvPatchField<Type>::"
194
"patchInternalField()"
195
);
196
197
return
Field<Type>::null
();
198
}
199
200
201
template
<
class
Type>
202
tmp<Field<Type>
>
slicedFvPatchField<Type>::patchNeighbourField
203
(
204
const
Field<Type>
& iField
205
)
const
206
{
207
notImplemented
208
(
209
"slicedFvPatchField<Type>::"
210
"patchNeighbourField(const DimensionedField<Type, volMesh>& iField)"
211
);
212
213
return
Field<Type>::null
();
214
}
215
216
217
template
<
class
Type>
218
tmp<Field<Type>
>
slicedFvPatchField<Type>::patchNeighbourField
()
const
219
{
220
notImplemented
221
(
222
"slicedFvPatchField<Type>::"
223
"patchNeighbourField()"
224
);
225
226
return
Field<Type>::null
();
227
}
228
229
230
template
<
class
Type>
231
tmp<Field<Type>
>
slicedFvPatchField<Type>::valueInternalCoeffs
232
(
233
const
tmp<scalarField>
&
234
)
const
235
{
236
notImplemented
237
(
238
"slicedFvPatchField<Type>::"
239
"valueInternalCoeffs(const tmp<scalarField>&)"
240
);
241
242
return
Field<Type>::null
();
243
}
244
245
246
template
<
class
Type>
247
tmp<Field<Type>
>
slicedFvPatchField<Type>::valueBoundaryCoeffs
248
(
249
const
tmp<scalarField>
&
250
)
const
251
{
252
notImplemented
253
(
254
"slicedFvPatchField<Type>::"
255
"valueBoundaryCoeffs(const tmp<scalarField>&)"
256
);
257
258
return
Field<Type>::null
();
259
}
260
261
262
template
<
class
Type>
263
tmp<Field<Type>
>
slicedFvPatchField<Type>::gradientInternalCoeffs
()
const
264
{
265
notImplemented
266
(
267
"slicedFvPatchField<Type>::"
268
"gradientInternalCoeffs()"
269
);
270
271
return
Field<Type>::null
();
272
}
273
274
275
template
<
class
Type>
276
tmp<Field<Type>
>
slicedFvPatchField<Type>::gradientBoundaryCoeffs
()
const
277
{
278
notImplemented
279
(
280
"slicedFvPatchField<Type>::"
281
"gradientBoundaryCoeffs()"
282
);
283
284
return
Field<Type>::null
();
285
}
286
287
288
template
<
class
Type>
289
void
slicedFvPatchField<Type>::write
(
Ostream
& os)
const
290
{
291
fvPatchField<Type>::write
(os);
292
this->writeEntry(
"value"
, os);
293
}
294
295
296
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297
298
}
// End namespace Foam
299
300
// ************************ vim: set sw=4 sts=4 et: ************************ //