16 #include <deal.II/base/point.h> 17 #include <deal.II/base/function_derivative.h> 18 #include <deal.II/lac/vector.h> 22 DEAL_II_NAMESPACE_OPEN
49 for (
unsigned int i=0; i<
incr.size (); ++i)
70 "match any known difference formula."));
82 for (
unsigned int i=0; i<incr.size (); ++i)
92 const unsigned int component)
const 95 ExcMessage (
"FunctionDerivative was not initialized for constant direction"));
100 return (f.value(p+incr[0], component)-f.value(p-incr[0], component))/(2*h);
102 return (f.value(p, component)-f.value(p-incr[0], component))/h;
104 return (-f.value(p+2*incr[0], component) + 8*f.value(p+incr[0], component)
105 -8*f.value(p-incr[0], component) + f.value(p-2*incr[0], component))/(12*h);
121 ExcMessage (
"FunctionDerivative was not initialized for constant direction"));
131 f.vector_value(p-incr[0], aux);
132 result.
sadd(1./(2*h), -1./(2*h), aux);
136 f.vector_value(p-incr[0], aux);
137 result.
sadd(1./h, -1./h, aux);
141 f.vector_value(p+2*incr[0], aux);
142 result.
add(-1., aux);
143 f.vector_value(p-incr[0], aux);
144 result.
add(-8., aux);
145 f.vector_value(p+incr[0], aux);
159 std::vector<double> &values,
160 const unsigned int component)
const 162 const unsigned int n = points.size();
163 const bool variable_direction = (incr.size() == 1) ?
false :
true;
164 if (variable_direction)
165 Assert (incr.size() == points.size(),
169 std::vector<double> aux(n);
171 std::vector<Point<dim> > paux(n);
178 for (
unsigned int i=0; i<n; ++i)
179 paux[i] = points[i]+incr[(variable_direction) ? i : 0U];
180 f.value_list(paux, values, component);
181 for (
unsigned int i=0; i<n; ++i)
182 paux[i] = points[i]-incr[(variable_direction) ? i : 0U];
183 f.value_list(paux, aux, component);
184 for (
unsigned int i=0; i<n; ++i)
185 values[i] = (values[i]-aux[i])/(2*h);
189 for (
unsigned int i=0; i<n; ++i)
190 paux[i] = points[i]-incr[(variable_direction) ? i : 0U];
191 f.value_list(paux, aux, component);
192 for (
unsigned int i=0; i<n; ++i)
193 values[i] = (values[i]-aux[i])/h;
196 for (
unsigned int i=0; i<n; ++i)
197 paux[i] = points[i]-2*incr[(variable_direction) ? i : 0U];
198 f.value_list(paux, values, component);
199 for (
unsigned int i=0; i<n; ++i)
200 paux[i] = points[i]+2*incr[(variable_direction) ? i : 0U];
201 f.value_list(paux, aux, component);
202 for (
unsigned int i=0; i<n; ++i)
204 for (
unsigned int i=0; i<n; ++i)
205 paux[i] = points[i]+incr[(variable_direction) ? i : 0U];
206 f.value_list(paux, aux, component);
207 for (
unsigned int i=0; i<n; ++i)
208 values[i] += 8.*aux[i];
209 for (
unsigned int i=0; i<n; ++i)
210 paux[i] = points[i]-incr[(variable_direction) ? i : 0U];
211 f.value_list(paux, aux, component);
212 for (
unsigned int i=0; i<n; ++i)
213 values[i] = (values[i] - 8.*aux[i])/(12*h);
228 return sizeof (*this);
238 DEAL_II_NAMESPACE_CLOSE
void sadd(const Number s, const Vector< Number > &V)
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
void set_formula(typename AutoDerivativeFunction< dim >::DifferenceFormula formula=AutoDerivativeFunction< dim >::Euler)
void add(const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::size_t memory_consumption() const
FunctionDerivative(const Function< dim > &f, const Point< dim > &direction, const double h=1.e-6)
virtual void vector_value(const Point< dim > &p, Vector< double > &value) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< Number > &values, const unsigned int component=0) const
static ::ExceptionBase & ExcNotImplemented()
std::vector< Tensor< 1, dim > > incr
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual void vector_value(const Point< dim > &p, Vector< Number > &values) const
void set_h(const double h)