From: Wolfgang Bangerth Date: Sun, 9 Jan 2022 17:35:47 +0000 (-0700) Subject: Make CellData look like our usual style. X-Git-Tag: v9.4.0-rc1~636^2~2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=5c9f3a4a2e1707ba0674d8163326ba0961a7ac63;p=dealii.git Make CellData look like our usual style. --- diff --git a/include/deal.II/fe/fe_tools_extrapolate.templates.h b/include/deal.II/fe/fe_tools_extrapolate.templates.h index 18617b8fac..d6b7fb68d7 100644 --- a/include/deal.II/fe/fe_tools_extrapolate.templates.h +++ b/include/deal.II/fe/fe_tools_extrapolate.templates.h @@ -99,92 +99,54 @@ namespace FETools }; - // A structure holding all data - // of cells needed from other processes - // for the extrapolate algorithm. + /** + * A structure holding all data of cells needed from other processes + * for the extrapolate algorithm. + */ struct CellData { + /** + * Default constructor. + */ CellData(); + /** + * Constructor setting the `dof_values` member to the right size. + */ CellData(const unsigned int dofs_per_cell); - Vector dof_values; - - unsigned int tree_index; - - typename dealii::internal::p4est::types::quadrant quadrant; - - int receiver; - + /** + * Comparison operator. + */ bool - operator<(const CellData &rhs) const - { - if (dealii::internal::p4est::functions::quadrant_compare( - &quadrant, &rhs.quadrant) < 0) - return true; - - return false; - } - - unsigned int - bytes_for_buffer() const - { - return (sizeof(unsigned int) + // dofs_per_cell - dof_values.size() * sizeof(value_type) + // dof_values - sizeof(unsigned int) + // tree_index - sizeof(typename dealii::internal::p4est::types< - dim>::quadrant)); // quadrant - } + operator<(const CellData &rhs) const; + /** + * Pack the data of this object into a char[] buffer. + */ void - pack_data(std::vector &buffer) const - { - buffer.resize(bytes_for_buffer()); - - char *ptr = buffer.data(); - - unsigned int n_dofs = dof_values.size(); - std::memcpy(ptr, &n_dofs, sizeof(unsigned int)); - ptr += sizeof(unsigned int); - - std::memcpy(ptr, dof_values.begin(), n_dofs * sizeof(value_type)); - ptr += n_dofs * sizeof(value_type); - - std::memcpy(ptr, &tree_index, sizeof(unsigned int)); - ptr += sizeof(unsigned int); - - std::memcpy( - ptr, - &quadrant, - sizeof(typename dealii::internal::p4est::types::quadrant)); - ptr += sizeof(typename dealii::internal::p4est::types::quadrant); - - Assert(ptr == buffer.data() + buffer.size(), ExcInternalError()); - } + pack_data(std::vector &buffer) const; + /** + * Unpack the data of the given buffer into the members of this object. + */ void - unpack_data(const std::vector &buffer) - { - const char * ptr = buffer.data(); - unsigned int n_dofs; - memcpy(&n_dofs, ptr, sizeof(unsigned int)); - ptr += sizeof(unsigned int); + unpack_data(const std::vector &buffer); - dof_values.reinit(n_dofs); - std::memcpy(dof_values.begin(), ptr, n_dofs * sizeof(value_type)); - ptr += n_dofs * sizeof(value_type); - std::memcpy(&tree_index, ptr, sizeof(unsigned int)); - ptr += sizeof(unsigned int); + /** + * The values of degrees of freedom associated with the current cell. + */ + Vector dof_values; - std::memcpy( - &quadrant, - ptr, - sizeof(typename dealii::internal::p4est::types::quadrant)); - ptr += sizeof(typename dealii::internal::p4est::types::quadrant); + /** + * The tree within the forest (i.e., the coarse cell) and which of its + * descendents we are currently working on. + */ + unsigned int tree_index; + typename dealii::internal::p4est::types::quadrant quadrant; - Assert(ptr == buffer.data() + buffer.size(), ExcInternalError()); - } + int receiver; }; // Problem: The function extrapolates a polynomial @@ -400,6 +362,8 @@ namespace FETools unsigned int round; }; + + template class ExtrapolateImplementation<1, 1, OutVector> { @@ -472,10 +436,90 @@ namespace FETools template ExtrapolateImplementation::CellData::CellData( const unsigned int dofs_per_cell) - : tree_index(0) + : dof_values(dofs_per_cell) + , tree_index(0) , receiver(0) + {} + + + + template + bool + ExtrapolateImplementation::CellData::operator<( + const CellData &rhs) const + { + if (dealii::internal::p4est::functions::quadrant_compare( + &quadrant, &rhs.quadrant) < 0) + return true; + + return false; + } + + + + template + void + ExtrapolateImplementation::CellData::pack_data( + std::vector &buffer) const + { + // Compute how much memory we need to pack the data of this + // structure into a char[] buffer. + const unsigned int bytes_for_buffer = + (sizeof(unsigned int) + // dofs_per_cell + dof_values.size() * sizeof(value_type) + // dof_values + sizeof(unsigned int) + // tree_index + sizeof( + typename dealii::internal::p4est::types::quadrant)); // quadrant + + buffer.resize(bytes_for_buffer); + + char *ptr = buffer.data(); + + unsigned int n_dofs = dof_values.size(); + std::memcpy(ptr, &n_dofs, sizeof(unsigned int)); + ptr += sizeof(unsigned int); + + std::memcpy(ptr, dof_values.begin(), n_dofs * sizeof(value_type)); + ptr += n_dofs * sizeof(value_type); + + std::memcpy(ptr, &tree_index, sizeof(unsigned int)); + ptr += sizeof(unsigned int); + + std::memcpy(ptr, + &quadrant, + sizeof( + typename dealii::internal::p4est::types::quadrant)); + ptr += sizeof(typename dealii::internal::p4est::types::quadrant); + + Assert(ptr == buffer.data() + buffer.size(), ExcInternalError()); + } + + + + template + void + ExtrapolateImplementation::CellData::unpack_data( + const std::vector &buffer) { - dof_values.reinit(dofs_per_cell); + const char * ptr = buffer.data(); + unsigned int n_dofs; + memcpy(&n_dofs, ptr, sizeof(unsigned int)); + ptr += sizeof(unsigned int); + + dof_values.reinit(n_dofs); + std::memcpy(dof_values.begin(), ptr, n_dofs * sizeof(value_type)); + ptr += n_dofs * sizeof(value_type); + + std::memcpy(&tree_index, ptr, sizeof(unsigned int)); + ptr += sizeof(unsigned int); + + std::memcpy(&quadrant, + ptr, + sizeof( + typename dealii::internal::p4est::types::quadrant)); + ptr += sizeof(typename dealii::internal::p4est::types::quadrant); + + Assert(ptr == buffer.data() + buffer.size(), ExcInternalError()); }