// ---------------------------------------------------------------------
//
-// Copyright (C) 2009 - 2014 by the deal.II authors
+// Copyright (C) 2009 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
* @ingroup MeshWorker
* @author Guido Kanschat, 2009
*/
- template <int dim, int spacedim = dim>
+ template <int dim, int spacedim = dim, typename Number=double>
class VectorDataBase :
public VectorSelector
{
* Virtual, but empty destructor.
*/
virtual ~VectorDataBase();
+
/**
* The only function added to VectorSelector is an abstract virtual
* function implemented in the derived class template and called by
* base element.
*/
virtual void fill(
- std::vector<std::vector<std::vector<double> > > &values,
- std::vector<std::vector<std::vector<Tensor<1,dim> > > > &gradients,
- std::vector<std::vector<std::vector<Tensor<2,dim> > > > &hessians,
+ std::vector<std::vector<std::vector<Number> > > &values,
+ std::vector<std::vector<std::vector<Tensor<1,dim,Number> > > > &gradients,
+ std::vector<std::vector<std::vector<Tensor<2,dim,Number> > > > &hessians,
const FEValuesBase<dim,spacedim> &fe,
const std::vector<types::global_dof_index> &index,
const unsigned int component,
* the active cells.
*/
virtual void mg_fill(
- std::vector<std::vector<std::vector<double> > > &values,
- std::vector<std::vector<std::vector<Tensor<1,dim> > > > &gradients,
- std::vector<std::vector<std::vector<Tensor<2,dim> > > > &hessians,
+ std::vector<std::vector<std::vector<Number> > > &values,
+ std::vector<std::vector<std::vector<Tensor<1,dim,Number> > > > &gradients,
+ std::vector<std::vector<std::vector<Tensor<2,dim,Number> > > > &hessians,
const FEValuesBase<dim,spacedim> &fe,
const unsigned int level,
const std::vector<types::global_dof_index> &index,
*/
template <class VECTOR, int dim, int spacedim = dim>
class VectorData :
- public VectorDataBase<dim, spacedim>
+ public VectorDataBase<dim, spacedim,typename VECTOR::value_type>
{
public:
/**
* Constructor.
*/
MGVectorData();
+
/**
* Constructor using a prefilled VectorSelector
*/
// ---------------------------------------------------------------------
//
-// Copyright (C) 2009 - 2014 by the deal.II authors
+// Copyright (C) 2009 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
namespace MeshWorker
{
- template <int dim, int spacedim>
- VectorDataBase<dim, spacedim>::~VectorDataBase()
+ template <int dim, int spacedim, typename Number>
+ VectorDataBase<dim, spacedim, Number>::~VectorDataBase()
{}
- template <int dim, int spacedim>
- VectorDataBase<dim, spacedim>::VectorDataBase(const VectorSelector &v)
+ template <int dim, int spacedim, typename Number>
+ VectorDataBase<dim, spacedim, Number>::VectorDataBase(const VectorSelector &v)
:
VectorSelector(v)
{}
- template <int dim, int spacedim>
- VectorDataBase<dim, spacedim>::VectorDataBase()
+ template <int dim, int spacedim, typename Number>
+ VectorDataBase<dim, spacedim, Number>::VectorDataBase()
{}
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename Number>
void
- VectorDataBase<dim, spacedim>::initialize(const AnyData &d)
+ VectorDataBase<dim, spacedim, Number>::initialize(const AnyData &d)
{
this->data = d;
VectorSelector::initialize(d);
}
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename Number>
void
- VectorDataBase<dim, spacedim>::fill(
- std::vector<std::vector<std::vector<double> > > &,
- std::vector<std::vector<std::vector<Tensor<1,dim> > > > &,
- std::vector<std::vector<std::vector<Tensor<2,dim> > > > &,
+ VectorDataBase<dim, spacedim, Number>::fill(
+ std::vector<std::vector<std::vector<Number> > > &,
+ std::vector<std::vector<std::vector<Tensor<1,dim,Number> > > > &,
+ std::vector<std::vector<std::vector<Tensor<2,dim,Number> > > > &,
const FEValuesBase<dim,spacedim> &,
const std::vector<types::global_dof_index> &,
const unsigned int,
}
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename Number>
void
- VectorDataBase<dim, spacedim>::mg_fill(
- std::vector<std::vector<std::vector<double> > > &,
- std::vector<std::vector<std::vector<Tensor<1,dim> > > > &,
- std::vector<std::vector<std::vector<Tensor<2,dim> > > > &,
+ VectorDataBase<dim, spacedim, Number>::mg_fill(
+ std::vector<std::vector<std::vector<Number> > > &,
+ std::vector<std::vector<std::vector<Tensor<1,dim,Number> > > > &,
+ std::vector<std::vector<std::vector<Tensor<2,dim,Number> > > > &,
const FEValuesBase<dim,spacedim> &,
const unsigned int,
const std::vector<types::global_dof_index> &,
template <class VECTOR, int dim, int spacedim>
VectorData<VECTOR, dim, spacedim>::VectorData(const VectorSelector &s)
:
- VectorDataBase<dim, spacedim>(s)
+ VectorDataBase<dim, spacedim, typename VECTOR::value_type>(s)
{}