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