/**
* Create a right hand side
- * vector for a point source at point @p p
- * for vector-valued finite elements.
- * Prior content of the
- * given @p rhs_vector vector is
+ * vector for a point source at point @p p. This variation of the function
+ * is meant for vector-valued problems with exactly dim components (it will
+ * also work for problems with more than dim components, and in this case
+ * simply consider only the first dim components of the shape functions).
+ * It computes a right hand side that corresponds to a forcing function that
+ * is equal to a delta function times a given direction.
+ * In other words, it creates a vector $F$ so that
+ * $F_i = \int_\Omega [\mathbf d \delta(x-p)] \cdot \phi_i(x) dx$.
+ * Note here that $\phi_i$ is a vector-valued function. $\mathbf d$ is
+ * the given direction of the source term $\mathbf d \delta(x-p)$ and
+ * corresponds to the @p direction argument to be passed to this function.
+ *
+ * Prior content of the given @p rhs_vector vector is
* deleted.
*
* See the general documentation of this
void create_point_source_vector(const Mapping<dim,spacedim> &mapping,
const DoFHandler<dim,spacedim> &dof,
const Point<spacedim> &p,
- const Point<dim> &orientation,
+ const Point<dim> &direction,
Vector<double> &rhs_vector);
/**
template <int dim, int spacedim>
void create_point_source_vector(const DoFHandler<dim,spacedim> &dof,
const Point<spacedim> &p,
- const Point<dim> &orientation,
+ const Point<dim> &direction,
Vector<double> &rhs_vector);
/**
void create_point_source_vector(const hp::MappingCollection<dim,spacedim> &mapping,
const hp::DoFHandler<dim,spacedim> &dof,
const Point<spacedim> &p,
- const Point<dim> &orientation,
+ const Point<dim> &direction,
Vector<double> &rhs_vector);
/**
template <int dim, int spacedim>
void create_point_source_vector(const hp::DoFHandler<dim,spacedim> &dof,
const Point<spacedim> &p,
- const Point<dim> &orientation,
+ const Point<dim> &direction,
Vector<double> &rhs_vector);
/**
rhs_vector = 0;
- std::pair<typename DoFHandler<dim,spacedim>::active_cell_iterator, Point<spacedim> >
+ const std::pair<typename DoFHandler<dim,spacedim>::active_cell_iterator, Point<spacedim> >
cell_point =
GridTools::find_active_cell_around_point (mapping, dof_handler, p);
- Quadrature<dim> q(GeometryInfo<dim>::project_to_unit_cell(cell_point.second));
+ const Quadrature<dim> q(GeometryInfo<dim>::project_to_unit_cell(cell_point.second));
const FEValuesExtractors::Vector vec (0);
FEValues<dim,spacedim> fe_values(mapping, dof_handler.get_fe(),