meshToMesh0.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) 2011-2015 OpenFOAM Foundation
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 "meshToMesh0.H"
27 #include "processorFvPatch.H"
28 #include "demandDrivenData.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 defineTypeNameAndDebug(meshToMesh0, 0);
35 }
36 
37 const Foam::scalar Foam::meshToMesh0::directHitTol = 1e-5;
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const fvMesh& meshFrom,
45  const fvMesh& meshTo,
46  const HashTable<word>& patchMap,
47  const wordList& cuttingPatchNames
48 )
49 :
50  fromMesh_(meshFrom),
51  toMesh_(meshTo),
52  patchMap_(patchMap),
53  cellAddressing_(toMesh_.nCells()),
54  boundaryAddressing_(toMesh_.boundaryMesh().size()),
55  inverseDistanceWeightsPtr_(NULL),
56  inverseVolumeWeightsPtr_(NULL),
57  cellToCellAddressingPtr_(NULL),
58  V_(0.0)
59 {
60  forAll(fromMesh_.boundaryMesh(), patchi)
61  {
62  fromMeshPatches_.insert
63  (
64  fromMesh_.boundaryMesh()[patchi].name(),
65  patchi
66  );
67  }
68 
69  forAll(toMesh_.boundaryMesh(), patchi)
70  {
71  toMeshPatches_.insert
72  (
73  toMesh_.boundaryMesh()[patchi].name(),
74  patchi
75  );
76  }
77 
78  forAll(cuttingPatchNames, i)
79  {
80  if (toMeshPatches_.found(cuttingPatchNames[i]))
81  {
82  cuttingPatches_.insert
83  (
84  cuttingPatchNames[i],
85  toMeshPatches_.find(cuttingPatchNames[i])()
86  );
87  }
88  else
89  {
90  WarningIn
91  (
92  "meshToMesh0::meshToMesh0"
93  "(const fvMesh& meshFrom, const fvMesh& meshTo,"
94  "const HashTable<word>& patchMap,"
95  "const wordList& cuttingPatchNames)"
96  ) << "Cannot find cutting-patch " << cuttingPatchNames[i]
97  << " in destination mesh" << endl;
98  }
99  }
100 
101  forAll(toMesh_.boundaryMesh(), patchi)
102  {
103  // Add the processor patches in the toMesh to the cuttingPatches list
104  if (isA<processorPolyPatch>(toMesh_.boundaryMesh()[patchi]))
105  {
106  cuttingPatches_.insert
107  (
108  toMesh_.boundaryMesh()[patchi].name(),
109  patchi
110  );
111  }
112  }
113 
114  calcAddressing();
115 }
116 
117 
119 (
120  const fvMesh& meshFrom,
121  const fvMesh& meshTo
122 )
123 :
124  fromMesh_(meshFrom),
125  toMesh_(meshTo),
126  cellAddressing_(toMesh_.nCells()),
127  boundaryAddressing_(toMesh_.boundaryMesh().size()),
128  inverseDistanceWeightsPtr_(NULL),
129  inverseVolumeWeightsPtr_(NULL),
130  cellToCellAddressingPtr_(NULL),
131  V_(0.0)
132 {
133  // check whether both meshes have got the same number
134  // of boundary patches
135  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
136  {
138  (
139  "meshToMesh0::meshToMesh0"
140  "(const fvMesh& meshFrom, const fvMesh& meshTo)"
141  ) << "Incompatible meshes: different number of patches, "
142  << "fromMesh = " << fromMesh_.boundary().size()
143  << ", toMesh = " << toMesh_.boundary().size()
144  << exit(FatalError);
145  }
146 
147  forAll(fromMesh_.boundaryMesh(), patchi)
148  {
149  if
150  (
151  fromMesh_.boundaryMesh()[patchi].name()
152  != toMesh_.boundaryMesh()[patchi].name()
153  )
154  {
156  (
157  "meshToMesh0::meshToMesh0"
158  "(const fvMesh& meshFrom, const fvMesh& meshTo)"
159  ) << "Incompatible meshes: different patch names for patch "
160  << patchi
161  << ", fromMesh = " << fromMesh_.boundary()[patchi].name()
162  << ", toMesh = " << toMesh_.boundary()[patchi].name()
163  << exit(FatalError);
164  }
165 
166  if
167  (
168  fromMesh_.boundaryMesh()[patchi].type()
169  != toMesh_.boundaryMesh()[patchi].type()
170  )
171  {
173  (
174  "meshToMesh0::meshToMesh0"
175  "(const fvMesh& meshFrom, const fvMesh& meshTo)"
176  ) << "Incompatible meshes: different patch types for patch "
177  << patchi
178  << ", fromMesh = " << fromMesh_.boundary()[patchi].type()
179  << ", toMesh = " << toMesh_.boundary()[patchi].type()
180  << exit(FatalError);
181  }
182 
183  fromMeshPatches_.insert
184  (
185  fromMesh_.boundaryMesh()[patchi].name(),
186  patchi
187  );
188 
189  toMeshPatches_.insert
190  (
191  toMesh_.boundaryMesh()[patchi].name(),
192  patchi
193  );
194 
195  patchMap_.insert
196  (
197  toMesh_.boundaryMesh()[patchi].name(),
198  fromMesh_.boundaryMesh()[patchi].name()
199  );
200  }
201 
202  calcAddressing();
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
207 
209 {
210  deleteDemandDrivenData(inverseDistanceWeightsPtr_);
211  deleteDemandDrivenData(inverseVolumeWeightsPtr_);
212  deleteDemandDrivenData(cellToCellAddressingPtr_);
213 }
214 
215 
216 // ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An STL-conforming hash table.
Definition: HashTable.H:61
void deleteDemandDrivenData(DataPtr &dataPtr)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Namespace for OpenFOAM.
meshToMesh0(const fvMesh &fromMesh, const fvMesh &toMesh, const HashTable< word > &patchMap, const wordList &cuttingPatchNames)
Construct from the two meshes, the patch name map for the patches.
Definition: meshToMesh0.C:43
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define WarningIn(functionName)
Report a warning using Foam::Warning.
#define forAll(list, i)
Definition: UList.H:421
label patchi
Template functions to aid in the implementation of demand driven data.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
~meshToMesh0()
Destructor.
Definition: meshToMesh0.C:208
error FatalError
defineTypeNameAndDebug(combustionModel, 0)