Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mg_transfer_global_coarsening.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2020 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_mg_transfer_global_coarsening_h
17#define dealii_mg_transfer_global_coarsening_h
18
20
22
25
27
29
30// Forward declarations
31#ifndef DOXYGEN
32namespace internal
33{
34 class MGTwoLevelTransferImplementation;
35}
36#endif
37
38
39
44{
53 {
58 bisect,
69 };
70
75 unsigned int
77 const unsigned int degree,
78 const PolynomialCoarseningSequenceType &p_sequence);
79
85 std::vector<unsigned int>
87 const unsigned int max_degree,
88 const PolynomialCoarseningSequenceType &p_sequence);
89
98 template <int dim, int spacedim>
99 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
101 const Triangulation<dim, spacedim> &tria);
102
103} // namespace MGTransferGlobalCoarseningTools
104
105
106
110template <int dim, typename VectorType>
112{
113public:
117 void
118 prolongate(VectorType &dst, const VectorType &src) const;
119
123 void
124 restrict_and_add(VectorType &dst, const VectorType &src) const;
125
132 void
133 interpolate(VectorType &dst, const VectorType &src) const;
134};
135
136
137
142template <int dim, typename Number>
144{
145public:
150 void
152 const DoFHandler<dim> & dof_handler_coarse,
153 const AffineConstraints<Number> &constraint_fine,
154 const AffineConstraints<Number> &constraint_coarse);
155
163 void
165 const DoFHandler<dim> & dof_handler_fine,
166 const DoFHandler<dim> & dof_handler_coarse,
167 const AffineConstraints<Number> &constraint_fine,
168 const AffineConstraints<Number> &constraint_coarse);
169
178 static bool
179 fast_polynomial_transfer_supported(const unsigned int fe_degree_fine,
180 const unsigned int fe_degree_coarse);
181
185 void
188
192 void
195
202 void
205
206private:
214 struct MGTransferScheme
215 {
219 unsigned int n_coarse_cells;
220
225
229 unsigned int dofs_per_cell_fine;
230
234 unsigned int degree_coarse;
235
239 unsigned int degree_fine;
240
244 std::vector<Number> weights;
245
250
255
260
265
270 std::vector<unsigned int> level_dof_indices_coarse;
271
276 std::vector<unsigned int> level_dof_indices_fine;
277 };
278
282 std::vector<MGTransferScheme> schemes;
283
290
294 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_fine;
295
299 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_coarse;
300
309
314
320
326 std::vector<unsigned int> constraint_coarse_distribute_indices;
327
334
339 std::vector<unsigned int> constraint_coarse_distribute_ptr;
340
345 std::vector<unsigned int> distribute_local_to_global_indices;
346
352
357 std::vector<unsigned int> distribute_local_to_global_ptr;
358
362 unsigned int n_components;
363
364 friend class internal::MGTwoLevelTransferImplementation;
365};
366
367
368
379template <int dim, typename VectorType>
381{
382public:
386 using Number = typename VectorType::value_type;
387
397 const std::function<void(const unsigned int, VectorType &)>
399
403 void
404 prolongate(const unsigned int to_level,
405 VectorType & dst,
406 const VectorType & src) const override;
407
411 virtual void
412 restrict_and_add(const unsigned int from_level,
413 VectorType & dst,
414 const VectorType & src) const override;
415
422 template <class InVector, int spacedim>
423 void
426 const InVector & src) const;
427
434 template <class OutVector, int spacedim>
435 void
437 OutVector & dst,
438 const MGLevelObject<VectorType> &src) const;
439
454 template <class InVector, int spacedim>
455 void
458 const InVector & src) const;
459
460private:
465
469 const std::function<void(const unsigned int, VectorType &)>
471};
472
473
474
475#ifndef DOXYGEN
476
477/* ----------------------- Inline functions --------------------------------- */
478
479
480
481template <int dim, typename VectorType>
484 const std::function<void(const unsigned int, VectorType &)>
485 &initialize_dof_vector)
486 : transfer(transfer)
487 , initialize_dof_vector(initialize_dof_vector)
488{}
489
490
491
492template <int dim, typename VectorType>
493void
495 const unsigned int to_level,
496 VectorType & dst,
497 const VectorType & src) const
498{
499 this->transfer[to_level].prolongate(dst, src);
500}
501
502
503
504template <int dim, typename VectorType>
505void
507 const unsigned int from_level,
508 VectorType & dst,
509 const VectorType & src) const
510{
511 this->transfer[from_level].restrict_and_add(dst, src);
512}
513
514
515
516template <int dim, typename VectorType>
517template <class InVector, int spacedim>
518void
520 const DoFHandler<dim, spacedim> &dof_handler,
522 const InVector & src) const
523{
524 (void)dof_handler;
525
526 Assert(
527 initialize_dof_vector,
529 "To be able to use this function, a function to initialize an internal "
530 "DoF vector has to be provided in the constructor of "
531 "MGTransferGlobalCoarsening."));
532
533 for (unsigned int level = dst.min_level(); level <= dst.max_level(); ++level)
534 initialize_dof_vector(level, dst[level]);
535
536 dst[dst.max_level()].copy_locally_owned_data_from(src);
537}
538
539
540
541template <int dim, typename VectorType>
542template <class OutVector, int spacedim>
543void
545 const DoFHandler<dim, spacedim> &dof_handler,
546 OutVector & dst,
547 const MGLevelObject<VectorType> &src) const
548{
549 (void)dof_handler;
550
551 dst.copy_locally_owned_data_from(src[src.max_level()]);
552}
553
554
555
556template <int dim, typename VectorType>
557template <class InVector, int spacedim>
558void
560 const DoFHandler<dim, spacedim> &dof_handler,
562 const InVector & src) const
563{
564 (void)dof_handler;
565
566 Assert(
567 initialize_dof_vector,
569 "To be able to use this function, a function to initialize an internal "
570 "DoF vector has to be provided in the constructor of "
571 "MGTransferGlobalCoarsening."));
572
573 const unsigned int min_level = transfer.min_level();
574 const unsigned int max_level = transfer.max_level();
575
576 AssertDimension(min_level, dst.min_level());
577 AssertDimension(max_level, dst.max_level());
578
579 for (unsigned int level = min_level; level <= max_level; ++level)
580 initialize_dof_vector(level, dst[level]);
581
582 dst[transfer.max_level()].copy_locally_owned_data_from(src);
583
584 for (unsigned int l = max_level; l > min_level; --l)
585 this->transfer[l].interpolate(dst[l - 1], dst[l]);
586}
587
588#endif
589
591
592#endif
unsigned int max_level() const
unsigned int min_level() const
const MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > & transfer
MGTransferGlobalCoarsening(const MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > &transfer, const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector={})
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
typename VectorType::value_type Number
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
const std::function< void(const unsigned int, VectorType &)> initialize_dof_vector
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
void reinit_geometric_transfer(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine, const AffineConstraints< Number > &constraint_coarse)
void prolongate(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
static bool fast_polynomial_transfer_supported(const unsigned int fe_degree_fine, const unsigned int fe_degree_coarse)
void reinit_polynomial_transfer(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine, const AffineConstraints< Number > &constraint_coarse)
void interpolate(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void restrict_and_add(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void restrict_and_add(VectorType &dst, const VectorType &src) const
void interpolate(VectorType &dst, const VectorType &src) const
void prolongate(VectorType &dst, const VectorType &src) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
unsigned int level
Definition: grid_out.cc:4590
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
static ::ExceptionBase & ExcMessage(std::string arg1)
void interpolate(const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, OutVector &u2)
std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > create_geometric_coarsening_sequence(const Triangulation< dim, spacedim > &tria)
unsigned int create_next_polynomial_coarsening_degree(const unsigned int degree, const PolynomialCoarseningSequenceType &p_sequence)
std::vector< unsigned int > create_polynomial_coarsening_sequence(const unsigned int max_degree, const PolynomialCoarseningSequenceType &p_sequence)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)