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 | Copyright (C) 2011-2012 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 (
44 ) const
45 {
46  const fvMesh& mesh = this->mesh();
47 
49  (
51  );
52 
53  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
55  (
57  (
58  IOobject
59  (
60  "snGradCorr("+vf.name()+')',
61  vf.instance(),
62  mesh,
63  IOobject::NO_READ,
64  IOobject::NO_WRITE
65  ),
66  mesh,
67  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
68  )
69  );
70 
71  Field<Type>& sfCorr = tsfCorr().internalField();
72 
73  const pointField& points = mesh.points();
74  const faceList& faces = mesh.faces();
75  const vectorField& Sf = mesh.Sf().internalField();
76  const vectorField& C = mesh.C().internalField();
77  const scalarField& magSf = mesh.magSf().internalField();
78  const labelList& owner = mesh.owner();
79  const labelList& neighbour = mesh.neighbour();
80 
81  forAll(sfCorr, facei)
82  {
84  (
86  );
87 
88  const face& fi = faces[facei];
89 
90  vector nf(Sf[facei]/magSf[facei]);
91 
92  for (label pi=0; pi<fi.size(); pi++)
93  {
94  // Next point index
95  label pj = (pi+1)%fi.size();
96 
97  // Edge normal in plane of face
98  vector edgen(nf^(points[fi[pj]] - points[fi[pi]]));
99 
100  // Edge centre field value
101  Type pvfe(0.5*(pvf[fi[pj]] + pvf[fi[pi]]));
102 
103  // Integrate face gradient
104  fgrad += edgen*pvfe;
105  }
106 
107  // Finalize face-gradient by dividing by face area
108  fgrad /= magSf[facei];
109 
110  // Calculate correction vector
111  vector dCorr(C[neighbour[facei]] - C[owner[facei]]);
112  dCorr /= (nf & dCorr);
113 
114  // if (mag(dCorr) > 2) dCorr *= 2/mag(dCorr);
115 
116  sfCorr[facei] = dCorr&fgrad;
117  }
118 
119  tsfCorr().boundaryField() = pTraits<Type>::zero;
120 
121  return tsfCorr;
122 }
123 
124 
125 template<class Type>
128 (
130 ) const
131 {
132  const fvMesh& mesh = this->mesh();
133 
134  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
136  (
138  (
139  IOobject
140  (
141  "snGradCorr("+vf.name()+')',
142  vf.instance(),
143  mesh,
144  IOobject::NO_READ,
145  IOobject::NO_WRITE
146  ),
147  mesh,
148  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
149  )
150  );
152 
153  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
154  {
155  ssf.replace
156  (
157  cmpt,
159  .fullGradCorrection(vf.component(cmpt))
160  );
161  }
162 
163  return tssf;
164 }
165 
166 
167 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
triangle.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
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:42
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, int(pTraits< arg1 >::rank)+int(pTraits< arg2 >::rank) >::type type
Definition: products.H:72
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
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:128
Foam::fv::faceCorrectedSnGrad::~faceCorrectedSnGrad
virtual ~faceCorrectedSnGrad()
Destructor.
Definition: faceCorrectedSnGrad.C:33
Foam::label
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
Foam::Field< Type >
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volPointInterpolation.H
faceCorrectedSnGrad.H
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::GeometricField::replace
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fv::faceCorrectedSnGrad
Simple central-difference snGrad scheme with non-orthogonal correction.
Definition: faceCorrectedSnGrad.H:54
Foam::outerProduct
Definition: products.H:64