volSurfaceMapping.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "volSurfaceMapping.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
37 ) const
38 {
39  // Grab labels for all faces in faMesh
40  const labelList& faceLabels = mesh_.faceLabels();
41 
42  auto tresult = tmp<Field<Type>>::New(faceLabels.size(), Zero);
43  auto& result = tresult.ref();
44 
45  // Get reference to volume mesh
46  const polyMesh& pMesh = mesh_();
47  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
48 
49  label patchID, faceID;
50 
51  // Grab droplet cloud source by identifying patch and face
52  forAll(faceLabels, i)
53  {
54  // Escape if face is beyond active faces, eg belongs to a face zone
55  if (faceLabels[i] < pMesh.nFaces())
56  {
57  patchID = bm.whichPatch(faceLabels[i]);
58  faceID = bm[patchID].whichFace(faceLabels[i]);
59 
60  result[i] = df[patchID][faceID];
61  }
62  }
63 
64  return tresult;
65 }
66 
67 
68 template<class Type>
70 (
72 ) const
73 {
74  // Grab labels for all faces in faMesh
75  const labelList& faceLabels = mesh_.faceLabels();
76 
77  auto tresult = tmp<Field<Type>>::New(faceLabels.size(), Zero);
78  auto& result = tresult.ref();
79 
80  // Get reference to volume mesh
81  const polyMesh& pMesh = mesh_();
82  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
83 
84  label patchID, faceID;
85 
86  // Grab droplet cloud source by identifying patch and face
87  forAll(faceLabels, i)
88  {
89  // Escape if face is beyond active faces, eg belongs to a face zone
90  if (faceLabels[i] < pMesh.nFaces())
91  {
92  patchID = bm.whichPatch(faceLabels[i]);
93  faceID = bm[patchID].whichFace(faceLabels[i]);
94 
95  result[i] = df[patchID].patchInternalField()()[faceID];
96  }
97  }
98 
99  return tresult;
100 }
101 
102 
103 template<class Type>
105 (
108 ) const
109 {
110  // Grab labels for all faces in faMesh
111  const labelList& faceLabels = mesh_.faceLabels();
112 
113  // Get reference to volume mesh
114  const polyMesh& pMesh = mesh_();
115  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
116 
117  label patchID, faceID;
118 
119  const Field<Type>& afi = af.internalField();
120 
121  forAll(faceLabels, i)
122  {
123  // Escape if face is beyond active faces, eg belongs to a face zone
124  if (faceLabels[i] < pMesh.nFaces())
125  {
126  patchID = bm.whichPatch(faceLabels[i]);
127  faceID = bm[patchID].whichFace(faceLabels[i]);
128 
129  bf[patchID][faceID] = afi[i];
130  }
131  }
132 }
133 
134 
135 template<class Type>
137 (
140 ) const
141 {
142  mapToVolume(taf(), bf);
143 
144  taf.clear();
145 }
146 
147 
148 template<class Type>
150 (
152  Field<Type>& f
153 ) const
154 {
155  const labelList& faceLabels = mesh_.faceLabels();
156 
157  const polyMesh& pMesh = mesh_();
158  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
159  label patchID, faceID;
160 
161  const Field<Type>& afi = af.internalField();
162 
163  forAll(faceLabels, i)
164  {
165  if (faceLabels[i] < pMesh.nFaces())
166  {
167  patchID = bm.whichPatch(faceLabels[i]);
168  faceID = bm[patchID].whichFace(faceLabels[i]);
169  f[faceID] = afi[i];
170  }
171  }
172 }
173 
174 
175 // ************************************************************************* //
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::volSurfaceMapping::mapToSurface
tmp< Field< Type > > mapToSurface(const typename GeometricField< Type, fvPatchField, volMesh >::Boundary &df) const
Map droplet cloud sources to surface.
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::GeometricField::internalField
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Definition: GeometricFieldI.H:43
Foam::volSurfaceMapping::mapToField
void mapToField(const GeometricField< Type, faPatchField, areaMesh > &af, Field< Type > &f) const
Definition: volSurfaceMapping.C:150
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Field
Generic templated field type.
Definition: Field.H:63
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:810
Foam::volSurfaceMapping::mapToVolume
void mapToVolume(const GeometricField< Type, faPatchField, areaMesh > &af, typename GeometricField< Type, fvPatchField, volMesh >::Boundary &bf) const
Map surface field to volume boundary field.
Definition: volSurfaceMapping.C:105
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
volSurfaceMapping.H
f
labelList f(nPoints)
Foam::volSurfaceMapping::mapInternalToSurface
tmp< Field< Type > > mapInternalToSurface(const typename GeometricField< Type, fvPatchField, volMesh >::Boundary &df) const
Map patch internal field to surface.
mapToVolume
vsm mapToVolume(Cs, Cvf.boundaryFieldRef())
Foam::List< label >
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53