* Constructor with constraints. Equivalent to the default constructor
* followed by initialize_constraints().
*/
+ MGTransferPrebuilt (const MGConstrainedDoFs &mg_constrained_dofs);
+
+ /**
+ * Constructor with constraints. Equivalent to the default constructor
+ * followed by initialize_constraints().
+ *
+ * @deprecated @p constraints is unused.
+ */
MGTransferPrebuilt (const ConstraintMatrix &constraints,
- const MGConstrainedDoFs &mg_constrained_dofs);
+ const MGConstrainedDoFs &mg_constrained_dofs) DEAL_II_DEPRECATED;
/**
* Destructor.
/**
* Initialize the constraints to be used in build_matrices().
*/
+ void initialize_constraints (const MGConstrainedDoFs &mg_constrained_dofs);
+
+ /**
+ * Initialize the constraints to be used in build_matrices().
+ *
+ * @deprecated @p constraints is unused.
+ */
void initialize_constraints (const ConstraintMatrix &constraints,
- const MGConstrainedDoFs &mg_constrained_dofs);
+ const MGConstrainedDoFs &mg_constrained_dofs) DEAL_II_DEPRECATED;
/**
* Reset the object to the state it had right after the default constructor.
* boundary.
*/
std::vector<std::vector<bool> > interface_dofs;
-
- /**
- * The constraints of the global system.
- */
- SmartPointer<const ConstraintMatrix, MGTransferPrebuilt<VectorType> > constraints;
};
// ---------------------------------------------------------------------
//
-// Copyright (C) 2001 - 2015 by the deal.II authors
+// Copyright (C) 2001 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
* discontinuous finite elements or with no local refinement.
*/
MGTransferBlockBase ();
+
/**
* Constructor with constraint matrices as well as mg_constrained_dofs.
*/
+ MGTransferBlockBase (const MGConstrainedDoFs &mg_constrained_dofs);
+
+ /**
+ * Constructor with constraint matrices as well as mg_constrained_dofs.
+ *
+ * @deprecated @p constraints is unused.
+ */
MGTransferBlockBase (const ConstraintMatrix &constraints,
- const MGConstrainedDoFs &mg_constrained_dofs);
+ const MGConstrainedDoFs &mg_constrained_dofs) DEAL_II_DEPRECATED;
+
/**
* Memory used by this object.
*/
*/
std::vector<std::vector<std::vector<std::pair<unsigned int, unsigned int> > > >
copy_indices;
- /**
- * The constraints of the global system.
- */
- SmartPointer<const ConstraintMatrix, MGTransferBlockBase> constraints;
+
/**
* The mg_constrained_dofs of the level systems.
*/
* discontinuous finite elements or with no local refinement.
*/
MGTransferBlockSelect ();
+
+ /**
+ * Constructor with constraint matrices as well as mg_constrained_dofs.
+ */
+ MGTransferBlockSelect (const MGConstrainedDoFs &mg_constrained_dofs);
+
/**
* Constructor with constraint matrices as well as mg_constrained_dofs.
+ *
+ * @deprecated @p constraints is unused.
*/
MGTransferBlockSelect (const ConstraintMatrix &constraints,
- const MGConstrainedDoFs &mg_constrained_dofs);
+ const MGConstrainedDoFs &mg_constrained_dofs) DEAL_II_DEPRECATED;
+
/**
* Destructor.
*/
template<typename VectorType>
-MGTransferPrebuilt<VectorType>::MGTransferPrebuilt (const ConstraintMatrix &c, const MGConstrainedDoFs &mg_c)
- :
- constraints(&c)
+MGTransferPrebuilt<VectorType>::MGTransferPrebuilt (const MGConstrainedDoFs &mg_c)
+{
+ this->mg_constrained_dofs = &mg_c;
+}
+
+
+
+template<typename VectorType>
+MGTransferPrebuilt<VectorType>::MGTransferPrebuilt (const ConstraintMatrix &/*c*/, const MGConstrainedDoFs &mg_c)
{
this->mg_constrained_dofs = &mg_c;
}
template <typename VectorType>
void MGTransferPrebuilt<VectorType>::initialize_constraints
-(const ConstraintMatrix &c, const MGConstrainedDoFs &mg_c)
+(const MGConstrainedDoFs &mg_c)
{
- constraints = &c;
this->mg_constrained_dofs = &mg_c;
}
+template <typename VectorType>
+void MGTransferPrebuilt<VectorType>::initialize_constraints
+(const ConstraintMatrix &/*c*/, const MGConstrainedDoFs &mg_c)
+{
+ initialize_constraints(mg_c);
+}
+
+
+
template <typename VectorType>
void MGTransferPrebuilt<VectorType>::clear ()
{
prolongation_matrices.resize(0);
prolongation_sparsities.resize(0);
interface_dofs.resize(0);
- constraints = 0;
}
// ---------------------------------------------------------------------
//
-// Copyright (C) 2000 - 2014 by the deal.II authors
+// Copyright (C) 2000 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
{}
+
+MGTransferBlockBase::MGTransferBlockBase (
+ const MGConstrainedDoFs &mg_c)
+ :
+ mg_constrained_dofs(&mg_c)
+{}
+
+
+
MGTransferBlockBase::MGTransferBlockBase (
- const ConstraintMatrix &c, const MGConstrainedDoFs &mg_c)
+ const ConstraintMatrix &/*c*/, const MGConstrainedDoFs &mg_c)
:
- constraints(&c),
mg_constrained_dofs(&mg_c)
{}
}
+
template <typename number>
void MGTransferSelect<number>::restrict_and_add (
const unsigned int from_level,
{}
+
template <typename number>
MGTransferBlockSelect<number>::MGTransferBlockSelect (
- const ConstraintMatrix &c, const MGConstrainedDoFs &mg_c)
- : MGTransferBlockBase(c, mg_c)
+ const MGConstrainedDoFs &mg_c)
+ : MGTransferBlockBase(mg_c)
{}
+
+
+template <typename number>
+MGTransferBlockSelect<number>::MGTransferBlockSelect (
+ const ConstraintMatrix &/*c*/, const MGConstrainedDoFs &mg_c)
+ : MGTransferBlockBase(mg_c)
+{}
+
+
+
template <typename number>
MGTransferBlockSelect<number>::~MGTransferBlockSelect ()
{}
+
template <typename number>
void MGTransferBlockSelect<number>::prolongate (
const unsigned int to_level,