17 #include <deal.II/base/logstream.h> 18 #include <deal.II/base/function.h> 20 #include <deal.II/lac/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/fe/fe_tools.h> 26 #include <deal.II/dofs/dof_accessor.h> 27 #include <deal.II/multigrid/mg_tools.h> 28 #include <deal.II/multigrid/mg_transfer_matrix_free.h> 30 #include <deal.II/matrix_free/shape_info.h> 31 #include <deal.II/matrix_free/fe_evaluation.h> 35 DEAL_II_NAMESPACE_OPEN
38 template<
int dim,
typename Number>
42 element_is_continuous(false),
49 template<
int dim,
typename Number>
62 template <
int dim,
typename Number>
68 template <
int dim,
typename Number>
77 template <
int dim,
typename Number>
101 compute_shift_within_children(
const unsigned int child,
102 const unsigned int fe_shift_1d,
107 unsigned int c_tensor_index[dim];
108 unsigned int tmp = child;
109 for (
unsigned int d=0; d<dim; ++d)
111 c_tensor_index[d] = tmp % 2;
114 const unsigned int n_child_dofs_1d = fe_degree + 1 + fe_shift_1d;
115 unsigned int factor = 1;
116 unsigned int shift = fe_shift_1d * c_tensor_index[0];
117 for (
unsigned int d=1; d<dim; ++d)
119 factor *= n_child_dofs_1d;
120 shift = shift + factor * fe_shift_1d * c_tensor_index[d];
130 void add_child_indices(
const unsigned int child,
131 const unsigned int fe_shift_1d,
132 const unsigned int fe_degree,
133 const std::vector<unsigned int> &lexicographic_numbering,
134 const std::vector<types::global_dof_index> &local_dof_indices,
137 const unsigned int n_child_dofs_1d = fe_degree + 1 + fe_shift_1d;
138 const unsigned int shift =
139 compute_shift_within_children<dim>(child, fe_shift_1d,
fe_degree);
141 local_dof_indices.size()/Utilities::fixed_power<dim>(fe_degree+1);
143 const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
145 for (
unsigned int k=0; k<(dim>2 ? (fe_degree+1) : 1); ++k)
146 for (
unsigned int j=0; j<(dim>1 ? (fe_degree+1) : 1); ++j)
147 for (
unsigned int i=0; i<(fe_degree+1); ++i, ++m)
149 const unsigned int index = c*n_scalar_cell_dofs+k*n_child_dofs_1d*
150 n_child_dofs_1d+j*n_child_dofs_1d+i;
152 indices[index] == local_dof_indices[lexicographic_numbering[m]],
154 indices[index] = local_dof_indices[lexicographic_numbering[m]];
162 template <
typename Number>
164 reinit_ghosted_vector(
const IndexSet &locally_owned,
165 std::vector<types::global_dof_index> &ghosted_level_dofs,
166 const MPI_Comm &communicator,
170 std::sort(ghosted_level_dofs.begin(), ghosted_level_dofs.end());
172 ghosted_dofs.
add_indices(ghosted_level_dofs.begin(),
173 std::unique(ghosted_level_dofs.begin(),
174 ghosted_level_dofs.end()));
175 ghosted_dofs.compress();
178 if (ghosted_level_vector.
size() == locally_owned.
size())
182 const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> part
184 ghosted_dofs.add_indices(part->ghost_indices());
190 ghosted_level_vector.
reinit(locally_owned, ghosted_dofs, communicator);
196 const std::vector<types::global_dof_index> &mine,
197 const std::vector<types::global_dof_index> &remote,
198 std::vector<unsigned int> &localized_indices)
200 localized_indices.resize(mine.size()+remote.size(),
202 for (
unsigned int i=0; i<mine.size(); ++i)
206 for (
unsigned int i=0; i<remote.size(); ++i)
214 template <
int dim,
typename Number>
230 std::string fe_name = mg_dof.
get_fe().base_element(0).get_name();
232 const std::size_t template_starts = fe_name.find_first_of(
'<');
233 Assert (fe_name[template_starts+1] == (dim==1?
'1':(dim==2?
'2':
'3')),
235 fe_name[template_starts+1] =
'1';
237 std_cxx11::shared_ptr<FiniteElement<1> > fe_1d
238 (FETools::get_fe_from_name<1>(fe_name));
246 mg_dof.
get_fe().dofs_per_cell);
270 ExcNotImplemented());
275 for (
unsigned int c=0; c<GeometryInfo<1>::max_children_per_cell; ++c)
276 for (
unsigned int j=0; j<basic_support_points.size(); ++j)
277 points_refined[shift*c+j][0] =
278 c*0.5 + 0.5 * basic_support_points[renumbering[j]][0];
280 n_child_dofs_1d = points_refined.size();
287 for (
unsigned int c=0; c<GeometryInfo<1>::max_children_per_cell; ++c)
292 < std::max(2.*(
double)std::numeric_limits<Number>::epsilon(),1e-12),
301 std::vector<std::vector<unsigned int> > coarse_level_indices(n_levels-1);
302 for (
unsigned int level=0; level<std::min(tria.
n_levels(),n_levels-1); ++level)
303 coarse_level_indices[level].resize(tria.
n_raw_cells(level),
305 std::vector<types::global_dof_index> local_dof_indices(mg_dof.
get_fe().dofs_per_cell);
314 for (
unsigned int level=n_levels-1; level > 0; --level)
317 std::vector<types::global_dof_index> global_level_dof_indices;
318 std::vector<types::global_dof_index> global_level_dof_indices_remote;
319 std::vector<types::global_dof_index> ghosted_level_dofs;
320 std::vector<types::global_dof_index> global_level_dof_indices_l0;
321 std::vector<types::global_dof_index> ghosted_level_dofs_l0;
325 cell != mg_dof.
end(level-1); ++cell)
328 if (!cell->has_children())
331 bool consider_cell =
false;
335 consider_cell =
true;
340 bool cell_is_remote = !consider_cell;
341 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
344 consider_cell =
true;
358 std::vector<types::global_dof_index> &next_indices =
359 cell_is_remote ? global_level_dof_indices_remote : global_level_dof_indices;
360 const std::size_t start_index = next_indices.size();
363 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
365 if (cell_is_remote && cell->child(c)->level_subdomain_id() !=
368 cell->child(c)->get_mg_dof_indices(local_dof_indices);
371 for (
unsigned int i=0; i<local_dof_indices.size(); ++i)
372 if (!owned_level_dofs.
is_element(local_dof_indices[i]))
373 ghosted_level_dofs.push_back(local_dof_indices[i]);
378 &next_indices[start_index]);
381 if (cell->child(c)->has_children() &&
386 const unsigned int child_index = coarse_level_indices[level][cell->child(c)->index()];
388 unsigned int parent_index =
counter;
398 std::make_pair(parent_index, c);
400 static_cast<unsigned short>(-1));
405 for (
unsigned int i=0; i<mg_dof.
get_fe().dofs_per_cell; ++i)
413 coarse_level_indices[level-1].size());
414 coarse_level_indices[level-1][cell->index()] = counter++;
421 if (level == 1 && !cell_is_remote)
423 cell->get_mg_dof_indices(local_dof_indices);
426 for (
unsigned int i=0; i<local_dof_indices.size(); ++i)
427 if (!owned_level_dofs_l0.
is_element(local_dof_indices[i]))
428 ghosted_level_dofs_l0.push_back(local_dof_indices[i]);
430 const std::size_t start_index = global_level_dof_indices_l0.size();
436 &global_level_dof_indices_l0[start_index]);
440 for (
unsigned int i=0; i<mg_dof.
get_fe().dofs_per_cell; ++i)
459 if (level < n_levels-1)
460 for (std::vector<std::pair<unsigned int,unsigned int> >::iterator
462 if (i->first >= tria.
n_cells(level))
464 i->first -= tria.
n_cells(level);
471 const MPI_Comm communicator =
475 ghosted_level_dofs, communicator,
479 copy_indices_to_mpi_local_numbers(*this->ghosted_level_vector[level].get_partitioner(),
480 global_level_dof_indices,
481 global_level_dof_indices_remote,
491 ghosted_level_dofs_l0, communicator,
492 this->ghosted_level_vector[0],
495 copy_indices_to_mpi_local_numbers(*this->ghosted_level_vector[0].get_partitioner(),
496 global_level_dof_indices_l0,
497 std::vector<types::global_dof_index>(),
507 for (
unsigned int level = 1; level<n_levels; ++level)
517 std::vector<unsigned int> degree_to_3 (n_child_dofs_1d);
519 for (
unsigned int i=1; i<n_child_dofs_1d-1; ++i)
521 degree_to_3.back() = 2;
525 weights_on_refined[level].resize(((n_owned_level_cells[level-1]+vec_size-1)/vec_size)*Utilities::fixed_power<dim>(3));
526 for (
unsigned int c=0; c<n_owned_level_cells[level-1]; ++c)
528 const unsigned int comp = c/vec_size;
529 const unsigned int v = c%vec_size;
531 for (
unsigned int k=0, m=0; k<(dim>2 ? n_child_dofs_1d : 1); ++k)
532 for (
unsigned int j=0; j<(dim>1 ? n_child_dofs_1d : 1); ++j)
534 unsigned int shift = 9*degree_to_3[k] + 3*degree_to_3[j];
535 for (
unsigned int i=0; i<n_child_dofs_1d; ++i, ++m)
536 weights_on_refined[level][comp*Utilities::fixed_power<dim>(3)+shift+degree_to_3[i]][v] = Number(1.)/
547 template <
int dim,
typename Number>
563 this->ghosted_level_vector[to_level] = 0.;
569 do_prolongate_add<0>(to_level, this->ghosted_level_vector[to_level],
570 this->ghosted_level_vector[to_level-1]);
572 do_prolongate_add<1>(to_level, this->ghosted_level_vector[to_level],
573 this->ghosted_level_vector[to_level-1]);
575 do_prolongate_add<2>(to_level, this->ghosted_level_vector[to_level],
576 this->ghosted_level_vector[to_level-1]);
578 do_prolongate_add<3>(to_level, this->ghosted_level_vector[to_level],
579 this->ghosted_level_vector[to_level-1]);
581 do_prolongate_add<4>(to_level, this->ghosted_level_vector[to_level],
582 this->ghosted_level_vector[to_level-1]);
584 do_prolongate_add<5>(to_level, this->ghosted_level_vector[to_level],
585 this->ghosted_level_vector[to_level-1]);
587 do_prolongate_add<6>(to_level, this->ghosted_level_vector[to_level],
588 this->ghosted_level_vector[to_level-1]);
590 do_prolongate_add<7>(to_level, this->ghosted_level_vector[to_level],
591 this->ghosted_level_vector[to_level-1]);
593 do_prolongate_add<8>(to_level, this->ghosted_level_vector[to_level],
594 this->ghosted_level_vector[to_level-1]);
596 do_prolongate_add<9>(to_level, this->ghosted_level_vector[to_level],
597 this->ghosted_level_vector[to_level-1]);
599 do_prolongate_add<10>(to_level, this->ghosted_level_vector[to_level],
600 this->ghosted_level_vector[to_level-1]);
602 AssertThrow(
false, ExcNotImplemented(
"Only degrees 0 up to 10 implemented."));
604 this->ghosted_level_vector[to_level].compress(VectorOperation::add);
605 dst = this->ghosted_level_vector[to_level];
610 template <
int dim,
typename Number>
626 this->ghosted_level_vector[from_level-1] = 0.;
629 do_restrict_add<0>(from_level, this->ghosted_level_vector[from_level-1],
630 this->ghosted_level_vector[from_level]);
632 do_restrict_add<1>(from_level, this->ghosted_level_vector[from_level-1],
633 this->ghosted_level_vector[from_level]);
635 do_restrict_add<2>(from_level, this->ghosted_level_vector[from_level-1],
636 this->ghosted_level_vector[from_level]);
638 do_restrict_add<3>(from_level, this->ghosted_level_vector[from_level-1],
639 this->ghosted_level_vector[from_level]);
641 do_restrict_add<4>(from_level, this->ghosted_level_vector[from_level-1],
642 this->ghosted_level_vector[from_level]);
644 do_restrict_add<5>(from_level, this->ghosted_level_vector[from_level-1],
645 this->ghosted_level_vector[from_level]);
647 do_restrict_add<6>(from_level, this->ghosted_level_vector[from_level-1],
648 this->ghosted_level_vector[from_level]);
650 do_restrict_add<7>(from_level, this->ghosted_level_vector[from_level-1],
651 this->ghosted_level_vector[from_level]);
653 do_restrict_add<8>(from_level, this->ghosted_level_vector[from_level-1],
654 this->ghosted_level_vector[from_level]);
656 do_restrict_add<9>(from_level, this->ghosted_level_vector[from_level-1],
657 this->ghosted_level_vector[from_level]);
659 do_restrict_add<10>(from_level, this->ghosted_level_vector[from_level-1],
660 this->ghosted_level_vector[from_level]);
662 AssertThrow(
false, ExcNotImplemented(
"Only degrees 0 up to 10 implemented."));
664 this->ghosted_level_vector[from_level-1].compress(VectorOperation::add);
665 dst += this->ghosted_level_vector[from_level-1];
672 template <
int dim,
typename Eval,
typename Number,
bool pro
longate>
674 perform_tensorized_op(
const Eval &evaluator,
689 evaluator.template values<0,prolongate,false>(t0, t2);
692 evaluator.template values<0,prolongate,false>(t0, t1);
693 evaluator.template values<1,prolongate,false>(t1, t2);
697 evaluator.template values<0,prolongate,false>(t0, t2);
698 evaluator.template values<1,prolongate,false>(t2, t1);
699 evaluator.template values<2,prolongate,false>(t1, t2);
702 Assert(
false, ExcNotImplemented());
705 t0 += Eval::dofs_per_cell;
706 t2 += Eval::n_q_points;
710 t0 += Eval::n_q_points;
711 t2 += Eval::dofs_per_cell;
716 template <
int dim,
int degree,
typename Number>
718 const unsigned int n_components,
721 Assert(degree > 0, ExcNotImplemented());
722 const int loop_length = 2*degree+1;
723 unsigned int degree_to_3 [loop_length];
725 for (
int i=1; i<loop_length-1; ++i)
727 degree_to_3[loop_length-1] = 2;
729 for (
int k=0; k<(dim>2 ? loop_length : 1); ++k)
730 for (
int j=0; j<(dim>1 ? loop_length : 1); ++j)
732 const unsigned int shift = 9*degree_to_3[k] + 3*degree_to_3[j];
733 data[0] *= weights[shift];
736 for (
int i=1; i<loop_length-1; ++i)
737 data[i] *= weights[shift+1];
738 data[loop_length-1] *= weights[shift+2];
746 template <
int dim,
typename Number>
747 template <
int degree>
755 const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
761 const unsigned int n_chunks = cell+vec_size > n_owned_level_cells[to_level-1] ?
762 n_owned_level_cells[to_level-1] - cell : vec_size;
765 for (
unsigned int v=0; v<n_chunks; ++v)
767 const unsigned int shift = compute_shift_within_children<dim>
773 for (
unsigned int k=0; k<(dim>2 ? (degree+1) : 1); ++k)
774 for (
unsigned int j=0; j<(dim>1 ? (degree+1) : 1); ++j)
775 for (
unsigned int i=0; i<(degree+1); ++i, ++m)
778 k*n_child_dofs_1d*n_child_dofs_1d+
779 j*n_child_dofs_1d+i]);
789 internal::MatrixFreeFunctions::tensor_symmetric, ExcNotImplemented());
793 (degree+1)*(degree+1));
794 typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+1,VectorizedArray<Number> > Evaluator;
795 Evaluator evaluator(
shape_info.shape_val_evenodd,
798 perform_tensorized_op<dim,Evaluator,Number,true>(evaluator,
802 weight_dofs_on_child<dim,degree,Number>(&
weights_on_refined[to_level][(cell/vec_size)*three_to_dim],
809 (degree+1)*(degree+1));
810 typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+2,VectorizedArray<Number> > Evaluator;
811 Evaluator evaluator(
shape_info.shape_val_evenodd,
814 perform_tensorized_op<dim,Evaluator,Number,true>(evaluator,
823 for (
unsigned int v=0; v<n_chunks; ++v)
834 template <
int dim,
typename Number>
835 template <
int degree>
843 const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
849 const unsigned int n_chunks = cell+vec_size > n_owned_level_cells[from_level-1] ?
850 n_owned_level_cells[from_level-1] - cell : vec_size;
856 for (
unsigned int v=0; v<n_chunks; ++v)
866 internal::MatrixFreeFunctions::tensor_symmetric, ExcNotImplemented());
870 (degree+1)*(degree+1));
871 typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+1,VectorizedArray<Number> > Evaluator;
872 Evaluator evaluator(
shape_info.shape_val_evenodd,
875 weight_dofs_on_child<dim,degree,Number>(&
weights_on_refined[from_level][(cell/vec_size)*three_to_dim],
878 perform_tensorized_op<dim,Evaluator,Number,false>(evaluator,
886 (degree+1)*(degree+1));
887 typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+2,VectorizedArray<Number> > Evaluator;
888 Evaluator evaluator(
shape_info.shape_val_evenodd,
891 perform_tensorized_op<dim,Evaluator,Number,false>(evaluator,
898 for (
unsigned int v=0; v<n_chunks; ++v)
900 const unsigned int shift = compute_shift_within_children<dim>
913 for (
unsigned int k=0; k<(dim>2 ? (degree+1) : 1); ++k)
914 for (
unsigned int j=0; j<(dim>1 ? (degree+1) : 1); ++j)
915 for (
unsigned int i=0; i<(degree+1); ++i, ++m)
917 k*n_child_dofs_1d*n_child_dofs_1d+
918 j*n_child_dofs_1d+i])
927 template <
int dim,
typename Number>
944 #include "mg_transfer_matrix_free.inst" 947 DEAL_II_NAMESPACE_CLOSE
unsigned int max_level() const
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
const types::subdomain_id invalid_subdomain_id
std::vector< unsigned int > n_owned_level_cells
cell_iterator begin(const unsigned int level=0) const
AlignedVector< VectorizedArray< Number > > evaluation_data
unsigned int n_raw_cells(const unsigned int level) const
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices_global_mine
std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > get_partitioner() const
unsigned int n_cells() const
unsigned int n_components
std::size_t memory_consumption() const
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
DEAL_VOLATILE unsigned int counter
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
cell_iterator end() const
#define AssertIndexRange(index, range)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
const unsigned int degree
void update_ghost_values() const
#define AssertThrow(cond, exc)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
unsigned int global_to_local(const types::global_dof_index global_index) const
const FiniteElement< dim, spacedim > & get_fe() const
const unsigned int dofs_per_line
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< parallel::distributed::Vector< Number > > > mg_constrained_dofs
unsigned int n_levels() const
virtual void restrict_and_add(const unsigned int from_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
Number local_element(const size_type local_index) const
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &mg_dof)
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
unsigned int global_dof_index
#define Assert(cond, exc)
virtual MPI_Comm get_communicator() const
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
void resize(const size_type size_in, const T &init=T())
MGLevelObject< parallel::distributed::Vector< Number > > ghosted_level_vector
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
bool element_is_continuous
virtual ~MGTransferMatrixFree()
std::vector< std::vector< unsigned int > > level_dof_indices
const unsigned int dofs_per_cell
void reinit(const size_type size, const bool omit_zeroing_entries=false)
const std::vector< Point< dim > > & get_unit_support_points() const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual void prolongate(const unsigned int to_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
virtual unsigned int n_global_levels() const
void do_prolongate_add(const unsigned int to_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
unsigned int n_child_cell_dofs
const Triangulation< dim, spacedim > & get_triangulation() const
internal::MatrixFreeFunctions::ShapeInfo< Number > shape_info
const unsigned int dofs_per_vertex
bool is_element(const size_type index) const
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
const types::global_dof_index invalid_dof_index
void do_restrict_add(const unsigned int from_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
size_type local_size() const
size_type n_elements() const
void build(const DoFHandler< dim, dim > &mg_dof)
virtual types::subdomain_id locally_owned_subdomain() const