pointLinear.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-2016 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 "pointLinear.H"
27 #include "fvMesh.H"
28 #include "volPointInterpolation.H"
29 #include "triangle.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class Type>
37 (
39 ) const
40 {
41  const fvMesh& mesh = this->mesh();
42 
44  (
46  );
47 
50 
51  Field<Type>& sfCorr = tsfCorr.ref().primitiveFieldRef();
52 
53  const pointField& points = mesh.points();
54  const pointField& C = mesh.C();
55  const faceList& faces = mesh.faces();
56  const scalarField& w = mesh.weights();
57  const labelList& owner = mesh.owner();
58  const labelList& neighbour = mesh.neighbour();
59 
60  forAll(sfCorr, facei)
61  {
62  point pi =
63  w[owner[facei]]*C[owner[facei]]
64  + (1.0 - w[owner[facei]])*C[neighbour[facei]];
65 
66  const face& f = faces[facei];
67 
69  (
70  pi,
71  points[f[0]],
72  points[f[f.size()-1]]
73  ).mag();
74 
75  scalar sumAt = at;
76  Type sumPsip = at*(1.0/3.0)*
77  (
78  sfCorr[facei]
79  + pvf[f[0]]
80  + pvf[f[f.size()-1]]
81  );
82 
83  for (label pointi=1; pointi<f.size(); pointi++)
84  {
86  (
87  pi,
88  points[f[pointi]],
89  points[f[pointi-1]]
90  ).mag();
91 
92  sumAt += at;
93  sumPsip += at*(1.0/3.0)*
94  (
95  sfCorr[facei]
96  + pvf[f[pointi]]
97  + pvf[f[pointi-1]]
98  );
99 
100  }
101 
102  sfCorr[facei] = sumPsip/sumAt - sfCorr[facei];
103  }
104 
105 
107  Boundary& bSfCorr = tsfCorr.ref().boundaryFieldRef();
108 
109  forAll(bSfCorr, patchi)
110  {
111  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
112 
113  if (pSfCorr.coupled())
114  {
115  const fvPatch& fvp = mesh.boundary()[patchi];
116  const scalarField& pWghts = mesh.weights().boundaryField()[patchi];
117  const labelUList& pOwner = fvp.faceCells();
118  const vectorField& pNbrC = mesh.C().boundaryField()[patchi];
119 
120  forAll(pOwner, facei)
121  {
122  label own = pOwner[facei];
123 
124  point pi =
125  pWghts[facei]*C[own]
126  + (1.0 - pWghts[facei])*pNbrC[facei];
127 
128  const face& f = faces[facei+fvp.start()];
129 
131  (
132  pi,
133  points[f[0]],
134  points[f[f.size()-1]]
135  ).mag();
136 
137  scalar sumAt = at;
138  Type sumPsip = at*(1.0/3.0)*
139  (
140  pSfCorr[facei]
141  + pvf[f[0]]
142  + pvf[f[f.size()-1]]
143  );
144 
145  for (label pointi=1; pointi<f.size(); pointi++)
146  {
148  (
149  pi,
150  points[f[pointi]],
151  points[f[pointi-1]]
152  ).mag();
153 
154  sumAt += at;
155  sumPsip += at*(1.0/3.0)*
156  (
157  pSfCorr[facei]
158  + pvf[f[pointi]]
159  + pvf[f[pointi-1]]
160  );
161 
162  }
163 
164  pSfCorr[facei] = sumPsip/sumAt - pSfCorr[facei];
165  }
166  }
167  else
168  {
169  pSfCorr = Zero;
170  }
171  }
172 
173  return tsfCorr;
174 }
175 
176 
177 namespace Foam
178 {
180 }
181 
182 // ************************************************************************* //
Graphite solid properties.
Definition: C.H:57
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: pointLinear.C:37
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:57
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Generic GeometricField class.
const volVectorField & C() const
Return cell centres as volVectorField.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
virtual bool coupled() const
Return true if this patch field is coupled.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dynamicFvMesh & mesh
const pointField & points
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
Pre-declare SubField and related Field type.
Definition: Field.H:57
makeSurfaceInterpolationScheme(cellCoBlended)
static const zero Zero
Definition: zero.H:91
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
labelList f(nPoints)
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensioned< scalar > mag(const dimensioned< Type > &)
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Face-point interpolation scheme class derived from linear and returns linear weighting factors but al...
Definition: pointLinear.H:52
Namespace for OpenFOAM.