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 | Copyright (C) 2011-2013 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "leastSquaresVectors.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(leastSquaresVectors, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 :
41  MeshObject<fvMesh, Foam::MoveableMeshObject, leastSquaresVectors>(mesh),
42  pVectors_
43  (
44  IOobject
45  (
46  "LeastSquaresP",
47  mesh_.pointsInstance(),
48  mesh_,
49  IOobject::NO_READ,
50  IOobject::NO_WRITE,
51  false
52  ),
53  mesh_,
55  ),
56  nVectors_
57  (
58  IOobject
59  (
60  "LeastSquaresN",
61  mesh_.pointsInstance(),
62  mesh_,
63  IOobject::NO_READ,
64  IOobject::NO_WRITE,
65  false
66  ),
67  mesh_,
69  )
70 {
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
84 {
85  if (debug)
86  {
87  Info<< "leastSquaresVectors::calcLeastSquaresVectors() :"
88  << "Calculating least square gradient vectors"
89  << 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(), symmTensor::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::GeometricBoundaryField& blsP =
120  pVectors_.boundaryField();
121 
122  forAll(blsP, 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(blsP, patchi)
174  {
175  fvsPatchVectorField& patchLsP = blsP[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  {
212  Info<< "leastSquaresVectors::calcLeastSquaresVectors() :"
213  << "Finished calculating least square gradient vectors"
214  << endl;
215  }
216 }
217 
218 
220 {
221  calcLeastSquaresVectors();
222  return true;
223 }
224 
225 
226 // ************************************************************************* //
volFields.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
p
p
Definition: pEqn.H:62
Foam::leastSquaresVectors::~leastSquaresVectors
virtual ~leastSquaresVectors()
Destructor.
Definition: invDistLeastSquaresVectors.C:77
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::fvsPatchScalarField
fvsPatchField< scalar > fvsPatchScalarField
Definition: fvsPatchFieldsFwd.H:43
Foam::fv::LeastSquaresVectors::calcLeastSquaresVectors
void calcLeastSquaresVectors()
Calculate Least-squares gradient vectors.
Definition: LeastSquaresVectors.C:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
leastSquaresVectors.H
C
volScalarField & C
Definition: readThermalProperties.H:104
Foam::leastSquaresVectors::movePoints
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
Definition: invDistLeastSquaresVectors.C:180
Foam::leastSquaresVectors::calcLeastSquaresVectors
void calcLeastSquaresVectors()
Construct Least-squares gradient vectors.
Definition: invDistLeastSquaresVectors.C:83
Foam::symmTensorField
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
Definition: primitiveFieldsFwd.H:51
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:49
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
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:59
Foam::fvsPatchVectorField
fvsPatchField< vector > fvsPatchVectorField
Definition: fvsPatchFieldsFwd.H:46
Foam::leastSquaresVectors::leastSquaresVectors
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
Definition: invDistLeastSquaresVectors.C:39
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:71
Foam::Info
messageStream Info
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::SymmTensor< scalar >::zero
static const SymmTensor zero
Definition: SymmTensor.H:77
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::labelUList
UList< label > labelUList
Definition: UList.H:63
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)