fourthGrad.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-2019 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 "fourthGrad.H"
27 #include "leastSquaresGrad.H"
28 #include "gaussGrad.H"
29 #include "fvMesh.H"
30 #include "volMesh.H"
31 #include "surfaceMesh.H"
32 #include "GeometricField.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 <
41  <
45  >
46 >
48 (
49  const GeometricField<Type, fvPatchField, volMesh>& vsf,
50  const word& name
51 ) const
52 {
53  // The fourth-order gradient is calculated in two passes. First,
54  // the standard least-square gradient is assembled. Then, the
55  // fourth-order correction is added to the second-order accurate
56  // gradient to complete the accuracy.
57 
58  typedef typename outerProduct<vector, Type>::type GradType;
59 
60  const fvMesh& mesh = vsf.mesh();
61 
62  // Assemble the second-order least-square gradient
63  // Calculate the second-order least-square gradient
64  tmp<GeometricField<GradType, fvPatchField, volMesh>> tsecondfGrad
65  = leastSquaresGrad<Type>(mesh).grad
66  (
67  vsf,
68  "leastSquaresGrad(" + vsf.name() + ")"
69  );
70  const GeometricField<GradType, fvPatchField, volMesh>& secondfGrad =
71  tsecondfGrad();
72 
73  tmp<GeometricField<GradType, fvPatchField, volMesh>> tfGrad
74  (
76  (
77  name,
78  secondfGrad
79  )
80  );
81  GeometricField<GradType, fvPatchField, volMesh>& fGrad = tfGrad.ref();
82 
83  const vectorField& C = mesh.C();
84 
85  const surfaceScalarField& lambda = mesh.weights();
86 
87  // Get reference to least square vectors
88  const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
89  const surfaceVectorField& ownLs = lsv.pVectors();
90  const surfaceVectorField& neiLs = lsv.nVectors();
91 
92  // owner/neighbour addressing
93  const labelUList& own = mesh.owner();
94  const labelUList& nei = mesh.neighbour();
95 
96  // Assemble the fourth-order gradient
97 
98  // Internal faces
99  forAll(own, facei)
100  {
101  Type dDotGradDelta = 0.5*
102  (
103  (C[nei[facei]] - C[own[facei]])
104  & (secondfGrad[nei[facei]] - secondfGrad[own[facei]])
105  );
106 
107  fGrad[own[facei]] -= lambda[facei]*ownLs[facei]*dDotGradDelta;
108  fGrad[nei[facei]] -= (1.0 - lambda[facei])*neiLs[facei]*dDotGradDelta;
109  }
110 
111  // Boundary faces
112  forAll(vsf.boundaryField(), patchi)
113  {
114  if (secondfGrad.boundaryField()[patchi].coupled())
115  {
116  const fvsPatchVectorField& patchOwnLs =
117  ownLs.boundaryField()[patchi];
118 
119  const scalarField& lambdap = lambda.boundaryField()[patchi];
120 
121  const fvPatch& p = vsf.boundaryField()[patchi].patch();
122 
123  const labelUList& faceCells = p.faceCells();
124 
125  // Build the d-vectors
126  vectorField pd(p.delta());
127 
128  const Field<GradType> neighbourSecondfGrad
129  (
130  secondfGrad.boundaryField()[patchi].patchNeighbourField()
131  );
132 
133  forAll(faceCells, patchFacei)
134  {
135  fGrad[faceCells[patchFacei]] -=
136  0.5*lambdap[patchFacei]*patchOwnLs[patchFacei]
137  *(
138  pd[patchFacei]
139  & (
140  neighbourSecondfGrad[patchFacei]
141  - secondfGrad[faceCells[patchFacei]]
142  )
143  );
144  }
145  }
146  }
147 
148  fGrad.correctBoundaryConditions();
150 
151  return tfGrad;
152 }
153 
154 
155 // ************************************************************************* //
fvsPatchField< vector > fvsPatchVectorField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
rDeltaT correctBoundaryConditions()
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
volVectorField vectorField(fieldObject, mesh)
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
Generic GeometricField class.
UList< label > labelUList
Definition: UList.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
volScalarField scalarField(fieldObject, mesh)
label patchi
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField