template <int dim> class Function;
template <int dim> class FunctionMap;
template <int dim> class Quadrature;
-template <int dim> class QGauss2;
+template <int dim> class QGauss;
template <typename number> class Vector;
template <typename number> class FullMatrix;
* See the general documentation of this
* class for further information.
*/
- template <int dim>
+ template <int dim, class VECTOR>
static void project (const Mapping<dim> &mapping,
const DoFHandler<dim> &dof,
const ConstraintMatrix &constraints,
const Quadrature<dim> &quadrature,
const Function<dim> &function,
- Vector<double> &vec,
+ VECTOR &vec,
const bool enforce_zero_boundary = false,
- const Quadrature<dim-1> &q_boundary = QGauss2<dim-1>(),
+ const Quadrature<dim-1> &q_boundary = QGauss<dim-1>(2),
const bool project_to_boundary_first = false);
/**
* quadrature formula is an invalid
* object since it makes no sense in 1d.
*/
+ template <class VECTOR>
static void project (const Mapping<1> &mapping,
const DoFHandler<1> &dof,
const ConstraintMatrix &constraints,
const Quadrature<1> &quadrature,
const Function<1> &function,
- Vector<double> &vec,
+ VECTOR &vec,
const bool enforce_zero_boundary = false,
const Quadrature<0> &q_boundary = *invalid_face_quadrature,
const bool project_to_boundary_first = false);
* function, see above, with
* <tt>mapping=MappingQ1@<dim@>()</tt>.
*/
- template <int dim>
+ template <int dim, class VECTOR>
static void project (const DoFHandler<dim> &dof,
const ConstraintMatrix &constraints,
const Quadrature<dim> &quadrature,
const Function<dim> &function,
- Vector<double> &vec,
+ VECTOR &vec,
const bool enforce_zero_boundary = false,
- const Quadrature<dim-1> &q_boundary = QGauss2<dim-1>(),
+ const Quadrature<dim-1> &q_boundary = QGauss<dim-1>(2),
const bool project_to_boundary_first = false);
/**
#if deal_II_dimension == 1
+template <class VECTOR>
void VectorTools::project (const Mapping<1> &,
const DoFHandler<1> &,
const ConstraintMatrix &,
const Quadrature<1> &,
const Function<1> &,
- Vector<double> &,
+ VECTOR &,
const bool ,
const Quadrature<0> &,
const bool )
{
- // this function should easily be implemented
- // using the template below. However some
- // changes have to be made since faces don't
- // exist in 1D. Maybe integrate the creation of
- // zero boundary values into the
- // project_boundary_values function?
+ // this function should easily be
+ // implemented using the template
+ // below. However some changes have
+ // to be made since faces don't
+ // exist in 1D. Maybe integrate the
+ // creation of zero boundary values
+ // into the project_boundary_values
+ // function?
Assert (false, ExcNotImplemented());
}
#endif
-template <int dim>
+template <int dim, class VECTOR>
void VectorTools::project (const Mapping<dim> &mapping,
const DoFHandler<dim> &dof,
const ConstraintMatrix &constraints,
const Quadrature<dim> &quadrature,
const Function<dim> &function,
- Vector<double> &vec,
+ VECTOR &vec_result,
const bool enforce_zero_boundary,
const Quadrature<dim-1> &q_boundary,
const bool project_to_boundary_first)
{
Assert (dof.get_fe().n_components() == function.n_components,
ExcInvalidFE());
+
+ Assert (vec_result.size() == dof.n_dofs(),
+ ExcDimensionMismatch (vec_result.size(), dof.n_dofs()));
const FiniteElement<dim> &fe = dof.get_fe();
// set up mass matrix and right hand side
- vec.reinit (dof.n_dofs());
+ Vector<double> vec (dof.n_dofs());
SparsityPattern sparsity(dof.n_dofs(),
dof.n_dofs(),
dof.max_couplings_between_dofs());
// distribute solution
constraints.distribute (vec);
+
+ // copy vec into vec_result. we
+ // can't use ve_result itself
+ // above, since it may be of
+ // another type than Vector<double>
+ // and that wouldn't necessarily go
+ // together with the matrix and
+ // other functions
+ std::copy (vec.begin(), vec.end(), vec_result.begin());
}
-template <int dim>
+template <int dim, class VECTOR>
void VectorTools::project (const DoFHandler<dim> &dof,
const ConstraintMatrix &constraints,
const Quadrature<dim> &quadrature,
const Function<dim> &function,
- Vector<double> &vec,
+ VECTOR &vec,
const bool enforce_zero_boundary,
const Quadrature<dim-1> &q_boundary,
const bool project_to_boundary_first)
// $Id$
// Version: $Name$
//
-// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005 by the deal.II authors
+// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
#include "vectors.instance.h"
#undef VEC
-template
-void VectorTools::project<deal_II_dimension>
-(const DoFHandler<deal_II_dimension> &,
- const ConstraintMatrix &,
- const Quadrature<deal_II_dimension> &,
- const Function<deal_II_dimension> &,
- Vector<double> &,
- const bool,
- const Quadrature<deal_II_dimension-1> &,
- const bool);
template
void VectorTools::create_right_hand_side<deal_II_dimension>
std::map<unsigned int,double> &,
const std::vector<bool> &);
-template
-void VectorTools::project<deal_II_dimension>
-(const Mapping<deal_II_dimension> &,
- const DoFHandler<deal_II_dimension> &,
- const ConstraintMatrix &,
- const Quadrature<deal_II_dimension> &,
- const Function<deal_II_dimension> &,
- Vector<double> &,
- const bool,
- const Quadrature<deal_II_dimension-1> &,
- const bool);
#endif
const Quadrature<deal_II_dimension>&,
const VEC&,
const unsigned int);
+
+template
+void VectorTools::project
+(const Mapping<deal_II_dimension> &,
+ const DoFHandler<deal_II_dimension> &,
+ const ConstraintMatrix &,
+ const Quadrature<deal_II_dimension> &,
+ const Function<deal_II_dimension> &,
+ VEC &,
+ const bool,
+ const Quadrature<deal_II_dimension-1> &,
+ const bool);
+
+template
+void VectorTools::project
+(const DoFHandler<deal_II_dimension> &,
+ const ConstraintMatrix &,
+ const Quadrature<deal_II_dimension> &,
+ const Function<deal_II_dimension> &,
+ VEC &,
+ const bool,
+ const Quadrature<deal_II_dimension-1> &,
+ const bool);
<h3>deal.II</h3>
<ol>
+ <li> <p>
+ New: The function <code>VectorTools::project</code> functions can now
+ also be used for vector arguments of type other than
+ <code>Vector<double></code>.
+ <br>
+ (WB 2006/04/17)
+ </p>
+
<li> <p>
New: The function <code>GridTools::get_finest_common_cells</code> can be
used to find those cells that are shared by two triangulations that are
subcelldata \
interpolate_boundary_values_* \
have_same_coarse_mesh_* \
- get_finest_common_cells_*
+ get_finest_common_cells_* \
+ project_*
# from above list of regular expressions, generate the real set of
# tests