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);
41 }
42 
43 Foam::scalar Foam::meshToMeshMethod::tolerance_ = 1e-6;
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  boundBox intersectBb(src_.bounds());
50  intersectBb &= tgt_.bounds();
51  intersectBb.inflate(0.01);
52 
54 
55  for (label srcCelli = 0; srcCelli < src_.nCells(); ++srcCelli)
56  {
57  boundBox cellBb(src_.cellBb(srcCelli));
58  if (intersectBb.overlaps(cellBb))
59  {
60  cells.append(srcCelli);
61  }
62  }
63 
64  if (debug)
65  {
66  Pout<< "participating source mesh cells: " << src_.nCells() << endl;
67  }
68 
69  return cells;
70 }
71 
72 
74 (
75  const label srcCelli,
76  const label tgtCelli
77 ) const
78 {
79  scalar threshold = tolerance_*src_.cellVolumes()[srcCelli];
80 
81  tetOverlapVolume overlapEngine;
82 
83  treeBoundBox bbTgtCell(tgt_.cellBb(tgtCelli));
84 
85  return overlapEngine.cellCellOverlapMinDecomp
86  (
87  src_,
88  srcCelli,
89  tgt_,
90  tgtCelli,
91  bbTgtCell,
92  threshold
93  );
94 }
95 
96 
98 (
99  const label srcCelli,
100  const label tgtCelli
101 ) const
102 {
103  tetOverlapVolume overlapEngine;
104 
105  treeBoundBox bbTgtCell(tgt_.cellBb(tgtCelli));
106 
107  scalar vol = overlapEngine.cellCellOverlapVolumeMinDecomp
108  (
109  src_,
110  srcCelli,
111  tgt_,
112  tgtCelli,
113  bbTgtCell
114  );
116  return vol;
117 }
118 
119 
122 (
123  const label srcCelli,
124  const label tgtCelli
125 )
126 {
127  tetOverlapVolume overlapEngine;
128 
129  treeBoundBox bbTgtCell(tgt_.cellBb(tgtCelli));
130 
131  Tuple2<scalar, point> volAndInertia =
132  overlapEngine.cellCellOverlapMomentMinDecomp
133  (
134  src_,
135  srcCelli,
136  tgt_,
137  tgtCelli,
138  bbTgtCell
139  );
140 
141  // Convert from inertia to centroid
142  if (volAndInertia.first() <= ROOTVSMALL)
143  {
144  volAndInertia.first() = 0.0;
145  volAndInertia.second() = Zero;
146  }
147  else
148  {
149  volAndInertia.second() /= volAndInertia.first();
150  }
151 
152  return volAndInertia;
153 }
154 
155 
157 (
158  const label celli,
159  const polyMesh& mesh,
160  const DynamicList<label>& visitedCells,
161  DynamicList<label>& nbrCellIDs
162 ) const
163 {
164  const labelList& nbrCells = mesh.cellCells()[celli];
165 
166  // filter out cells already visited from cell neighbours
167  for (const label nbrCelli : nbrCells)
168  {
169  if (!visitedCells.contains(nbrCelli))
170  {
171  nbrCellIDs.push_uniq(nbrCelli);
172  }
173  }
174 }
175 
176 
178 (
179  labelListList& srcToTgtAddr,
180  scalarListList& srcToTgtWght,
181  labelListList& tgtToSrcAddr,
182  scalarListList& tgtToSrcWght
183 ) const
184 {
185  srcToTgtAddr.setSize(src_.nCells());
186  srcToTgtWght.setSize(src_.nCells());
187  tgtToSrcAddr.setSize(tgt_.nCells());
188  tgtToSrcWght.setSize(tgt_.nCells());
189 
190  if (!src_.nCells())
191  {
192  return false;
193  }
194  else if (!tgt_.nCells())
195  {
196  if (debug)
197  {
198  Pout<< "mesh interpolation: have " << src_.nCells() << " source "
199  << " cells but no target cells" << endl;
200  }
201 
202  return false;
203  }
204 
205  return true;
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
210 
211 Foam::meshToMeshMethod::meshToMeshMethod
212 (
213  const polyMesh& src,
214  const polyMesh& tgt
215 )
216 :
217  src_(src),
218  tgt_(tgt),
219  V_(0.0)
220 {
221  if (!src_.nCells() || !tgt_.nCells())
222  {
223  if (debug)
224  {
225  Pout<< "mesh interpolation: cells not on processor: Source cells = "
226  << src_.nCells() << ", target cells = " << tgt_.nCells()
227  << endl;
228  }
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
234 
236 {}
237 
238 
239 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240 
242 (
243  const polyMesh& mesh1,
244  const polyMesh& mesh2,
245  const labelListList& mesh1ToMesh2Addr
246 ) const
247 {
248  Pout<< "Source size = " << mesh1.nCells() << endl;
249  Pout<< "Target size = " << mesh2.nCells() << endl;
250 
251  word fName("addressing_" + mesh1.name() + "_to_" + mesh2.name());
252 
253  if (Pstream::parRun())
254  {
255  fName = fName + "_proc" + Foam::name(Pstream::myProcNo());
256  }
257 
258  OFstream os(src_.time().path()/fName + ".obj");
259 
260  label vertI = 0;
261  forAll(mesh1ToMesh2Addr, i)
262  {
263  const labelList& addr = mesh1ToMesh2Addr[i];
264  forAll(addr, j)
265  {
266  label celli = addr[j];
267  const vector& c0 = mesh1.cellCentres()[i];
268 
269  const cell& c = mesh2.cells()[celli];
270  const pointField pts(c.points(mesh2.faces(), mesh2.points()));
271  forAll(pts, j)
272  {
273  const point& p = pts[j];
274  os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
275  vertI++;
276  os << "v " << c0.x() << ' ' << c0.y() << ' ' << c0.z()
277  << nl;
278  vertI++;
279  os << "l " << vertI - 1 << ' ' << vertI << nl;
280  }
281  }
282  }
283 }
284 
285 
286 // ************************************************************************* //
const polyMesh & src_
Reference to the source mesh.
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition: boundBoxI.H:381
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
virtual ~meshToMeshMethod()
Destructor.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
const cellList & cells() const
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
label push_uniq(const T &val)
Append an element if not already in the list.
Definition: ListI.H:285
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
void writeConnectivity(const polyMesh &mesh1, const polyMesh &mesh2, const labelListList &mesh1ToMesh2Addr) const
Write the connectivity (debugging)
virtual scalar interVol(const label srcCelli, const label tgtCelli) const
Return the intersection volume between two cells.
static scalar tolerance_
Tolerance used in volume overlap calculations.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
A class for handling words, derived from Foam::string.
Definition: word.H:63
Calculates the overlap volume of two cells using tetrahedral decomposition.
labelList maskCells() const
Return src cell IDs for the overlap region.
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
const vectorField & cellCentres() const
const polyMesh & tgt_
Reference to the target mesh.
virtual bool initialise(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght) const
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
Base class for mesh-to-mesh calculation methods.
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.
vector point
Point is a vector.
Definition: point.H:37
virtual Tuple2< scalar, point > interVolAndCentroid(const label srcCellI, const label tgtCellI)
Return the intersection volume and centroid between two cells.
label nCells() const noexcept
Number of mesh cells.
const dimensionedScalar c
Speed of light in a vacuum.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
const T1 & first() const noexcept
Access the first element.
Definition: Tuple2.H:132
virtual bool intersect(const label srcCelli, const label tgtCelli) const
Return the true if cells intersect.
boundBox cellBb(const label celli) const
The bounding box for given cell index.
const labelListList & cellCells() const
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
scalar cellCellOverlapVolumeMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB) const
Calculates the overlap volume.
Namespace for OpenFOAM.
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.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127