fvPatchTemplates.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 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fvPatch.H"
30
31// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const UList<Type>& f
37) const
38{
39 return patchInternalField(f, this->faceCells());
40}
41
42
43template<class Type>
45(
46 const UList<Type>& f,
48) const
49{
50 auto tpif = tmp<Field<Type>>::New(size());
51 auto& pif = tpif.ref();
52
53 forAll(pif, facei)
54 {
55 pif[facei] = f[faceCells[facei]];
56 }
57
58 return tpif;
59}
60
61
62template<class Type>
64(
65 const UList<Type>& f,
66 Field<Type>& pif
67) const
68{
69 pif.resize(size());
70
71 const labelUList& faceCells = this->faceCells();
72
73 forAll(pif, facei)
74 {
75 pif[facei] = f[faceCells[facei]];
76 }
77}
78
79
80template<class GeometricField, class Type>
82(
83 const GeometricField& gf
84) const
85{
86 return gf.boundaryField()[index()];
87}
88
89
90// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
PatchField< Type > Patch
The patch field type for the GeometricBoundaryField.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: faPatchField.C:175
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
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:113
A class for managing temporary objects.
Definition: tmp.H:65
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333