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