<ol>
+ <li> Fixed: parallel::distributed::Vector now detects of the size of MPI
+ messages exceeds 2GB or if the local range exceeds the size of 32-bit
+ integers and throws an exception about unsupported range of operation.
+ <br>
+ (Martin Kronbichler, 2015/11/26)
+ </li>
+
<li> Fixed: GridGenerator::extract_boundary_mesh() in 3d could generate
surface cells that did not uniformly had a right- or left-handed coordinate
system associated with them when viewed from one side of the surface. This
{
AssertIndexRange (local_index, local_size() + n_ghost_indices_data);
if (local_index < local_size())
- return local_range_data.first + local_index;
+ return local_range_data.first + types::global_dof_index(local_index);
else
return ghost_indices_data.nth_index_in_set (local_index-local_size());
}
* multiple threads. This may or may not be desired when working also with
* MPI.
*
+ * <h4>Limitations regarding the vector size</h4>
+ *
+ * This vector class is based on two different number types for indexing.
+ * The so-called global index type encodes the overall size of the vector.
+ * Its type is <code>types::global_dof_index</type>. The largest possible
+ * value is <code>2^32-1</code> or approximately four billion in case 64
+ * bit integers are disabled at configuration of deal.II (default case) or
+ * <code>2^64-1</code> or approximately <code>10^19</code> if 64 bit
+ * integers are enabled (see the glossary entry on
+ * @ref GlobalDoFIndex for further information).
+ *
+ * The second relevant index type is the local index used within one MPI
+ * rank. As opposed to the global index, the implementation assumes
+ * 32-bit unsigned integers unconditionally. In other words, to actually
+ * use a vector with more than four billion entries, you need to use MPI
+ * with more than one rank (which in general is a safe assumption since
+ * four billion entries consume at least 16 GB of memory for floats or 32
+ * GB of memory for doubles). If more than 4 billion local elements are
+ * present, the implementation tries to detect that, which triggers an
+ * exception will be thrown and aborts the code. Note, however, that the
+ * detection of overflow is tricky and the detection mechanism might fail
+ * in some circumstances. Therefore, it is strongly recommended to not
+ * rely on this class to detect the case if local elements are more than
+ * four billion.
+ *
* @author Katharina Kormann, Martin Kronbichler, 2010, 2011
*/
template <typename Number>
import_data = new Number[part.n_import_indices()];
for (unsigned int i=0; i<n_import_targets; i++)
{
+ AssertThrow (static_cast<size_type>(part.import_targets()[i].second)*
+ sizeof(Number) <
+ static_cast<size_type>(std::numeric_limits<int>::max()),
+ ExcMessage("Index overflow: Maximum message size in MPI is 2GB. "
+ "The number of ghost entries times the size of 'Number' "
+ "exceeds this value. This is not supported."));
MPI_Recv_init (&import_data[current_index_start],
part.import_targets()[i].second*sizeof(Number),
MPI_BYTE,
current_index_start = part.local_size();
for (unsigned int i=0; i<n_ghost_targets; i++)
{
+ AssertThrow (static_cast<size_type>(part.ghost_targets()[i].second)*
+ sizeof(Number) <
+ static_cast<size_type>(std::numeric_limits<int>::max()),
+ ExcMessage("Index overflow: Maximum message size in MPI is 2GB. "
+ "The number of ghost entries times the size of 'Number' "
+ "exceeds this value. This is not supported."));
MPI_Send_init (&this->val[current_index_start],
part.ghost_targets()[i].second*sizeof(Number),
MPI_BYTE,
(locally_owned_indices.nth_index_in_set(0),
locally_owned_indices.nth_index_in_set(0) +
locally_owned_indices.n_elements());
+ AssertThrow (local_range_data.second-local_range_data.first <
+ static_cast<types::global_dof_index>(std::numeric_limits<unsigned int>::max()),
+ ExcMessage("Index overflow: This class supports at most 2^32-1 locally owned vector entries"));
locally_owned_range_data.set_size (locally_owned_indices.size());
locally_owned_range_data.add_range (local_range_data.first, local_range_data.second);
locally_owned_range_data.compress();
ghost_indices_data.set_size(locally_owned_range_data.size());
ghost_indices_data.subtract_set (locally_owned_range_data);
ghost_indices_data.compress();
+ AssertThrow (ghost_indices_data.n_elements() <
+ static_cast<types::global_dof_index>(std::numeric_limits<unsigned int>::max()),
+ ExcMessage("Index overflow: This class supports at most 2^32-1 ghost elements"));
n_ghost_indices_data = ghost_indices_data.n_elements();
have_ghost_indices =