From: Bruno Turcksin Date: Tue, 7 Jul 2020 00:40:01 +0000 (+0000) Subject: Add get_vector_partitioner and get_dof_handler functions X-Git-Tag: v9.3.0-rc1~1314^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=3e130500b26a708d8544c34eac38454f6ebf64cb;p=dealii.git Add get_vector_partitioner and get_dof_handler functions --- diff --git a/include/deal.II/matrix_free/cuda_matrix_free.h b/include/deal.II/matrix_free/cuda_matrix_free.h index e621bc1033..71e0aad607 100644 --- a/include/deal.II/matrix_free/cuda_matrix_free.h +++ b/include/deal.II/matrix_free/cuda_matrix_free.h @@ -339,12 +339,31 @@ namespace CUDAWrappers initialize_dof_vector( LinearAlgebra::distributed::Vector &vec) const; + /** + * Return the partitioner that represents the locally owned data and the + * ghost indices where access is needed to for the cell loop. The + * partitioner is constructed from the locally owned dofs and ghost dofs + * given by the respective fields. If you want to have specific information + * about these objects, you can query them with the respective access + * functions. If you just want to initialize a (parallel) vector, you should + * usually prefer this data structure as the data exchange information can + * be reused from one vector to another. + */ + const std::shared_ptr & + get_vector_partitioner() const; + /** * Free all the memory allocated. */ void free(); + /** + * Return the DoFHandler. + */ + const DoFHandler & + get_dof_handler() const; + /** * Return an approximation of the memory consumption of this class in bytes. */ @@ -592,6 +611,11 @@ namespace CUDAWrappers */ std::vector row_start; + /** + * Pointer to the DoFHandler associated with the object. + */ + const DoFHandler *dof_handler; + friend class internal::ReinitHelper; }; @@ -654,7 +678,7 @@ namespace CUDAWrappers } - + /*----------------------- Helper functions ---------------------------------*/ /** * Compute the quadrature point index in the local cell of a given thread. * @@ -709,6 +733,30 @@ namespace CUDAWrappers return *(data->q_points + data->padding_length * cell + q_point_id_in_cell(n_q_points_1d)); } + + + /*----------------------- Inline functions ---------------------------------*/ + +# ifndef DOXYGEN + + template + const std::shared_ptr & + MatrixFree::get_vector_partitioner() const + { + return partitioner; + } + + template + const DoFHandler & + MatrixFree::get_dof_handler() const + { + Assert(dof_handler != nullptr, ExcNotInitialized()); + + return *dof_handler; + } + +# endif + } // namespace CUDAWrappers DEAL_II_NAMESPACE_CLOSE diff --git a/include/deal.II/matrix_free/cuda_matrix_free.templates.h b/include/deal.II/matrix_free/cuda_matrix_free.templates.h index a4ce00591c..d488792a67 100644 --- a/include/deal.II/matrix_free/cuda_matrix_free.templates.h +++ b/include/deal.II/matrix_free/cuda_matrix_free.templates.h @@ -612,6 +612,7 @@ namespace CUDAWrappers , constrained_dofs(nullptr) , padding_length(0) , my_id(-1) + , dof_handler(nullptr) {} @@ -845,12 +846,14 @@ namespace CUDAWrappers void MatrixFree::internal_reinit( const Mapping & mapping, - const DoFHandler & dof_handler, + const DoFHandler & dof_handler_, const AffineConstraints &constraints, const Quadrature<1> & quad, std::shared_ptr comm, const AdditionalData additional_data) { + dof_handler = &dof_handler_; + if (typeid(Number) == typeid(double)) cudaDeviceSetSharedMemConfig(cudaSharedMemBankSizeEightByte); @@ -868,9 +871,9 @@ namespace CUDAWrappers // TODO: only free if we actually need arrays of different length free(); - n_dofs = dof_handler.n_dofs(); + n_dofs = dof_handler->n_dofs(); - const FiniteElement &fe = dof_handler.get_fe(); + const FiniteElement &fe = dof_handler->get_fe(); fe_degree = fe.degree; // TODO this should be a templated parameter @@ -946,14 +949,14 @@ namespace CUDAWrappers cells_per_block = cells_per_block_shmem(dim, fe_degree); internal::ReinitHelper helper( - this, mapping, fe, quad, shape_info, dof_handler, update_flags); + this, mapping, fe, quad, shape_info, *dof_handler, update_flags); // Create a graph coloring using CellFilter = FilteredIterator::active_cell_iterator>; CellFilter begin(IteratorFilters::LocallyOwnedCell(), - dof_handler.begin_active()); - CellFilter end(IteratorFilters::LocallyOwnedCell(), dof_handler.end()); + dof_handler->begin_active()); + CellFilter end(IteratorFilters::LocallyOwnedCell(), dof_handler->end()); std::vector> graph; if (begin != end) @@ -976,10 +979,10 @@ namespace CUDAWrappers graph.resize(3, std::vector()); std::vector ghost_vertices( - dof_handler.get_triangulation().n_vertices(), false); + dof_handler->get_triangulation().n_vertices(), false); for (const auto cell : - dof_handler.get_triangulation().active_cell_iterators()) + dof_handler->get_triangulation().active_cell_iterators()) if (cell->is_ghost()) for (unsigned int i = 0; i < GeometryInfo::vertices_per_cell; @@ -1032,10 +1035,10 @@ namespace CUDAWrappers IndexSet locally_relevant_dofs; if (comm) { - DoFTools::extract_locally_relevant_dofs(dof_handler, + DoFTools::extract_locally_relevant_dofs(*dof_handler, locally_relevant_dofs); partitioner = std::make_shared( - dof_handler.locally_owned_dofs(), locally_relevant_dofs, *comm); + dof_handler->locally_owned_dofs(), locally_relevant_dofs, *comm); } for (unsigned int i = 0; i < n_colors; ++i) { @@ -1094,7 +1097,7 @@ namespace CUDAWrappers } else { - const unsigned int n_local_dofs = dof_handler.n_dofs(); + const unsigned int n_local_dofs = dof_handler->n_dofs(); unsigned int i_constraint = 0; for (unsigned int i = 0; i < n_local_dofs; ++i) {