interpolationCellPointI.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-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29
30template<class Type>
32(
33 const cellPointWeight& cpw
34) const
35{
36 const barycentric& weights = cpw.weights();
37 const triFace& faceVertices = cpw.faceVertices();
38
39 Type t = this->psi_[cpw.cell()]*weights[0];
40 t += psip_[faceVertices[0]]*weights[1];
41 t += psip_[faceVertices[1]]*weights[2];
42 t += psip_[faceVertices[2]]*weights[3];
43
44 return t;
45}
46
47
48template<class Type>
50(
51 const vector& position,
52 const label celli,
53 const label facei
54) const
55{
56 return interpolate(cellPointWeight(this->pMesh_, position, celli, facei));
57}
58
59
60template<class Type>
62(
64 const tetIndices& tetIs,
65 const label facei
66) const
67{
68 // Assumes that the position is consistent with the supplied
69 // tetIndices. Does not pay attention to whether or not facei is
70 // supplied or not - the result will be essentially the same.
71 // Performs a consistency check, however.
72
73 if (facei >= 0)
74 {
75 if (facei != tetIs.face())
76 {
78 << "specified face " << facei << " inconsistent with the face "
79 << "stored by tetIndices: " << tetIs.face()
80 << exit(FatalError);
81 }
82 }
84 const triFace faceVertices(tetIs.faceTriIs(this->pMesh_));
85
86 return
87 (
88 this->psi_[tetIs.cell()]*coordinates[0]
89 + psip_[faceVertices[0]]*coordinates[1]
90 + psip_[faceVertices[1]]*coordinates[2]
91 + psip_[faceVertices[2]]*coordinates[3]
92 );
93}
94
95
96// ************************************************************************* //
Foam::cellPointWeight.
label cell() const noexcept
Cell index.
const barycentric & weights() const noexcept
Interpolation weights.
const triFace & faceVertices() const noexcept
Interpolation addressing for points on face.
bool interpolate() const noexcept
Same as isPointData()
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
label face() const noexcept
Return the face index.
Definition: tetIndices.H:139
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:55
label cell() const noexcept
Return the cell index.
Definition: tetIndices.H:133
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:72
mesh interpolate(rAU)
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130