*
* @author Wolfgang Bangerth, Guido Kanschat, 1999, 2000, 2001, 2002
*/
-template<typename number>
+template <typename number>
class MGTransferPrebuilt : public MGTransferBase<Vector<number> >
{
public:
* on each of the levels
* separately, i.a. an @p{MGVector}.
*/
- template<int dim, class InVector>
+ template <int dim, class InVector>
void
copy_to_mg (const MGDoFHandler<dim> &mg_dof,
MGLevelObject<Vector<number> > &dst,
* constrained degrees of freedom
* are set to zero.
*/
- template<int dim, class OutVector>
+ template <int dim, class OutVector>
void
copy_from_mg (const MGDoFHandler<dim> &mg_dof,
OutVector &dst,
* function, but probably not for
* continuous elements.
*/
- template<int dim, class OutVector>
+ template <int dim, class OutVector>
void
copy_from_mg_add (const MGDoFHandler<dim> &mg_dof,
OutVector &dst,
*/
std::vector<boost::shared_ptr<SparseMatrix<double> > > prolongation_matrices;
#endif
+
+ /**
+ * Structure that is used to
+ * disambiguate calls to
+ * @p{copy_to_mg} for 1d and
+ * non-1d. We provide two
+ * functions of @p{copy_to_mg},
+ * where the 1d function takes an
+ * argument of type
+ * @p{is_1d<true>} and the other
+ * one of type @p{is_1d<false>}.
+ */
+ template <bool> struct is_1d;
+
+ /**
+ * Implementation of the
+ * @p{copy_to_mg} function for
+ * 1d. We have to resort to some
+ * template trickery because we
+ * can't specialize template
+ * functions on the (outer)
+ * template of the class, without
+ * also fully specializing the
+ * inner (member function)
+ * template parameters. However,
+ * it can be done by adding the
+ * additional argument that
+ * converts template
+ * specialization into function
+ * overloading.
+ */
+ template <int dim, class InVector>
+ void
+ copy_to_mg (const MGDoFHandler<dim> &mg_dof,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src,
+ const is_1d<true> &) const;
+
+ /**
+ * Same for all other space
+ * dimensions.
+ */
+ template <int dim, class InVector>
+ void
+ copy_to_mg (const MGDoFHandler<dim> &mg_dof,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src,
+ const is_1d<false> &) const;
};
* on each of the levels
* separately, i.a. an @p{MGVector}.
*/
- template<int dim, class InVector>
+ template <int dim, class InVector>
void
copy_to_mg (const MGDoFHandler<dim> &mg_dof,
MGLevelObject<BlockVector<number> > &dst,
* constrained degrees of freedom
* are set to zero.
*/
- template<int dim, class OutVector>
+ template <int dim, class OutVector>
void
copy_from_mg (const MGDoFHandler<dim> &mg_dof,
OutVector &dst,
* function, but probably not for
* continuous elements.
*/
- template<int dim, class OutVector>
+ template <int dim, class OutVector>
void
copy_from_mg_add (const MGDoFHandler<dim> &mg_dof,
OutVector &dst,
const MGLevelObject<BlockVector<number> > &src) const;
+
+ private:
+ /**
+ * Structure that is used to
+ * disambiguate calls to
+ * @p{copy_to_mg} for 1d and
+ * non-1d. We provide two
+ * functions of @p{copy_to_mg},
+ * where the 1d function takes an
+ * argument of type
+ * @p{is_1d<true>} and the other
+ * one of type @p{is_1d<false>}.
+ */
+ template <bool> struct is_1d;
+
+ /**
+ * Implementation of the
+ * @p{copy_to_mg} function for
+ * 1d. We have to resort to some
+ * template trickery because we
+ * can't specialize template
+ * functions on the (outer)
+ * template of the class, without
+ * also fully specializing the
+ * inner (member function)
+ * template parameters. However,
+ * it can be done by adding the
+ * additional argument that
+ * converts template
+ * specialization into function
+ * overloading.
+ */
+ template <int dim, class InVector>
+ void
+ copy_to_mg (const MGDoFHandler<dim> &mg_dof,
+ MGLevelObject<BlockVector<number> > &dst,
+ const InVector &src,
+ const is_1d<true> &) const;
+
+ /**
+ * Same for all other space
+ * dimensions.
+ */
+ template <int dim, class InVector>
+ void
+ copy_to_mg (const MGDoFHandler<dim> &mg_dof,
+ MGLevelObject<BlockVector<number> > &dst,
+ const InVector &src,
+ const is_1d<false> &) const;
};
Vector<number> &dst,
const Vector<number> &src) const;
+ /**
+ * Structure that is used to
+ * disambiguate calls to
+ * @p{copy_to_mg} for 1d and
+ * non-1d. We provide two
+ * functions of @p{copy_to_mg},
+ * where the 1d function takes an
+ * argument of type
+ * @p{is_1d<true>} and the other
+ * one of type @p{is_1d<false>}.
+ */
+ template <bool> struct is_1d;
+
/**
* Transfer from a vector on the
* global grid to vectors defined
* on each of the levels
* separately, i.a. an @p{MGVector}.
*/
- template<int dim, class InVector>
+ template <int dim, class InVector>
void
copy_to_mg (const MGDoFHandler<dim> &mg_dof,
MGLevelObject<Vector<number> > &dst,
* constrained degrees of freedom
* are set to zero.
*/
- template<int dim, class OutVector>
+ template <int dim, class OutVector>
void
copy_from_mg (const MGDoFHandler<dim> &mg_dof,
OutVector &dst,
* function, but probably not for
* continuous elements.
*/
- template<int dim, class OutVector>
+ template <int dim, class OutVector>
void
copy_from_mg_add (const MGDoFHandler<dim> &mg_dof,
OutVector &dst,
* Selected component.
*/
unsigned int selected;
+
+
+ /**
+ * Implementation of the
+ * @p{copy_to_mg} function for
+ * 1d. We have to resort to some
+ * template trickery because we
+ * can't specialize template
+ * functions on the (outer)
+ * template of the class, without
+ * also fully specializing the
+ * inner (member function)
+ * template parameters. However,
+ * it can be done by adding the
+ * additional argument that
+ * converts template
+ * specialization into function
+ * overloading.
+ */
+ template <int dim, class InVector>
+ void
+ copy_to_mg (const MGDoFHandler<dim> &mg_dof,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src,
+ const is_1d<true> &) const;
+
+ /**
+ * Same for all other space
+ * dimensions.
+ */
+ template <int dim, class InVector>
+ void
+ copy_to_mg (const MGDoFHandler<dim> &mg_dof,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src,
+ const is_1d<false> &) const;
};
//
//----------------------------------------------------------------------------
-#ifndef __deal2__mg_tools_templates_h
-#define __deal2__mg_tools_templates_h
+#ifndef __deal2__mg_transfer_templates_h
+#define __deal2__mg_transfer_templates_h
#include <lac/sparse_matrix.h>
#include <dofs/dof_constraints.h>
/* --------------------- MGTransferPrebuilt -------------- */
-//TODO[GK]: this file must really be changed: it contains #if's for deal_II_dimension, but we can't use this in headers, since application programs might not use this way to select the dimension. the only way in header files is to have proper template specializations
+template <typename number>
+template <int dim, class InVector>
+void
+MGTransferPrebuilt<number>::copy_to_mg (const MGDoFHandler<dim> &mg_dof_handler,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src) const
+{
+ // forward to the correct
+ // specialization
+ copy_to_mg (mg_dof_handler, dst, src, is_1d<dim==1>());
+}
-//TODO:[?] This function needs to be specially implemented, since in 2d mode we use faces
-#if deal_II_dimension == 1
template <typename number>
template <int dim, class InVector>
void
-MGTransferPrebuilt<number>::copy_to_mg (
- const MGDoFHandler<dim>&,
- MGLevelObject<Vector<number> >&,
- const InVector&) const
+MGTransferPrebuilt<number>::copy_to_mg (const MGDoFHandler<dim> &,
+ MGLevelObject<Vector<number> > &,
+ const InVector &,
+ const is_1d<true> &) const
{
Assert(false, ExcNotImplemented());
}
-#else
-
template <typename number>
template <int dim, class InVector>
void
-MGTransferPrebuilt<number>::copy_to_mg (
- const MGDoFHandler<dim>& mg_dof_handler,
- MGLevelObject<Vector<number> >& dst,
- const InVector& src) const
+MGTransferPrebuilt<number>::copy_to_mg (const MGDoFHandler<dim> &mg_dof_handler,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src,
+ const is_1d<false> &) const
{
// Make src a real finite element function
// InVector src = osrc;
};
}
-#endif
template <typename number>
/* --------------------- MGTransferSelect -------------- */
-//TODO:[?] This function needs to be specially implemented, since in 2d mode we use faces
-#if deal_II_dimension == 1
+
+template <typename number>
+template <int dim, class InVector>
+void
+MGTransferSelect<number>::
+copy_to_mg (const MGDoFHandler<dim> &mg_dof_handler,
+ MGLevelObject<Vector<number> > &dst,
+ const InVector &src) const
+{
+ // forward to the correct
+ // specialization
+ copy_to_mg (mg_dof_handler, dst, src, is_1d<dim==1>());
+}
+
+
template <typename number>
template <int dim, class InVector>
MGTransferSelect<number>::
copy_to_mg (const MGDoFHandler<dim> &,
MGLevelObject<Vector<number> > &,
- const InVector &) const
+ const InVector &,
+ const is_1d<true> &) const
{
Assert(false, ExcNotImplemented());
}
-#else
template <typename number>
MGTransferSelect<number>::
copy_to_mg (const MGDoFHandler<dim> &mg_dof_handler,
MGLevelObject<Vector<number> > &dst,
- const InVector &osrc) const
+ const InVector &osrc,
+ const is_1d<false> &) const
{
// Make src a real finite element function
InVector src = osrc;
};
}
-#endif
template <typename number>
/* --------------------- MGTransferBlock -------------- */
-#if deal_II_dimension == 1
+
+template <typename number>
+template <int dim, class InVector>
+void
+MGTransferBlock<number>::
+copy_to_mg (const MGDoFHandler<dim> &mg_dof_handler,
+ MGLevelObject<BlockVector<number> > &dst,
+ const InVector &src) const
+{
+ // forward to the correct
+ // specialization
+ copy_to_mg (mg_dof_handler, dst, src, is_1d<dim==1>());
+}
+
+
template <typename number>
template <int dim, class InVector>
MGTransferBlock<number>::
copy_to_mg (const MGDoFHandler<dim> &,
MGLevelObject<BlockVector<number> > &,
- const InVector &) const
+ const InVector &,
+ const is_1d<true> &) const
{
Assert(false, ExcNotImplemented());
}
-#else
template <typename number>
MGTransferBlock<number>::
copy_to_mg (const MGDoFHandler<dim> &mg_dof_handler,
MGLevelObject<BlockVector<number> > &dst,
- const InVector &src) const
+ const InVector &src,
+ const is_1d<false> &) const
{
// Make src a real finite element
// function
};
}
-#endif
+
template <typename number>