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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::cubic
28 
29 Group
30  grpFvSurfaceInterpolationSchemes
31 
32 Description
33  Cubic interpolation scheme class derived from linear and returns
34  linear weighting factors but also applies an explicit correction.
35 
36 SourceFiles
37  cubic.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef cubic_H
42 #define cubic_H
43 
44 #include "linear.H"
45 #include "gaussGrad.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class cubic Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class Type>
57 class cubic
58 :
59  public linear<Type>
60 {
61  // Private Member Functions
62 
63  //- No copy construct
64  cubic(const cubic&) = delete;
65 
66  //- No copy assignment
67  void operator=(const cubic&) = delete;
68 
69 
70 public:
71 
72  //- Runtime type information
73  TypeName("cubic");
74 
75 
76  // Constructors
77 
78  //- Construct from mesh
79  cubic(const fvMesh& mesh)
80  :
81  linear<Type>(mesh)
82  {}
83 
84  //- Construct from mesh and Istream
86  (
87  const fvMesh& mesh,
88  Istream&
89  )
90  :
91  linear<Type>(mesh)
92  {}
93 
94  //- Construct from mesh, faceFlux and Istream
96  (
97  const fvMesh& mesh,
98  const surfaceScalarField&,
99  Istream&
100  )
101  :
102  linear<Type>(mesh)
103  {}
104 
105 
106  // Member Functions
107 
108  //- Return true if this scheme uses an explicit correction
109  virtual bool corrected() const
110  {
111  return true;
112  }
113 
114  //- Return the explicit correction to the face-interpolate
117  (
119  ) const
120  {
121  const fvMesh& mesh = this->mesh();
122 
123  // calculate the appropriate interpolation factors
125 
126  const surfaceScalarField kSc
127  (
128  lambda*(scalar(1) - lambda*(scalar(3) - scalar(2)*lambda))
129  );
130 
131  const surfaceScalarField kVecP(sqr(scalar(1) - lambda)*lambda);
132  const surfaceScalarField kVecN(sqr(lambda)*(lambda - scalar(1)));
133 
135  (
137  (
138  IOobject
139  (
140  "cubic::correction(" + vf.name() +')',
141  mesh.time().timeName(),
142  mesh,
145  false
146  ),
148  )
149  );
151  tsfCorr.ref();
152 
153  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
154  {
155  sfCorr.replace
156  (
157  cmpt,
158  sfCorr.component(cmpt)
159  + (
161  <
162  typename outerProduct
163  <
164  vector,
165  typename pTraits<Type>::cmptType
166  >::type
168  (
170  <typename pTraits<Type>::cmptType>(mesh)
171  .grad(vf.component(cmpt)),
172  kVecP,
173  kVecN
174  ) & mesh.Sf()
175  )/mesh.magSf()/mesh.surfaceInterpolation::deltaCoeffs()
176  );
177  }
178 
180  Boundary& sfCorrbf = sfCorr.boundaryFieldRef();
181 
182  forAll(sfCorrbf, pi)
183  {
184  if (!sfCorrbf[pi].coupled())
185  {
186  sfCorrbf[pi] = Zero;
187  }
188  }
189 
190  return tsfCorr;
191  }
192 };
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 } // End namespace Foam
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 #endif
202 
203 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::GeometricField::replace
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::surfaceInterpolationScheme::interpolate
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.
Definition: surfaceInterpolationScheme.C:126
Foam::cubic::TypeName
TypeName("cubic")
Runtime type information.
Foam::fv::gaussGrad
Basic second-order gradient scheme using face-interpolation and Gauss' theorem.
Definition: gaussGrad.H:63
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
gaussGrad.H
Foam::cubic::correction
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:116
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::cubic::corrected
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: cubic.H:108
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::surfaceInterpolation::weights
virtual const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
Definition: surfaceInterpolation.C:102
Foam::cubic::cubic
cubic(const fvMesh &mesh)
Construct from mesh.
Definition: cubic.H:78
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::surfaceInterpolationScheme
Abstract base class for surface interpolation schemes.
Definition: surfaceInterpolationScheme.H:57
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::surfaceInterpolationScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: surfaceInterpolationScheme.H:144
Foam::linear
Central-differencing interpolation scheme class.
Definition: linear.H:55
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::outerProduct
Definition: products.H:106
Foam::cubic
Cubic interpolation scheme class derived from linear and returns linear weighting factors but also ap...
Definition: cubic.H:56
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:319