interpolationCellPatchConstrained.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-------------------------------------------------------------------------------
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
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
32
33template<class Type>
35(
37)
38:
39 interpolation<Type>(psi)
40{}
41
42
43// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44
45template<class Type>
47(
48 const vector& pt,
49 const label celli,
50 const label facei
51) const
52{
53 if (facei >= 0 && facei >= this->psi_.mesh().nInternalFaces())
54 {
55 // Use boundary value
56 const polyBoundaryMesh& pbm = this->psi_.mesh().boundaryMesh();
57 label patchi = pbm.patchID()[facei-this->psi_.mesh().nInternalFaces()];
58 label patchFacei = pbm[patchi].whichFace(facei);
59
60 return this->psi_.boundaryField()[patchi][patchFacei];
61 }
62 else
63 {
64 return this->psi_[celli];
65 }
66}
67
68
69// ************************************************************************* //
Generic GeometricField class.
Uses the cell value for any point in the cell apart from a boundary face where it uses the boundary v...
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const labelList & patchID() const
Per boundary face label the patch index.
const polyMesh & mesh() const noexcept
Return the mesh reference.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
bool interpolate() const noexcept
Same as isPointData()
const volScalarField & psi