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