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 | Copyright (C) 2013-2014 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 "meshToMeshMethod.H"
27 #include "tetOverlapVolume.H"
28 #include "OFstream.H"
29 #include "Time.H"
30 #include "treeBoundBox.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(meshToMeshMethod, 0);
37  defineRunTimeSelectionTable(meshToMeshMethod, components);
38 }
39 
40 Foam::scalar Foam::meshToMeshMethod::tolerance_ = 1e-6;
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 {
46  boundBox intersectBb
47  (
48  max(src_.bounds().min(), tgt_.bounds().min()),
49  min(src_.bounds().max(), tgt_.bounds().max())
50  );
51 
52  intersectBb.inflate(0.01);
53 
54  const cellList& srcCells = src_.cells();
55  const faceList& srcFaces = src_.faces();
56  const pointField& srcPts = src_.points();
57 
59  forAll(srcCells, srcI)
60  {
61  boundBox cellBb(srcCells[srcI].points(srcFaces, srcPts), false);
62  if (intersectBb.overlaps(cellBb))
63  {
64  cells.append(srcI);
65  }
66  }
67 
68  if (debug)
69  {
70  Pout<< "participating source mesh cells: " << cells.size() << endl;
71  }
72 
73  return cells;
74 }
75 
76 
78 (
79  const label srcCellI,
80  const label tgtCellI
81 ) const
82 {
83  scalar threshold = tolerance_*src_.cellVolumes()[srcCellI];
84 
85  tetOverlapVolume overlapEngine;
86 
87  treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCellI]);
88 
89  return overlapEngine.cellCellOverlapMinDecomp
90  (
91  src_,
92  srcCellI,
93  tgt_,
94  tgtCellI,
95  bbTgtCell,
96  threshold
97  );
98 }
99 
100 
102 (
103  const label srcCellI,
104  const label tgtCellI
105 ) const
106 {
107  tetOverlapVolume overlapEngine;
108 
109  treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCellI]);
110 
111  scalar vol = overlapEngine.cellCellOverlapVolumeMinDecomp
112  (
113  src_,
114  srcCellI,
115  tgt_,
116  tgtCellI,
117  bbTgtCell
118  );
119 
120  return vol;
121 }
122 
123 
126 (
127  const label srcCellI,
128  const label tgtCellI
129 )
130 {
131  tetOverlapVolume overlapEngine;
132 
133  treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCellI]);
134 
135  Tuple2<scalar, point> volAndInertia =
136  overlapEngine.cellCellOverlapMomentMinDecomp
137  (
138  src_,
139  srcCellI,
140  tgt_,
141  tgtCellI,
142  bbTgtCell
143  );
144 
145  // Convert from inertia to centroid
146  if (volAndInertia.first() <= ROOTVSMALL)
147  {
148  volAndInertia.first() = 0.0;
149  volAndInertia.second() = vector::zero;
150  }
151  else
152  {
153  volAndInertia.second() /= volAndInertia.first();
154  }
155 
156  return volAndInertia;
157 }
158 
159 
161 (
162  const label cellI,
163  const polyMesh& mesh,
164  const DynamicList<label>& visitedCells,
165  DynamicList<label>& nbrCellIDs
166 ) const
167 {
168  const labelList& nbrCells = mesh.cellCells()[cellI];
169 
170  // filter out cells already visited from cell neighbours
171  forAll(nbrCells, i)
172  {
173  label nbrCellI = nbrCells[i];
174 
175  if
176  (
177  (findIndex(visitedCells, nbrCellI) == -1)
178  && (findIndex(nbrCellIDs, nbrCellI) == -1)
179  )
180  {
181  nbrCellIDs.append(nbrCellI);
182  }
183  }
184 }
185 
186 
188 (
189  labelListList& srcToTgtAddr,
190  scalarListList& srcToTgtWght,
191  labelListList& tgtToSrcAddr,
192  scalarListList& tgtToSrcWght
193 ) const
194 {
195  srcToTgtAddr.setSize(src_.nCells());
196  srcToTgtWght.setSize(src_.nCells());
197  tgtToSrcAddr.setSize(tgt_.nCells());
198  tgtToSrcWght.setSize(tgt_.nCells());
199 
200  if (!src_.nCells())
201  {
202  return false;
203  }
204  else if (!tgt_.nCells())
205  {
206  if (debug)
207  {
208  Pout<< "mesh interpolation: hhave " << src_.nCells() << " source "
209  << " cells but no target cells" << endl;
210  }
211 
212  return false;
213  }
214 
215  return true;
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
220 
222 (
223  const polyMesh& src,
224  const polyMesh& tgt
225 )
226 :
227  src_(src),
228  tgt_(tgt),
229  V_(0.0)
230 {
231  if (!src_.nCells() || !tgt_.nCells())
232  {
233  if (debug)
234  {
235  Pout<< "mesh interpolation: cells not on processor: Source cells = "
236  << src_.nCells() << ", target cells = " << tgt_.nCells()
237  << endl;
238  }
239  }
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
244 
246 {}
247 
248 
249 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
250 
252 (
253  const polyMesh& mesh1,
254  const polyMesh& mesh2,
255  const labelListList& mesh1ToMesh2Addr
256 ) const
257 {
258  Pout<< "Source size = " << mesh1.nCells() << endl;
259  Pout<< "Target size = " << mesh2.nCells() << endl;
260 
261  word fName("addressing_" + mesh1.name() + "_to_" + mesh2.name());
262 
263  if (Pstream::parRun())
264  {
265  fName = fName + "_proc" + Foam::name(Pstream::myProcNo());
266  }
267 
268  OFstream os(src_.time().path()/fName + ".obj");
269 
270  label vertI = 0;
271  forAll(mesh1ToMesh2Addr, i)
272  {
273  const labelList& addr = mesh1ToMesh2Addr[i];
274  forAll(addr, j)
275  {
276  label cellI = addr[j];
277  const vector& c0 = mesh1.cellCentres()[i];
278 
279  const cell& c = mesh2.cells()[cellI];
280  const pointField pts(c.points(mesh2.faces(), mesh2.points()));
281  forAll(pts, j)
282  {
283  const point& p = pts[j];
284  os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
285  vertI++;
286  os << "v " << c0.x() << ' ' << c0.y() << ' ' << c0.z()
287  << nl;
288  vertI++;
289  os << "l " << vertI - 1 << ' ' << vertI << nl;
290  }
291  }
292  }
293 }
294 
295 
296 // ************************************************************************* //
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::meshToMeshMethod::~meshToMeshMethod
virtual ~meshToMeshMethod()
Destructor.
Definition: meshToMeshMethod.C:245
Foam::DynamicList< label >
Foam::Tuple2::second
const Type2 & second() const
Return second.
Definition: Tuple2.H:106
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:103
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:210
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
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:126
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::Tuple2::first
const Type1 & first() const
Return first.
Definition: Tuple2.H:94
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::meshToMeshMethod::interVol
virtual scalar interVol(const label srcCellI, const label tgtCellI) const
Return the intersection volume between two cells.
Definition: meshToMeshMethod.C:102
Foam::meshToMeshMethod::tolerance_
static scalar tolerance_
Tolerance used in volume overlap calculations.
Definition: meshToMeshMethod.H:65
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
OFstream.H
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::meshToMeshMethod::writeConnectivity
void writeConnectivity(const polyMesh &mesh1, const polyMesh &mesh2, const labelListList &mesh1ToMesh2Addr) const
Write the connectivity (debugging)
Definition: meshToMeshMethod.C:252
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
treeBoundBox.H
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::boundBox::overlaps
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
Foam::meshToMeshMethod::tgt_
const polyMesh & tgt_
Reference to the target mesh.
Definition: meshToMeshMethod.H:59
Foam::primitiveMesh::cellCells
const labelListList & cellCells() const
Definition: primitiveMeshCellCells.C:100
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::tetOverlapVolume::cellCellOverlapMinDecomp
static void cellCellOverlapMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB, tetsOp &combineTetsOp)
Cell overlap calculation.
Definition: tetOverlapVolumeTemplates.C:99
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::meshToMeshMethod::src_
const polyMesh & src_
Reference to the source mesh.
Definition: meshToMeshMethod.H:56
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:129
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::meshToMeshMethod::appendNbrCells
virtual void appendNbrCells(const label tgtCellI, const polyMesh &mesh, const DynamicList< label > &visitedTgtCells, DynamicList< label > &nbrTgtCellIDs) const
Append target cell neihgbour cells to cellIDs list.
Definition: meshToMeshMethod.C:161
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::meshToMeshMethod::maskCells
labelList maskCells() const
Return src cell IDs for the overlap region.
Definition: meshToMeshMethod.C:44
meshToMeshMethod.H
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
tetOverlapVolume.H
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::meshToMeshMethod::initialise
virtual bool initialise(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght) const
Definition: meshToMeshMethod.C:188
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::meshToMeshMethod::meshToMeshMethod
meshToMeshMethod(const meshToMeshMethod &)
Disallow default bitwise copy construct.
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::Tuple2
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
Foam::tetOverlapVolume
Calculates the overlap volume of two cells using tetrahedral decomposition.
Definition: tetOverlapVolume.H:54
Foam::meshToMeshMethod::intersect
virtual bool intersect(const label srcCellI, const label tgtCellI) const
Return the true if cells intersect.
Definition: meshToMeshMethod.C:78
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47