17 #include <deal.II/base/logstream.h> 18 #include <deal.II/base/function.h> 19 #include <deal.II/base/vectorization.h> 20 #include <deal.II/lac/la_parallel_vector.h> 21 #include <deal.II/grid/tria.h> 22 #include <deal.II/grid/tria_iterator.h> 23 #include <deal.II/dofs/dof_tools.h> 24 #include <deal.II/fe/fe.h> 25 #include <deal.II/dofs/dof_accessor.h> 26 #include <deal.II/multigrid/mg_tools.h> 27 #include <deal.II/multigrid/mg_transfer_matrix_free.h> 28 #include <deal.II/multigrid/mg_transfer_internal.h> 30 #include <deal.II/matrix_free/evaluation_kernels.h> 34 DEAL_II_NAMESPACE_OPEN
37 template <
int dim,
typename Number>
41 element_is_continuous(false),
48 template <
int dim,
typename Number>
52 element_is_continuous(false),
61 template <
int dim,
typename Number>
65 this->mg_constrained_dofs = &mg_c;
70 template <
int dim,
typename Number>
75 element_is_continuous =
false;
77 n_child_cell_dofs = 0;
78 level_dof_indices.clear();
79 parent_child_connect.clear();
80 dirichlet_indices.clear();
81 n_owned_level_cells.clear();
82 prolongation_matrix_1d.clear();
83 evaluation_data.clear();
84 weights_on_refined.clear();
88 template <
int dim,
typename Number>
92 this->fill_and_communicate_copy_indices(mg_dof);
94 std::vector<std::vector<Number> > weights_unvectorized;
98 internal::MGTransfer::setup_transfer<dim,Number>(mg_dof,
99 this->mg_constrained_dofs,
102 parent_child_connect,
105 weights_unvectorized,
106 this->copy_indices_global_mine,
107 this->ghosted_level_vector);
110 element_is_continuous = elem_info.element_is_continuous;
112 n_child_cell_dofs = elem_info.n_child_cell_dofs;
115 prolongation_matrix_1d.resize(elem_info.prolongation_matrix_1d.size());
116 for (
unsigned int i=0; i<elem_info.prolongation_matrix_1d.size(); i++)
117 prolongation_matrix_1d[i] = elem_info.prolongation_matrix_1d[i];
121 const unsigned int n_levels = mg_dof.get_triangulation().n_global_levels();
123 const unsigned int n_weights_per_cell = Utilities::fixed_power<dim>(3);
124 weights_on_refined.resize(n_levels-1);
125 for (
unsigned int level = 1; level<n_levels; ++level)
127 weights_on_refined[level-1].resize(((n_owned_level_cells[level-1]+vec_size-1)/vec_size)*n_weights_per_cell);
129 for (
unsigned int c=0; c<n_owned_level_cells[level-1]; ++c)
131 const unsigned int comp = c/vec_size;
132 const unsigned int v = c%vec_size;
133 for (
unsigned int i = 0; i<n_weights_per_cell; ++i)
136 weights_on_refined[level-1][comp*n_weights_per_cell+i][v] = weights_unvectorized[level-1][c*n_weights_per_cell+i];
141 evaluation_data.resize(n_child_cell_dofs);
146 template <
int dim,
typename Number>
152 Assert ((to_level >= 1) && (to_level<=level_dof_indices.size()),
160 this->ghosted_level_vector[to_level-1].copy_locally_owned_data_from(src);
161 this->ghosted_level_vector[to_level-1].update_ghost_values();
162 this->ghosted_level_vector[to_level] = 0.;
168 do_prolongate_add<0>(to_level, this->ghosted_level_vector[to_level],
169 this->ghosted_level_vector[to_level-1]);
170 else if (fe_degree == 1)
171 do_prolongate_add<1>(to_level, this->ghosted_level_vector[to_level],
172 this->ghosted_level_vector[to_level-1]);
173 else if (fe_degree == 2)
174 do_prolongate_add<2>(to_level, this->ghosted_level_vector[to_level],
175 this->ghosted_level_vector[to_level-1]);
176 else if (fe_degree == 3)
177 do_prolongate_add<3>(to_level, this->ghosted_level_vector[to_level],
178 this->ghosted_level_vector[to_level-1]);
179 else if (fe_degree == 4)
180 do_prolongate_add<4>(to_level, this->ghosted_level_vector[to_level],
181 this->ghosted_level_vector[to_level-1]);
182 else if (fe_degree == 5)
183 do_prolongate_add<5>(to_level, this->ghosted_level_vector[to_level],
184 this->ghosted_level_vector[to_level-1]);
185 else if (fe_degree == 6)
186 do_prolongate_add<6>(to_level, this->ghosted_level_vector[to_level],
187 this->ghosted_level_vector[to_level-1]);
188 else if (fe_degree == 7)
189 do_prolongate_add<7>(to_level, this->ghosted_level_vector[to_level],
190 this->ghosted_level_vector[to_level-1]);
191 else if (fe_degree == 8)
192 do_prolongate_add<8>(to_level, this->ghosted_level_vector[to_level],
193 this->ghosted_level_vector[to_level-1]);
194 else if (fe_degree == 9)
195 do_prolongate_add<9>(to_level, this->ghosted_level_vector[to_level],
196 this->ghosted_level_vector[to_level-1]);
197 else if (fe_degree == 10)
198 do_prolongate_add<10>(to_level, this->ghosted_level_vector[to_level],
199 this->ghosted_level_vector[to_level-1]);
201 do_prolongate_add<-1>(to_level, this->ghosted_level_vector[to_level],
202 this->ghosted_level_vector[to_level-1]);
210 template <
int dim,
typename Number>
216 Assert ((from_level >= 1) && (from_level<=level_dof_indices.size()),
224 this->ghosted_level_vector[from_level].copy_locally_owned_data_from(src);
225 this->ghosted_level_vector[from_level].update_ghost_values();
226 this->ghosted_level_vector[from_level-1] = 0.;
229 do_restrict_add<0>(from_level, this->ghosted_level_vector[from_level-1],
230 this->ghosted_level_vector[from_level]);
231 else if (fe_degree == 1)
232 do_restrict_add<1>(from_level, this->ghosted_level_vector[from_level-1],
233 this->ghosted_level_vector[from_level]);
234 else if (fe_degree == 2)
235 do_restrict_add<2>(from_level, this->ghosted_level_vector[from_level-1],
236 this->ghosted_level_vector[from_level]);
237 else if (fe_degree == 3)
238 do_restrict_add<3>(from_level, this->ghosted_level_vector[from_level-1],
239 this->ghosted_level_vector[from_level]);
240 else if (fe_degree == 4)
241 do_restrict_add<4>(from_level, this->ghosted_level_vector[from_level-1],
242 this->ghosted_level_vector[from_level]);
243 else if (fe_degree == 5)
244 do_restrict_add<5>(from_level, this->ghosted_level_vector[from_level-1],
245 this->ghosted_level_vector[from_level]);
246 else if (fe_degree == 6)
247 do_restrict_add<6>(from_level, this->ghosted_level_vector[from_level-1],
248 this->ghosted_level_vector[from_level]);
249 else if (fe_degree == 7)
250 do_restrict_add<7>(from_level, this->ghosted_level_vector[from_level-1],
251 this->ghosted_level_vector[from_level]);
252 else if (fe_degree == 8)
253 do_restrict_add<8>(from_level, this->ghosted_level_vector[from_level-1],
254 this->ghosted_level_vector[from_level]);
255 else if (fe_degree == 9)
256 do_restrict_add<9>(from_level, this->ghosted_level_vector[from_level-1],
257 this->ghosted_level_vector[from_level]);
258 else if (fe_degree == 10)
259 do_restrict_add<10>(from_level, this->ghosted_level_vector[from_level-1],
260 this->ghosted_level_vector[from_level]);
263 do_restrict_add<-1>(from_level, this->ghosted_level_vector[from_level-1],
264 this->ghosted_level_vector[from_level]);
267 dst += this->ghosted_level_vector[from_level-1];
274 template <
int dim,
int degree,
typename Number>
277 const unsigned int fe_degree,
282 const int loop_length = degree != -1 ? 2*degree+1 : 2*fe_degree+1;
283 unsigned int degree_to_3 [100];
285 for (
int i=1; i<loop_length-1; ++i)
287 degree_to_3[loop_length-1] = 2;
289 for (
int k=0; k<(dim>2 ? loop_length : 1); ++k)
290 for (
int j=0; j<(dim>1 ? loop_length : 1); ++j)
292 const unsigned int shift = 9*degree_to_3[k] + 3*degree_to_3[j];
293 data[0] *= weights[shift];
296 for (
int i=1; i<loop_length-1; ++i)
297 data[i] *= weights[shift+1];
298 data[loop_length-1] *= weights[shift+2];
306 template <
int dim,
typename Number>
307 template <
int degree>
314 const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
315 const unsigned int n_child_dofs_1d = 2*degree_size - element_is_continuous;
316 const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
319 for (
unsigned int cell=0; cell < n_owned_level_cells[to_level-1];
322 const unsigned int n_lanes = cell+vec_size > n_owned_level_cells[to_level-1] ?
323 n_owned_level_cells[to_level-1] - cell : vec_size;
326 for (
unsigned int v=0; v<n_lanes; ++v)
328 const unsigned int shift = internal::MGTransfer::compute_shift_within_children<dim>
329 (parent_child_connect[to_level-1][cell+v].second,
330 fe_degree+1-element_is_continuous, fe_degree);
331 const unsigned int *indices = &level_dof_indices[to_level-1][parent_child_connect[to_level-1][cell+v].first*n_child_cell_dofs+shift];
334 for (
unsigned int k=0; k<(dim>2 ? degree_size : 1); ++k)
335 for (
unsigned int j=0; j<(dim>1 ? degree_size : 1); ++j)
336 for (
unsigned int i=0; i<degree_size; ++i, ++m)
337 evaluation_data[m][v] =
339 k*n_child_dofs_1d*n_child_dofs_1d+
340 j*n_child_dofs_1d+i]);
343 for (std::vector<unsigned short>::const_iterator i=dirichlet_indices[to_level-1][cell+v].begin(); i!=dirichlet_indices[to_level-1][cell+v].end(); ++i)
344 evaluation_data[*i][v] = 0.;
349 degree_size * n_child_dofs_1d);
351 if (element_is_continuous)
358 ::do_forward(prolongation_matrix_1d,
359 evaluation_data.begin()+ c*Utilities::fixed_power<dim>(degree_size),
360 evaluation_data.begin()+ c*n_scalar_cell_dofs,
361 fe_degree+1, 2*fe_degree+1);
362 weight_dofs_on_child<dim,degree,Number>(&weights_on_refined[to_level-1][(cell/vec_size)*three_to_dim],
364 evaluation_data.begin());
371 ::do_forward(prolongation_matrix_1d,
372 evaluation_data.begin() + c*Utilities::fixed_power<dim>(degree_size),
373 evaluation_data.begin() + c*n_scalar_cell_dofs,
374 fe_degree+1, 2*fe_degree+2);
378 const unsigned int *indices = &level_dof_indices[to_level][cell*
380 for (
unsigned int v=0; v<n_lanes; ++v)
382 for (
unsigned int i=0; i<n_child_cell_dofs; ++i)
384 indices += n_child_cell_dofs;
391 template <
int dim,
typename Number>
392 template <
int degree>
399 const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
400 const unsigned int n_child_dofs_1d = 2*degree_size - element_is_continuous;
401 const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
404 for (
unsigned int cell=0; cell < n_owned_level_cells[from_level-1];
407 const unsigned int n_lanes = cell+vec_size > n_owned_level_cells[from_level-1] ?
408 n_owned_level_cells[from_level-1] - cell : vec_size;
412 const unsigned int *indices = &level_dof_indices[from_level][cell*
414 for (
unsigned int v=0; v<n_lanes; ++v)
416 for (
unsigned int i=0; i<n_child_cell_dofs; ++i)
418 indices += n_child_cell_dofs;
423 degree_size * n_child_dofs_1d);
425 if (element_is_continuous)
427 weight_dofs_on_child<dim,degree,Number>(&weights_on_refined[from_level-1][(cell/vec_size)*three_to_dim],
429 &evaluation_data[0]);
433 ::do_backward(prolongation_matrix_1d,
false,
434 evaluation_data.begin() + c*n_scalar_cell_dofs,
435 evaluation_data.begin() + c*Utilities::fixed_power<dim>(degree_size),
436 fe_degree+1, 2*fe_degree+1);
443 ::do_backward(prolongation_matrix_1d,
false,
444 evaluation_data.begin() + c*n_scalar_cell_dofs,
445 evaluation_data.begin() + c*Utilities::fixed_power<dim>(degree_size),
446 fe_degree+1, 2*fe_degree+2);
450 for (
unsigned int v=0; v<n_lanes; ++v)
452 const unsigned int shift = internal::MGTransfer::compute_shift_within_children<dim>
453 (parent_child_connect[from_level-1][cell+v].second,
454 fe_degree+1-element_is_continuous, fe_degree);
456 n_child_cell_dofs+n_child_cell_dofs-1,
457 level_dof_indices[from_level-1].size());
458 const unsigned int *indices = &level_dof_indices[from_level-1][parent_child_connect[from_level-1][cell+v].first*n_child_cell_dofs+shift];
462 for (std::vector<unsigned short>::const_iterator i=dirichlet_indices[from_level-1][cell+v].begin(); i!=dirichlet_indices[from_level-1][cell+v].end(); ++i)
463 evaluation_data[*i][v] = 0.;
465 for (
unsigned int k=0; k<(dim>2 ? degree_size : 1); ++k)
466 for (
unsigned int j=0; j<(dim>1 ? degree_size : 1); ++j)
467 for (
unsigned int i=0; i<degree_size; ++i, ++m)
469 k*n_child_dofs_1d*n_child_dofs_1d+
470 j*n_child_dofs_1d+i])
471 += evaluation_data[m][v];
479 template <
int dim,
typename Number>
496 template <
int dim,
typename Number>
498 : same_for_all (true)
505 template <
int dim,
typename Number>
507 (
const std::vector<MGConstrainedDoFs> &mg_c)
508 : same_for_all (
false)
510 for (
unsigned int i=0; i<mg_c.size(); ++i)
511 matrix_free_transfer_vector.emplace_back(mg_c[i]);
516 template <
int dim,
typename Number>
521 ExcMessage(
"This object was initialized with support for usage with " 522 "one DoFHandler for each block, but this method assumes " 523 "that the same DoFHandler is used for all the blocks!"));
526 matrix_free_transfer_vector[0].initialize_constraints(mg_c);
531 template <
int dim,
typename Number>
533 (
const std::vector<MGConstrainedDoFs> &mg_c)
536 ExcMessage(
"This object was initialized with support for using " 537 "the same DoFHandler for all the blocks, but this " 538 "method assumes that there is a separate DoFHandler " 542 for (
unsigned int i = 0; i<mg_c.size(); ++i)
543 matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
548 template <
int dim,
typename Number>
551 matrix_free_transfer_vector.clear();
556 template <
int dim,
typename Number>
561 matrix_free_transfer_vector[0].build(mg_dof);
566 template <
int dim,
typename Number>
571 for (
unsigned int i=0; i<mg_dof.size(); ++i)
572 matrix_free_transfer_vector[i].build(*mg_dof[i]);
577 template <
int dim,
typename Number>
583 const unsigned int n_blocks = src.
n_blocks();
589 for (
unsigned int b = 0; b < n_blocks; ++b)
591 const unsigned int data_block = same_for_all ? 0 : b;
592 matrix_free_transfer_vector[data_block].prolongate(to_level, dst.
block(b), src.
block(b));
598 template <
int dim,
typename Number>
604 const unsigned int n_blocks = src.
n_blocks();
610 for (
unsigned int b = 0; b < n_blocks; ++b)
612 const unsigned int data_block = same_for_all ? 0 : b;
613 matrix_free_transfer_vector[data_block].restrict_and_add(from_level, dst.
block(b), src.
block(b));
619 template <
int dim,
typename Number>
623 std::size_t total_memory_consumption= 0;
624 for (
const auto &el: matrix_free_transfer_vector)
625 total_memory_consumption += el.memory_consumption();
626 return total_memory_consumption;
632 #include "mg_transfer_matrix_free.inst" 635 DEAL_II_NAMESPACE_CLOSE
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
#define AssertDimension(dim1, dim2)
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void copy_locally_owned_data_from(const Vector< Number2 > &src)
std::size_t memory_consumption() const
void build(const DoFHandler< dim, dim > &mg_dof)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::size_t memory_consumption() const
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
size_type local_size() const
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
unsigned int n_blocks() const
MGTransferBlockMatrixFree()=default
Number local_element(const size_type local_index) const
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const
static ::ExceptionBase & ExcNotImplemented()
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
BlockType & block(const unsigned int i)
void build(const DoFHandler< dim, dim > &mg_dof)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)