faceCorrectedSnGrad.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-2022 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 "faceCorrectedSnGrad.H"
27 #include "volPointInterpolation.H"
28 #include "triangle.H"
29 
30 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
31 
32 template<class Type>
34 {}
35 
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
39 template<class Type>
42 (
43  const VolField<Type>& vf
44 ) const
45 {
46  const fvMesh& mesh = this->mesh();
47 
49  (
51  );
52 
53  // construct SurfaceField<Type>
54  tmp<SurfaceField<Type>> tsfCorr
55  (
57  (
58  "snGradCorr("+vf.name()+')',
59  mesh,
61  )
62  );
63 
64  Field<Type>& sfCorr = tsfCorr.ref().primitiveFieldRef();
65 
66  const pointField& points = mesh.points();
67  const faceList& faces = mesh.faces();
68  const vectorField& Sf = mesh.Sf();
69  const vectorField& C = mesh.C();
70  const scalarField& magSf = mesh.magSf();
71  const labelList& owner = mesh.owner();
72  const labelList& neighbour = mesh.neighbour();
73 
74  forAll(sfCorr, facei)
75  {
77  (
79  );
80 
81  const face& fi = faces[facei];
82 
83  vector nf(Sf[facei]/magSf[facei]);
84 
85  for (label pi=0; pi<fi.size(); pi++)
86  {
87  // Next point index
88  label pj = (pi+1)%fi.size();
89 
90  // Edge normal in plane of face
91  vector edgen(nf^(points[fi[pj]] - points[fi[pi]]));
92 
93  // Edge centre field value
94  Type pvfe(0.5*(pvf[fi[pj]] + pvf[fi[pi]]));
95 
96  // Integrate face gradient
97  fgrad += edgen*pvfe;
98  }
99 
100  // Finalise face-gradient by dividing by face area
101  fgrad /= magSf[facei];
102 
103  // Calculate correction vector
104  vector dCorr(C[neighbour[facei]] - C[owner[facei]]);
105  dCorr /= (nf & dCorr);
106 
107  // if (mag(dCorr) > 2) dCorr *= 2/mag(dCorr);
108 
109  sfCorr[facei] = dCorr&fgrad;
110  }
111 
112  tsfCorr.ref().boundaryFieldRef() = Zero;
113 
114  return tsfCorr;
115 }
116 
117 
118 template<class Type>
121 (
122  const VolField<Type>& vf
123 ) const
124 {
125  const fvMesh& mesh = this->mesh();
126 
127  // construct SurfaceField<Type>
129  (
131  (
132  "snGradCorr("+vf.name()+')',
133  mesh,
135  )
136  );
137  SurfaceField<Type>& ssf = tssf.ref();
138 
139  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
140  {
141  ssf.replace
142  (
143  cmpt,
145  .fullGradCorrection(vf.component(cmpt))
146  );
147  }
148 
149  return tssf;
150 }
151 
152 
153 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Graphite solid properties.
Definition: C.H:51
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricField class.
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
const word & name() const
Return name.
Definition: IOobject.H:310
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const volVectorField & C() const
Return cell centres.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:451
const surfaceVectorField & Sf() const
Return cell face area vectors.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:457
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Simple central-difference snGrad scheme with non-orthogonal correction.
tmp< SurfaceField< Type > > fullGradCorrection(const VolField< Type > &) const
Return the explicit correction to the faceCorrectedSnGrad.
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &) const
Return the explicit correction to the faceCorrectedSnGrad.
virtual ~faceCorrectedSnGrad()
Destructor.
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
Traits class for primitives.
Definition: pTraits.H:53
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1361
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const pointField & points
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
uint8_t direction
Definition: direction.H:45