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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "slicedFvPatchField.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38  const Field<Type>& completeField
39 )
40 :
42 {
43  // Set the fvPatchField to a slice of the given complete field
44  UList<Type>::shallowCopy(p.patchSlice(completeField));
45 }
46 
47 
48 template<class Type>
50 (
51  const fvPatch& p,
53 )
54 :
56 {}
57 
58 
59 template<class Type>
61 (
62  const fvPatch& p,
64  const dictionary& dict
65 )
66 :
67  fvPatchField<Type>(p, iF, dict, false)
68 {
70 }
71 
72 
73 template<class Type>
75 (
76  const slicedFvPatchField<Type>& ptf,
77  const fvPatch& p,
79  const fvPatchFieldMapper& mapper
80 )
81 :
82  fvPatchField<Type>(ptf, p, iF, mapper)
83 {
85 }
86 
87 
88 template<class Type>
90 (
91  const slicedFvPatchField<Type>& ptf,
93 )
94 :
96 {
97  // Transfer the slice from the argument
99 }
100 
101 
102 template<class Type>
105 {
106  return tmp<fvPatchField<Type>>
107  (
108  new slicedFvPatchField<Type>(*this)
109  );
110 }
111 
112 
113 template<class Type>
115 (
116  const slicedFvPatchField<Type>& ptf
117 )
118 :
120  (
121  ptf.patch(),
122  ptf.internalField(),
123  Field<Type>()
124  )
125 {
126  // Transfer the slice from the argument
128 }
129 
130 
131 template<class Type>
134 (
136 ) const
137 {
138  return tmp<fvPatchField<Type>>
139  (
140  new slicedFvPatchField<Type>(*this, iF)
141  );
142 }
143 
144 
145 template<class Type>
147 {
148  // Set the fvPatchField storage pointer to nullptr before its destruction
149  // to protect the field it a slice of.
151 }
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
156 template<class Type>
158 {
160 
161  return Field<Type>::null();
162 }
163 
164 
165 template<class Type>
168 {
170 
171  return Field<Type>::null();
172 }
173 
174 
175 template<class Type>
177 {
179 }
180 
181 
182 template<class Type>
185 (
186  const Field<Type>& iField
187 ) const
188 {
190 
191  return Field<Type>::null();
192 }
193 
194 
195 template<class Type>
198 {
200 
201  return Field<Type>::null();
202 }
203 
204 
205 template<class Type>
208 (
209  const tmp<scalarField>&
210 ) const
211 {
213 
214  return Field<Type>::null();
215 }
216 
217 
218 template<class Type>
221 (
222  const tmp<scalarField>&
223 ) const
224 {
226 
227  return Field<Type>::null();
228 }
229 
230 
231 template<class Type>
234 {
236 
237  return Field<Type>::null();
238 }
239 
240 
241 template<class Type>
244 {
246 
247  return Field<Type>::null();
248 }
249 
250 
251 template<class Type>
253 {
255  this->writeEntry("value", os);
256 }
257 
258 
259 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::fvPatchField::internalField
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:363
Foam::slicedFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: slicedFvPatchField.C:157
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::slicedFvPatchField::gradientInternalCoeffs
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: slicedFvPatchField.C:233
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::slicedFvPatchField::valueBoundaryCoeffs
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: slicedFvPatchField.C:221
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Foam::Field
Generic templated field type.
Definition: Field.H:63
slicedFvPatchField.H
Foam::slicedFvPatchField::slicedFvPatchField
slicedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &)
Construct from patch, internal field and field to slice.
Definition: slicedFvPatchField.C:35
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::slicedFvPatchField::clone
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
Definition: slicedFvPatchField.C:104
Foam::slicedFvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField of the values on the patch or on the.
Definition: slicedFvPatchField.C:197
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::slicedFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: slicedFvPatchField.C:208
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::slicedFvPatchField
Specialization of fvPatchField which creates the underlying fvPatchField as a slice of the given comp...
Definition: slicedFvPatchField.H:64
Foam::slicedFvPatchField::~slicedFvPatchField
virtual ~slicedFvPatchField()
Destructor.
Definition: slicedFvPatchField.C:146
Foam::slicedFvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: slicedFvPatchField.C:243
Foam::UList
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:103
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::slicedFvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: slicedFvPatchField.C:167
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::slicedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: slicedFvPatchField.C:252
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54