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-2016 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const fvPatch& p,
35  const Field<Type>& completeField
36 )
37 :
39 {
40  // Set the fvPatchField to a slice of the given complete field
41  UList<Type>::shallowCopy(p.patchSlice(completeField));
42 }
43 
44 
45 template<class Type>
47 (
48  const fvPatch& p,
50 )
51 :
53 {}
54 
55 
56 template<class Type>
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
65 {
67 }
68 
69 
70 template<class Type>
72 (
73  const slicedFvPatchField<Type>& ptf,
74  const fvPatch& p,
76  const fvPatchFieldMapper& mapper
77 )
78 :
79  fvPatchField<Type>(ptf, p, iF, mapper)
80 {
82 }
83 
84 
85 template<class Type>
87 (
88  const slicedFvPatchField<Type>& ptf,
90 )
91 :
93 {
94  // Transfer the slice from the argument
96 }
97 
98 
99 template<class Type>
102 {
103  return tmp<fvPatchField<Type>>
104  (
105  new slicedFvPatchField<Type>(*this)
106  );
107 }
108 
109 
110 template<class Type>
112 (
113  const slicedFvPatchField<Type>& ptf
114 )
115 :
117  (
118  ptf.patch(),
119  ptf.internalField(),
120  Field<Type>()
121  )
122 {
123  // Transfer the slice from the argument
125 }
126 
127 
128 template<class Type>
131 (
133 ) const
134 {
135  return tmp<fvPatchField<Type>>
136  (
137  new slicedFvPatchField<Type>(*this, iF)
138  );
139 }
140 
141 
142 template<class Type>
144 {
145  // Set the fvPatchField storage pointer to NULL before its destruction
146  // to protect the field it a slice of.
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
153 template<class Type>
155 {
157 
158  return Field<Type>::null();
159 }
160 
161 
162 template<class Type>
164 {
166 }
167 
168 
169 template<class Type>
172 {
174 
175  return Field<Type>::null();
176 }
177 
178 
179 template<class Type>
181 {
183 }
184 
185 
186 template<class Type>
189 (
190  const Field<Type>& iField
191 ) const
192 {
194 
195  return Field<Type>::null();
196 }
197 
198 
199 template<class Type>
202 {
204 
205  return Field<Type>::null();
206 }
207 
208 
209 template<class Type>
212 (
213  const tmp<scalarField>&
214 ) const
215 {
217 
218  return Field<Type>::null();
219 }
220 
221 
222 template<class Type>
225 (
226  const tmp<scalarField>&
227 ) const
228 {
230 
231  return Field<Type>::null();
232 }
233 
234 
235 template<class Type>
238 {
240 
241  return Field<Type>::null();
242 }
243 
244 
245 template<class Type>
248 {
250 
251  return Field<Type>::null();
252 }
253 
254 
255 template<class Type>
257 {
259  this->writeEntry("value", os);
260 }
261 
262 
263 // ************************************************************************* //
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:339
virtual ~slicedFvPatchField()
Destructor.
Specialization of fvPatchField which creates the underlying fvPatchField as a slice of the given comp...
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
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 > > patchInternalField() const
Return internal field next to patch as patch field.
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
slicedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &)
Construct from patch, internal field and field to slice.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:54
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
runTime write()
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField of the values on the patch or on the.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:346