1
+ // Copyright (C) 2024 Paul T. 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>
164
+ void assemble_transfer_matrix (la::MatSet<T> auto mat_set,
165
+ const dolfinx::fem::FunctionSpace<T>& V_from,
166
+ const dolfinx::fem::FunctionSpace<T>& V_to,
167
+ const std::vector<std::int64_t >& inclusion_map,
168
+ std::function<T(std::int32_t )> weight)
169
+ {
170
+ // TODO: P1 and which elements supported?
171
+ auto mesh_from = V_from.mesh ();
172
+ auto mesh_to = V_to.mesh ();
173
+ assert (mesh_from);
174
+ assert (mesh_to);
175
+
176
+ MPI_Comm comm = mesh_from->comm ();
177
+ {
178
+ // Check comms equal
179
+ int result;
180
+ MPI_Comm_compare (comm, mesh_to->comm (), &result);
181
+ assert (result == MPI_CONGRUENT);
182
+ }
183
+ assert (mesh_from->topology ()->dim () == mesh_to->topology ()->dim ());
184
+
185
+ auto to_v_to_f = mesh_to->topology ()->connectivity (0 , 1 );
186
+ auto to_f_to_v = mesh_to->topology ()->connectivity (1 , 0 );
187
+ assert (to_v_to_f);
188
+ assert (to_f_to_v);
189
+
190
+ auto dofmap_from = V_from.dofmap ();
191
+ auto dofmap_to = V_to.dofmap ();
192
+ assert (dofmap_from);
193
+ assert (dofmap_to);
194
+
195
+ assert (mesh_to->topology ()->index_map (0 ));
196
+ assert (mesh_from->topology ()->index_map (0 ));
197
+ const common::IndexMap& im_to = *mesh_to->topology ()->index_map (0 );
198
+ const common::IndexMap& im_from = *mesh_from->topology ()->index_map (0 );
199
+
200
+ assert (inclusion_map.size () == im_from.size_global ());
201
+
202
+ for (int dof_from_global = 0 ; dof_from_global < im_from.size_global ();
203
+ dof_from_global++)
204
+ {
205
+ std::int64_t dof_to_global = inclusion_map[dof_from_global];
206
+
207
+ std::vector<std::int32_t > local_dof_to_v{0 };
208
+ im_to.global_to_local (std::vector<std::int64_t >{dof_to_global},
209
+ local_dof_to_v);
210
+
211
+ auto local_dof_to = local_dof_to_v[0 ];
212
+
213
+ bool is_remote = (local_dof_to == -1 );
214
+ bool is_ghost = local_dof_to >= im_to.size_local ();
215
+ if (is_remote || is_ghost)
216
+ continue ;
217
+
218
+ std::vector<std::int32_t > dof_from_v{0 };
219
+ im_from.global_to_local (std::vector<std::int64_t >{dof_from_global},
220
+ dof_from_v);
221
+
222
+ for (auto e : to_v_to_f->links (local_dof_to))
223
+ {
224
+ for (auto n : to_f_to_v->links (e))
225
+ {
226
+ // For now we only support distance <= 1 -> this should be somehow
227
+ // asserted.
228
+ std::int32_t distance = (n == local_dof_to) ? 0 : 1 ;
229
+ mat_set (std::vector<int32_t >{dof_from_v[0 ]}, std::vector<int32_t >{n},
230
+ std::vector<double >{weight (distance)});
231
+ }
232
+ }
233
+ }
234
+ }
235
+
236
+ } // namespace dolfinx::transfer
0 commit comments