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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "meshToMesh0.H"
30 #include "processorFvPatch.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 defineTypeNameAndDebug(meshToMesh0, 0);
37 }
38 
39 const Foam::scalar Foam::meshToMesh0::directHitTol = 1e-5;
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const fvMesh& meshFrom,
47  const fvMesh& meshTo,
48  const HashTable<word>& patchMap,
49  const wordList& cuttingPatchNames
50 )
51 :
52  fromMesh_(meshFrom),
53  toMesh_(meshTo),
54  patchMap_(patchMap),
55  cellAddressing_(toMesh_.nCells()),
56  boundaryAddressing_(toMesh_.boundaryMesh().size()),
57  inverseDistanceWeightsPtr_(nullptr),
58  inverseVolumeWeightsPtr_(nullptr),
59  cellToCellAddressingPtr_(nullptr),
60  V_(0)
61 {
62  forAll(fromMesh_.boundaryMesh(), patchi)
63  {
64  fromMeshPatches_.insert
65  (
66  fromMesh_.boundaryMesh()[patchi].name(),
67  patchi
68  );
69  }
70 
71  forAll(toMesh_.boundaryMesh(), patchi)
72  {
73  toMeshPatches_.insert
74  (
75  toMesh_.boundaryMesh()[patchi].name(),
76  patchi
77  );
78  }
79 
80  forAll(cuttingPatchNames, i)
81  {
82  if (toMeshPatches_.found(cuttingPatchNames[i]))
83  {
84  cuttingPatches_.insert
85  (
86  cuttingPatchNames[i],
87  toMeshPatches_.find(cuttingPatchNames[i])()
88  );
89  }
90  else
91  {
93  << "Cannot find cutting-patch " << cuttingPatchNames[i]
94  << " in destination mesh" << endl;
95  }
96  }
97 
98  forAll(toMesh_.boundaryMesh(), patchi)
99  {
100  // Add the processor patches in the toMesh to the cuttingPatches list
101  if (isA<processorPolyPatch>(toMesh_.boundaryMesh()[patchi]))
102  {
103  cuttingPatches_.insert
104  (
105  toMesh_.boundaryMesh()[patchi].name(),
106  patchi
107  );
108  }
109  }
110 
111  calcAddressing();
112 }
113 
114 
116 (
117  const fvMesh& meshFrom,
118  const fvMesh& meshTo
119 )
120 :
121  fromMesh_(meshFrom),
122  toMesh_(meshTo),
123  cellAddressing_(toMesh_.nCells()),
124  boundaryAddressing_(toMesh_.boundaryMesh().size()),
125  inverseDistanceWeightsPtr_(nullptr),
126  inverseVolumeWeightsPtr_(nullptr),
127  cellToCellAddressingPtr_(nullptr),
128  V_(0)
129 {
130  // check whether both meshes have got the same number
131  // of boundary patches
132  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
133  {
135  << "Incompatible meshes: different number of patches, "
136  << "fromMesh = " << fromMesh_.boundary().size()
137  << ", toMesh = " << toMesh_.boundary().size()
138  << exit(FatalError);
139  }
140 
141  forAll(fromMesh_.boundaryMesh(), patchi)
142  {
143  if
144  (
145  fromMesh_.boundaryMesh()[patchi].name()
146  != toMesh_.boundaryMesh()[patchi].name()
147  )
148  {
150  << "Incompatible meshes: different patch names for patch "
151  << patchi
152  << ", fromMesh = " << fromMesh_.boundary()[patchi].name()
153  << ", toMesh = " << toMesh_.boundary()[patchi].name()
154  << exit(FatalError);
155  }
156 
157  if
158  (
159  fromMesh_.boundaryMesh()[patchi].type()
160  != toMesh_.boundaryMesh()[patchi].type()
161  )
162  {
164  << "Incompatible meshes: different patch types for patch "
165  << patchi
166  << ", fromMesh = " << fromMesh_.boundary()[patchi].type()
167  << ", toMesh = " << toMesh_.boundary()[patchi].type()
168  << exit(FatalError);
169  }
170 
171  fromMeshPatches_.insert
172  (
173  fromMesh_.boundaryMesh()[patchi].name(),
174  patchi
175  );
176 
177  toMeshPatches_.insert
178  (
179  toMesh_.boundaryMesh()[patchi].name(),
180  patchi
181  );
182 
183  patchMap_.insert
184  (
185  toMesh_.boundaryMesh()[patchi].name(),
186  fromMesh_.boundaryMesh()[patchi].name()
187  );
188  }
189 
190  calcAddressing();
191 }
192 
193 
194 // ************************************************************************* //
processorFvPatch.H
meshToMesh0.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::meshToMesh0::meshToMesh0
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:45
Foam::FatalError
error FatalError
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::List< word >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328