BUG: cloudSet.C: force early construction of tetBasePtIs to avoid demand-driven comms
[OpenFOAM-2.0.x.git] / src / sampling / meshToMeshInterpolation / meshToMesh / meshToMesh.C
blobc5f95879501941ce7bfb29ccd481ef35cb7c19d8
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "meshToMesh.H"
27 #include "processorFvPatch.H"
28 #include "demandDrivenData.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 defineTypeNameAndDebug(Foam::meshToMesh, 0);
34 const Foam::scalar Foam::meshToMesh::directHitTol = 1e-5;
36 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
38 Foam::meshToMesh::meshToMesh
40     const fvMesh& meshFrom,
41     const fvMesh& meshTo,
42     const HashTable<word>& patchMap,
43     const wordList& cuttingPatchNames
46     fromMesh_(meshFrom),
47     toMesh_(meshTo),
48     patchMap_(patchMap),
49     cellAddressing_(toMesh_.nCells()),
50     boundaryAddressing_(toMesh_.boundaryMesh().size()),
51     inverseDistanceWeightsPtr_(NULL)
53     forAll(fromMesh_.boundaryMesh(), patchi)
54     {
55         fromMeshPatches_.insert
56         (
57             fromMesh_.boundaryMesh()[patchi].name(),
58             patchi
59         );
60     }
62     forAll(toMesh_.boundaryMesh(), patchi)
63     {
64         toMeshPatches_.insert
65         (
66             toMesh_.boundaryMesh()[patchi].name(),
67             patchi
68         );
69     }
71     forAll(cuttingPatchNames, i)
72     {
73         if (toMeshPatches_.found(cuttingPatchNames[i]))
74         {
75             cuttingPatches_.insert
76             (
77                 cuttingPatchNames[i],
78                 toMeshPatches_.find(cuttingPatchNames[i])()
79             );
80         }
81         else
82         {
83             WarningIn
84             (
85                 "meshToMesh::meshToMesh"
86                 "(const fvMesh& meshFrom, const fvMesh& meshTo,"
87                 "const HashTable<word>& patchMap,"
88                 "const wordList& cuttingPatchNames)"
89             )   << "Cannot find cutting-patch " << cuttingPatchNames[i]
90                 << " in destination mesh" << endl;
91         }
92     }
94     forAll(toMesh_.boundaryMesh(), patchi)
95     {
96         // Add the processor patches in the toMesh to the cuttingPatches list
97         if (isA<processorPolyPatch>(toMesh_.boundaryMesh()[patchi]))
98         {
99             cuttingPatches_.insert
100             (
101                 toMesh_.boundaryMesh()[patchi].name(),
102                 patchi
103             );
104         }
105     }
107     calcAddressing();
111 Foam::meshToMesh::meshToMesh
113     const fvMesh& meshFrom,
114     const fvMesh& meshTo
117     fromMesh_(meshFrom),
118     toMesh_(meshTo),
119     cellAddressing_(toMesh_.nCells()),
120     boundaryAddressing_(toMesh_.boundaryMesh().size()),
121     inverseDistanceWeightsPtr_(NULL)
123     // check whether both meshes have got the same number
124     // of boundary patches
125     if (fromMesh_.boundary().size() != toMesh_.boundary().size())
126     {
127         FatalErrorIn
128         (
129             "meshToMesh::meshToMesh"
130             "(const fvMesh& meshFrom, const fvMesh& meshTo)"
131         )   << "Incompatible meshes: different number of patches, "
132             << "fromMesh = " << fromMesh_.boundary().size()
133             << ", toMesh = " << toMesh_.boundary().size()
134             << exit(FatalError);
135     }
137     forAll(fromMesh_.boundaryMesh(), patchi)
138     {
139         if
140         (
141             fromMesh_.boundaryMesh()[patchi].name()
142          != toMesh_.boundaryMesh()[patchi].name()
143         )
144         {
145             FatalErrorIn
146             (
147                 "meshToMesh::meshToMesh"
148                 "(const fvMesh& meshFrom, const fvMesh& meshTo)"
149             )   << "Incompatible meshes: different patch names for patch "
150                 << patchi
151                 << ", fromMesh = " << fromMesh_.boundary()[patchi].name()
152                 << ", toMesh = " << toMesh_.boundary()[patchi].name()
153                 << exit(FatalError);
154         }
156         if
157         (
158             fromMesh_.boundaryMesh()[patchi].type()
159          != toMesh_.boundaryMesh()[patchi].type()
160         )
161         {
162             FatalErrorIn
163             (
164                 "meshToMesh::meshToMesh"
165                 "(const fvMesh& meshFrom, const fvMesh& meshTo)"
166             )   << "Incompatible meshes: different patch types for patch "
167                 << patchi
168                 << ", fromMesh = " << fromMesh_.boundary()[patchi].type()
169                 << ", toMesh = " << toMesh_.boundary()[patchi].type()
170                 << exit(FatalError);
171         }
173         fromMeshPatches_.insert
174         (
175             fromMesh_.boundaryMesh()[patchi].name(),
176             patchi
177         );
179         toMeshPatches_.insert
180         (
181             toMesh_.boundaryMesh()[patchi].name(),
182             patchi
183         );
185         patchMap_.insert
186         (
187             toMesh_.boundaryMesh()[patchi].name(),
188             fromMesh_.boundaryMesh()[patchi].name()
189         );
190     }
192     calcAddressing();
196 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
198 Foam::meshToMesh::~meshToMesh()
200     deleteDemandDrivenData(inverseDistanceWeightsPtr_);
204 // ************************************************************************* //