fvPatchTemplates.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-2018 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 "fvPatch.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const UList<Type>& f
34 ) const
35 {
36  tmp<Field<Type>> tpif(new Field<Type>(size()));
37  Field<Type>& pif = tpif.ref();
38 
39  const labelUList& faceCells = this->faceCells();
40 
41  forAll(pif, facei)
42  {
43  pif[facei] = f[faceCells[facei]];
44  }
45 
46  return tpif;
47 }
48 
49 
50 template<class Type>
52 (
53  const UList<Type>& f,
54  Field<Type>& pif
55 ) const
56 {
57  pif.setSize(size());
58 
59  const labelUList& faceCells = this->faceCells();
60 
61  forAll(pif, facei)
62  {
63  pif[facei] = f[faceCells[facei]];
64  }
65 }
66 
67 
68 template<class GeometricField, class Type>
70 (
71  const GeometricField& gf
72 ) const
73 {
74  return gf.boundaryField()[index()];
75 }
76 
77 
78 template<class GeometricField, class Type>
80 (
81  GeometricField& gf
82 ) const
83 {
84  return gf.boundaryFieldRef()[index()];
85 }
86 
87 
88 // ************************************************************************* //
const GeometricField::Patch & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Generic GeometricField class.
Pre-declare SubField and related Field type.
Definition: Field.H:57
tmp< Field< Type > > patchInternalField(const UList< Type > &) const
Return given internal field next to patch as patch field.
void setSize(const label)
Reset size of List.
Definition: List.C:281
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
PatchField< Type > Patch
Type of the patch field of which the.
A class for managing temporary objects.
Definition: PtrList.H:53