All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
leastSquaresStencil.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  grpleastSquares
26 Class
27  Foam::fvsc::leastSquares::leastSquaresStencil
28 Description
29  Realisation least squares method for calculationg of differential operators
30 \*---------------------------------------------------------------------------*/
31 
32 #include "leastSquaresStencil.H"
33 #include "polyMesh.H"
34 #include "fvMesh.H"
35 #include "word.H"
36 #include "IOstream.H"
37 #include "Ostream.H"
38 #include <HashTable.H>
39 #include "addToRunTimeSelectionTable.H"
40 #include "faceSet.H"
41 
42 namespace Foam
43 {
44 namespace fvsc
45 {
46  defineTypeNameAndDebug(leastSquares,0);
48  (
49  fvscStencil,
50  leastSquares,
51  components
52  );
53 }
54 }
55 
56 
57 // constructors
59 :
60  fvscStencil(io),
61  leastSquaresBase(mesh_)
62 {
63  faceSet degenerateFacesSet
64  (
65  //faceSetHeader
66  mesh_,
67  "degenerateStencilFaces",
68  IOobject::READ_IF_PRESENT,
69  IOobject::NO_WRITE
70  );
71 
72  if (degenerateFacesSet.size() > 0) //.headerOk())
73  {
74  Info << "Found set with faces for reduced approximation QGD terms" << endl;
75 
76  //read list of degenerated faces
77 
78 
79  const labelList degenerateFaces = degenerateFacesSet.toc();
80 
81  labelHashSet intDegFaces;
82  List<labelHashSet> procDegFaces (procDegFaces_.size());
83 
84  forAll(degenerateFaces, iDegFace)
85  {
86  label faceId = degenerateFaces[iDegFace];
87 
88  if (mesh_.isInternalFace(faceId))
89  {
90  intDegFaces.insert(faceId);
91  }
92  else
93  {
94  label patchId = mesh_.boundaryMesh().whichPatch(faceId);
95  label iPatch = -1;
96  forAll(procPairs_, iPP)
97  {
98  if (procPairs_[iPP] == patchId)
99  {
100  iPatch = iPP;
101  break;
102  }
103  }
104 
105  if (iPatch > -1)
106  {
107  procDegFaces[iPatch].insert
108  (
109  faceId
110  -
111  mesh_.boundaryMesh()[patchId].start()
112  );
113  }
114  }
115  }
116 
117  internalDegFaces_.append(intDegFaces.toc());
118  forAll(procDegFaces, iProcPair)
119  {
120  if (procDegFaces[iProcPair].toc().size() > 0)
121  {
122  procDegFaces_[iProcPair].append
123  (
124  procDegFaces[iProcPair].toc()
125  );
126  }
127  }
128  }
129  else
130  {
131  Info << "Set \"degenerateStencilFaces\" with faces for reduced approximation QGD terms not found" << endl;
132  }
133 }
134 
136 {
137 }
138 
139 //- Calculate gradient of volume vector field on the faces.
140 //
141 // \param iVF Internal vector field.
142 // Allowable values: constant reference to the volVectorField.
143 //
144 // \return Gradient of iVF (tensor field) which was computed on the faces of mesh.
145 Foam::tmp<Foam::surfaceTensorField> Foam::fvsc::leastSquares::Grad(const volVectorField& iVF)
146 {
147  surfaceVectorField gradComp0col = Grad(iVF.component(0));
148  surfaceVectorField gradComp1col = Grad(iVF.component(1));
149  surfaceVectorField gradComp2col = Grad(iVF.component(2));
150 
151  tmp<surfaceTensorField> tgradIVF(0* nf_* fvc::snGrad(iVF));
152  surfaceTensorField& gradIVF = tgradIVF.ref();
153 
154  //set internal field
155  gradIVF.primitiveFieldRef().replace(0, gradComp0col.primitiveField().component(0));
156  gradIVF.primitiveFieldRef().replace(1, gradComp1col.primitiveField().component(0));
157  gradIVF.primitiveFieldRef().replace(2, gradComp2col.primitiveField().component(0));
158 
159  gradIVF.primitiveFieldRef().replace(3, gradComp0col.primitiveField().component(1));
160  gradIVF.primitiveFieldRef().replace(4, gradComp1col.primitiveField().component(1));
161  gradIVF.primitiveFieldRef().replace(5, gradComp2col.primitiveField().component(1));
162 
163  gradIVF.primitiveFieldRef().replace(6, gradComp0col.primitiveField().component(2));
164  gradIVF.primitiveFieldRef().replace(7, gradComp1col.primitiveField().component(2));
165  gradIVF.primitiveFieldRef().replace(8, gradComp2col.primitiveField().component(2));
166 
167  //set external fields
168  forAll(mesh_.boundaryMesh(), patchi)
169  {
170  forAll(mesh_.boundary()[patchi], facei)
171  {
172  gradIVF.boundaryFieldRef()[patchi][facei].component(0) =
173  gradComp0col.boundaryField()[patchi][facei].component(0);
174  gradIVF.boundaryFieldRef()[patchi][facei].component(1) =
175  gradComp1col.boundaryField()[patchi][facei].component(0);
176  gradIVF.boundaryFieldRef()[patchi][facei].component(2) =
177  gradComp2col.boundaryField()[patchi][facei].component(0);
178 
179  gradIVF.boundaryFieldRef()[patchi][facei].component(3) =
180  gradComp0col.boundaryField()[patchi][facei].component(1);
181  gradIVF.boundaryFieldRef()[patchi][facei].component(4) =
182  gradComp1col.boundaryField()[patchi][facei].component(1);
183  gradIVF.boundaryFieldRef()[patchi][facei].component(5) =
184  gradComp2col.boundaryField()[patchi][facei].component(1);
185 
186  gradIVF.boundaryFieldRef()[patchi][facei].component(6) =
187  gradComp0col.boundaryField()[patchi][facei].component(2);
188  gradIVF.boundaryFieldRef()[patchi][facei].component(7) =
189  gradComp1col.boundaryField()[patchi][facei].component(2);
190  gradIVF.boundaryFieldRef()[patchi][facei].component(8) =
191  gradComp2col.boundaryField()[patchi][facei].component(2);
192  }
193  }
194 
195  return tgradIVF;
196 };
197 
198 //- Calculate divergence of volume vector field on the faces.
199 //
200 // \param iVF Internal vector field.
201 // Allowable values: constant reference to the volVectorField.
202 //
203 // \return Divergence of iVF (scalar field) which was computed on the faces of mesh.
204 Foam::tmp<Foam::surfaceScalarField> Foam::fvsc::leastSquares::Div(const volVectorField& iVF)
205 {
206  surfaceVectorField gradComp0 = Grad(iVF.component(0));
207  surfaceVectorField gradComp1 = Grad(iVF.component(1));
208  surfaceVectorField gradComp2 = Grad(iVF.component(2));
209 
210  tmp<surfaceScalarField> tdivIVF(0 * (nf_ & fvc::snGrad(iVF)));
211  surfaceScalarField& divIVF = tdivIVF.ref();
212 
213  divIVF.primitiveFieldRef() = gradComp0.primitiveField().component(0)
214  + gradComp1.primitiveField().component(1)
215  + gradComp2.primitiveField().component(2);
216 
217  forAll(mesh_.boundary(), patchi)
218  {
219  divIVF.boundaryFieldRef()[patchi] =
220  gradComp0.boundaryField()[patchi].component(0)
221  +
222  gradComp1.boundaryField()[patchi].component(1)
223  +
224  gradComp2.boundaryField()[patchi].component(2);
225  }
226 
227  return tdivIVF;
228 };
229 
230 //- Calculate divergence of volume tensor field on the faces.
231 //
232 // \param iTF Internal tensor field.
233 // Allowable values: constant reference to the volTensorField.
234 //
235 // \return Divergence of iTF (vector field) which was computed on the faces of mesh.
236 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::leastSquares::Div(const volTensorField& iTF)
237 {
238  tmp<surfaceVectorField> gradComp0 (Grad(iTF.component(0)));
239  tmp<surfaceVectorField> gradComp1 (Grad(iTF.component(1)));
240  tmp<surfaceVectorField> gradComp2 (Grad(iTF.component(2)));
241 
242  tmp<surfaceVectorField> gradComp3 (Grad(iTF.component(3)));
243  tmp<surfaceVectorField> gradComp4 (Grad(iTF.component(4)));
244  tmp<surfaceVectorField> gradComp5 (Grad(iTF.component(5)));
245 
246  tmp<surfaceVectorField> gradComp6 (Grad(iTF.component(6)));
247  tmp<surfaceVectorField> gradComp7 (Grad(iTF.component(7)));
248  tmp<surfaceVectorField> gradComp8 (Grad(iTF.component(8)));
249 
250  tmp<surfaceScalarField> divComp0 (gradComp0().component(0) + gradComp3().component(1) + gradComp6().component(2));
251  tmp<surfaceScalarField> divComp1 (gradComp1().component(0) + gradComp4().component(1) + gradComp7().component(2));
252  tmp<surfaceScalarField> divComp2 (gradComp2().component(0) + gradComp5().component(1) + gradComp8().component(2));
253 
254  tmp<surfaceVectorField> tdivITF(0*nf_*fvc::snGrad(iTF.component(0)));
255  surfaceVectorField& divITF = tdivITF.ref();
256 
257  divITF.primitiveFieldRef().replace(0, divComp0().primitiveField());
258  divITF.primitiveFieldRef().replace(1, divComp1().primitiveField());
259  divITF.primitiveFieldRef().replace(2, divComp2().primitiveField());
260 
261  forAll(mesh_.boundary(), patchi)
262  {
263  forAll(mesh_.boundary()[patchi], facei)
264  {
265  divITF.boundaryFieldRef()[patchi][facei].component(0) =
266  divComp0().boundaryField()[patchi][facei];
267  divITF.boundaryFieldRef()[patchi][facei].component(1) =
268  divComp1().boundaryField()[patchi][facei];
269  divITF.boundaryFieldRef()[patchi][facei].component(2) =
270  divComp2().boundaryField()[patchi][facei];
271  }
272  }
273 
274  return tdivITF;
275 }
276 
277 //
278 //END-OF-FILE
279 //
280 
281 
List< DynamicList< label > > procDegFaces_
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.
tmp< surfaceScalarField > Div(const volVectorField &iVF)
Calculate divergence of volume vector field on the faces.
addToRunTimeSelectionTable(fvscStencil, GaussVolPoint, components)
const fvMesh & mesh_
Definition: fvscStencil.H:46
HashSet< label, Hash< label > > labelHashSet
leastSquares(const IOobject &)
Construct from IOobject. Optional flag for if IOobject is the.
forAll(Y, i)
Definition: QGDYEqn.H:36
defineTypeNameAndDebug(fvscStencil, 0)
DynamicList< label > internalDegFaces_