filmGaussGrad.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) 2023 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 "filmGaussGrad.H"
27 #include "filmFvPatch.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 <
35 >
37 (
38  const VolField<Type>& vsf,
39  const word& name
40 ) const
41 {
42  typedef typename outerProduct<vector, Type>::type GradType;
43 
44  const fvMesh& mesh = vsf.mesh();
45 
46  tmp<VolField<GradType>> tfGrad
47  (
49  (
50  name,
51  mesh,
52  dimensioned<GradType>
53  (
54  "zero",
55  vsf.dimensions()/dimLength,
56  Zero
57  ),
58  extrapolatedCalculatedFvPatchField<GradType>::typeName
59  )
60  );
61  VolField<GradType>& fGrad = tfGrad.ref();
62 
63  SurfaceField<Type> ssf(this->tinterpScheme_().interpolate(vsf));
64 
65  const labelUList& owner = mesh.owner();
66  const labelUList& neighbour = mesh.neighbour();
67  const vectorField& Sf = mesh.Sf();
68  const scalarField& V = mesh.V();
69 
70  Field<GradType>& ifGrad = fGrad;
71  const Field<Type>& issf = ssf;
72 
73  forAll(owner, facei)
74  {
75  GradType Sfssf = Sf[facei]*issf[facei];
76 
77  ifGrad[owner[facei]] += Sfssf;
78  ifGrad[neighbour[facei]] -= Sfssf;
79  }
80 
81  forAll(mesh.boundary(), patchi)
82  {
83  const fvPatch& p = mesh.boundary()[patchi];
84  const labelUList& pFaceCells = p.faceCells();
85  const vectorField& pSf = mesh.Sf().boundaryField()[patchi];
86  const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
87 
88  if (isA<filmFvPatch>(p))
89  {
90  // Scale the film patch contribution by V/delta
91 
92  const scalarField& deltaCoeffs = p.deltaCoeffs();
93 
94  forAll(p, facei)
95  {
96  ifGrad[pFaceCells[facei]] +=
97  0.5*V[pFaceCells[facei]]*deltaCoeffs[facei]
98  *(pSf[facei]/mag(pSf[facei]))*pssf[facei];
99  }
100  }
101  else
102  {
103  forAll(p, facei)
104  {
105  ifGrad[pFaceCells[facei]] += pSf[facei]*pssf[facei];
106  }
107  }
108  }
109 
110  ifGrad /= mesh.V();
111 
112  fGrad.correctBoundaryConditions();
114 
115  return tfGrad;
116 }
117 
118 
119 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
virtual tmp< VolField< typename outerProduct< vector, Type >::type > > calcGrad(const VolField< Type > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
A class for managing temporary objects.
Definition: tmp.H:55
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
label patchi
U correctBoundaryConditions()
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
dimensioned< scalar > mag(const dimensioned< Type > &)
UList< label > labelUList
Definition: UList.H:65
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField & p