calculateMeshToMesh0Weights.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) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "meshToMesh0.H"
29 #include "tetOverlapVolume.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::meshToMesh0::calculateInverseDistanceWeights() const
34 {
35  if (debug)
36  {
38  << "Calculating inverse distance weighting factors" << endl;
39  }
40 
41  if (inverseDistanceWeightsPtr_)
42  {
44  << "weighting factors already calculated"
45  << exit(FatalError);
46  }
47 
48  //- Initialise overlap volume to zero
49  V_ = 0.0;
50 
51  inverseDistanceWeightsPtr_ = new scalarListList(toMesh_.nCells());
52  scalarListList& invDistCoeffs = *inverseDistanceWeightsPtr_;
53 
54  // get reference to source mesh data
55  const labelListList& cc = fromMesh_.cellCells();
56  const vectorField& centreFrom = fromMesh_.C();
57  const vectorField& centreTo = toMesh_.C();
58 
59  forAll(cellAddressing_, celli)
60  {
61  if (cellAddressing_[celli] != -1)
62  {
63  const vector& target = centreTo[celli];
64  scalar m = mag(target - centreFrom[cellAddressing_[celli]]);
65 
66  const labelList& neighbours = cc[cellAddressing_[celli]];
67 
68  // if the nearest cell is a boundary cell or there is a direct hit,
69  // pick up the value
70  label directCelli = -1;
71  if (m < directHitTol || neighbours.empty())
72  {
73  directCelli = celli;
74  }
75  else
76  {
77  forAll(neighbours, ni)
78  {
79  scalar nm = mag(target - centreFrom[neighbours[ni]]);
80  if (nm < directHitTol)
81  {
82  directCelli = neighbours[ni];
83  break;
84  }
85  }
86  }
87 
88 
89  if (directCelli != -1)
90  {
91  // Direct hit
92  invDistCoeffs[directCelli].setSize(1);
93  invDistCoeffs[directCelli][0] = 1.0;
94  V_ += fromMesh_.V()[cellAddressing_[directCelli]];
95  }
96  else
97  {
98  invDistCoeffs[celli].setSize(neighbours.size() + 1);
99 
100  // The first coefficient corresponds to the centre cell.
101  // The rest is ordered in the same way as the cellCells list.
102  scalar invDist = 1.0/m;
103  invDistCoeffs[celli][0] = invDist;
104  scalar sumInvDist = invDist;
105 
106  // now add the neighbours
107  forAll(neighbours, ni)
108  {
109  invDist = 1.0/mag(target - centreFrom[neighbours[ni]]);
110  invDistCoeffs[celli][ni + 1] = invDist;
111  sumInvDist += invDist;
112  }
113 
114  // divide by the total inverse-distance
115  forAll(invDistCoeffs[celli], i)
116  {
117  invDistCoeffs[celli][i] /= sumInvDist;
118  }
119 
120 
121  V_ +=
122  invDistCoeffs[celli][0]
123  *fromMesh_.V()[cellAddressing_[celli]];
124  for (label i = 1; i < invDistCoeffs[celli].size(); i++)
125  {
126  V_ +=
127  invDistCoeffs[celli][i]*fromMesh_.V()[neighbours[i-1]];
128  }
129  }
130  }
131  }
132 }
133 
134 
135 void Foam::meshToMesh0::calculateInverseVolumeWeights() const
136 {
137  if (debug)
138  {
140  << "Calculating inverse volume weighting factors" << endl;
141  }
142 
143  if (inverseVolumeWeightsPtr_)
144  {
146  << "weighting factors already calculated"
147  << exit(FatalError);
148  }
149 
150  //- Initialise overlap volume to zero
151  V_ = 0.0;
152 
153  inverseVolumeWeightsPtr_ = new scalarListList(toMesh_.nCells());
154  scalarListList& invVolCoeffs = *inverseVolumeWeightsPtr_;
155 
156  const labelListList& cellToCell = cellToCellAddressing();
157 
158  tetOverlapVolume overlapEngine;
159 
160  forAll(cellToCell, celli)
161  {
162  const labelList& overlapCells = cellToCell[celli];
163 
164  if (overlapCells.size() > 0)
165  {
166  invVolCoeffs[celli].setSize(overlapCells.size());
167 
168  forAll(overlapCells, j)
169  {
170  label cellFrom = overlapCells[j];
171  treeBoundBox bbFromMesh
172  (
173  pointField
174  (
175  fromMesh_.points(),
176  fromMesh_.cellPoints()[cellFrom]
177  )
178  );
179 
180  scalar v = overlapEngine.cellCellOverlapVolumeMinDecomp
181  (
182  toMesh_,
183  celli,
184 
185  fromMesh_,
186  cellFrom,
187  bbFromMesh
188  );
189  invVolCoeffs[celli][j] = v/toMesh_.V()[celli];
190 
191  V_ += v;
192  }
193  }
194  }
195 }
196 
197 
198 void Foam::meshToMesh0::calculateCellToCellAddressing() const
199 {
200  if (debug)
201  {
203  << "Calculating cell to cell addressing" << endl;
204  }
205 
206  if (cellToCellAddressingPtr_)
207  {
209  << "addressing already calculated"
210  << exit(FatalError);
211  }
212 
213  //- Initialise overlap volume to zero
214  V_ = 0.0;
215 
216  tetOverlapVolume overlapEngine;
217 
218  cellToCellAddressingPtr_ = new labelListList(toMesh_.nCells());
219  labelListList& cellToCell = *cellToCellAddressingPtr_;
220 
221 
222  forAll(cellToCell, iTo)
223  {
224  const labelList overLapCells =
225  overlapEngine.overlappingCells(fromMesh_, toMesh_, iTo);
226  if (overLapCells.size() > 0)
227  {
228  cellToCell[iTo].setSize(overLapCells.size());
229  forAll(overLapCells, j)
230  {
231  cellToCell[iTo][j] = overLapCells[j];
232  V_ += fromMesh_.V()[overLapCells[j]];
233  }
234  }
235  }
236 }
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 const Foam::scalarListList& Foam::meshToMesh0::inverseDistanceWeights() const
241 {
242  if (!inverseDistanceWeightsPtr_)
243  {
244  calculateInverseDistanceWeights();
245  }
246 
247  return *inverseDistanceWeightsPtr_;
248 }
249 
250 
251 const Foam::scalarListList& Foam::meshToMesh0::inverseVolumeWeights() const
252 {
253  if (!inverseVolumeWeightsPtr_)
254  {
255  calculateInverseVolumeWeights();
256  }
257 
258  return *inverseVolumeWeightsPtr_;
259 }
260 
261 
262 const Foam::labelListList& Foam::meshToMesh0::cellToCellAddressing() const
263 {
264  if (!cellToCellAddressingPtr_)
265  {
266  calculateCellToCellAddressing();
267  }
268 
269  return *cellToCellAddressingPtr_;
270 }
271 
272 
273 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::scalarListList
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
meshToMesh0.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
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:62
Foam::primitiveMesh::cellCells
const labelListList & cellCells() const
Definition: primitiveMeshCellCells.C:102
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:358
Foam::FatalError
error FatalError
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
tetOverlapVolume.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::List< scalarList >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:190