meshToMesh0.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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_(nullptr),
56  inverseVolumeWeightsPtr_(nullptr),
57  cellToCellAddressingPtr_(nullptr),
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  {
91  << "Cannot find cutting-patch " << cuttingPatchNames[i]
92  << " in destination mesh" << endl;
93  }
94  }
95 
96  forAll(toMesh_.boundaryMesh(), patchi)
97  {
98  // Add the processor patches in the toMesh to the cuttingPatches list
99  if (isA<processorPolyPatch>(toMesh_.boundaryMesh()[patchi]))
100  {
101  cuttingPatches_.insert
102  (
103  toMesh_.boundaryMesh()[patchi].name(),
104  patchi
105  );
106  }
107  }
108 
109  calcAddressing();
110 }
111 
112 
114 (
115  const fvMesh& meshFrom,
116  const fvMesh& meshTo
117 )
118 :
119  fromMesh_(meshFrom),
120  toMesh_(meshTo),
121  cellAddressing_(toMesh_.nCells()),
122  boundaryAddressing_(toMesh_.boundaryMesh().size()),
123  inverseDistanceWeightsPtr_(nullptr),
124  inverseVolumeWeightsPtr_(nullptr),
125  cellToCellAddressingPtr_(nullptr),
126  V_(0.0)
127 {
128  // check whether both meshes have got the same number
129  // of boundary patches
130  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
131  {
133  << "Incompatible meshes: different number of patches, "
134  << "fromMesh = " << fromMesh_.boundary().size()
135  << ", toMesh = " << toMesh_.boundary().size()
136  << exit(FatalError);
137  }
138 
139  forAll(fromMesh_.boundaryMesh(), patchi)
140  {
141  if
142  (
143  fromMesh_.boundaryMesh()[patchi].name()
144  != toMesh_.boundaryMesh()[patchi].name()
145  )
146  {
148  << "Incompatible meshes: different patch names for patch "
149  << patchi
150  << ", fromMesh = " << fromMesh_.boundary()[patchi].name()
151  << ", toMesh = " << toMesh_.boundary()[patchi].name()
152  << exit(FatalError);
153  }
154 
155  if
156  (
157  fromMesh_.boundaryMesh()[patchi].type()
158  != toMesh_.boundaryMesh()[patchi].type()
159  )
160  {
162  << "Incompatible meshes: different patch types for patch "
163  << patchi
164  << ", fromMesh = " << fromMesh_.boundary()[patchi].type()
165  << ", toMesh = " << toMesh_.boundary()[patchi].type()
166  << exit(FatalError);
167  }
168 
169  fromMeshPatches_.insert
170  (
171  fromMesh_.boundaryMesh()[patchi].name(),
172  patchi
173  );
174 
175  toMeshPatches_.insert
176  (
177  toMesh_.boundaryMesh()[patchi].name(),
178  patchi
179  );
180 
181  patchMap_.insert
182  (
183  toMesh_.boundaryMesh()[patchi].name(),
184  fromMesh_.boundaryMesh()[patchi].name()
185  );
186  }
187 
188  calcAddressing();
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
193 
195 {
196  deleteDemandDrivenData(inverseDistanceWeightsPtr_);
197  deleteDemandDrivenData(inverseVolumeWeightsPtr_);
198  deleteDemandDrivenData(cellToCellAddressingPtr_);
199 }
200 
201 
202 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
An STL-conforming hash table.
Definition: HashTable.H:61
defineTypeNameAndDebug(combustionModel, 0)
~meshToMesh0()
Destructor.
Definition: meshToMesh0.C:194
Template functions to aid in the implementation of demand driven data.
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
void deleteDemandDrivenData(DataPtr &dataPtr)
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
Namespace for OpenFOAM.