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  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faceCorrectedSnGrad.H"
30 #include "volPointInterpolation.H"
31 #include "triangle.H"
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
38 (
40 ) const
41 {
42  const fvMesh& mesh = this->mesh();
43 
45  (
47  );
48 
49  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
51  (
53  (
54  IOobject
55  (
56  "snGradCorr("+vf.name()+')',
57  vf.instance(),
58  mesh,
59  IOobject::NO_READ,
60  IOobject::NO_WRITE
61  ),
62  mesh,
63  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
64  )
65  );
66 
67  Field<Type>& sfCorr = tsfCorr.ref().primitiveFieldRef();
68 
69  const pointField& points = mesh.points();
70  const faceList& faces = mesh.faces();
71  const vectorField& Sf = mesh.Sf();
72  const vectorField& C = mesh.C();
73  const scalarField& magSf = mesh.magSf();
74  const labelList& owner = mesh.owner();
75  const labelList& neighbour = mesh.neighbour();
76 
77  forAll(sfCorr, facei)
78  {
80  (
82  );
83 
84  const face& fi = faces[facei];
85 
86  const vector nf(Sf[facei]/magSf[facei]);
87 
88  for (label pi = 0; pi < fi.size(); ++pi)
89  {
90  // Next point index
91  const label pj = fi.fcIndex(pi);
92 
93  // Edge normal in plane of face
94  const vector edgen(nf^(points[fi[pj]] - points[fi[pi]]));
95 
96  // Edge centre field value
97  const Type pvfe(0.5*(pvf[fi[pj]] + pvf[fi[pi]]));
98 
99  // Integrate face gradient
100  fgrad += edgen*pvfe;
101  }
102 
103  // Finalize face-gradient by dividing by face area
104  fgrad /= magSf[facei];
105 
106  // Calculate correction vector
107  vector dCorr(C[neighbour[facei]] - C[owner[facei]]);
108  dCorr /= (nf & dCorr);
109 
110  sfCorr[facei] = dCorr&fgrad;
111  }
112 
113  tsfCorr.ref().boundaryFieldRef() = Zero;
114 
115  return tsfCorr;
116 }
117 
118 
119 template<class Type>
122 (
124 ) const
125 {
126  const fvMesh& mesh = this->mesh();
127 
128  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
130  (
132  (
133  IOobject
134  (
135  "snGradCorr("+vf.name()+')',
136  vf.instance(),
137  mesh,
138  IOobject::NO_READ,
139  IOobject::NO_WRITE
140  ),
141  mesh,
142  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
143  )
144  );
146 
147  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
148  {
149  ssf.replace
150  (
151  cmpt,
153  .fullGradCorrection(vf.component(cmpt))
154  );
155  }
156 
157  return tssf;
158 }
159 
160 
161 // ************************************************************************* //
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.
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
Definition: faceCorrectedSnGrad.C:38
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
Definition: faceCorrectedSnGrad.C:122
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
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:85
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
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
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
Surface gradient scheme with full explicit non-orthogonal correction.
Definition: faceCorrectedSnGrad.H:69
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