interpolationCellPatchConstrained.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const VolField<Type>& psi
35 )
36 :
38 {}
39 
40 
41 template<class Type>
43 (
45 )
46 :
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
53 template<class Type>
55 (
56  const vector& pt,
57  const label celli,
58  const label facei
59 ) const
60 {
61  if (facei >= 0 && facei >= this->psi_.mesh().nInternalFaces())
62  {
63  // Use boundary value
64  const polyBoundaryMesh& pbm = this->psi_.mesh().boundaryMesh();
65  label patchi = pbm.patchID()[facei-this->psi_.mesh().nInternalFaces()];
66  label patchFacei = pbm[patchi].whichFace(facei);
67 
68  return this->psi_.boundaryField()[patchi][patchFacei];
69  }
70  else
71  {
72  return this->psi_[celli];
73  }
74 }
75 
76 
77 // ************************************************************************* //
Generic GeometricField class.
Uses the cell value for any point in the cell apart from a boundary face where it uses the boundary v...
interpolationCellPatchConstrained(const VolField< Type > &psi)
Construct from components.
Type interpolate(const vector &position, const label celli, const label facei=-1) const
Interpolate field to the given point in the given cell.
Foam::polyBoundaryMesh.
const labelList & patchID() const
Per boundary face label the patch index.
const polyMesh & mesh() const
Return the mesh reference.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
label patchi
const volScalarField & psi
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