All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
fvscStencil.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
25  grpfvscStencil
26 Class
27  Foam::fvsc::fvscStencil
28 Description
29  This is a method for calculation the differential operators without
30  tangential derivatives. They are further used in the calculation of
31  the QGD terms.
32 \*---------------------------------------------------------------------------*/
33 
34 #include "fvscStencil.H"
35 #include "volFields.H"
36 #include "fvMesh.H"
37 #include "Time.H"
38 #include "addToRunTimeSelectionTable.H"
39 #include "coupledFvsPatchFields.H"
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace fvsc
46 {
47  defineTypeNameAndDebug(fvscStencil, 0);
48  defineRunTimeSelectionTable(fvscStencil, components);
49 }
50 }
51 
52 namespace Foam
53 {
54 namespace fvsc
55 {
56 
57 PtrList<fvscStencil> fvscStencil::stencils_(0);
58 
59 autoPtr<fvscStencil> fvscStencil::New
60 (
61  const word& fvscType,
62  const fvMesh& mesh
63 )
64 {
65  Info<< "Selecting finite volume surface calculus stencil type " << fvscType << endl;
66 
67  componentsConstructorTable::iterator cstrIter =
68  componentsConstructorTablePtr_->find(fvscType);
69 
70  if (cstrIter == componentsConstructorTablePtr_->end())
71  {
72  FatalErrorIn
73  (
74  "fvscStencil::New(const word&, const fvMesh&)"
75  ) << "Unknown Model type " << fvscType << nl << nl
76  << "Valid model types are:" << nl
77  << componentsConstructorTablePtr_->sortedToc()
78  << exit(FatalError);
79  }
80 
81  return autoPtr<fvscStencil>
82  (
83  cstrIter()
84  (
85  IOobject
86  (
87  fvscType,
88  mesh.time().timeName(),
89  mesh,
90  IOobject::NO_READ,
91  IOobject::NO_WRITE
92  )
93  )
94  );
95 }
96 
97 //tmp<fvscStencil> fvscStencil::lookupOrNew
99 (
100  const word& name,
101  const fvMesh& mesh
102 )
103 {
104  if (!mesh.thisDb().foundObject<fvscStencil>(name))
105  {
106  stencils_.append
107  (
108  fvscStencil::New(name,mesh)
109  );
110  stencils_.last().checkIn();
111  }
112 
113  return
114  const_cast<fvscStencil&>
115  (
116  mesh.thisDb().lookupObject<fvscStencil>(name)
117  );
118 }
119 //
120 fvscStencil::fvscStencil(const IOobject& io)
121 :
122  regIOobject(io, false),
123  refCount(),
124  mesh_(refCast<const fvMesh>(io.db())),
125  runTime_(mesh_.time()),
126  nf_
127  (
128  mesh_.Sf() / mesh_.magSf()
129  )
130 {
131  nf_.rename("nf");
132 }
133 
135 {
136 }
137 
138 }; //namespace fvsc
139 
140 }; //namespace Foam
141 
142 
143 //END-OF-FILE
144 
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
defineRunTimeSelectionTable(fvscStencil, components)
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
defineTypeNameAndDebug(fvscStencil, 0)