leastSquaresVectors.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
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 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "leastSquaresVectors.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(leastSquaresVectors, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 :
43  MeshObject<fvMesh, Foam::MoveableMeshObject, leastSquaresVectors>(mesh),
44  pVectors_
45  (
46  IOobject
47  (
48  "LeastSquaresP",
49  mesh_.pointsInstance(),
50  mesh_,
51  IOobject::NO_READ,
52  IOobject::NO_WRITE,
53  false
54  ),
55  mesh_,
57  ),
58  nVectors_
59  (
60  IOobject
61  (
62  "LeastSquaresN",
63  mesh_.pointsInstance(),
64  mesh_,
65  IOobject::NO_READ,
66  IOobject::NO_WRITE,
67  false
68  ),
69  mesh_,
71  )
72 {
73  calcLeastSquaresVectors();
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
78 
80 {}
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 void Foam::leastSquaresVectors::calcLeastSquaresVectors()
86 {
87  if (debug)
88  {
89  InfoInFunction << "Calculating least square gradient vectors" << endl;
90  }
91 
92  const fvMesh& mesh = mesh_;
93 
94  // Set local references to mesh data
95  const labelUList& owner = mesh_.owner();
96  const labelUList& neighbour = mesh_.neighbour();
97 
98  const volVectorField& C = mesh.C();
99  const surfaceScalarField& w = mesh.weights();
100  const surfaceScalarField& magSf = mesh.magSf();
101 
102 
103  // Set up temporary storage for the dd tensor (before inversion)
104  symmTensorField dd(mesh_.nCells(), Zero);
105 
106  forAll(owner, facei)
107  {
108  label own = owner[facei];
109  label nei = neighbour[facei];
110 
111  vector d = C[nei] - C[own];
112  symmTensor wdd = (magSf[facei]/magSqr(d))*sqr(d);
113 
114  dd[own] += (1 - w[facei])*wdd;
115  dd[nei] += w[facei]*wdd;
116  }
117 
118 
119  surfaceVectorField::Boundary& pVectorsBf =
120  pVectors_.boundaryFieldRef();
121 
122  forAll(pVectorsBf, patchi)
123  {
124  const fvsPatchScalarField& pw = w.boundaryField()[patchi];
125  const fvsPatchScalarField& pMagSf = magSf.boundaryField()[patchi];
126 
127  const fvPatch& p = pw.patch();
128  const labelUList& faceCells = p.patch().faceCells();
129 
130  // Build the d-vectors
131  vectorField pd(p.delta());
132 
133  if (pw.coupled())
134  {
135  forAll(pd, patchFacei)
136  {
137  const vector& d = pd[patchFacei];
138 
139  dd[faceCells[patchFacei]] +=
140  ((1 - pw[patchFacei])*pMagSf[patchFacei]/magSqr(d))*sqr(d);
141  }
142  }
143  else
144  {
145  forAll(pd, patchFacei)
146  {
147  const vector& d = pd[patchFacei];
148 
149  dd[faceCells[patchFacei]] +=
150  (pMagSf[patchFacei]/magSqr(d))*sqr(d);
151  }
152  }
153  }
154 
155 
156  // Invert the dd tensor
157  const symmTensorField invDd(inv(dd));
158 
159 
160  // Revisit all faces and calculate the pVectors_ and nVectors_ vectors
161  forAll(owner, facei)
162  {
163  label own = owner[facei];
164  label nei = neighbour[facei];
165 
166  vector d = C[nei] - C[own];
167  scalar magSfByMagSqrd = magSf[facei]/magSqr(d);
168 
169  pVectors_[facei] = (1 - w[facei])*magSfByMagSqrd*(invDd[own] & d);
170  nVectors_[facei] = -w[facei]*magSfByMagSqrd*(invDd[nei] & d);
171  }
172 
173  forAll(pVectorsBf, patchi)
174  {
175  fvsPatchVectorField& patchLsP = pVectorsBf[patchi];
176 
177  const fvsPatchScalarField& pw = w.boundaryField()[patchi];
178  const fvsPatchScalarField& pMagSf = magSf.boundaryField()[patchi];
179 
180  const fvPatch& p = pw.patch();
181  const labelUList& faceCells = p.faceCells();
182 
183  // Build the d-vectors
184  vectorField pd(p.delta());
185 
186  if (pw.coupled())
187  {
188  forAll(pd, patchFacei)
189  {
190  const vector& d = pd[patchFacei];
191 
192  patchLsP[patchFacei] =
193  ((1 - pw[patchFacei])*pMagSf[patchFacei]/magSqr(d))
194  *(invDd[faceCells[patchFacei]] & d);
195  }
196  }
197  else
198  {
199  forAll(pd, patchFacei)
200  {
201  const vector& d = pd[patchFacei];
202 
203  patchLsP[patchFacei] =
204  pMagSf[patchFacei]*(1.0/magSqr(d))
205  *(invDd[faceCells[patchFacei]] & d);
206  }
207  }
208  }
209 
210  if (debug)
211  {
213  << "Finished calculating least square gradient vectors" << endl;
214  }
215 }
216 
217 
219 {
220  calcLeastSquaresVectors();
221  return true;
222 }
223 
224 
225 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
volFields.H
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::leastSquaresVectors::~leastSquaresVectors
virtual ~leastSquaresVectors()
Destructor.
Definition: invDistLeastSquaresVectors.C:79
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::fvsPatchScalarField
fvsPatchField< scalar > fvsPatchScalarField
Definition: fvsPatchFieldsFwd.H:45
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
leastSquaresVectors.H
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::leastSquaresVectors::movePoints
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
Definition: invDistLeastSquaresVectors.C:179
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::symmTensorField
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
Definition: primitiveFieldsFwd.H:56
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::fvsPatchVectorField
fvsPatchField< vector > fvsPatchVectorField
Definition: fvsPatchFieldsFwd.H:48
Foam::leastSquaresVectors::leastSquaresVectors
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
Definition: invDistLeastSquaresVectors.C:41
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:50
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:79
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)