Reference documentation for deal.II version 8.5.1
auto_derivative_function.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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/auto_derivative_function.h>
18 #include <deal.II/lac/vector.h>
19 
20 #include <cmath>
21 
22 DEAL_II_NAMESPACE_OPEN
23 
24 template <int dim>
26 AutoDerivativeFunction (const double hh,
27  const unsigned int n_components,
28  const double initial_time)
29  :
30  Function<dim>(n_components, initial_time),
31  h(1),
32  ht(dim),
33  formula(Euler)
34 {
35  set_h(hh);
36  set_formula();
37 }
38 
39 
40 template <int dim>
42 {}
43 
44 
45 
46 template <int dim>
47 void
49 {
50  // go through all known formulas, reject ones we don't know about
51  // and don't handle in the member functions of this class
52  switch (form)
53  {
54  case Euler:
55  case UpwindEuler:
56  case FourthOrder:
57  break;
58  default:
59  Assert(false, ExcMessage("The argument passed to this function does not "
60  "match any known difference formula."));
61  }
62 
63  formula = form;
64 }
65 
66 
67 template <int dim>
68 void
70 {
71  h=hh;
72  for (unsigned int i=0; i<dim; ++i)
73  ht[i][i]=h;
74 }
75 
76 
77 template <int dim>
80  const unsigned int comp) const
81 {
82  Tensor<1,dim> grad;
83  switch (formula)
84  {
85  case UpwindEuler:
86  {
87  Point<dim> q1;
88  for (unsigned int i=0; i<dim; ++i)
89  {
90  q1=p-ht[i];
91  grad[i]=(this->value(p, comp)-this->value(q1, comp))/h;
92  }
93  break;
94  }
95  case Euler:
96  {
97  Point<dim> q1, q2;
98  for (unsigned int i=0; i<dim; ++i)
99  {
100  q1=p+ht[i];
101  q2=p-ht[i];
102  grad[i]=(this->value(q1, comp)-this->value(q2, comp))/(2*h);
103  }
104  break;
105  }
106  case FourthOrder:
107  {
108  Point<dim> q1, q2, q3, q4;
109  for (unsigned int i=0; i<dim; ++i)
110  {
111  q2=p+ht[i];
112  q1=q2+ht[i];
113  q3=p-ht[i];
114  q4=q3-ht[i];
115  grad[i]=(- this->value(q1, comp)
116  +8*this->value(q2, comp)
117  -8*this->value(q3, comp)
118  + this->value(q4, comp))/(12*h);
119 
120  }
121  break;
122  }
123  default:
124  Assert(false, ExcNotImplemented());
125  }
126  return grad;
127 }
128 
129 
130 template <int dim>
131 void
134  std::vector<Tensor<1,dim> > &gradients) const
135 {
136  Assert (gradients.size() == this->n_components,
137  ExcDimensionMismatch(gradients.size(), this->n_components));
138 
139  switch (formula)
140  {
141  case UpwindEuler:
142  {
143  Point<dim> q1;
144  Vector<double> v(this->n_components), v1(this->n_components);
145  const double h_inv=1./h;
146  for (unsigned int i=0; i<dim; ++i)
147  {
148  q1=p-ht[i];
149  this->vector_value(p, v);
150  this->vector_value(q1, v1);
151 
152  for (unsigned int comp=0; comp<this->n_components; ++comp)
153  gradients[comp][i]=(v(comp)-v1(comp))*h_inv;
154  }
155  break;
156  }
157 
158  case Euler:
159  {
160  Point<dim> q1, q2;
161  Vector<double> v1(this->n_components), v2(this->n_components);
162  const double h_inv_2=1./(2*h);
163  for (unsigned int i=0; i<dim; ++i)
164  {
165  q1=p+ht[i];
166  q2=p-ht[i];
167  this->vector_value(q1, v1);
168  this->vector_value(q2, v2);
169 
170  for (unsigned int comp=0; comp<this->n_components; ++comp)
171  gradients[comp][i]=(v1(comp)-v2(comp))*h_inv_2;
172  }
173  break;
174  }
175 
176  case FourthOrder:
177  {
178  Point<dim> q1, q2, q3, q4;
180  v1(this->n_components), v2(this->n_components),
181  v3(this->n_components), v4(this->n_components);
182  const double h_inv_12=1./(12*h);
183  for (unsigned int i=0; i<dim; ++i)
184  {
185  q2=p+ht[i];
186  q1=q2+ht[i];
187  q3=p-ht[i];
188  q4=q3-ht[i];
189  this->vector_value(q1, v1);
190  this->vector_value(q2, v2);
191  this->vector_value(q3, v3);
192  this->vector_value(q4, v4);
193 
194  for (unsigned int comp=0; comp<this->n_components; ++comp)
195  gradients[comp][i]=(-v1(comp)+8*v2(comp)-8*v3(comp)+v4(comp))*h_inv_12;
196  }
197  break;
198  }
199 
200  default:
201  Assert(false, ExcNotImplemented());
202  }
203 }
204 
205 
206 template <int dim>
207 void
209 gradient_list (const std::vector<Point<dim> > &points,
210  std::vector<Tensor<1,dim> > &gradients,
211  const unsigned int comp) const
212 {
213  Assert (gradients.size() == points.size(),
214  ExcDimensionMismatch(gradients.size(), points.size()));
215 
216  switch (formula)
217  {
218  case UpwindEuler:
219  {
220  Point<dim> q1;
221  for (unsigned int p=0; p<points.size(); ++p)
222  for (unsigned int i=0; i<dim; ++i)
223  {
224  q1=points[p]-ht[i];
225  gradients[p][i]=(this->value(points[p], comp)-this->value(q1, comp))/h;
226  }
227  break;
228  }
229 
230  case Euler:
231  {
232  Point<dim> q1, q2;
233  for (unsigned int p=0; p<points.size(); ++p)
234  for (unsigned int i=0; i<dim; ++i)
235  {
236  q1=points[p]+ht[i];
237  q2=points[p]-ht[i];
238  gradients[p][i]=(this->value(q1, comp)-this->value(q2, comp))/(2*h);
239  }
240  break;
241  }
242 
243  case FourthOrder:
244  {
245  Point<dim> q1, q2, q3, q4;
246  for (unsigned int p=0; p<points.size(); ++p)
247  for (unsigned int i=0; i<dim; ++i)
248  {
249  q2=points[p]+ht[i];
250  q1=q2+ht[i];
251  q3=points[p]-ht[i];
252  q4=q3-ht[i];
253  gradients[p][i]=(- this->value(q1, comp)
254  +8*this->value(q2, comp)
255  -8*this->value(q3, comp)
256  + this->value(q4, comp))/(12*h);
257  }
258  break;
259  }
260 
261  default:
262  Assert(false, ExcNotImplemented());
263  }
264 }
265 
266 
267 
268 template <int dim>
269 void
271 vector_gradient_list (const std::vector<Point<dim> > &points,
272  std::vector<std::vector<Tensor<1,dim> > > &gradients) const
273 {
274  Assert (gradients.size() == points.size(),
275  ExcDimensionMismatch(gradients.size(), points.size()));
276  for (unsigned int p=0; p<points.size(); ++p)
277  Assert (gradients[p].size() == this->n_components,
278  ExcDimensionMismatch(gradients.size(), this->n_components));
279 
280  switch (formula)
281  {
282  case UpwindEuler:
283  {
284  Point<dim> q1;
285  for (unsigned int p=0; p<points.size(); ++p)
286  for (unsigned int i=0; i<dim; ++i)
287  {
288  q1=points[p]-ht[i];
289  for (unsigned int comp=0; comp<this->n_components; ++comp)
290  gradients[p][comp][i]=(this->value(points[p], comp)-this->value(q1, comp))/h;
291  }
292  break;
293  }
294 
295  case Euler:
296  {
297  Point<dim> q1, q2;
298  for (unsigned int p=0; p<points.size(); ++p)
299  for (unsigned int i=0; i<dim; ++i)
300  {
301  q1=points[p]+ht[i];
302  q2=points[p]-ht[i];
303  for (unsigned int comp=0; comp<this->n_components; ++comp)
304  gradients[p][comp][i]=(this->value(q1, comp) -
305  this->value(q2, comp))/(2*h);
306  }
307  break;
308  }
309 
310  case FourthOrder:
311  {
312  Point<dim> q1, q2, q3, q4;
313  for (unsigned int p=0; p<points.size(); ++p)
314  for (unsigned int i=0; i<dim; ++i)
315  {
316  q2=points[p]+ht[i];
317  q1=q2+ht[i];
318  q3=points[p]-ht[i];
319  q4=q3-ht[i];
320  for (unsigned int comp=0; comp<this->n_components; ++comp)
321  gradients[p][comp][i]=(- this->value(q1, comp)
322  +8*this->value(q2, comp)
323  -8*this->value(q3, comp)
324  + this->value(q4, comp))/(12*h);
325  }
326  break;
327  }
328 
329  default:
330  Assert(false, ExcNotImplemented());
331  }
332 }
333 
334 
335 template <int dim>
338 {
339  switch (ord)
340  {
341  case 0:
342  case 1:
343  return UpwindEuler;
344  case 2:
345  return Euler;
346  case 3:
347  case 4:
348  return FourthOrder;
349  default:
350  Assert(false, ExcNotImplemented());
351  }
352  return Euler;
353 }
354 
355 
356 template class AutoDerivativeFunction<1>;
357 template class AutoDerivativeFunction<2>;
358 template class AutoDerivativeFunction<3>;
359 
360 DEAL_II_NAMESPACE_CLOSE
virtual void vector_gradient(const Point< dim > &p, std::vector< Tensor< 1, dim > > &gradients) const
static DifferenceFormula get_formula_of_order(const unsigned int ord)
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const
void set_formula(const DifferenceFormula formula=Euler)
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)
virtual void vector_gradient_list(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const
AutoDerivativeFunction(const double h, const unsigned int n_components=1, const double initial_time=0.0)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)