FreeFOAM The Cross-Platform CFD Toolkit
calculateMeshToMeshWeights.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) 1991-2010 OpenCFD Ltd.
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 <sampling/meshToMesh.H>
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void meshToMesh::calculateInverseDistanceWeights() const
36 {
37  if (debug)
38  {
39  Info<< "meshToMesh::calculateInverseDistanceWeights() : "
40  << "calculating inverse distance weighting factors" << endl;
41  }
42 
43  if (inverseDistanceWeightsPtr_)
44  {
45  FatalErrorIn("meshToMesh::calculateInverseDistanceWeights()")
46  << "weighting factors already calculated"
47  << exit(FatalError);
48  }
49 
50  inverseDistanceWeightsPtr_ = new scalarListList(toMesh_.nCells());
51  scalarListList& invDistCoeffs = *inverseDistanceWeightsPtr_;
52 
53  // get reference to source mesh data
54  const labelListList& cc = fromMesh_.cellCells();
55  const vectorField& centreFrom = fromMesh_.C().internalField();
56  const vectorField& centreTo = toMesh_.C().internalField();
57 
58  forAll (cellAddressing_, celli)
59  {
60  if (cellAddressing_[celli] != -1)
61  {
62  const vector& target = centreTo[celli];
63  scalar m = mag(target - centreFrom[cellAddressing_[celli]]);
64 
65  const labelList& neighbours = cc[cellAddressing_[celli]];
66 
67  // if the nearest cell is a boundary cell or there is a direct hit,
68  // pick up the value
69  if
70  (
71  m < directHitTol // Direct hit
72  || neighbours.empty()
73  )
74  {
75  invDistCoeffs[celli].setSize(1);
76  invDistCoeffs[celli][0] = 1.0;
77  }
78  else
79  {
80  invDistCoeffs[celli].setSize(neighbours.size() + 1);
81 
82  // The first coefficient corresponds to the centre cell.
83  // The rest is ordered in the same way as the cellCells list.
84  scalar invDist = 1.0/m;
85  invDistCoeffs[celli][0] = invDist;
86  scalar sumInvDist = invDist;
87 
88  // now add the neighbours
89  forAll (neighbours, ni)
90  {
91  invDist = 1.0/mag(target - centreFrom[neighbours[ni]]);
92  invDistCoeffs[celli][ni + 1] = invDist;
93  sumInvDist += invDist;
94  }
95 
96  // divide by the total inverse-distance
97  forAll (invDistCoeffs[celli], i)
98  {
99  invDistCoeffs[celli][i] /= sumInvDist;
100  }
101  }
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 const scalarListList& meshToMesh::inverseDistanceWeights() const
110 {
111  if (!inverseDistanceWeightsPtr_)
112  {
113  calculateInverseDistanceWeights();
114  }
115 
116  return *inverseDistanceWeightsPtr_;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace Foam
123 
124 // ************************ vim: set sw=4 sts=4 et: ************************ //