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-2021 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 (
71  const Field<Type>& f
72 ) const
73 {
74  const labelList& faceLabels = mesh_.faceLabels();
75 
76  auto tresult = tmp<Field<Type>>::New(faceLabels.size(), Zero);
77  auto& result = tresult.ref();
78 
79  const polyMesh& pMesh = mesh_();
80  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
81  label patchID, faceID;
82 
83  forAll(faceLabels, i)
84  {
85  if (faceLabels[i] < pMesh.nFaces())
86  {
87  patchID = bm.whichPatch(faceLabels[i]);
88  faceID = bm[patchID].whichFace(faceLabels[i]);
89 
90  result[i] = f[faceID];
91  }
92  }
93 
94  return tresult;
95 }
96 
97 
98 template<class Type>
100 (
102 ) const
103 {
104  // Grab labels for all faces in faMesh
105  const labelList& faceLabels = mesh_.faceLabels();
106 
107  auto tresult = tmp<Field<Type>>::New(faceLabels.size(), Zero);
108  auto& result = tresult.ref();
109 
110  // Get reference to volume mesh
111  const polyMesh& pMesh = mesh_();
112  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
113 
114  label patchID, faceID;
115 
116  // Grab droplet cloud source by identifying patch and face
117  forAll(faceLabels, i)
118  {
119  // Escape if face is beyond active faces, eg belongs to a face zone
120  if (faceLabels[i] < pMesh.nFaces())
121  {
122  patchID = bm.whichPatch(faceLabels[i]);
123  faceID = bm[patchID].whichFace(faceLabels[i]);
124 
125  result[i] = df[patchID].patchInternalField()()[faceID];
126  }
127  }
128 
129  return tresult;
130 }
131 
132 
133 template<class Type>
135 (
138 ) const
139 {
140  // Grab labels for all faces in faMesh
141  const labelList& faceLabels = mesh_.faceLabels();
142 
143  // Get reference to volume mesh
144  const polyMesh& pMesh = mesh_();
145  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
146 
147  label patchID, faceID;
148 
149  const Field<Type>& afi = af.internalField();
150 
151  forAll(faceLabels, i)
152  {
153  // Escape if face is beyond active faces, eg belongs to a face zone
154  if (faceLabels[i] < pMesh.nFaces())
155  {
156  patchID = bm.whichPatch(faceLabels[i]);
157  faceID = bm[patchID].whichFace(faceLabels[i]);
158 
159  bf[patchID][faceID] = afi[i];
160  }
161  }
162 }
163 
164 
165 template<class Type>
167 (
170 ) const
171 {
172  mapToVolume(taf(), bf);
173 
174  taf.clear();
175 }
176 
177 
178 template<class Type>
180 (
182  Field<Type>& f
183 ) const
184 {
185  const Field<Type>& afi = af.internalField();
186 
187  mapToField(afi, f);
188 }
189 
190 
191 template<class Type>
193 (
194  const Field<Type>& af,
195  Field<Type>& f
196 ) const
197 {
198  const labelList& faceLabels = mesh_.faceLabels();
199 
200  const polyMesh& pMesh = mesh_();
201  const polyBoundaryMesh& bm = pMesh.boundaryMesh();
202  label patchID, faceID;
203 
204  forAll(faceLabels, i)
205  {
206  if (faceLabels[i] < pMesh.nFaces())
207  {
208  patchID = bm.whichPatch(faceLabels[i]);
209  faceID = bm[patchID].whichFace(faceLabels[i]);
210  f[faceID] = af[i];
211  }
212  }
213 }
214 
215 
216 // ************************************************************************* //
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
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::volSurfaceMapping::mapToSurface
tmp< Field< Type > > mapToSurface(const typename GeometricField< Type, fvPatchField, volMesh >::Boundary &df) const
Map Boundary field 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
Map surface field to field.
Definition: volSurfaceMapping.C:180
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:812
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:135
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::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53