* ExcShapeFunctionNotPrimitive. In that case, use the
* shape_value_component() function.
*
- * @param function_no Number of the shape function to be evaluated. Note
+ * @param i Number of the shape function $\varphi_i$ to be evaluated. Note
* that this number runs from zero to dofs_per_cell, even in the case of an
* FEFaceValues or FESubfaceValues object.
*
- * @param point_no Number of the quadrature point at which function is to be
+ * @param q_point Number of the quadrature point at which function is to be
* evaluated
*
* @dealiiRequiresUpdateFlags{update_values}
*/
const double &
- shape_value(const unsigned int function_no,
- const unsigned int point_no) const;
+ shape_value(const unsigned int i, const unsigned int q_point) const;
/**
* Compute one vector component of the value of a shape function at a
* shape function is not primitive, but then it is necessary since the other
* function cannot be used.
*
- * @param function_no Number of the shape function to be evaluated.
+ * @param i Number of the shape function $\varphi_i$ to be evaluated.
*
- * @param point_no Number of the quadrature point at which function is to be
+ * @param q_point Number of the quadrature point at which function is to be
* evaluated.
*
* @param component vector component to be evaluated.
* @dealiiRequiresUpdateFlags{update_values}
*/
double
- shape_value_component(const unsigned int function_no,
- const unsigned int point_no,
+ shape_value_component(const unsigned int i,
+ const unsigned int q_point,
const unsigned int component) const;
/**
- * Compute the gradient of the <tt>function_no</tt>th shape function at the
+ * Compute the gradient of the <tt>i</tt>th shape function at the
* <tt>quadrature_point</tt>th quadrature point with respect to real cell
* coordinates. If you want to get the derivative in one of the coordinate
* directions, use the appropriate function of the Tensor class to extract
* The same holds for the arguments of this function as for the
* shape_value() function.
*
- * @param function_no Number of the shape function to be evaluated.
+ * @param i Number of the shape function $\varphi_i$ to be evaluated.
*
- * @param quadrature_point Number of the quadrature point at which function
+ * @param q_point Number of the quadrature point at which function
* is to be evaluated.
*
* @dealiiRequiresUpdateFlags{update_gradients}
*/
const Tensor<1, spacedim> &
- shape_grad(const unsigned int function_no,
- const unsigned int quadrature_point) const;
+ shape_grad(const unsigned int i, const unsigned int q_point) const;
/**
* Return one vector component of the gradient of a shape function at a
* @dealiiRequiresUpdateFlags{update_gradients}
*/
Tensor<1, spacedim>
- shape_grad_component(const unsigned int function_no,
- const unsigned int point_no,
+ shape_grad_component(const unsigned int i,
+ const unsigned int q_point,
const unsigned int component) const;
/**
- * Second derivatives of the <tt>function_no</tt>th shape function at the
- * <tt>point_no</tt>th quadrature point with respect to real cell
+ * Second derivatives of the <tt>i</tt>th shape function at the
+ * <tt>q_point</tt>th quadrature point with respect to real cell
* coordinates. If you want to get the derivatives in one of the coordinate
* directions, use the appropriate function of the Tensor class to extract
* one component. Since only a reference to the hessian values is returned,
* @dealiiRequiresUpdateFlags{update_hessians}
*/
const Tensor<2, spacedim> &
- shape_hessian(const unsigned int function_no,
- const unsigned int point_no) const;
+ shape_hessian(const unsigned int i, const unsigned int q_point) const;
/**
* Return one vector component of the hessian of a shape function at a
* @dealiiRequiresUpdateFlags{update_hessians}
*/
Tensor<2, spacedim>
- shape_hessian_component(const unsigned int function_no,
- const unsigned int point_no,
+ shape_hessian_component(const unsigned int i,
+ const unsigned int q_point,
const unsigned int component) const;
/**
- * Third derivatives of the <tt>function_no</tt>th shape function at the
- * <tt>point_no</tt>th quadrature point with respect to real cell
+ * Third derivatives of the <tt>i</tt>th shape function at the
+ * <tt>q_point</tt>th quadrature point with respect to real cell
* coordinates. If you want to get the 3rd derivatives in one of the
* coordinate directions, use the appropriate function of the Tensor class
* to extract one component. Since only a reference to the 3rd derivative
* @dealiiRequiresUpdateFlags{update_3rd_derivatives}
*/
const Tensor<3, spacedim> &
- shape_3rd_derivative(const unsigned int function_no,
- const unsigned int point_no) const;
+ shape_3rd_derivative(const unsigned int i, const unsigned int q_point) const;
/**
* Return one vector component of the third derivative of a shape function
* @dealiiRequiresUpdateFlags{update_3rd_derivatives}
*/
Tensor<3, spacedim>
- shape_3rd_derivative_component(const unsigned int function_no,
- const unsigned int point_no,
+ shape_3rd_derivative_component(const unsigned int i,
+ const unsigned int q_point,
const unsigned int component) const;
/** @} */
quadrature_point_indices() const;
/**
- * Position of the <tt>q</tt>th quadrature point in real space.
+ * Return the location of the <tt>q_point</tt>th quadrature point in
+ * real space.
*
* @dealiiRequiresUpdateFlags{update_quadrature_points}
*/
const Point<spacedim> &
- quadrature_point(const unsigned int q) const;
+ quadrature_point(const unsigned int q_point) const;
/**
* Return a reference to the vector of quadrature points in real space.
/**
* Mapped quadrature weight. If this object refers to a volume evaluation
* (i.e. the derived class is of type FEValues), then this is the Jacobi
- * determinant times the weight of the *<tt>i</tt>th unit quadrature point.
+ * determinant times the weight of the <tt>q_point</tt>th unit quadrature
+ * point.
*
* For surface evaluations (i.e. classes FEFaceValues or FESubfaceValues),
* it is the mapped surface element times the weight of the quadrature
* @dealiiRequiresUpdateFlags{update_JxW_values}
*/
double
- JxW(const unsigned int quadrature_point) const;
+ JxW(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by JxW().
* @dealiiRequiresUpdateFlags{update_jacobians}
*/
const DerivativeForm<1, dim, spacedim> &
- jacobian(const unsigned int quadrature_point) const;
+ jacobian(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_jacobian_grads}
*/
const DerivativeForm<2, dim, spacedim> &
- jacobian_grad(const unsigned int quadrature_point) const;
+ jacobian_grad(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_jacobian_pushed_forward_grads}
*/
const Tensor<3, spacedim> &
- jacobian_pushed_forward_grad(const unsigned int quadrature_point) const;
+ jacobian_pushed_forward_grad(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_jacobian_2nd_derivatives}
*/
const DerivativeForm<3, dim, spacedim> &
- jacobian_2nd_derivative(const unsigned int quadrature_point) const;
+ jacobian_2nd_derivative(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_jacobian_pushed_forward_2nd_derivatives}
*/
const Tensor<4, spacedim> &
- jacobian_pushed_forward_2nd_derivative(
- const unsigned int quadrature_point) const;
+ jacobian_pushed_forward_2nd_derivative(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_jacobian_3rd_derivatives}
*/
const DerivativeForm<4, dim, spacedim> &
- jacobian_3rd_derivative(const unsigned int quadrature_point) const;
+ jacobian_3rd_derivative(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_jacobian_pushed_forward_3rd_derivatives}
*/
const Tensor<5, spacedim> &
- jacobian_pushed_forward_3rd_derivative(
- const unsigned int quadrature_point) const;
+ jacobian_pushed_forward_3rd_derivative(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* @dealiiRequiresUpdateFlags{update_inverse_jacobians}
*/
const DerivativeForm<1, spacedim, dim> &
- inverse_jacobian(const unsigned int quadrature_point) const;
+ inverse_jacobian(const unsigned int q_point) const;
/**
* Return a reference to the array holding the values returned by
* Return the normal vector at a quadrature point. If you call this
* function for a face (i.e., when using a FEFaceValues or FESubfaceValues
* object), then this function returns the outward normal vector to
- * the cell at the <tt>i</tt>th quadrature point of the face.
+ * the cell at the <tt>q_point</tt>th quadrature point of the face.
*
* In contrast, if you call this function for a cell of codimension one
* (i.e., when using a `FEValues<dim,spacedim>` object with
* @dealiiRequiresUpdateFlags{update_normal_vectors}
*/
const Tensor<1, spacedim> &
- normal_vector(const unsigned int i) const;
+ normal_vector(const unsigned int q_point) const;
/**
* Return the normal vectors at all quadrature points represented by
const hp::QCollection<dim - 1> & quadrature);
/**
- * Boundary form of the transformation of the cell at the <tt>i</tt>th
+ * Boundary form of the transformation of the cell at the <tt>q_point</tt>th
* quadrature point. See
* @ref GlossBoundaryForm.
*
* @dealiiRequiresUpdateFlags{update_boundary_forms}
*/
const Tensor<1, spacedim> &
- boundary_form(const unsigned int i) const;
+ boundary_form(const unsigned int q_point) const;
/**
* Return the list of outward normal vectors times the Jacobian of the
template <int dim, int spacedim>
inline const double &
FEValuesBase<dim, spacedim>::shape_value(const unsigned int i,
- const unsigned int j) const
+ const unsigned int q_point) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
Assert(this->update_flags & update_values,
// if the entire FE is primitive,
// then we can take a short-cut:
if (fe->is_primitive())
- return this->finite_element_output.shape_values(i, j);
+ return this->finite_element_output.shape_values(i, q_point);
else
{
// otherwise, use the mapping
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() +
fe->system_to_component_index(i).first];
- return this->finite_element_output.shape_values(row, j);
+ return this->finite_element_output.shape_values(row, q_point);
}
}
inline double
FEValuesBase<dim, spacedim>::shape_value_component(
const unsigned int i,
- const unsigned int j,
+ const unsigned int q_point,
const unsigned int component) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
const unsigned int row =
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() + component];
- return this->finite_element_output.shape_values(row, j);
+ return this->finite_element_output.shape_values(row, q_point);
}
template <int dim, int spacedim>
inline const Tensor<1, spacedim> &
FEValuesBase<dim, spacedim>::shape_grad(const unsigned int i,
- const unsigned int j) const
+ const unsigned int q_point) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
Assert(this->update_flags & update_gradients,
// if the entire FE is primitive,
// then we can take a short-cut:
if (fe->is_primitive())
- return this->finite_element_output.shape_gradients[i][j];
+ return this->finite_element_output.shape_gradients[i][q_point];
else
{
// otherwise, use the mapping
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() +
fe->system_to_component_index(i).first];
- return this->finite_element_output.shape_gradients[row][j];
+ return this->finite_element_output.shape_gradients[row][q_point];
}
}
inline Tensor<1, spacedim>
FEValuesBase<dim, spacedim>::shape_grad_component(
const unsigned int i,
- const unsigned int j,
+ const unsigned int q_point,
const unsigned int component) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
const unsigned int row =
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() + component];
- return this->finite_element_output.shape_gradients[row][j];
+ return this->finite_element_output.shape_gradients[row][q_point];
}
template <int dim, int spacedim>
inline const Tensor<2, spacedim> &
FEValuesBase<dim, spacedim>::shape_hessian(const unsigned int i,
- const unsigned int j) const
+ const unsigned int q_point) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
Assert(this->update_flags & update_hessians,
// if the entire FE is primitive,
// then we can take a short-cut:
if (fe->is_primitive())
- return this->finite_element_output.shape_hessians[i][j];
+ return this->finite_element_output.shape_hessians[i][q_point];
else
{
// otherwise, use the mapping
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() +
fe->system_to_component_index(i).first];
- return this->finite_element_output.shape_hessians[row][j];
+ return this->finite_element_output.shape_hessians[row][q_point];
}
}
inline Tensor<2, spacedim>
FEValuesBase<dim, spacedim>::shape_hessian_component(
const unsigned int i,
- const unsigned int j,
+ const unsigned int q_point,
const unsigned int component) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
const unsigned int row =
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() + component];
- return this->finite_element_output.shape_hessians[row][j];
+ return this->finite_element_output.shape_hessians[row][q_point];
}
template <int dim, int spacedim>
inline const Tensor<3, spacedim> &
-FEValuesBase<dim, spacedim>::shape_3rd_derivative(const unsigned int i,
- const unsigned int j) const
+FEValuesBase<dim, spacedim>::shape_3rd_derivative(
+ const unsigned int i,
+ const unsigned int q_point) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
Assert(this->update_flags & update_3rd_derivatives,
// if the entire FE is primitive,
// then we can take a short-cut:
if (fe->is_primitive())
- return this->finite_element_output.shape_3rd_derivatives[i][j];
+ return this->finite_element_output.shape_3rd_derivatives[i][q_point];
else
{
// otherwise, use the mapping
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() +
fe->system_to_component_index(i).first];
- return this->finite_element_output.shape_3rd_derivatives[row][j];
+ return this->finite_element_output.shape_3rd_derivatives[row][q_point];
}
}
inline Tensor<3, spacedim>
FEValuesBase<dim, spacedim>::shape_3rd_derivative_component(
const unsigned int i,
- const unsigned int j,
+ const unsigned int q_point,
const unsigned int component) const
{
AssertIndexRange(i, fe->n_dofs_per_cell());
const unsigned int row =
this->finite_element_output
.shape_function_to_row_table[i * fe->n_components() + component];
- return this->finite_element_output.shape_3rd_derivatives[row][j];
+ return this->finite_element_output.shape_3rd_derivatives[row][q_point];
}
template <int dim, int spacedim>
inline const Tensor<3, spacedim> &
FEValuesBase<dim, spacedim>::jacobian_pushed_forward_grad(
- const unsigned int i) const
+ const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobian_pushed_forward_grads,
ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobian_pushed_forward_grads[i];
+ return this->mapping_output.jacobian_pushed_forward_grads[q_point];
}
template <int dim, int spacedim>
inline const DerivativeForm<3, dim, spacedim> &
-FEValuesBase<dim, spacedim>::jacobian_2nd_derivative(const unsigned int i) const
+FEValuesBase<dim, spacedim>::jacobian_2nd_derivative(
+ const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobian_2nd_derivatives,
ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobian_2nd_derivatives[i];
+ return this->mapping_output.jacobian_2nd_derivatives[q_point];
}
template <int dim, int spacedim>
inline const Tensor<4, spacedim> &
FEValuesBase<dim, spacedim>::jacobian_pushed_forward_2nd_derivative(
- const unsigned int i) const
+ const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobian_pushed_forward_2nd_derivatives,
ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_2nd_derivatives"));
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
+ return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[q_point];
}
template <int dim, int spacedim>
inline const DerivativeForm<4, dim, spacedim> &
-FEValuesBase<dim, spacedim>::jacobian_3rd_derivative(const unsigned int i) const
+FEValuesBase<dim, spacedim>::jacobian_3rd_derivative(
+ const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobian_3rd_derivatives,
ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobian_3rd_derivatives[i];
+ return this->mapping_output.jacobian_3rd_derivatives[q_point];
}
template <int dim, int spacedim>
inline const Tensor<5, spacedim> &
FEValuesBase<dim, spacedim>::jacobian_pushed_forward_3rd_derivative(
- const unsigned int i) const
+ const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobian_pushed_forward_3rd_derivatives,
ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_3rd_derivatives"));
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
+ return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[q_point];
}
template <int dim, int spacedim>
inline const Point<spacedim> &
-FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int i) const
+FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int q_point) const
{
Assert(this->update_flags & update_quadrature_points,
ExcAccessToUninitializedField("update_quadrature_points"));
- AssertIndexRange(i, this->mapping_output.quadrature_points.size());
+ AssertIndexRange(q_point, this->mapping_output.quadrature_points.size());
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.quadrature_points[i];
+ return this->mapping_output.quadrature_points[q_point];
}
template <int dim, int spacedim>
inline double
-FEValuesBase<dim, spacedim>::JxW(const unsigned int i) const
+FEValuesBase<dim, spacedim>::JxW(const unsigned int q_point) const
{
Assert(this->update_flags & update_JxW_values,
ExcAccessToUninitializedField("update_JxW_values"));
- AssertIndexRange(i, this->mapping_output.JxW_values.size());
+ AssertIndexRange(q_point, this->mapping_output.JxW_values.size());
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.JxW_values[i];
+ return this->mapping_output.JxW_values[q_point];
}
template <int dim, int spacedim>
inline const DerivativeForm<1, dim, spacedim> &
-FEValuesBase<dim, spacedim>::jacobian(const unsigned int i) const
+FEValuesBase<dim, spacedim>::jacobian(const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobians,
ExcAccessToUninitializedField("update_jacobians"));
- AssertIndexRange(i, this->mapping_output.jacobians.size());
+ AssertIndexRange(q_point, this->mapping_output.jacobians.size());
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobians[i];
+ return this->mapping_output.jacobians[q_point];
}
template <int dim, int spacedim>
inline const DerivativeForm<2, dim, spacedim> &
-FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int i) const
+FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int q_point) const
{
Assert(this->update_flags & update_jacobian_grads,
ExcAccessToUninitializedField("update_jacobians_grads"));
- AssertIndexRange(i, this->mapping_output.jacobian_grads.size());
+ AssertIndexRange(q_point, this->mapping_output.jacobian_grads.size());
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.jacobian_grads[i];
+ return this->mapping_output.jacobian_grads[q_point];
}
template <int dim, int spacedim>
inline const DerivativeForm<1, spacedim, dim> &
-FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int i) const
+FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int q_point) const
{
Assert(this->update_flags & update_inverse_jacobians,
ExcAccessToUninitializedField("update_inverse_jacobians"));
- AssertIndexRange(i, this->mapping_output.inverse_jacobians.size());
+ AssertIndexRange(q_point, this->mapping_output.inverse_jacobians.size());
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.inverse_jacobians[i];
+ return this->mapping_output.inverse_jacobians[q_point];
}
template <int dim, int spacedim>
inline const Tensor<1, spacedim> &
-FEValuesBase<dim, spacedim>::normal_vector(const unsigned int i) const
+FEValuesBase<dim, spacedim>::normal_vector(const unsigned int q_point) const
{
Assert(this->update_flags & update_normal_vectors,
(typename FEValuesBase<dim, spacedim>::ExcAccessToUninitializedField(
"update_normal_vectors")));
- AssertIndexRange(i, this->mapping_output.normal_vectors.size());
+ AssertIndexRange(q_point, this->mapping_output.normal_vectors.size());
Assert(present_cell.is_initialized(), ExcNotReinited());
- return this->mapping_output.normal_vectors[i];
+ return this->mapping_output.normal_vectors[q_point];
}
template <int dim, int spacedim>
inline const Tensor<1, spacedim> &
-FEFaceValuesBase<dim, spacedim>::boundary_form(const unsigned int i) const
+FEFaceValuesBase<dim, spacedim>::boundary_form(const unsigned int q_point) const
{
- AssertIndexRange(i, this->mapping_output.boundary_forms.size());
+ AssertIndexRange(q_point, this->mapping_output.boundary_forms.size());
Assert(this->update_flags & update_boundary_forms,
(typename FEValuesBase<dim, spacedim>::ExcAccessToUninitializedField(
"update_boundary_forms")));
- return this->mapping_output.boundary_forms[i];
+ return this->mapping_output.boundary_forms[q_point];
}
#endif // DOXYGEN