All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
fvscStencil.H
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
25  grpfvscStencil
26 Class
27  Foam::fvsc::fvscStencil
28 Source file
29  fvscStencil.C
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvMesh.H"
33 #include "regIOobject.H"
34 #include "runTimeSelectionTables.H"
35 #include "surfaceFields.H"
36 
37 #ifndef fvscStencil_H
38 #define fvscStencil_H
39 
40 namespace Foam
41 {
42 
43 namespace fvsc
44 {
45 
46 class fvscStencil : public regIOobject, public refCount
47 {
48 
49 protected:
50 
51  //-
52  const fvMesh& mesh_;
53 
54  //-
55  const Time& runTime_;
56 
57  //-
58  surfaceVectorField nf_;
59 
60  //-
61  static PtrList<fvscStencil> stencils_;
62 
63 public:
64 
65  //-
66  TypeName("fvscStencil");
67 
68  //-
70  (
71  autoPtr,
73  components,
74  (
75  const IOobject& io
76  ),
77  (io)
78  );
79 
80  //- Construct from components
82  (
83  const IOobject& io
84  );
85 
86  //- Return a reference to the selected fvscStencil model
87  static autoPtr<fvscStencil> New
88  (
89  const word& name,
90  const fvMesh& mesh
91  );
92 
93  //-
94  //static tmp<fvscStencil> lookupOrNew
95  static fvscStencil& lookupOrNew
96  (
97  const word& nname,
98  const fvMesh& mesh
99  );
100 
101  //-
102  virtual ~fvscStencil();
103 
104  //-
105  virtual tmp<surfaceVectorField> Grad(const volScalarField& vF)
106  {
107  notImplemented("tmp<surfaceVectorField> Grad(const volScalarField& vF)");
108  return tmp<surfaceVectorField>(nullptr);
109  }
110 
111  //-
112  virtual tmp<surfaceTensorField> Grad(const volVectorField& iVF)
113  {
114  notImplemented("tmp<surfaceTensorField> Grad(const volVectorField& vF)");
115  return tmp<surfaceTensorField>(nullptr);
116  }
117 
118  //-
119  virtual tmp<surfaceScalarField> Div(const volVectorField& iVF)
120  {
121  notImplemented("tmp<surfaceScalarField> Grad(const volVectorField& vF)");
122  return tmp<surfaceScalarField>(nullptr);
123  }
124 
125  //-
126  virtual tmp<surfaceVectorField> Div(const volTensorField& iTF)
127  {
128  notImplemented("tmp<surfaceVectorField> Grad(const volTensorField& vF)");
129  return tmp<surfaceVectorField>(nullptr);
130  }
131 
132 
133  virtual bool writeData(Ostream&) const
134  {
135  return true;
136  }
137 };
139 }
140 
141 }
142 
143 #endif
144 
145 //END-OF-FILE
146 
This is a method for calculation the differential operators without tangential derivatives. They are further used in the calculation of the QGD terms.
Definition: fvscStencil.H:38
TypeName("fvscStencil")
virtual bool writeData(Ostream &) const
Definition: fvscStencil.H:154
virtual tmp< surfaceVectorField > Grad(const volScalarField &vF)
Definition: fvscStencil.H:120
const Time & runTime_
Definition: fvscStencil.H:51
virtual tmp< surfaceScalarField > Div(const volVectorField &iVF)
Definition: fvscStencil.H:138
const fvMesh & mesh_
Definition: fvscStencil.H:46
surfaceVectorField nf_
Definition: fvscStencil.H:56
Methods calculating of differential operators.
static fvscStencil & lookupOrNew(const word &nname, const fvMesh &mesh)
static tmp&lt;fvscStencil&gt; lookupOrNew
Definition: fvscStencil.C:92
fvscStencil(const IOobject &io)
Construct from components.
Definition: fvscStencil.C:113
static autoPtr< fvscStencil > New(const word &name, const fvMesh &mesh)
Return a reference to the selected fvscStencil model.
Definition: fvscStencil.C:53
static PtrList< fvscStencil > stencils_
Definition: fvscStencil.H:61
declareRunTimeSelectionTable(autoPtr, fvscStencil, components,(const IOobject &io),(io))