cubic.H
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 Class
25  Foam::cubic
26 
27 Description
28  Cubic interpolation scheme class derived from linear and returns
29  linear weighting factors but also applies an explicit correction.
30 
31 SourceFiles
32  cubic.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef cubic_H
37 #define cubic_H
38 
39 #include "linear.H"
40 #include "gaussGrad.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class cubic Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 template<class Type>
52 class cubic
53 :
54  public linear<Type>
55 {
56  // Private Member Functions
57 
58  //- Disallow default bitwise copy construct
59  cubic(const cubic&);
60 
61  //- Disallow default bitwise assignment
62  void operator=(const cubic&);
63 
64 
65 public:
66 
67  //- Runtime type information
68  TypeName("cubic");
69 
70 
71  // Constructors
72 
73  //- Construct from mesh
74  cubic(const fvMesh& mesh)
75  :
76  linear<Type>(mesh)
77  {}
78 
79  //- Construct from mesh and Istream
81  (
82  const fvMesh& mesh,
83  Istream&
84  )
85  :
86  linear<Type>(mesh)
87  {}
88 
89  //- Construct from mesh, faceFlux and Istream
91  (
92  const fvMesh& mesh,
93  const surfaceScalarField&,
94  Istream&
95  )
96  :
97  linear<Type>(mesh)
98  {}
99 
100 
101  // Member Functions
102 
103  //- Return true if this scheme uses an explicit correction
104  virtual bool corrected() const
105  {
106  return true;
107  }
108 
109  //- Return the explicit correction to the face-interpolate
112  (
114  ) const
115  {
116  const fvMesh& mesh = this->mesh();
117 
118  // calculate the appropriate interpolation factors
119  const surfaceScalarField& lambda = mesh.weights();
120 
121  const surfaceScalarField kSc
122  (
123  lambda*(scalar(1) - lambda*(scalar(3) - scalar(2)*lambda))
124  );
125 
126  const surfaceScalarField kVecP(sqr(scalar(1) - lambda)*lambda);
127  const surfaceScalarField kVecN(sqr(lambda)*(lambda - scalar(1)));
128 
130  (
132  (
133  IOobject
134  (
135  "cubic::correction(" + vf.name() +')',
136  mesh.time().timeName(),
137  mesh,
140  false
141  ),
143  )
144  );
146  tsfCorr.ref();
147 
148  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
149  {
150  sfCorr.replace
151  (
152  cmpt,
153  sfCorr.component(cmpt)
154  + (
156  <
157  typename outerProduct
158  <
159  vector,
160  typename pTraits<Type>::cmptType
161  >::type
163  (
165  <typename pTraits<Type>::cmptType>(mesh)
166  .grad(vf.component(cmpt)),
167  kVecP,
168  kVecN
169  ) & mesh.Sf()
170  )/mesh.magSf()/mesh.surfaceInterpolation::deltaCoeffs()
171  );
172  }
173 
175  Boundary& sfCorrbf = sfCorr.boundaryFieldRef();
176 
177  forAll(sfCorrbf, pi)
178  {
179  if (!sfCorrbf[pi].coupled())
180  {
181  sfCorrbf[pi] = Zero;
182  }
183  }
184 
185  return tsfCorr;
186  }
187 };
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace Foam
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 #endif
197 
198 // ************************************************************************* //
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
uint8_t direction
Definition: direction.H:46
const surfaceVectorField & Sf() const
Return cell face area vectors.
const fvMesh & mesh() const
Return mesh reference.
Central-differencing interpolation scheme class.
Definition: linear.H:50
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Traits class for primitives.
Definition: pTraits.H:50
Basic second-order gradient scheme using face-interpolation and Gauss&#39; theorem.
Definition: gaussGrad.H:57
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Cubic interpolation scheme class derived from linear and returns linear weighting factors but also ap...
Definition: cubic.H:51
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: cubic.H:103
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
static const zero Zero
Definition: zero.H:91
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: cubic.H:111
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
TypeName("cubic")
Runtime type information.
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Abstract base class for surface interpolation schemes.
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243