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 Group grpfvsc
25  This group contains common part of QGD solvers.
26 Class
27  Foam::fvsc
28 Description
29  Methods calculating of differential operators
30 \*---------------------------------------------------------------------------*/
31 
32 
33 #include "fvsc.H"
34 #include "fvscStencil.H"
35 #include "volFields.H"
36 #include "surfaceFields.H"
37 #include "prismMatcher.H"
38 
39 namespace Foam
40 {
41 namespace fvsc
42 {
43  word fvscOpName(const Foam::fvMesh& mesh, Foam::word termName);
44 }
45 }
46 
47 Foam::word
48 Foam::fvsc::fvscOpName(const Foam::fvMesh& mesh, Foam::word termName)
49 {
50  word opname = "none";
51  if (mesh.schemesDict().subDict("fvsc").found(termName))
52  {
53  mesh.schemesDict().subDict("fvsc").lookup(termName) >> opname;
54  }
55  else
56  {
57  mesh.schemesDict().subDict("fvsc").lookup("default") >> opname;
58  }
59 
60  if (((opname == "leastSquares") or (opname == "leastSquaresOpt")) and (mesh.nGeometricD() == 3))
61  {
62  FatalErrorIn("Foam::fvsc::fvscOpName") << "Can't use leastSquares or leastSquaresOpt in 3D case." << nl << exit(FatalError);
63  }
64 
65  if (opname == "GaussVolPoint")
66  {
67  const labelList wedgeBC = mesh.boundaryMesh().findIndices("wedge", true);
68  if (wedgeBC.size() > 0)
69  {
70  prismMatcher prism;
71  forAll (mesh.cells(), cellI)
72  {
73  if (prism.isA(mesh, cellI))
74  {
75  FatalErrorInFunction
76  << "GaussVolPoint scheme does not support solving axisymmetric cases with wedge BC and prism cells." << endl
77  << "Try to set leastSquares scheme."
78  << exit(FatalError);
79  }
80  }
81  }
82  }
83 
84  return opname;
85 }
86 
87 Foam::tmp<Foam::surfaceVectorField>
88 Foam::fvsc::grad(const volScalarField& vf)
89 {
90  word tname = "grad(" + vf.name() + ")";
91 
92  fvscStencil& Stencil = fvscStencil::lookupOrNew
93  (
94  fvscOpName(vf.mesh(),tname),
95  vf.mesh()
96  );
97 
98  tmp<surfaceVectorField> tGrad(Stencil.Grad(vf));
99 
100  return tGrad;
101 }
102 
103 Foam::tmp<Foam::surfaceVectorField>
104 Foam::fvsc::grad(const tmp<volScalarField>& tvf)
105 {
106  return Foam::fvsc::grad(tvf());
107 }
108 
109 Foam::tmp<Foam::surfaceTensorField>
110 Foam::fvsc::grad(const volVectorField& vf)
111 {
112  word tname = "grad(" + vf.name() + ")";
113 
114  fvscStencil& Stencil = fvscStencil::lookupOrNew
115  (
116  fvscOpName(vf.mesh(),tname),
117  vf.mesh()
118  );
119 
120  return Stencil.Grad(vf);
121 }
122 
123 Foam::tmp<Foam::surfaceTensorField>
124 Foam::fvsc::grad(const tmp<volVectorField>& tvf)
125 {
126  return Foam::fvsc::grad(tvf());
127 }
128 
129 Foam::tmp<Foam::surfaceScalarField>
130 Foam::fvsc::div(const volVectorField& vf)
131 {
132  word tname = "div(" + vf.name() + ")";
133 
134  fvscStencil& Stencil = fvscStencil::lookupOrNew
135  (
136  fvscOpName(vf.mesh(),tname),
137  vf.mesh()
138  );
139 
140  return Stencil.Div(vf);
141 }
142 
143 Foam::tmp<Foam::surfaceScalarField>
144 Foam::fvsc::div(const tmp<volVectorField>& tvf)
145 {
146  return Foam::fvsc::div(tvf());
147 }
148 
149 Foam::tmp<Foam::surfaceVectorField>
150 Foam::fvsc::div(const volTensorField& vf)
151 {
152  word tname = "div(" + vf.name() + ")";
153 
154  fvscStencil& Stencil = fvscStencil::lookupOrNew
155  (
156  fvscOpName(vf.mesh(),tname),
157  vf.mesh()
158  );
159 
160  return Stencil.Div(vf);
161 }
162 
163 Foam::tmp<Foam::surfaceVectorField>
164 Foam::fvsc::div(const tmp<volTensorField>& tvf)
165 {
166  return Foam::fvsc::div(tvf());
167 }
168 
169 
170 //END-OF-FILE
171 
172 
word fvscOpName(const Foam::fvMesh &mesh, Foam::word termName)
forAll(Y, i)
Definition: QGDYEqn.H:36
Methods calculating of differential operators.
tmp< surfaceScalarField > div(const volVectorField &vF)
tmp< surfaceVectorField > grad(const volScalarField &vF)