slicedFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const fvPatch& p,
35  const Field<Type>& completeField
36 )
37 :
38  fvPatchField<Type>(p, iF, Field<Type>())
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  const fvPatchField<Type>& pf
51 )
52 :
53  fvPatchField<Type>(p, iF, Field<Type>())
54 {
55  // Set the fvPatchField values to the given fvPatchField
57 }
58 
59 
60 template<class Type>
62 (
63  const slicedFvPatchField<Type>& ptf,
65 )
66 :
67  fvPatchField<Type>(ptf.patch(), iF, Field<Type>())
68 {
69  // Transfer the slice from the argument
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
76 template<class Type>
78 {
79  // Set the fvPatchField storage pointer to nullptr before its destruction
80  // to protect the field it a slice of.
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class Type>
89 {
91 
92  return Field<Type>::null();
93 }
94 
95 
96 template<class Type>
98 {
100 }
101 
102 
103 template<class Type>
106 (
107  const tmp<scalarField>&
108 ) const
109 {
111 
112  return Field<Type>::null();
113 }
114 
115 
116 template<class Type>
119 (
120  const tmp<scalarField>&
121 ) const
122 {
124 
125  return Field<Type>::null();
126 }
127 
128 
129 template<class Type>
132 {
134 
135  return Field<Type>::null();
136 }
137 
138 
139 template<class Type>
142 {
144 
145  return Field<Type>::null();
146 }
147 
148 
149 template<class Type>
151 {
153  writeEntry(os, "value", *this);
154 }
155 
156 
157 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:82
static const Field< Type > & null()
Return a null field.
Definition: Field.H:111
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
void shallowCopy(const UList< T > &)
Copy the pointer held by the given UList.
Definition: UListI.H:156
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Specialisation of fvPatchField which creates the underlying fvPatchField as a slice of the given comp...
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
virtual ~slicedFvPatchField()
Destructor.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
slicedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &)
Construct from patch, internal field and field to slice.
A class for managing temporary objects.
Definition: tmp.H:55
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
volScalarField & p