Reference documentation for deal.II version 8.5.1
function_derivative.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/point.h>
17 #include <deal.II/base/function_derivative.h>
18 #include <deal.II/lac/vector.h>
19 
20 #include <cmath>
21 
22 DEAL_II_NAMESPACE_OPEN
23 
24 template <int dim>
26  const Point<dim> &dir,
27  const double h)
28  :
29  AutoDerivativeFunction<dim> (h, f.n_components, f.get_time()),
30  f(f),
31  h(h),
32  incr(1, h *dir)
33 {
34  set_formula();
35 }
36 
37 
38 
39 template <int dim>
41  const std::vector<Point<dim> > &dir,
42  const double h)
43  :
44  AutoDerivativeFunction<dim> (h, f.n_components, f.get_time()),
45  f(f),
46  h(h),
47  incr(dir.size())
48 {
49  for (unsigned int i=0; i<incr.size (); ++i)
50  incr[i] = h*dir[i];
51  set_formula();
52 }
53 
54 
55 
56 template <int dim>
57 void
59 {
60  // go through all known formulas, reject ones we don't know about
61  // and don't handle in the member functions of this class
62  switch (form)
63  {
67  break;
68  default:
69  Assert(false, ExcMessage("The argument passed to this function does not "
70  "match any known difference formula."));
71  }
72 
73  formula = form;
74 }
75 
76 
77 
78 template <int dim>
79 void
80 FunctionDerivative<dim>::set_h (const double new_h)
81 {
82  for (unsigned int i=0; i<incr.size (); ++i)
83  incr[i] *= new_h/h;
84  h = new_h;
85 }
86 
87 
88 
89 template <int dim>
90 double
92  const unsigned int component) const
93 {
94  Assert (incr.size() == 1,
95  ExcMessage ("FunctionDerivative was not initialized for constant direction"));
96 
97  switch (formula)
98  {
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);
106  default:
107  Assert(false, ExcNotImplemented());
108  }
109  return 0.;
110 }
111 
112 
113 
114 template <int dim>
115 void
117  const Point<dim> &p,
118  Vector<double> &result) const
119 {
120  Assert (incr.size() == 1,
121  ExcMessage ("FunctionDerivative was not initialized for constant direction"));
122  Vector<double> aux(result.size());
123 
124  // Formulas are the same as in
125  // value, but here we have to use
126  // Vector arithmetic
127  switch (formula)
128  {
130  f.vector_value(p+incr[0], result);
131  f.vector_value(p-incr[0], aux);
132  result.sadd(1./(2*h), -1./(2*h), aux);
133  return;
135  f.vector_value(p, result);
136  f.vector_value(p-incr[0], aux);
137  result.sadd(1./h, -1./h, aux);
138  return;
140  f.vector_value(p-2*incr[0], result);
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);
146  result.add(8., aux);
147  result/=(12.*h);
148  return;
149  default:
150  Assert(false, ExcNotImplemented());
151  }
152 }
153 
154 
155 
156 template <int dim>
157 void
158 FunctionDerivative<dim>::value_list (const std::vector<Point<dim> > &points,
159  std::vector<double> &values,
160  const unsigned int component) const
161 {
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(),
166  ExcDimensionMismatch(incr.size(), points.size()));
167 
168  // Vector of auxiliary values
169  std::vector<double> aux(n);
170  // Vector of auxiliary points
171  std::vector<Point<dim> > paux(n);
172  // Use the same formulas as in
173  // value, but with vector
174  // arithmetic
175  switch (formula)
176  {
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);
186  return;
188  f.value_list(points, values, component);
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;
194  return;
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)
203  values[i] -= aux[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);
214  return;
215  default:
216  Assert(false, ExcNotImplemented());
217  }
218 }
219 
220 
221 
222 template <int dim>
223 std::size_t
225 {
226  // only simple data elements, so
227  // use sizeof operator
228  return sizeof (*this);
229 }
230 
231 
232 
233 // explicit instantiations
234 template class FunctionDerivative<1>;
235 template class FunctionDerivative<2>;
236 template class FunctionDerivative<3>;
237 
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)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::size_t size() const
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)