// $Id$
// Version: $Name: $
//
-// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006 by the deal.II authors
+// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
* Likewise for second order
* derivatives.
*/
- typedef std::vector<std::vector<Tensor<2,dim> > > GradGradVector;
+ typedef std::vector<std::vector<Tensor<2,dim> > > HessianVector;
/**
* Store the values of the shape
* for the layout of the data in
* this field.
*/
- GradGradVector shape_2nd_derivatives;
+ HessianVector shape_hessians;
/**
* Store an array of weights
* given shape function occupies
* in the #shape_values,
* #shape_gradients and
- * #shape_2nd_derivatives
+ * #shape_hessians
* arrays. If all shape functions
* are primitive, then this is
* the identity mapping. If, on
* and you have to walk over all (or only the non-zero) components of
* the shape function using this set of functions.
*
- * <li> get_function_values(), get_function_grads(), etc.: Compute a
+ * <li> get_function_values(), get_function_gradients(), etc.: Compute a
* finite element function or its derivative in quadrature points.
*
* <li> reinit: initialize the FEValues object for a certain cell.
* function.
*/
const Tensor<2,dim> &
+ shape_hessian (const unsigned int function_no,
+ const unsigned int point_no) const;
+
+ /**
+ * @deprecated Wrapper for shape_hessian()
+ */
+ const Tensor<2,dim> &
shape_2nd_derivative (const unsigned int function_no,
const unsigned int point_no) const;
* is scalar, then only component
* zero is allowed and the return
* value equals that of the
- * shape_2nd_derivative()
+ * shape_hessian()
* function. If the finite
* element is vector valued but
* all shape functions are
* non-zero in only one
* component), then the value
* returned by
- * shape_2nd_derivative()
+ * shape_hessian()
* equals that of this function
* for exactly one
* component. This function is
* function cannot be used.
*/
Tensor<2,dim>
+ shape_hessian_component (const unsigned int function_no,
+ const unsigned int point_no,
+ const unsigned int component) const;
+
+ /**
+ * @deprecated Wrapper for shape_hessian_component()
+ */
+ Tensor<2,dim>
shape_2nd_derivative_component (const unsigned int function_no,
const unsigned int point_no,
const unsigned int component) const;
* finite element in use is a scalar one,
* i.e. has only one vector component. If
* it is a vector-valued one, then use
- * the other get_function_grads()
+ * the other get_function_gradients()
* function.
*
* The function assumes that the
* on the unit cell).
*/
template <class InputVector>
+ void get_function_gradients (const InputVector &fe_function,
+ std::vector<Tensor<1,dim> > &gradients) const;
+
+ /**
+ * @deprecated Use
+ * get_function_gradients() instead.
+ */
+ template <class InputVector>
void get_function_grads (const InputVector &fe_function,
std::vector<Tensor<1,dim> > &gradients) const;
* on the unit cell).
*/
template <class InputVector>
+ void get_function_gradients (const InputVector &fe_function,
+ std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
+
+
+ /**
+ * @deprecated Use
+ * get_function_gradients() instead.
+ */
+ template <class InputVector>
void get_function_grads (const InputVector &fe_function,
std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
* corresponding arguments.
*/
template <class InputVector>
+ void get_function_gradients (const InputVector& fe_function,
+ const VectorSlice<const std::vector<unsigned int> >& indices,
+ std::vector<Tensor<1,dim> >& gradients) const;
+
+ /**
+ * @deprecated Use
+ * get_function_gradients() instead.
+ */
+ template <class InputVector>
void get_function_grads (const InputVector& fe_function,
const VectorSlice<const std::vector<unsigned int> >& indices,
std::vector<Tensor<1,dim> >& gradients) const;
* corresponding arguments.
*/
template <class InputVector>
+ void get_function_gradients (const InputVector& fe_function,
+ const VectorSlice<const std::vector<unsigned int> >& indices,
+ std::vector<std::vector<Tensor<1,dim> > >& gradients,
+ bool quadrature_points_fastest = false) const;
+
+ /**
+ * @deprecated Use
+ * get_function_gradients() instead.
+ */
+ template <class InputVector>
void get_function_grads (const InputVector& fe_function,
const VectorSlice<const std::vector<unsigned int> >& indices,
std::vector<std::vector<Tensor<1,dim> > >& gradients,
* points.
*
* The function assumes that the
- * @p second_derivatives object
+ * @p hessians object
* already has the correct size.
*
* This function may only be used if the
* i.e. has only one vector component. If
* it is a vector-valued one, then use
* the other
- * get_function_2nd_derivatives()
+ * get_function_hessians()
* function.
*
* The actual data type of the input
*/
template <class InputVector>
void
- get_function_2nd_derivatives (const InputVector& fe_function,
- std::vector<Tensor<2,dim> >& second_derivatives) const;
+ get_function_hessians (const InputVector& fe_function,
+ std::vector<Tensor<2,dim> >& hessians) const;
/**
* @p cell at the quadrature points.
*
* The function assumes that the
- * @p second_derivatives object already has
+ * @p hessians object already has
* the right size.
*
* This function does the same as
*/
template <class InputVector>
void
- get_function_2nd_derivatives (const InputVector &fe_function,
- std::vector<std::vector<Tensor<2,dim> > > &second_derivatives,
- bool quadrature_points_fastest = false) const;
+ get_function_hessians (const InputVector &fe_function,
+ std::vector<std::vector<Tensor<2,dim> > > &hessians,
+ bool quadrature_points_fastest = false) const;
+
+ /**
+ * @deprecated Wrapper for get_function_hessians()
+ */
+ template <class InputVector>
+ void
+ get_function_2nd_derivatives (const InputVector&,
+ std::vector<Tensor<2,dim> >&) const;
+
+ /**
+ * @deprecated Wrapper for get_function_hessians()
+ */
+ template <class InputVector>
+ void
+ get_function_2nd_derivatives (const InputVector&,
+ std::vector<std::vector<Tensor<2,dim> > >&,
+ bool = false) const;
+
//@}
/**
* about degrees of
* freedom. These functions are,
* above all, the
- * <tt>get_function_value/grads/2nd_derivatives</tt>
+ * <tt>get_function_value/gradients/hessians</tt>
* functions. If you want to call
* these functions, you have to
* call the @p reinit variants
* information about degrees of
* freedom. These functions are,
* above all, the
- * <tt>get_function_value/grads/2nd_derivatives</tt>
+ * <tt>get_function_value/gradients/hessians</tt>
* functions. If you want to call
* these functions, you have to
* call the @p reinit variants
* information about degrees of
* freedom. These functions are,
* above all, the
- * <tt>get_function_value/grads/2nd_derivatives</tt>
+ * <tt>get_function_value/gradients/hessians</tt>
* functions. If you want to call
* these functions, you have to
* call the @p reinit variants
template <int dim>
inline
const Tensor<2,dim> &
-FEValuesBase<dim>::shape_2nd_derivative (const unsigned int i,
+FEValuesBase<dim>::shape_hessian (const unsigned int i,
const unsigned int j) const
{
- Assert (this->update_flags & update_second_derivatives,
+ Assert (this->update_flags & update_hessians,
ExcAccessToUninitializedField());
Assert (fe->is_primitive (i),
ExcShapeFunctionNotPrimitive(i));
- Assert (i<this->shape_2nd_derivatives.size(),
- ExcIndexRange (i, 0, this->shape_2nd_derivatives.size()));
- Assert (j<this->shape_2nd_derivatives[0].size(),
- ExcIndexRange (j, 0, this->shape_2nd_derivatives[0].size()));
+ Assert (i<this->shape_hessians.size(),
+ ExcIndexRange (i, 0, this->shape_hessians.size()));
+ Assert (j<this->shape_hessians[0].size(),
+ ExcIndexRange (j, 0, this->shape_hessians[0].size()));
// if the entire FE is primitive,
// then we can take a short-cut:
if (fe->is_primitive())
- return this->shape_2nd_derivatives[i][j];
+ return this->shape_hessians[i][j];
else
// otherwise, use the mapping
// between shape function numbers
// question to which vector
// component the call of this
// function refers
- return this->shape_2nd_derivatives[this->shape_function_to_row_table[i]][j];
+ return this->shape_hessians[this->shape_function_to_row_table[i]][j];
}
+template <int dim>
+inline
+const Tensor<2,dim> &
+FEValuesBase<dim>::shape_2nd_derivative (const unsigned int i,
+ const unsigned int j) const
+{
+ return shape_hessian(i,j);
+}
+
+
template <int dim>
inline
Tensor<2,dim>
-FEValuesBase<dim>::shape_2nd_derivative_component (const unsigned int i,
+FEValuesBase<dim>::shape_hessian_component (const unsigned int i,
const unsigned int j,
const unsigned int component) const
{
- Assert (this->update_flags & update_second_derivatives,
+ Assert (this->update_flags & update_hessians,
ExcAccessToUninitializedField());
Assert (component < fe->n_components(),
ExcIndexRange(component, 0, fe->n_components()));
if (fe->is_primitive(i))
{
if (component == fe->system_to_component_index(i).first)
- return this->shape_2nd_derivatives[this->shape_function_to_row_table[i]][j];
+ return this->shape_hessians[this->shape_function_to_row_table[i]][j];
else
return Tensor<2,dim>();
}
std::count (fe->get_nonzero_components(i).begin(),
fe->get_nonzero_components(i).begin()+component,
true));
- return this->shape_2nd_derivatives[row][j];
+ return this->shape_hessians[row][j];
};
}
+template <int dim>
+inline
+Tensor<2,dim>
+FEValuesBase<dim>::shape_2nd_derivative_component (const unsigned int i,
+ const unsigned int j,
+ const unsigned int component) const
+{
+ return shape_hessian_component(i,j,component);
+}
+
template <int dim>
inline
}
+template <int dim>
+template <class InputVector>
+void
+FEValuesBase<dim>::get_function_grads (const InputVector &fe_function,
+ std::vector<Tensor<1,dim> > &gradients) const
+{
+ get_function_gradients(fe_function, gradients);
+}
+
+
+template <int dim>
+template <class InputVector>
+void FEValuesBase<dim>::get_function_grads (
+ const InputVector& fe_function,
+ const VectorSlice<const std::vector<unsigned int> >& indices,
+ std::vector<Tensor<1,dim> > &values) const
+{
+ get_function_gradients(fe_function, indices, values);
+}
+
+
+template <int dim>
+template <class InputVector>
+void
+FEValuesBase<dim>::
+get_function_grads (const InputVector &fe_function,
+ std::vector<std::vector<Tensor<1,dim> > > &gradients) const
+{
+ get_function_gradients(fe_function, gradients);
+}
+
+
+template <int dim>
+template <class InputVector>
+void FEValuesBase<dim>::get_function_grads (
+ const InputVector& fe_function,
+ const VectorSlice<const std::vector<unsigned int> >& indices,
+ std::vector<std::vector<Tensor<1,dim> > >& values,
+ bool q_points_fastest) const
+{
+ get_function_gradients(fe_function, indices, values, q_points_fastest);
+}
+
+
+
+template <int dim>
+template <class InputVector>
+inline void
+FEValuesBase<dim>::
+get_function_2nd_derivatives (const InputVector &fe_function,
+ std::vector<Tensor<2,dim> > &hessians) const
+{
+ get_function_hessians(fe_function, hessians);
+}
+
+
+
+template <int dim>
+template <class InputVector>
+void
+FEValuesBase<dim>::
+get_function_2nd_derivatives (const InputVector &fe_function,
+ std::vector<std::vector<Tensor<2,dim> > > &hessians,
+ bool quadrature_points_fastest) const
+{
+ get_function_hessians(fe_function, hessians, quadrature_points_fastest);
+}
+
+
+
/*------------------------ Inline functions: FEValues ----------------------------*/
= ("You have previously called the FEValues::reinit function with a\n"
"cell iterator of type Triangulation<dim>::cell_iterator. However,\n"
"when you do this, you cannot call some functions in the FEValues\n"
- "class, such as the get_function_values/grads/2nd_derivatives\n"
+ "class, such as the get_function_values/gradients/hessians\n"
"functions. If you need these functions, then you need to call\n"
"FEValues::reinit with an iterator type that allows to extract\n"
"degrees of freedom, such as DoFHandler<dim>::cell_iterator.");
this->shape_gradients.resize (n_nonzero_shape_components,
std::vector<Tensor<1,dim> > (n_quadrature_points));
- if (flags & update_second_derivatives)
- this->shape_2nd_derivatives.resize (n_nonzero_shape_components,
+ if (flags & update_hessians)
+ this->shape_hessians.resize (n_nonzero_shape_components,
std::vector<Tensor<2,dim> > (n_quadrature_points));
if (flags & update_q_points)
template <class InputVector>
void
FEValuesBase<dim>::
-get_function_grads (const InputVector &fe_function,
+get_function_gradients (const InputVector &fe_function,
std::vector<Tensor<1,dim> > &gradients) const
{
Assert (this->update_flags & update_gradients, ExcAccessToUninitializedField());
template <int dim>
template <class InputVector>
-void FEValuesBase<dim>::get_function_grads (
+void FEValuesBase<dim>::get_function_gradients (
const InputVector& fe_function,
const VectorSlice<const std::vector<unsigned int> >& indices,
std::vector<Tensor<1,dim> > &values) const
template <class InputVector>
void
FEValuesBase<dim>::
-get_function_grads (const InputVector &fe_function,
+get_function_gradients (const InputVector &fe_function,
std::vector<std::vector<Tensor<1,dim> > > &gradients) const
{
Assert (gradients.size() == n_quadrature_points,
template <int dim>
template <class InputVector>
-void FEValuesBase<dim>::get_function_grads (
+void FEValuesBase<dim>::get_function_gradients (
const InputVector& fe_function,
const VectorSlice<const std::vector<unsigned int> >& indices,
std::vector<std::vector<Tensor<1,dim> > >& values,
template <class InputVector>
void
FEValuesBase<dim>::
-get_function_2nd_derivatives (const InputVector &fe_function,
- std::vector<Tensor<2,dim> > &second_derivatives) const
+get_function_hessians (const InputVector &fe_function,
+ std::vector<Tensor<2,dim> > &hessians) const
{
Assert (fe->n_components() == 1,
ExcDimensionMismatch(fe->n_components(), 1));
- Assert (second_derivatives.size() == n_quadrature_points,
- ExcDimensionMismatch(second_derivatives.size(), n_quadrature_points));
- Assert (this->update_flags & update_second_derivatives, ExcAccessToUninitializedField());
+ Assert (hessians.size() == n_quadrature_points,
+ ExcDimensionMismatch(hessians.size(), n_quadrature_points));
+ Assert (this->update_flags & update_hessians, ExcAccessToUninitializedField());
Assert (present_cell.get() != 0,
ExcMessage ("FEValues object is not reinit'ed to any cell"));
Assert (fe_function.size() == present_cell->n_dofs_for_dof_handler(),
present_cell->get_interpolated_dof_values(fe_function, dof_values);
// initialize with zero
- std::fill_n (second_derivatives.begin(), n_quadrature_points, Tensor<2,dim>());
+ std::fill_n (hessians.begin(), n_quadrature_points, Tensor<2,dim>());
// add up contributions of trial
// functions. note that here we
for (unsigned int point=0; point<n_quadrature_points; ++point)
for (unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
{
- Tensor<2,dim> tmp = this->shape_2nd_derivative(shape_func,point);
+ Tensor<2,dim> tmp = this->shape_hessian(shape_func,point);
tmp *= dof_values(shape_func);
- second_derivatives[point] += tmp;
+ hessians[point] += tmp;
};
}
template <class InputVector>
void
FEValuesBase<dim>::
-get_function_2nd_derivatives (const InputVector &fe_function,
+get_function_hessians (const InputVector &fe_function,
std::vector<std::vector<Tensor<2,dim> > > &second_derivs,
bool quadrature_points_fastest) const
{
Assert (second_derivs[i].size() == n_components,
ExcDimensionMismatch(second_derivs[i].size(), n_components));
- Assert (this->update_flags & update_second_derivatives, ExcAccessToUninitializedField());
+ Assert (this->update_flags & update_hessians, ExcAccessToUninitializedField());
Assert (present_cell.get() != 0,
ExcMessage ("FEValues object is not reinit'ed to any cell"));
Assert (fe_function.size() == present_cell->n_dofs_for_dof_handler(),
for (unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
if (fe->is_primitive(shape_func))
{
- Tensor<2,dim> tmp(shape_2nd_derivative(shape_func,point));
+ Tensor<2,dim> tmp(shape_hessian(shape_func,point));
tmp *= dof_values(shape_func);
second_derivs[fe->system_to_component_index(shape_func).first][point]
+= tmp;
else
for (unsigned int c=0; c<n_components; ++c)
{
- Tensor<2,dim> tmp = this->shape_2nd_derivative_component(shape_func,point,c);
+ Tensor<2,dim> tmp = this->shape_hessian_component(shape_func,point,c);
tmp *= dof_values(shape_func);
second_derivs[c][point] += tmp;
}
for (unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
if (fe->is_primitive(shape_func))
{
- Tensor<2,dim> tmp(shape_2nd_derivative(shape_func,point));
+ Tensor<2,dim> tmp(shape_hessian(shape_func,point));
tmp *= dof_values(shape_func);
second_derivs[point][fe->system_to_component_index(shape_func).first]
+= tmp;
else
for (unsigned int c=0; c<n_components; ++c)
{
- Tensor<2,dim> tmp = this->shape_2nd_derivative_component(shape_func,point,c);
+ Tensor<2,dim> tmp = this->shape_hessian_component(shape_func,point,c);
tmp *= dof_values(shape_func);
second_derivs[point][c] += tmp;
}
{
return (MemoryConsumption::memory_consumption (this->shape_values) +
MemoryConsumption::memory_consumption (this->shape_gradients) +
- MemoryConsumption::memory_consumption (this->shape_2nd_derivatives) +
+ MemoryConsumption::memory_consumption (this->shape_hessians) +
MemoryConsumption::memory_consumption (this->JxW_values) +
MemoryConsumption::memory_consumption (this->quadrature_points) +
MemoryConsumption::memory_consumption (this->normal_vectors) +