correctedCellVolumeWeightMethod.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) 2015 OpenFOAM Foundation
9  Copyright (C) 2015 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 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(correctedCellVolumeWeightMethod, 0);
38  (
39  meshToMeshMethod,
40  correctedCellVolumeWeightMethod,
41  components
42  );
43 }
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 (
49  labelListList& srcToTgtCellAddr,
50  scalarListList& srcToTgtCellWght,
51  pointListList& srcToTgtCellVec,
52  labelListList& tgtToSrcCellAddr,
53  scalarListList& tgtToSrcCellWght,
54  pointListList& tgtToSrcCellVec,
55  const label srcSeedI,
56  const label tgtSeedI,
57  const labelList& srcCellIDs,
58  boolList& mapFlag,
59  label& startSeedI
60 )
61 {
62  label srcCellI = srcSeedI;
63  label tgtCellI = tgtSeedI;
64 
65  List<DynamicList<label>> srcToTgtAddr(src_.nCells());
66  List<DynamicList<scalar>> srcToTgtWght(src_.nCells());
67  List<DynamicList<point>> srcToTgtVec(src_.nCells());
68 
69  List<DynamicList<label>> tgtToSrcAddr(tgt_.nCells());
70  List<DynamicList<scalar>> tgtToSrcWght(tgt_.nCells());
71  List<DynamicList<point>> tgtToSrcVec(tgt_.nCells());
72 
73  // list of tgt cell neighbour cells
74  DynamicList<label> nbrTgtCells(10);
75 
76  // list of tgt cells currently visited for srcCellI to avoid multiple hits
77  DynamicList<label> visitedTgtCells(10);
78 
79  // list to keep track of tgt cells used to seed src cells
80  labelList seedCells(src_.nCells(), -1);
81  seedCells[srcCellI] = tgtCellI;
82 
83  const scalarField& srcVol = src_.cellVolumes();
84  const pointField& srcCc = src_.cellCentres();
85  const pointField& tgtCc = tgt_.cellCentres();
86 
87  do
88  {
89  nbrTgtCells.clear();
90  visitedTgtCells.clear();
91 
92  // append initial target cell and neighbours
93  nbrTgtCells.append(tgtCellI);
94  appendNbrCells(tgtCellI, tgt_, visitedTgtCells, nbrTgtCells);
95 
96  do
97  {
98  tgtCellI = nbrTgtCells.remove();
99  visitedTgtCells.append(tgtCellI);
100 
101  Tuple2<scalar, point> vol = interVolAndCentroid
102  (
103  srcCellI,
104  tgtCellI
105  );
106 
107  // accumulate addressing and weights for valid intersection
108  if (vol.first()/srcVol[srcCellI] > tolerance_)
109  {
110  // store src/tgt cell pair
111  srcToTgtAddr[srcCellI].append(tgtCellI);
112  srcToTgtWght[srcCellI].append(vol.first());
113  srcToTgtVec[srcCellI].append(vol.second()-tgtCc[tgtCellI]);
114 
115  tgtToSrcAddr[tgtCellI].append(srcCellI);
116  tgtToSrcWght[tgtCellI].append(vol.first());
117  tgtToSrcVec[tgtCellI].append(vol.second()-srcCc[srcCellI]);
118 
119  appendNbrCells(tgtCellI, tgt_, visitedTgtCells, nbrTgtCells);
120 
121  // accumulate intersection volume
122  V_ += vol.first();
123  }
124  }
125  while (!nbrTgtCells.empty());
126 
127  mapFlag[srcCellI] = false;
128 
129  // find new source seed cell
130  setNextCells
131  (
132  startSeedI,
133  srcCellI,
134  tgtCellI,
135  srcCellIDs,
136  mapFlag,
137  visitedTgtCells,
138  seedCells
139  );
140  }
141  while (srcCellI != -1);
142 
143  // transfer addressing into persistent storage
144  forAll(srcToTgtCellAddr, i)
145  {
146  srcToTgtCellAddr[i].transfer(srcToTgtAddr[i]);
147  srcToTgtCellWght[i].transfer(srcToTgtWght[i]);
148  srcToTgtCellVec[i].transfer(srcToTgtVec[i]);
149  }
150 
151  forAll(tgtToSrcCellAddr, i)
152  {
153  tgtToSrcCellAddr[i].transfer(tgtToSrcAddr[i]);
154  tgtToSrcCellWght[i].transfer(tgtToSrcWght[i]);
155  tgtToSrcCellVec[i].transfer(tgtToSrcVec[i]);
156  }
157 
158 
159  if (debug%2)
160  {
161  // At this point the overlaps are still in volume so we could
162  // get out the relative error
163  forAll(srcToTgtCellAddr, cellI)
164  {
165  scalar srcVol = src_.cellVolumes()[cellI];
166  scalar tgtVol = sum(srcToTgtCellWght[cellI]);
167 
168  if (mag(srcVol) > ROOTVSMALL && mag((tgtVol-srcVol)/srcVol) > 1e-6)
169  {
171  << "At cell " << cellI << " cc:"
172  << src_.cellCentres()[cellI]
173  << " vol:" << srcVol
174  << " total overlap volume:" << tgtVol
175  << endl;
176  }
177  }
178 
179  forAll(tgtToSrcCellAddr, cellI)
180  {
181  scalar tgtVol = tgt_.cellVolumes()[cellI];
182  scalar srcVol = sum(tgtToSrcCellWght[cellI]);
183 
184  if (mag(tgtVol) > ROOTVSMALL && mag((srcVol-tgtVol)/tgtVol) > 1e-6)
185  {
187  << "At cell " << cellI << " cc:"
188  << tgt_.cellCentres()[cellI]
189  << " vol:" << tgtVol
190  << " total overlap volume:" << srcVol
191  << endl;
192  }
193  }
194  }
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
199 
201 (
202  const polyMesh& src,
203  const polyMesh& tgt
204 )
205 :
206  cellVolumeWeightMethod(src, tgt)
207 {}
208 
209 
210 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
211 
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
219 (
220  labelListList& srcToTgtAddr,
221  scalarListList& srcToTgtWght,
222  pointListList& srcToTgtVec,
223 
224  labelListList& tgtToSrcAddr,
225  scalarListList& tgtToSrcWght,
226  pointListList& tgtToSrcVec
227 )
228 {
229  bool ok = initialise
230  (
231  srcToTgtAddr,
232  srcToTgtWght,
233  tgtToSrcAddr,
234  tgtToSrcWght
235  );
236 
237  if (!ok)
238  {
239  return;
240  }
241 
242  srcToTgtVec.setSize(srcToTgtAddr.size());
243  tgtToSrcVec.setSize(tgtToSrcAddr.size());
244 
245 
246  // (potentially) participating source mesh cells
247  const labelList srcCellIDs(maskCells());
248 
249  // list to keep track of whether src cell can be mapped
250  boolList mapFlag(src_.nCells(), false);
251  boolUIndList(mapFlag, srcCellIDs) = true;
252 
253  // find initial point in tgt mesh
254  label srcSeedI = -1;
255  label tgtSeedI = -1;
256  label startSeedI = 0;
257 
258  bool startWalk =
259  findInitialSeeds
260  (
261  srcCellIDs,
262  mapFlag,
263  startSeedI,
264  srcSeedI,
265  tgtSeedI
266  );
267 
268  if (startWalk)
269  {
270  calculateAddressing
271  (
272  srcToTgtAddr,
273  srcToTgtWght,
274  srcToTgtVec,
275  tgtToSrcAddr,
276  tgtToSrcWght,
277  tgtToSrcVec,
278  srcSeedI,
279  tgtSeedI,
280  srcCellIDs,
281  mapFlag,
282  startSeedI
283  );
284  }
285  else
286  {
287  // if meshes are collocated, after inflating the source mesh bounding
288  // box tgt mesh cells may be transferred, but may still not overlap
289  // with the source mesh
290  return;
291  }
292 }
293 
294 
295 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
correctedCellVolumeWeightMethod.H
Foam::DynamicList< label >
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::correctedCellVolumeWeightMethod::calculateAddressing
void calculateAddressing(labelListList &srcToTgtCellAddr, scalarListList &srcToTgtCellWght, pointListList &srcToTgtCellVec, labelListList &tgtToSrcCellAddr, scalarListList &tgtToSrcCellWght, pointListList &tgtToSrcCellVec, const label srcSeedI, const label tgtSeedI, const labelList &srcCellIDs, boolList &mapFlag, label &startSeedI)
Calculate the mesh-to-mesh addressing and weights.
Definition: correctedCellVolumeWeightMethod.C:48
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::DynamicList::clear
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
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::boolUIndList
UIndirectList< bool > boolUIndList
UIndirectList of bools.
Definition: UIndirectList.H:54
Foam::Field< scalar >
Foam::correctedCellVolumeWeightMethod::calculate
virtual void calculate(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, pointListList &srcToTgtVec, labelListList &tgtToSrcAddr, scalarListList &tgtToSrcWght, pointListList &tgtToSrcVec)
Calculate addressing and weights and optionally offset vectors.
Definition: correctedCellVolumeWeightMethod.C:219
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::correctedCellVolumeWeightMethod::~correctedCellVolumeWeightMethod
virtual ~correctedCellVolumeWeightMethod()
Destructor.
Definition: correctedCellVolumeWeightMethod.C:212
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::correctedCellVolumeWeightMethod::correctedCellVolumeWeightMethod
correctedCellVolumeWeightMethod(const correctedCellVolumeWeightMethod &)=delete
No copy construct.
Foam::List< labelList >
Foam::cellVolumeWeightMethod
Cell-volume-weighted mesh-to-mesh interpolation class.
Definition: cellVolumeWeightMethod.H:52
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::Tuple2::second
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
Foam::Tuple2< scalar, point >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::Tuple2::first
const T1 & first() const noexcept
Return first.
Definition: Tuple2.H:118
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328