faceCorrectedSnGrad.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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
27 
28 #include "faceCorrectedSnGrad.H"
29 #include "volPointInterpolation.H"
30 #include "triangle.H"
31 
32 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 {}
37 
38 
39 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40 
41 template<class Type>
44 (
46 ) const
47 {
48  const fvMesh& mesh = this->mesh();
49 
51  (
53  );
54 
55  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
57  (
59  (
60  IOobject
61  (
62  "snGradCorr("+vf.name()+')',
63  vf.instance(),
64  mesh,
65  IOobject::NO_READ,
66  IOobject::NO_WRITE
67  ),
68  mesh,
69  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
70  )
71  );
72 
73  Field<Type>& sfCorr = tsfCorr.ref().primitiveFieldRef();
74 
75  const pointField& points = mesh.points();
76  const faceList& faces = mesh.faces();
77  const vectorField& Sf = mesh.Sf();
78  const vectorField& C = mesh.C();
79  const scalarField& magSf = mesh.magSf();
80  const labelList& owner = mesh.owner();
81  const labelList& neighbour = mesh.neighbour();
82 
83  forAll(sfCorr, facei)
84  {
86  (
88  );
89 
90  const face& fi = faces[facei];
91 
92  vector nf(Sf[facei]/magSf[facei]);
93 
94  for (label pi=0; pi<fi.size(); pi++)
95  {
96  // Next point index
97  label pj = (pi+1)%fi.size();
98 
99  // Edge normal in plane of face
100  vector edgen(nf^(points[fi[pj]] - points[fi[pi]]));
101 
102  // Edge centre field value
103  Type pvfe(0.5*(pvf[fi[pj]] + pvf[fi[pi]]));
104 
105  // Integrate face gradient
106  fgrad += edgen*pvfe;
107  }
108 
109  // Finalize face-gradient by dividing by face area
110  fgrad /= magSf[facei];
111 
112  // Calculate correction vector
113  vector dCorr(C[neighbour[facei]] - C[owner[facei]]);
114  dCorr /= (nf & dCorr);
115 
116  // if (mag(dCorr) > 2) dCorr *= 2/mag(dCorr);
117 
118  sfCorr[facei] = dCorr&fgrad;
119  }
120 
121  tsfCorr.ref().boundaryFieldRef() = Zero;
122 
123  return tsfCorr;
124 }
125 
126 
127 template<class Type>
130 (
132 ) const
133 {
134  const fvMesh& mesh = this->mesh();
135 
136  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
138  (
140  (
141  IOobject
142  (
143  "snGradCorr("+vf.name()+')',
144  vf.instance(),
145  mesh,
146  IOobject::NO_READ,
147  IOobject::NO_WRITE
148  ),
149  mesh,
150  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
151  )
152  );
154 
155  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
156  {
157  ssf.replace
158  (
159  cmpt,
161  .fullGradCorrection(vf.component(cmpt))
162  );
163  }
164 
165  return tssf;
166 }
167 
168 
169 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
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.
triangle.H
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::fv::faceCorrectedSnGrad::fullGradCorrection
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fullGradCorrection(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the faceCorrectedSnGrad.
Definition: faceCorrectedSnGrad.C:44
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:114
Foam::fv::faceCorrectedSnGrad::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the faceCorrectedSnGrad.
Definition: faceCorrectedSnGrad.C:130
Foam::fv::faceCorrectedSnGrad::~faceCorrectedSnGrad
virtual ~faceCorrectedSnGrad()
Destructor.
Definition: faceCorrectedSnGrad.C:35
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Field
Generic templated field type.
Definition: Field.H:63
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
volPointInterpolation.H
faceCorrectedSnGrad.H
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::List< face >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fv::faceCorrectedSnGrad
Simple central-difference snGrad scheme with non-orthogonal correction.
Definition: faceCorrectedSnGrad.H:59
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::outerProduct
Definition: products.H:106