All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
reducedFaceNormalStencil.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) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
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  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22  You should have received a copy of the GNU General Public License
23  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 Class
25  Foam::fvsc::reduced::reducedFaceNormalStencil
26 Description
27  Methods calculating of differential operators without tangential direvatives
28 \*---------------------------------------------------------------------------*/
29 
30 
31 
32 #include "fvc.H"
34 #include "addToRunTimeSelectionTable.H"
35 
36 namespace Foam
37 {
38 namespace fvsc
39 {
40  defineTypeNameAndDebug(reduced,0);
42  (
45  components
46  );
47 }
48 }
49 
50 // constructors
51 Foam::fvsc::reduced::reduced(const IOobject& io)
52 :
53  fvscStencil(io)
54 {
55 }
56 
58 {
59 }
60 
61 //- Calculate gradient of volume scalar function on the faces
62 //
63 // \param iF Internal scalar field.
64 // Allowable values: constant reference to the volScalarField.
65 //
66 // \return Gradient of iF (vector field) which was computed on the faces of mesh.
67 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::reduced::Grad(const volScalarField& vF)
68 {
69  tmp<surfaceVectorField> tgradIF(nf_ * fvc::snGrad(vF));
70 
71  return tgradIF;
72 };
73 
74 //- Calculate gradient of volume vector field on the faces.
75 //
76 // \param iVF Internal vector field.
77 // Allowable values: constant reference to the volVectorField.
78 //
79 // \return Gradient of iVF (tensor field) which was computed on the faces of mesh.
80 Foam::tmp<Foam::surfaceTensorField> Foam::fvsc::reduced::Grad(const volVectorField& iVF)
81 {
82 
83  tmp<surfaceTensorField> tgradIVF(nf_ * fvc::snGrad(iVF));
84 
85  return tgradIVF;
86 };
87 
88 Foam::tmp<Foam::surfaceScalarField> Foam::fvsc::reduced::Div(const volVectorField& iVF)
89 {
90  tmp<surfaceScalarField> tdivIVF(nf_ & fvc::snGrad(iVF));
91 
92  return tdivIVF;
93 };
94 
95 //- Calculate divergence of volume tensor field on the faces.
96 //
97 // \param iTF Internal tensor field.
98 // Allowable values: constant reference to the volTensorField.
99 //
100 // \return Divergence of iTF (vector field) which was computed on the faces of mesh.
101 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::reduced::Div(const volTensorField& iTF)
102 {
103  tmp<surfaceVectorField> tdivITF(nf_ & fvc::snGrad(iTF));
104 
105  return tdivITF;
106 }
107 
108 
109 //
110 //END-OF-FILE
111 //
112 
113 
tmp< surfaceScalarField > Div(const volVectorField &iVF)
addToRunTimeSelectionTable(fvscStencil, GaussVolPoint, components)
reduced(const IOobject &)
Construct from IOobject.
surfaceVectorField nf_
Definition: fvscStencil.H:54
fvscStencil(const IOobject &io)
Construct from components.
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.
defineTypeNameAndDebug(fvscStencil, 0)