Fix MGTransfer for 64 bit indices (and computations beyond 32 bit).
BlockVector<float>;
BlockVector<long double>;
+ parallel::distributed::Vector<double>;
+
@DEAL_II_EXPAND_TRILINOS_VECTOR@;
@DEAL_II_EXPAND_TRILINOS_MPI_VECTOR@;
}
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/block_sparsity_pattern.h>
#include <deal.II/lac/trilinos_sparse_matrix.h>
+#include <deal.II/lac/parallel_vector.h>
#include <deal.II/lac/vector_memory.h>
};
#ifdef DEAL_II_WITH_TRILINOS
+ template <typename Number>
+ struct MatrixSelector<parallel::distributed::Vector<Number> >
+ {
+ typedef ::dealii::TrilinosWrappers::SparsityPattern Sparsity;
+ typedef ::dealii::TrilinosWrappers::SparseMatrix Matrix;
+
+ template <class DSP, class DH>
+ static void reinit(Matrix &matrix, Sparsity &, int level, const DSP &dsp, DH &dh)
+ {
+ matrix.reinit(dh.locally_owned_mg_dofs(level+1),
+ dh.locally_owned_mg_dofs(level),
+ dsp, MPI_COMM_WORLD, true);
+ }
+
+ };
template <>
struct MatrixSelector<dealii::TrilinosWrappers::MPI::Vector>
{
* The data is organized as follows: one vector per level. Each element of
* these vectors contains first the global index, then the level index.
*/
- std::vector<std::vector<std::pair<types::global_dof_index, unsigned int> > >
+ std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
copy_indices;
/**
*
* Organization of the data is like for @p copy_indices_mine.
*/
- std::vector<std::vector<std::pair<types::global_dof_index, unsigned int> > >
+ std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
copy_indices_to_me;
/**
*
* Organization of the data is like for @p copy_indices_mine.
*/
- std::vector<std::vector<std::pair<types::global_dof_index, unsigned int> > >
+ std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
copy_indices_from_me;
}
}
+ /**
+ * Adjust vectors on all levels to correct size. Here, we just count the
+ * numbers of degrees of freedom on each level and @p reinit each level
+ * vector to this length.
+ */
+ template <int dim, typename number, int spacedim>
+ void
+ reinit_vector (const dealii::DoFHandler<dim,spacedim> &mg_dof,
+ const std::vector<unsigned int> &,
+ MGLevelObject<parallel::distributed::Vector<number> > &v)
+ {
+ const parallel::Triangulation<dim,spacedim> *tria =
+ (dynamic_cast<const parallel::Triangulation<dim,spacedim>*>
+ (&mg_dof.get_tria()));
+
+ for (unsigned int level=v.min_level();
+ level<=v.max_level(); ++level)
+ {
+ const IndexSet vector_index_set = v[level].locally_owned_elements();
+ if (vector_index_set.size() != mg_dof.locally_owned_mg_dofs(level).size() ||
+ mg_dof.locally_owned_mg_dofs(level) != vector_index_set)
+ v[level].reinit(mg_dof.locally_owned_mg_dofs(level), tria != 0 ?
+ tria->get_communicator() : MPI_COMM_SELF);
+ else
+ v[level] = 0.;
+ }
+ }
+
+
#ifdef DEAL_II_WITH_TRILINOS
/**
* Adjust vectors on all levels to correct size. Here, we just count the
--level;
VECTOR &dst_level = dst[level];
- typedef std::vector<std::pair<types::global_dof_index, unsigned int> >::const_iterator IT;
+ typedef std::vector<std::pair<types::global_dof_index, types::global_dof_index> >::const_iterator IT;
for (IT i= copy_indices[level].begin();
i != copy_indices[level].end(); ++i)
dst_level(i->second) = src(i->first);
dst = 0;
for (unsigned int level=0; level<mg_dof_handler.get_tria().n_global_levels(); ++level)
{
- typedef std::vector<std::pair<types::global_dof_index, unsigned int> >::const_iterator IT;
+ typedef std::vector<std::pair<types::global_dof_index, types::global_dof_index> >::const_iterator IT;
// First copy all indices local to this process
if (constraints==0)
// functions
for (unsigned int level=0; level<mg_dof_handler.get_tria().n_global_levels(); ++level)
{
- typedef std::vector<std::pair<types::global_dof_index, unsigned int> >::const_iterator IT;
+ typedef std::vector<std::pair<types::global_dof_index, types::global_dof_index> >::const_iterator IT;
if (constraints==0)
for (IT i= copy_indices[level].begin();
i != copy_indices[level].end(); ++i)
// for a cell and one of its children
std::vector<types::global_dof_index> dof_indices_parent (dofs_per_cell);
std::vector<types::global_dof_index> dof_indices_child (dofs_per_cell);
+ std::vector<types::global_dof_index> entries (dofs_per_cell);
// for each level: first build the sparsity
// pattern of the matrices and then build the
// level which have children
for (unsigned int level=0; level<n_levels-1; ++level)
{
-
- // reset the dimension of the structure.
- // note that for the number of entries
- // per row, the number of parent dofs
- // coupling to a child dof is
- // necessary. this, of course, is the
- // number of degrees of freedom per
+ // find the locally relevant level dofs for setting the writable rows of
+ // the sparsity pattern
+ IndexSet level_p1_relevant_dofs = mg_dof.locally_owned_mg_dofs(level+1);
+ std::vector<types::global_dof_index> relevant_dofs;
+ for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
+ cell != mg_dof.end(level); ++cell)
+ if (cell->has_children() &&
+ ( mg_dof.get_tria().locally_owned_subdomain()==numbers::invalid_subdomain_id
+ || cell->level_subdomain_id()==mg_dof.get_tria().locally_owned_subdomain()
+ ))
+ for (unsigned int child=0; child<cell->n_children(); ++child)
+ {
+ cell->child(child)->get_mg_dof_indices (dof_indices_child);
+ for (unsigned int i=0; i<dofs_per_cell; ++i)
+ if (!level_p1_relevant_dofs.is_element(dof_indices_child[i]))
+ relevant_dofs.push_back(dof_indices_child[i]);
+ }
+ std::sort(relevant_dofs.begin(), relevant_dofs.end());
+ level_p1_relevant_dofs.add_indices(relevant_dofs.begin(),
+ std::unique(relevant_dofs.begin(),
+ relevant_dofs.end()));
+
+ // reset the dimension of the structure. note that for the number of
+ // entries per row, the number of parent dofs coupling to a child dof is
+ // necessary. this, of course, is the number of degrees of freedom per
// cell
- // increment dofs_per_cell
- // since a useless diagonal
- // element will be stored
+ //
+ // increment dofs_per_cell since a useless diagonal element will be
+ // stored
DynamicSparsityPattern dsp (sizes[level+1],
- sizes[level]);
- std::vector<types::global_dof_index> entries (dofs_per_cell);
+ sizes[level],
+ level_p1_relevant_dofs);
for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
cell != mg_dof.end(level); ++cell)
if (cell->has_children() &&
if (global_mine && level_mine)
copy_indices[level].push_back(
- std::pair<unsigned int, unsigned int> (global_dof_indices[i], level_dof_indices[i]));
+ std::make_pair (global_dof_indices[i], level_dof_indices[i]));
else if (level_mine)
copy_indices_from_me[level].push_back(
- std::pair<unsigned int, unsigned int> (global_dof_indices[i], level_dof_indices[i]));
+ std::make_pair (global_dof_indices[i], level_dof_indices[i]));
else if (global_mine)
copy_indices_to_me[level].push_back(
- std::pair<unsigned int, unsigned int> (global_dof_indices[i], level_dof_indices[i]));
+ std::make_pair (global_dof_indices[i], level_dof_indices[i]));
else
continue;
// If we are in debugging mode, we order the copy indices, so we get
// more reliable output for regression texts
#ifdef DEBUG
- std::less<std::pair<types::global_dof_index, unsigned int> > compare;
+ std::less<std::pair<types::global_dof_index, types::global_dof_index> > compare;
for (unsigned int level=0; level<copy_indices.size(); ++level)
std::sort(copy_indices[level].begin(), copy_indices[level].end(), compare);
for (unsigned int level=0; level<copy_indices_from_me.size(); ++level)
}
for (deal_II_dimension : DIMENSIONS; V1,V2 : DEAL_II_VEC_TEMPLATES; S1, S2 : REAL_SCALARS)
- {
+ {
template void
MGTransferPrebuilt<V1<S1> >::copy_to_mg (
const DoFHandler<deal_II_dimension>&, MGLevelObject<V1<S1> >&, const V2<S2>&) const;
const MGLevelObject<V1<S1> >&) const;
}
+for (deal_II_dimension : DIMENSIONS; S2 : REAL_SCALARS)
+ {
+ template void
+ MGTransferPrebuilt<parallel::distributed::Vector<double> >::copy_to_mg (
+ const DoFHandler<deal_II_dimension>&, MGLevelObject<parallel::distributed::Vector<double> >&, const parallel::distributed::Vector<S2>&) const;
+ template void
+ MGTransferPrebuilt<parallel::distributed::Vector<double> >::copy_from_mg (const DoFHandler<deal_II_dimension>&, parallel::distributed::Vector<S2>&,
+ const MGLevelObject<parallel::distributed::Vector<double> >&) const;
+ template void
+ MGTransferPrebuilt<parallel::distributed::Vector<double> >::copy_from_mg_add (const DoFHandler<deal_II_dimension>&, parallel::distributed::Vector<S2>&,
+ const MGLevelObject<parallel::distributed::Vector<double> >&) const;
+ }
+
for(deal_II_dimension : DIMENSIONS)
- {
+ {
#ifdef DEAL_II_WITH_TRILINOS
template void
const MGLevelObject<TrilinosWrappers::MPI::Vector>&) const;
#endif
}
-