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-2020 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 #include "objectRegistry.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const UList<Type>& f
35 ) const
36 {
37  tmp<Field<Type>> tpif(new Field<Type>(size()));
38  Field<Type>& pif = tpif.ref();
39 
40  const labelUList& faceCells = this->faceCells();
41 
42  forAll(pif, facei)
43  {
44  pif[facei] = f[faceCells[facei]];
45  }
46 
47  return tpif;
48 }
49 
50 
51 template<class Type>
53 (
54  const UList<Type>& f,
55  Field<Type>& pif
56 ) const
57 {
58  pif.setSize(size());
59 
60  const labelUList& faceCells = this->faceCells();
61 
62  forAll(pif, facei)
63  {
64  pif[facei] = f[faceCells[facei]];
65  }
66 }
67 
68 
69 template<class GeometricField, class Type>
71 (
72  const GeometricField& gf
73 ) const
74 {
75  return gf.boundaryField()[index()];
76 }
77 
78 
79 template<class GeometricField, class Type>
81 (
82  GeometricField& gf
83 ) const
84 {
85  return gf.boundaryFieldRef()[index()];
86 }
87 
88 
89 template<class GeometricField, class Type>
91 (
92  const word& name
93 ) const
94 {
95  return patchField<GeometricField, Type>
96  (
97  db().template lookupObject<GeometricField>(name)
98  );
99 }
100 
101 
102 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic GeometricField class.
PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual label size() const
Return size.
Definition: fvPatch.H:138
const GeometricField::Patch & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
tmp< Field< Type > > patchInternalField(const UList< Type > &) const
Return given internal field next to patch as patch field.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
labelList f(nPoints)