All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
GaussVolPointBase3D.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 Class
25  Foam::fvsc::GaussVolPoint::GaussVolPointBase3D
26 SourceFile
27  GaussVolPointBase3D.C
28 \*---------------------------------------------------------------------------*/
29 
30 #include "surfaceFieldsFwd.H"
31 #include "volFieldsFwd.H"
32 #include "pointFieldsFwd.H"
33 #include "tmp.H"
34 #include "List.H"
35 #include "DynamicList.H"
36 #include "label.H"
37 #include "scalar.H"
38 #include "autoPtr.H"
39 
40 #ifndef GaussVolPointBase3D_H
41 #define GaussVolPointBase3D_H
42 
43 namespace Foam
44 {
45 
46 class fvMesh;
47 class face;
48 class volPointInterpolation;
49 
50 typedef List<face> faceList;
51 
52 namespace fvsc
53 {
54 
56 {
57 
58  //-
59  const volPointInterpolation& volPoint_;
60 
61  //+
62  tmp<surfaceVectorField> nfRef_;
63 
64  //+
65  List<List<label> > bgfid_;
66 
67  //+
68  List<bool> processorPatch_;
69 
70  /* coefficients for quad faces */
71 
72  //+
73  List<label> qf_;
74 
75  //+
76  List<List<scalar> > aqx_;
77 
78  //+
79  List<List<scalar> > aqy_;
80 
81  //+
82  List<List<scalar> > aqz_;
83 
84  //+
85  List<scalar> vq_;
86 
87  //+
88  List<List<label> > bqf_;
89 
90  //+
91  List<List<List<scalar> > > baqx_;
92 
93  //+
94  List<List<List<scalar> > > baqy_;
95 
96  //+
97  List<List<List<scalar> > > baqz_;
98 
99  //+
100  List<List<scalar> > bvq_;
101 
102  /* coefficients for tri faces */
103  //
104  List<label> tf_;
105 
106  //
107  List<List<scalar> > atx_;
108 
109  //
110  List<List<scalar> > aty_;
111 
112  //
113  List<List<scalar> > atz_;
114 
115  //
116  List<scalar> vt_;
117 
118  //
119  List<List<label> > btf_;
120 
121  //
122  List<List<List<scalar> > > batx_;
123 
124  //
125  List<List<List<scalar> > > baty_;
126 
127  //
128  List<List<List<scalar> > > batz_;
129 
130  //-
131  List<List<scalar> > bvt_;
132 
133  /* coefficients for other faces */
134 
135  //-
136  List<List<scalar> > bmvON_;
137 
138  //+
139  DynamicList<label> of_;
140 
141  //+
142  List<List<label> > bof_;
143 
144 protected:
145 
146  //-
147  void calcGradfIF
148  (
149  const volScalarField& sf,
150  const pointScalarField& pf,
151  const faceList& faces,
152  surfaceVectorField& gradf,
153  const surfaceVectorField& dfdn
154  );
155 
156  //-
157  void calcGradfBF
158  (
159  const volScalarField& sf,
160  const pointScalarField& pf,
161  const faceList& faces,
162  surfaceVectorField& gradf,
163  const surfaceVectorField& dfdn
164  );
165 
166  //-
167  void calcGradfIF
168  (
169  const volVectorField& sf,
170  const pointVectorField& pf,
171  const faceList& faces,
172  surfaceTensorField& gradf,
173  const surfaceTensorField& dfdn
174  );
175 
176  //-
177  void calcGradfBF
178  (
179  const volVectorField& sf,
180  const pointVectorField& pf,
181  const faceList& faces,
182  surfaceTensorField& gradf,
183  const surfaceTensorField& dfdn
184  );
185 
186  //-
187  void calcDivfIF
188  (
189  const volVectorField& sf,
190  const pointVectorField& pf,
191  const faceList& faces,
192  surfaceScalarField& divf,
193  const surfaceScalarField& dfdn
194  );
195 
196  //-
197  void calcDivfBF
198  (
199  const volVectorField& sf,
200  const pointVectorField& pf,
201  const faceList& faces,
202  surfaceScalarField& divf,
203  const surfaceScalarField& dfdn
204  );
205 
206  //-
207  void calcDivfIF
208  (
209  const volTensorField& tf,
210  const pointTensorField& pf,
211  const faceList& faces,
212  surfaceVectorField& divf,
213  const surfaceVectorField& dfdn
214  );
215 
216  //-
217  void calcDivfBF
218  (
219  const volTensorField& sf,
220  const pointTensorField& pf,
221  const faceList& faces,
222  surfaceVectorField& divf,
223  const surfaceVectorField& dfdn
224  );
225 
226  //- Calculate weights for triangles
227  void triCalcWeights
228  (
229  const fvMesh& m
230  );
231 
232  //- Calcualte weights for quads
233  void quaCalcWeights
234  (
235  const fvMesh& m
236  );
237 
238 public:
239 
240  //-
241  GaussVolPointBase3D(const fvMesh& mesh);
242 
243  //-
245 
246  //-
247  void faceGrad(const volScalarField& f, surfaceVectorField& gradf);
248 
249  //-
250  void faceGrad(const volVectorField& f, surfaceTensorField& gradf);
251 
252  //-
253  void faceDiv(const volVectorField& f, surfaceScalarField& divf);
254 
255  //-
256  void faceDiv(const volTensorField& f, surfaceVectorField& divf);
257 
258 };
259 
260 }
261 
262 }
263 
264 #endif
265 
266 //
267 //END-OF-FILE
268 //
269 
void faceDiv(const volVectorField &f, surfaceScalarField &divf)
pf
Definition: updateFields.H:21
void faceGrad(const volScalarField &f, surfaceVectorField &gradf)
void triCalcWeights(const fvMesh &m)
Calculate weights for triangles.
List< face > faceList
void calcDivfBF(const volVectorField &sf, const pointVectorField &pf, const faceList &faces, surfaceScalarField &divf, const surfaceScalarField &dfdn)
void calcGradfBF(const volScalarField &sf, const pointScalarField &pf, const faceList &faces, surfaceVectorField &gradf, const surfaceVectorField &dfdn)
Methods calculating of differential operators.
GaussVolPointBase3D(const fvMesh &mesh)
void calcGradfIF(const volScalarField &sf, const pointScalarField &pf, const faceList &faces, surfaceVectorField &gradf, const surfaceVectorField &dfdn)
void quaCalcWeights(const fvMesh &m)
Calcualte weights for quads.
void calcDivfIF(const volVectorField &sf, const pointVectorField &pf, const faceList &faces, surfaceScalarField &divf, const surfaceScalarField &dfdn)