All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
fvsc.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
26 Description
27  Methods calculating of differential operators
28 \*---------------------------------------------------------------------------*/
29 
30 
31 #include "fvsc.H"
32 #include "fvscStencil.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 
36 namespace Foam
37 {
38 namespace fvsc
39 {
40  word fvscOpName(const Foam::fvMesh& mesh, Foam::word termName);
41 }
42 }
43 
44 Foam::word
45 Foam::fvsc::fvscOpName(const Foam::fvMesh& mesh, Foam::word termName)
46 {
47  word opname = "none";
48  if (mesh.schemesDict().subDict("fvsc").found(termName))
49  {
50  mesh.schemesDict().subDict("fvsc").lookup(termName) >> opname;
51  }
52  else
53  {
54  mesh.schemesDict().subDict("fvsc").lookup("default") >> opname;
55  }
56 
57  return opname;
58 }
59 
60 Foam::tmp<Foam::surfaceVectorField>
61 Foam::fvsc::grad(const volScalarField& vf)
62 {
63  word tname = "grad(" + vf.name() + ")";
64 
66  (
67  fvscOpName(vf.mesh(),tname),
68  vf.mesh()
69  );
70 
71  tmp<surfaceVectorField> tGrad(Stencil.Grad(vf));
72 
73  return tGrad;
74 }
75 
76 Foam::tmp<Foam::surfaceVectorField>
77 Foam::fvsc::grad(const tmp<volScalarField>& tvf)
78 {
79  return Foam::fvsc::grad(tvf());
80 }
81 
82 Foam::tmp<Foam::surfaceTensorField>
83 Foam::fvsc::grad(const volVectorField& vf)
84 {
85  word tname = "grad(" + vf.name() + ")";
86 
88  (
89  fvscOpName(vf.mesh(),tname),
90  vf.mesh()
91  );
92 
93  return Stencil.Grad(vf);
94 }
95 
96 Foam::tmp<Foam::surfaceTensorField>
97 Foam::fvsc::grad(const tmp<volVectorField>& tvf)
98 {
99  return Foam::fvsc::grad(tvf());
100 }
101 
102 Foam::tmp<Foam::surfaceScalarField>
103 Foam::fvsc::div(const volVectorField& vf)
104 {
105  word tname = "div(" + vf.name() + ")";
106 
108  (
109  fvscOpName(vf.mesh(),tname),
110  vf.mesh()
111  );
112 
113  return Stencil.Div(vf);
114 }
115 
116 Foam::tmp<Foam::surfaceScalarField>
117 Foam::fvsc::div(const tmp<volVectorField>& tvf)
118 {
119  return Foam::fvsc::div(tvf());
120 }
121 
122 Foam::tmp<Foam::surfaceVectorField>
123 Foam::fvsc::div(const volTensorField& vf)
124 {
125  word tname = "div(" + vf.name() + ")";
126 
128  (
129  fvscOpName(vf.mesh(),tname),
130  vf.mesh()
131  );
132 
133  return Stencil.Div(vf);
134 }
135 
136 Foam::tmp<Foam::surfaceVectorField>
137 Foam::fvsc::div(const tmp<volTensorField>& tvf)
138 {
139  return Foam::fvsc::div(tvf());
140 }
141 
142 
143 //END-OF-FILE
144 
145 
word fvscOpName(const Foam::fvMesh &mesh, Foam::word termName)
static fvscStencil & lookupOrNew(const word &nname, const fvMesh &mesh)
static tmp&lt;fvscStencil&gt; lookupOrNew
fvscStencil(const IOobject &io)
Construct from components.
Methods calculating of differential operators.
tmp< surfaceScalarField > div(const volVectorField &vF)
tmp< surfaceVectorField > grad(const volScalarField &vF)