1+ // Copyright (C) 2024 Paul Kühner
2+ //
3+ // This file is part of DOLFINX (https://www.fenicsproject.org)
4+ //
5+ // SPDX-License-Identifier: LGPL-3.0-or-later
6+
7+ #pragma once
8+
9+ #include < algorithm>
10+ #include < concepts>
11+ #include < cstdint>
12+ #include < iterator>
13+ #include < vector>
14+
15+ #include " dolfinx/common/IndexMap.h"
16+ #include " dolfinx/fem/FunctionSpace.h"
17+ #include " dolfinx/la/SparsityPattern.h"
18+ #include " dolfinx/la/utils.h"
19+ #include " dolfinx/mesh/Mesh.h"
20+
21+ namespace dolfinx ::transfer
22+ {
23+
24+ template <std::floating_point T>
25+ std::vector<std::int64_t >
26+ inclusion_mapping (const dolfinx::mesh::Mesh<T>& mesh_from,
27+ const dolfinx::mesh::Mesh<T>& mesh_to)
28+ {
29+
30+ const common::IndexMap& im_from = *mesh_from.topology ()->index_map (0 );
31+ const common::IndexMap& im_to = *mesh_to.topology ()->index_map (0 );
32+
33+ std::vector<std::int64_t > map (im_from.size_global (), -1 );
34+
35+ std::span<const T> x_from = mesh_from.geometry ().x ();
36+ std::span<const T> x_to = mesh_to.geometry ().x ();
37+
38+ auto build_global_to_local = [&](const auto & im)
39+ {
40+ return [&](std::int32_t idx)
41+ {
42+ std::array<std::int64_t , 1 > tmp;
43+ im.local_to_global (std::vector<std::int32_t >{idx}, tmp);
44+ return tmp[0 ];
45+ };
46+ };
47+
48+ auto to_global_to = build_global_to_local (im_to);
49+ auto to_global_from = build_global_to_local (im_from);
50+
51+ for (std::int32_t i = 0 ; i < im_from.size_local (); i++)
52+ {
53+ std::ranges::subrange vertex_from (std::next (x_from.begin (), 3 * i),
54+ std::next (x_from.begin (), 3 * (i + 1 )));
55+ for (std::int64_t j = 0 ; j < im_to.size_local () + im_to.num_ghosts (); j++)
56+ {
57+ std::ranges::subrange vertex_to (std::next (x_to.begin (), 3 * j),
58+ std::next (x_to.begin (), 3 * (j + 1 )));
59+
60+ if (std::ranges::equal (
61+ vertex_from, vertex_to, [](T a, T b)
62+ { return std::abs (a - b) <= std::numeric_limits<T>::epsilon (); }))
63+ {
64+ map[to_global_from (i)] = to_global_to (j);
65+ break ;
66+ }
67+ }
68+ }
69+
70+ if (dolfinx::MPI::size (mesh_to.comm ()) == 1 )
71+ {
72+ // no communication required
73+ assert (std::ranges::all_of (map, [](auto e) { return e >= 0 ; }));
74+ return map;
75+ }
76+
77+ // map holds at this point for every original local index the corresponding
78+ // mapped global index. All other entries are still -1, but are available on
79+ // other processes.
80+ std::vector<std::int64_t > result (map.size (), -1 );
81+ MPI_Allreduce (map.data (), result.data (), map.size (),
82+ dolfinx::MPI::mpi_type<std::int64_t >(), MPI_MAX,
83+ mesh_from.comm ());
84+
85+ assert (std::ranges::all_of (result, [](auto e) { return e >= 0 ; }));
86+ return result;
87+ }
88+
89+ template <std::floating_point T>
90+ la::SparsityPattern
91+ create_sparsity_pattern (const dolfinx::fem::FunctionSpace<T>& V_from,
92+ const dolfinx::fem::FunctionSpace<T>& V_to,
93+ const std::vector<std::int64_t >& inclusion_map)
94+ {
95+ // TODO: P1 and which elements supported?
96+ auto mesh_from = V_from.mesh ();
97+ auto mesh_to = V_to.mesh ();
98+ assert (mesh_from);
99+ assert (mesh_to);
100+
101+ MPI_Comm comm = mesh_from->comm ();
102+ {
103+ // Check comms equal
104+ int result;
105+ MPI_Comm_compare (comm, mesh_to->comm (), &result);
106+ assert (result == MPI_CONGRUENT);
107+ }
108+ assert (mesh_from->topology ()->dim () == mesh_to->topology ()->dim ());
109+
110+ auto to_v_to_f = mesh_to->topology ()->connectivity (0 , 1 );
111+ auto to_f_to_v = mesh_to->topology ()->connectivity (1 , 0 );
112+ assert (to_v_to_f);
113+ assert (to_f_to_v);
114+
115+ auto dofmap_from = V_from.dofmap ();
116+ auto dofmap_to = V_to.dofmap ();
117+ assert (dofmap_from);
118+ assert (dofmap_to);
119+
120+ assert (mesh_to->topology ()->index_map (0 ));
121+ assert (mesh_from->topology ()->index_map (0 ));
122+ const common::IndexMap& im_to = *mesh_to->topology ()->index_map (0 );
123+ const common::IndexMap& im_from = *mesh_from->topology ()->index_map (0 );
124+
125+ dolfinx::la::SparsityPattern sp (
126+ comm, {dofmap_from->index_map , dofmap_to->index_map },
127+ {dofmap_from->index_map_bs (), dofmap_to->index_map_bs ()});
128+
129+ assert (inclusion_map.size () == im_from.size_global ());
130+ for (int dof_from_global = 0 ; dof_from_global < im_from.size_global ();
131+ dof_from_global++)
132+ {
133+ std::int64_t dof_to_global = inclusion_map[dof_from_global];
134+
135+ std::vector<std::int32_t > local_dof_to_v{0 };
136+ im_to.global_to_local (std::vector<std::int64_t >{dof_to_global},
137+ local_dof_to_v);
138+
139+ auto local_dof_to = local_dof_to_v[0 ];
140+
141+ bool is_remote = (local_dof_to == -1 );
142+ bool is_ghost = local_dof_to >= im_to.size_local ();
143+ if (is_remote || is_ghost)
144+ continue ;
145+
146+ std::vector<std::int32_t > dof_from_v{0 };
147+ im_from.global_to_local (std::vector<std::int64_t >{dof_from_global},
148+ dof_from_v);
149+
150+ std::ranges::for_each (
151+ to_v_to_f->links (local_dof_to),
152+ [&](auto e)
153+ {
154+ sp.insert (
155+ std::vector<int32_t >(to_f_to_v->links (e).size (), dof_from_v[0 ]),
156+ to_f_to_v->links (e));
157+ });
158+ }
159+ sp.finalize ();
160+ return sp;
161+ }
162+
163+ template <std::floating_point T, typename F>
164+ requires std::invocable<F&, std::int32_t >
165+ void assemble_transfer_matrix (la::MatSet<T> auto mat_set,
166+ const dolfinx::fem::FunctionSpace<T>& V_from,
167+ const dolfinx::fem::FunctionSpace<T>& V_to,
168+ const std::vector<std::int64_t >& inclusion_map,
169+ F weight)
170+ {
171+ // TODO: P1 and which elements supported?
172+ auto mesh_from = V_from.mesh ();
173+ auto mesh_to = V_to.mesh ();
174+ assert (mesh_from);
175+ assert (mesh_to);
176+
177+ MPI_Comm comm = mesh_from->comm ();
178+ {
179+ // Check comms equal
180+ int result;
181+ MPI_Comm_compare (comm, mesh_to->comm (), &result);
182+ assert (result == MPI_CONGRUENT);
183+ }
184+ assert (mesh_from->topology ()->dim () == mesh_to->topology ()->dim ());
185+
186+ auto to_v_to_f = mesh_to->topology ()->connectivity (0 , 1 );
187+ auto to_f_to_v = mesh_to->topology ()->connectivity (1 , 0 );
188+ assert (to_v_to_f);
189+ assert (to_f_to_v);
190+
191+ auto dofmap_from = V_from.dofmap ();
192+ auto dofmap_to = V_to.dofmap ();
193+ assert (dofmap_from);
194+ assert (dofmap_to);
195+
196+ assert (mesh_to->topology ()->index_map (0 ));
197+ assert (mesh_from->topology ()->index_map (0 ));
198+ const common::IndexMap& im_to = *mesh_to->topology ()->index_map (0 );
199+ const common::IndexMap& im_from = *mesh_from->topology ()->index_map (0 );
200+
201+ assert (inclusion_map.size () == im_from.size_global ());
202+
203+ for (int dof_from_global = 0 ; dof_from_global < im_from.size_global ();
204+ dof_from_global++)
205+ {
206+ std::int64_t dof_to_global = inclusion_map[dof_from_global];
207+
208+ std::vector<std::int32_t > local_dof_to_v{0 };
209+ im_to.global_to_local (std::vector<std::int64_t >{dof_to_global},
210+ local_dof_to_v);
211+
212+ auto local_dof_to = local_dof_to_v[0 ];
213+
214+ bool is_remote = (local_dof_to == -1 );
215+ bool is_ghost = local_dof_to >= im_to.size_local ();
216+ if (is_remote || is_ghost)
217+ continue ;
218+
219+ std::vector<std::int32_t > dof_from_v{0 };
220+ im_from.global_to_local (std::vector<std::int64_t >{dof_from_global},
221+ dof_from_v);
222+
223+ for (auto e : to_v_to_f->links (local_dof_to))
224+ {
225+ for (auto n : to_f_to_v->links (e))
226+ {
227+ // For now we only support distance <= 1 -> this should be somehow
228+ // asserted.
229+ std::int32_t distance = (n == local_dof_to) ? 0 : 1 ;
230+ mat_set (std::vector<int32_t >{dof_from_v[0 ]}, std::vector<int32_t >{n},
231+ std::vector<double >{weight (distance)});
232+ }
233+ }
234+ }
235+ }
236+
237+ } // namespace dolfinx::transfer
0 commit comments