meshToMeshMethod.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) 2013-2014 OpenFOAM Foundation
9  Copyright (C) 2018-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 "meshToMeshMethod.H"
30 #include "tetOverlapVolume.H"
31 #include "OFstream.H"
32 #include "Time.H"
33 #include "treeBoundBox.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(meshToMeshMethod, 0);
40  defineRunTimeSelectionTable(meshToMeshMethod, components);
41 }
42 
43 Foam::scalar Foam::meshToMeshMethod::tolerance_ = 1e-6;
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  boundBox intersectBb
50  (
51  max(src_.bounds().min(), tgt_.bounds().min()),
52  min(src_.bounds().max(), tgt_.bounds().max())
53  );
54 
55  intersectBb.inflate(0.01);
56 
57  const cellList& srcCells = src_.cells();
58  const faceList& srcFaces = src_.faces();
59  const pointField& srcPts = src_.points();
60 
62  forAll(srcCells, srcI)
63  {
64  boundBox cellBb(srcCells[srcI].points(srcFaces, srcPts), false);
65  if (intersectBb.overlaps(cellBb))
66  {
67  cells.append(srcI);
68  }
69  }
70 
71  if (debug)
72  {
73  Pout<< "participating source mesh cells: " << cells.size() << endl;
74  }
75 
76  return cells;
77 }
78 
79 
81 (
82  const label srcCelli,
83  const label tgtCelli
84 ) const
85 {
86  scalar threshold = tolerance_*src_.cellVolumes()[srcCelli];
87 
88  tetOverlapVolume overlapEngine;
89 
90  // Note: avoid demand-driven construction of cellPoints
91  // treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCelli]);
92  const UList<label>& cellFaces = tgt_.cells()[tgtCelli];
93  treeBoundBox bbTgtCell(tgt_.points(), tgt_.faces()[cellFaces[0]]);
94  for (label i = 1; i < cellFaces.size(); ++i)
95  {
96  bbTgtCell.add(tgt_.points(), tgt_.faces()[cellFaces[i]]);
97  }
98 
99  return overlapEngine.cellCellOverlapMinDecomp
100  (
101  src_,
102  srcCelli,
103  tgt_,
104  tgtCelli,
105  bbTgtCell,
106  threshold
107  );
108 }
109 
110 
112 (
113  const label srcCelli,
114  const label tgtCelli
115 ) const
116 {
117  tetOverlapVolume overlapEngine;
118 
119  // Note: avoid demand-driven construction of cellPoints
120  // treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCelli]);
121  const UList<label>& cellFaces = tgt_.cells()[tgtCelli];
122  treeBoundBox bbTgtCell(tgt_.points(), tgt_.faces()[cellFaces[0]]);
123  for (label i = 1; i < cellFaces.size(); ++i)
124  {
125  bbTgtCell.add(tgt_.points(), tgt_.faces()[cellFaces[i]]);
126  }
127 
128  scalar vol = overlapEngine.cellCellOverlapVolumeMinDecomp
129  (
130  src_,
131  srcCelli,
132  tgt_,
133  tgtCelli,
134  bbTgtCell
135  );
136 
137  return vol;
138 }
139 
140 
143 (
144  const label srcCelli,
145  const label tgtCelli
146 )
147 {
148  tetOverlapVolume overlapEngine;
149 
150  // Note: avoid demand-driven construction of cellPoints
151  // treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCelli]);
152  const UList<label>& cellFaces = tgt_.cells()[tgtCelli];
153  treeBoundBox bbTgtCell(tgt_.points(), tgt_.faces()[cellFaces[0]]);
154  for (label i = 1; i < cellFaces.size(); ++i)
155  {
156  bbTgtCell.add(tgt_.points(), tgt_.faces()[cellFaces[i]]);
157  }
158 
159  Tuple2<scalar, point> volAndInertia =
160  overlapEngine.cellCellOverlapMomentMinDecomp
161  (
162  src_,
163  srcCelli,
164  tgt_,
165  tgtCelli,
166  bbTgtCell
167  );
168 
169  // Convert from inertia to centroid
170  if (volAndInertia.first() <= ROOTVSMALL)
171  {
172  volAndInertia.first() = 0.0;
173  volAndInertia.second() = Zero;
174  }
175  else
176  {
177  volAndInertia.second() /= volAndInertia.first();
178  }
179 
180  return volAndInertia;
181 }
182 
183 
185 (
186  const label celli,
187  const polyMesh& mesh,
188  const DynamicList<label>& visitedCells,
189  DynamicList<label>& nbrCellIDs
190 ) const
191 {
192  const labelList& nbrCells = mesh.cellCells()[celli];
193 
194  // filter out cells already visited from cell neighbours
195  for (const label nbrCelli : nbrCells)
196  {
197  if (!visitedCells.found(nbrCelli))
198  {
199  nbrCellIDs.appendUniq(nbrCelli);
200  }
201  }
202 }
203 
204 
206 (
207  labelListList& srcToTgtAddr,
208  scalarListList& srcToTgtWght,
209  labelListList& tgtToSrcAddr,
210  scalarListList& tgtToSrcWght
211 ) const
212 {
213  srcToTgtAddr.setSize(src_.nCells());
214  srcToTgtWght.setSize(src_.nCells());
215  tgtToSrcAddr.setSize(tgt_.nCells());
216  tgtToSrcWght.setSize(tgt_.nCells());
217 
218  if (!src_.nCells())
219  {
220  return false;
221  }
222  else if (!tgt_.nCells())
223  {
224  if (debug)
225  {
226  Pout<< "mesh interpolation: have " << src_.nCells() << " source "
227  << " cells but no target cells" << endl;
228  }
229 
230  return false;
231  }
232 
233  return true;
234 }
235 
236 
237 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
238 
239 Foam::meshToMeshMethod::meshToMeshMethod
240 (
241  const polyMesh& src,
242  const polyMesh& tgt
243 )
244 :
245  src_(src),
246  tgt_(tgt),
247  V_(0.0)
248 {
249  if (!src_.nCells() || !tgt_.nCells())
250  {
251  if (debug)
252  {
253  Pout<< "mesh interpolation: cells not on processor: Source cells = "
254  << src_.nCells() << ", target cells = " << tgt_.nCells()
255  << endl;
256  }
257  }
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
262 
264 {}
265 
266 
267 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
268 
270 (
271  const polyMesh& mesh1,
272  const polyMesh& mesh2,
273  const labelListList& mesh1ToMesh2Addr
274 ) const
275 {
276  Pout<< "Source size = " << mesh1.nCells() << endl;
277  Pout<< "Target size = " << mesh2.nCells() << endl;
278 
279  word fName("addressing_" + mesh1.name() + "_to_" + mesh2.name());
280 
281  if (Pstream::parRun())
282  {
283  fName = fName + "_proc" + Foam::name(Pstream::myProcNo());
284  }
285 
286  OFstream os(src_.time().path()/fName + ".obj");
287 
288  label vertI = 0;
289  forAll(mesh1ToMesh2Addr, i)
290  {
291  const labelList& addr = mesh1ToMesh2Addr[i];
292  forAll(addr, j)
293  {
294  label celli = addr[j];
295  const vector& c0 = mesh1.cellCentres()[i];
296 
297  const cell& c = mesh2.cells()[celli];
298  const pointField pts(c.points(mesh2.faces(), mesh2.points()));
299  forAll(pts, j)
300  {
301  const point& p = pts[j];
302  os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
303  vertI++;
304  os << "v " << c0.x() << ' ' << c0.y() << ' ' << c0.z()
305  << nl;
306  vertI++;
307  os << "l " << vertI - 1 << ' ' << vertI << nl;
308  }
309  }
310  }
311 }
312 
313 
314 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::meshToMeshMethod::~meshToMeshMethod
virtual ~meshToMeshMethod()
Destructor.
Definition: meshToMeshMethod.C:263
Foam::DynamicList< label >
Foam::tetOverlapVolume::cellCellOverlapVolumeMinDecomp
scalar cellCellOverlapVolumeMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB) const
Calculates the overlap volume.
Definition: tetOverlapVolume.C:96
Foam::meshToMeshMethod::appendNbrCells
virtual void appendNbrCells(const label tgtCelli, const polyMesh &mesh, const DynamicList< label > &visitedTgtCells, DynamicList< label > &nbrTgtCellIDs) const
Append target cell neighbour cells to cellIDs list.
Definition: meshToMeshMethod.C:185
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::boundBox::inflate
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:175
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::meshToMeshMethod::interVolAndCentroid
virtual Tuple2< scalar, point > interVolAndCentroid(const label srcCellI, const label tgtCellI)
Return the intersection volume and centroid between two cells.
Definition: meshToMeshMethod.C:143
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::meshToMeshMethod::interVol
virtual scalar interVol(const label srcCelli, const label tgtCelli) const
Return the intersection volume between two cells.
Definition: meshToMeshMethod.C:112
Foam::meshToMeshMethod::tolerance_
static scalar tolerance_
Tolerance used in volume overlap calculations.
Definition: meshToMeshMethod.H:68
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::meshToMeshMethod::writeConnectivity
void writeConnectivity(const polyMesh &mesh1, const polyMesh &mesh2, const labelListList &mesh1ToMesh2Addr) const
Write the connectivity (debugging)
Definition: meshToMeshMethod.C:270
Foam::meshToMeshMethod::intersect
virtual bool intersect(const label srcCelli, const label tgtCelli) const
Return the true if cells intersect.
Definition: meshToMeshMethod.C:81
Foam::Field< vector >
treeBoundBox.H
Foam::meshToMeshMethod::tgt_
const polyMesh & tgt_
Reference to the target mesh.
Definition: meshToMeshMethod.H:62
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::primitiveMesh::cellCells
const labelListList & cellCells() const
Definition: primitiveMeshCellCells.C:102
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::meshToMeshMethod::src_
const polyMesh & src_
Reference to the source mesh.
Definition: meshToMeshMethod.H:59
Foam::tetOverlapVolume::cellCellOverlapMomentMinDecomp
Tuple2< scalar, point > cellCellOverlapMomentMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB) const
Calculates the overlap volume and moment.
Definition: tetOverlapVolume.C:122
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::meshToMeshMethod::maskCells
labelList maskCells() const
Return src cell IDs for the overlap region.
Definition: meshToMeshMethod.C:47
meshToMeshMethod.H
tetOverlapVolume.H
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:450
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Time.H
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::meshToMeshMethod::initialise
virtual bool initialise(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght) const
Definition: meshToMeshMethod.C:206
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::Tuple2::second
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:60
Foam::tetOverlapVolume
Calculates the overlap volume of two cells using tetrahedral decomposition.
Definition: tetOverlapVolume.H:56
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::Tuple2::first
const T1 & first() const noexcept
Return first.
Definition: Tuple2.H:118
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::DynamicList::appendUniq
label appendUniq(const T &val)
Append an element if not already in the list.
Definition: DynamicListI.H:689
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191