Reference documentation for deal.II version 8.5.1
derivative_approximation.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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/quadrature_lib.h>
17 #include <deal.II/base/work_stream.h>
18 #include <deal.II/lac/vector.h>
19 #include <deal.II/lac/block_vector.h>
20 #include <deal.II/lac/la_vector.h>
21 #include <deal.II/lac/la_parallel_vector.h>
22 #include <deal.II/lac/la_parallel_block_vector.h>
23 #include <deal.II/lac/petsc_vector.h>
24 #include <deal.II/lac/petsc_block_vector.h>
25 #include <deal.II/lac/trilinos_vector.h>
26 #include <deal.II/lac/trilinos_block_vector.h>
27 #include <deal.II/grid/tria_iterator.h>
28 #include <deal.II/grid/grid_tools.h>
29 #include <deal.II/grid/filtered_iterator.h>
30 #include <deal.II/dofs/dof_accessor.h>
31 #include <deal.II/dofs/dof_handler.h>
32 #include <deal.II/fe/fe.h>
33 #include <deal.II/fe/fe_values.h>
34 #include <deal.II/hp/fe_values.h>
35 #include <deal.II/fe/mapping_q1.h>
36 #include <deal.II/hp/q_collection.h>
37 #include <deal.II/hp/fe_collection.h>
38 #include <deal.II/hp/mapping_collection.h>
39 #include <deal.II/numerics/derivative_approximation.h>
40 
41 #include <cmath>
42 
43 DEAL_II_NAMESPACE_OPEN
44 
45 
46 
47 namespace
48 {
49  template <typename T>
50  inline T sqr (const T t)
51  {
52  return t*t;
53  }
54 }
55 
56 // --------------- First the classes and functions that describe individual
57 // --------------- derivatives
58 
60 {
61  namespace internal
62  {
71  template <int dim>
72  class Gradient
73  {
74  public:
79  static const UpdateFlags update_flags;
80 
86 
92 
99  template <class InputVector, int spacedim>
100  static ProjectedDerivative
102  const InputVector &solution,
103  const unsigned int component);
104 
109  static double derivative_norm (const Derivative &d);
110 
118  static void symmetrize (Derivative &derivative_tensor);
119  };
120 
121  // static variables
122  template <int dim>
124 
125 
126  template <int dim>
127  template <class InputVector, int spacedim>
128  inline
132  const InputVector &solution,
133  const unsigned int component)
134  {
135  if (fe_values.get_fe().n_components() == 1)
136  {
137  std::vector<typename InputVector::value_type> values (1);
138  fe_values.get_function_values (solution, values);
139  return values[0];
140  }
141  else
142  {
143  std::vector<Vector<typename InputVector::value_type> > values
144  (1, Vector<typename InputVector::value_type>(fe_values.get_fe().n_components()));
145  fe_values.get_function_values (solution, values);
146  return values[0](component);
147  }
148  }
149 
150 
151 
152  template <int dim>
153  inline
154  double
156  {
157  double s = 0;
158  for (unsigned int i=0; i<dim; ++i)
159  s += d[i]*d[i];
160  return std::sqrt(s);
161  }
162 
163 
164 
165  template <int dim>
166  inline
167  void
169  {
170  // nothing to do here
171  }
172 
173 
174 
183  template <int dim>
185  {
186  public:
191  static const UpdateFlags update_flags;
192 
198 
204 
211  template <class InputVector, int spacedim>
212  static ProjectedDerivative
214  const InputVector &solution,
215  const unsigned int component);
216 
224  static double derivative_norm (const Derivative &d);
225 
235  static void symmetrize (Derivative &derivative_tensor);
236  };
237 
238  template <int dim>
240 
241 
242  template <int dim>
243  template <class InputVector, int spacedim>
244  inline
248  const InputVector &solution,
249  const unsigned int component)
250  {
251  if (fe_values.get_fe().n_components() == 1)
252  {
253  std::vector<Tensor<1,dim,typename InputVector::value_type> > values (1);
254  fe_values.get_function_gradients (solution, values);
255  return ProjectedDerivative(values[0]);
256  }
257  else
258  {
259  std::vector<std::vector<Tensor<1,dim,typename InputVector::value_type> > > values
260  (1, std::vector<Tensor<1,dim,typename InputVector::value_type> >(fe_values.get_fe().n_components()));
261  fe_values.get_function_gradients (solution, values);
262  return ProjectedDerivative(values[0][component]);
263  };
264  }
265 
266 
267 
268  template <>
269  inline
270  double
272  derivative_norm (const Derivative &d)
273  {
274  return std::fabs (d[0][0]);
275  }
276 
277 
278 
279  template <>
280  inline
281  double
283  derivative_norm (const Derivative &d)
284  {
285  // note that d should be a
286  // symmetric 2x2 tensor, so the
287  // eigenvalues are:
288  //
289  // 1/2(a+b\pm\sqrt((a-b)^2+4c^2))
290  //
291  // if the d_11=a, d_22=b,
292  // d_12=d_21=c
293  const double radicand = ::sqr(d[0][0] - d[1][1]) +
294  4*::sqr(d[0][1]);
295  const double eigenvalues[2]
296  = { 0.5*(d[0][0] + d[1][1] + std::sqrt(radicand)),
297  0.5*(d[0][0] + d[1][1] - std::sqrt(radicand))
298  };
299 
300  return std::max (std::fabs (eigenvalues[0]),
301  std::fabs (eigenvalues[1]));
302  }
303 
304 
305 
306  template <>
307  inline
308  double
310  derivative_norm (const Derivative &d)
311  {
312  /*
313  compute the three eigenvalues of the tensor @p{d} and take the
314  largest. one could use the following maple script to generate C
315  code:
316 
317  with(linalg);
318  readlib(C);
319  A:=matrix(3,3,[[a00,a01,a02],[a01,a11,a12],[a02,a12,a22]]);
320  E:=eigenvals(A);
321  EE:=vector(3,[E[1],E[2],E[3]]);
322  C(EE);
323 
324  Unfortunately, with both optimized and non-optimized output, at some
325  places the code `sqrt(-1.0)' is emitted, and I don't know what
326  Maple intends to do with it. This happens both with Maple4 and
327  Maple5.
328 
329  Fortunately, Roger Young provided the following Fortran code, which
330  is transcribed below to C. The code uses an algorithm that uses the
331  invariants of a symmetric matrix. (The translated algorithm is
332  augmented by a test for R>0, since R==0 indicates that all three
333  eigenvalues are equal.)
334 
335 
336  PROGRAM MAIN
337 
338  C FIND EIGENVALUES OF REAL SYMMETRIC MATRIX
339  C (ROGER YOUNG, 2001)
340 
341  IMPLICIT NONE
342 
343  REAL*8 A11, A12, A13, A22, A23, A33
344  REAL*8 I1, J2, J3, AM
345  REAL*8 S11, S12, S13, S22, S23, S33
346  REAL*8 SS12, SS23, SS13
347  REAL*8 R,R3, XX,YY, THETA
348  REAL*8 A1,A2,A3
349  REAL*8 PI
350  PARAMETER (PI=3.141592653587932384D0)
351  REAL*8 A,B,C, TOL
352  PARAMETER (TOL=1.D-14)
353 
354  C DEFINE A TEST MATRIX
355 
356  A11 = -1.D0
357  A12 = 5.D0
358  A13 = 3.D0
359  A22 = -2.D0
360  A23 = 0.5D0
361  A33 = 4.D0
362 
363 
364  I1 = A11 + A22 + A33
365  AM = I1/3.D0
366 
367  S11 = A11 - AM
368  S22 = A22 - AM
369  S33 = A33 - AM
370  S12 = A12
371  S13 = A13
372  S23 = A23
373 
374  SS12 = S12*S12
375  SS23 = S23*S23
376  SS13 = S13*S13
377 
378  J2 = S11*S11 + S22*S22 + S33*S33
379  J2 = J2 + 2.D0*(SS12 + SS23 + SS13)
380  J2 = J2/2.D0
381 
382  J3 = S11**3 + S22**3 + S33**3
383  J3 = J3 + 3.D0*S11*(SS12 + SS13)
384  J3 = J3 + 3.D0*S22*(SS12 + SS23)
385  J3 = J3 + 3.D0*S33*(SS13 + SS23)
386  J3 = J3 + 6.D0*S12*S23*S13
387  J3 = J3/3.D0
388 
389  R = SQRT(4.D0*J2/3.D0)
390  R3 = R*R*R
391  XX = 4.D0*J3/R3
392 
393  YY = 1.D0 - DABS(XX)
394  IF(YY.LE.0.D0)THEN
395  IF(YY.GT.(-TOL))THEN
396  WRITE(6,*)'Equal roots: XX= ',XX
397  A = -(XX/DABS(XX))*SQRT(J2/3.D0)
398  B = AM + A
399  C = AM - 2.D0*A
400  WRITE(6,*)B,' (twice) ',C
401  STOP
402  ELSE
403  WRITE(6,*)'Error: XX= ',XX
404  STOP
405  ENDIF
406  ENDIF
407 
408  THETA = (ACOS(XX))/3.D0
409 
410  A1 = AM + R*COS(THETA)
411  A2 = AM + R*COS(THETA + 2.D0*PI/3.D0)
412  A3 = AM + R*COS(THETA + 4.D0*PI/3.D0)
413 
414  WRITE(6,*)A1,A2,A3
415 
416  STOP
417  END
418 
419  */
420 
421  const double am = trace(d) / 3.;
422 
423  // s := d - trace(d) I
424  Tensor<2,3> s = d;
425  for (unsigned int i=0; i<3; ++i)
426  s[i][i] -= am;
427 
428  const double ss01 = s[0][1] * s[0][1],
429  ss12 = s[1][2] * s[1][2],
430  ss02 = s[0][2] * s[0][2];
431 
432  const double J2 = (s[0][0]*s[0][0] + s[1][1]*s[1][1] + s[2][2]*s[2][2]
433  + 2 * (ss01 + ss02 + ss12)) / 2.;
434  const double J3 = (std::pow(s[0][0],3) + std::pow(s[1][1],3) + std::pow(s[2][2],3)
435  + 3. * s[0][0] * (ss01 + ss02)
436  + 3. * s[1][1] * (ss01 + ss12)
437  + 3. * s[2][2] * (ss02 + ss12)
438  + 6. * s[0][1] * s[0][2] * s[1][2]) / 3.;
439 
440  const double R = std::sqrt (4. * J2 / 3.);
441 
442  double EE[3] = { 0, 0, 0 };
443  // the eigenvalues are away from
444  // @p{am} in the order of R. thus,
445  // if R<<AM, then we have the
446  // degenerate case with three
447  // identical eigenvalues. check
448  // this first
449  if (R <= 1e-14*std::fabs(am))
450  EE[0] = EE[1] = EE[2] = am;
451  else
452  {
453  // at least two eigenvalues are
454  // distinct
455  const double R3 = R*R*R;
456  const double XX = 4. * J3 / R3;
457  const double YY = 1. - std::fabs(XX);
458 
459  Assert (YY > -1e-14, ExcInternalError());
460 
461  if (YY < 0)
462  {
463  // two roots are equal
464  const double a = (XX>0 ? -1. : 1.) * R / 2;
465  EE[0] = EE[1] = am + a;
466  EE[2] = am - 2.*a;
467  }
468  else
469  {
470  const double theta = std::acos(XX) / 3.;
471  EE[0] = am + R*std::cos(theta);
472  EE[1] = am + R*std::cos(theta + 2./3.*numbers::PI);
473  EE[2] = am + R*std::cos(theta + 4./3.*numbers::PI);
474  };
475  };
476 
477  return std::max (std::fabs (EE[0]),
478  std::max (std::fabs (EE[1]),
479  std::fabs (EE[2])));
480  }
481 
482 
483 
484  template <int dim>
485  inline
486  double
489  {
490  // computing the spectral norm is
491  // not so simple in general. it is
492  // feasible for dim==3 as shown
493  // above, since then there are
494  // still closed form expressions of
495  // the roots of the characteristic
496  // polynomial, and they can easily
497  // be computed using
498  // maple. however, for higher
499  // dimensions, some other method
500  // needs to be employed. maybe some
501  // steps of the power method would
502  // suffice?
503  Assert (false, ExcNotImplemented());
504  return 0;
505  }
506 
507 
508 
509  template <int dim>
510  inline
511  void
513  {
514  // symmetrize non-diagonal entries
515  for (unsigned int i=0; i<dim; ++i)
516  for (unsigned int j=i+1; j<dim; ++j)
517  {
518  const double s = (d[i][j] + d[j][i]) / 2;
519  d[i][j] = d[j][i] = s;
520  };
521  }
522 
523 
524 
525  template <int dim>
526  class ThirdDerivative
527  {
528  public:
533  static const UpdateFlags update_flags;
534 
540  typedef Tensor<3,dim> Derivative;
541 
546  typedef Tensor<2,dim> ProjectedDerivative;
547 
554  template <class InputVector, int spacedim>
555  static ProjectedDerivative
556  get_projected_derivative (const FEValues<dim,spacedim> &fe_values,
557  const InputVector &solution,
558  const unsigned int component);
559 
567  static double derivative_norm (const Derivative &d);
568 
578  static void symmetrize (Derivative &derivative_tensor);
579  };
580 
581  template <int dim>
582  const UpdateFlags ThirdDerivative<dim>::update_flags = update_hessians;
583 
584 
585  template <int dim>
586  template <class InputVector, int spacedim>
587  inline
588  typename ThirdDerivative<dim>::ProjectedDerivative
589  ThirdDerivative<dim>::
590  get_projected_derivative (const FEValues<dim,spacedim> &fe_values,
591  const InputVector &solution,
592  const unsigned int component)
593  {
594  if (fe_values.get_fe().n_components() == 1)
595  {
596  std::vector<Tensor<2,dim,typename InputVector::value_type> > values (1);
597  fe_values.get_function_hessians (solution, values);
598  return ProjectedDerivative(values[0]);
599  }
600  else
601  {
602  std::vector<std::vector<Tensor<2,dim,typename InputVector::value_type> > > values
603  (1, std::vector<Tensor<2,dim,typename InputVector::value_type> >(fe_values.get_fe().n_components()));
604  fe_values.get_function_hessians (solution, values);
605  return ProjectedDerivative(values[0][component]);
606  };
607  }
608 
609 
610 
611  template <>
612  inline
613  double
614  ThirdDerivative<1>::
615  derivative_norm (const Derivative &d)
616  {
617  return std::fabs (d[0][0][0]);
618  }
619 
620 
621 
622  template <int dim>
623  inline
624  double
625  ThirdDerivative<dim>::
626  derivative_norm (const Derivative &d)
627  {
628  // return the Frobenius-norm. this is a
629  // member function of Tensor<rank_,dim>
630  return d.norm();
631  }
632 
633 
634  template <int dim>
635  inline
636  void
637  ThirdDerivative<dim>::symmetrize (Derivative &d)
638  {
639  // symmetrize non-diagonal entries
640 
641  // first do it in the case, that i,j,k are
642  // pairwise different (which can onlky happen
643  // in dim >= 3)
644  for (unsigned int i=0; i<dim; ++i)
645  for (unsigned int j=i+1; j<dim; ++j)
646  for (unsigned int k=j+1; k<dim; ++k)
647  {
648  const double s = (d[i][j][k] +
649  d[i][k][j] +
650  d[j][i][k] +
651  d[j][k][i] +
652  d[k][i][j] +
653  d[k][j][i]) / 6;
654  d[i][j][k]
655  = d[i][k][j]
656  = d[j][i][k]
657  = d[j][k][i]
658  = d[k][i][j]
659  = d[k][j][i]
660  = s;
661  }
662  // now do the case, where two indices are
663  // equal
664  for (unsigned int i=0; i<dim; ++i)
665  for (unsigned int j=i+1; j<dim; ++j)
666  {
667  // case 1: index i (lower one) is
668  // double
669  const double s = (d[i][i][j] +
670  d[i][j][i] +
671  d[j][i][i] ) / 3;
672  d[i][i][j]
673  = d[i][j][i]
674  = d[j][i][i]
675  = s;
676 
677  // case 2: index j (higher one) is
678  // double
679  const double t = (d[i][j][j] +
680  d[j][i][j] +
681  d[j][j][i] ) / 3;
682  d[i][j][j]
683  = d[j][i][j]
684  = d[j][j][i]
685  = t;
686  }
687  }
688 
689 
690  template <int order, int dim>
691  class DerivativeSelector
692  {
693  public:
699  typedef void DerivDescr;
700 
701  };
702 
703  template <int dim>
704  class DerivativeSelector<1,dim>
705  {
706  public:
707 
708  typedef Gradient<dim> DerivDescr;
709  };
710 
711  template <int dim>
712  class DerivativeSelector<2,dim>
713  {
714  public:
715 
716  typedef SecondDerivative<dim> DerivDescr;
717  };
718 
719  template <int dim>
720  class DerivativeSelector<3,dim>
721  {
722  public:
723 
724  typedef ThirdDerivative<dim> DerivDescr;
725  };
726  }
727 }
728 
729 // Dummy structures and dummy function used for WorkStream
730 namespace DerivativeApproximation
731 {
732  namespace internal
733  {
734  namespace Assembler
735  {
736  struct Scratch
737  {
738  Scratch() {}
739  };
740 
741  struct CopyData
742  {
743  CopyData() {}
744  };
745  }
746  }
747 }
748 
749 // ------------------------------- now for the functions that do the
750 // ------------------------------- actual work
751 
752 namespace DerivativeApproximation
753 {
754  namespace internal
755  {
760  template <class DerivativeDescription, int dim, template <int, int> class DoFHandlerType,
761  class InputVector, int spacedim>
762  void
763  approximate_cell (const Mapping<dim,spacedim> &mapping,
764  const DoFHandlerType<dim,spacedim> &dof_handler,
765  const InputVector &solution,
766  const unsigned int component,
767  const TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
768  false> > &cell,
769  typename DerivativeDescription::Derivative &derivative)
770  {
771  QMidpoint<dim> midpoint_rule;
772 
773  // create collection objects from
774  // single quadratures, mappings,
775  // and finite elements. if we have
776  // an hp DoFHandler,
777  // dof_handler.get_fe() returns a
778  // collection of which we do a
779  // shallow copy instead
780  const hp::QCollection<dim> q_collection (midpoint_rule);
781  const hp::FECollection<dim> fe_collection(dof_handler.get_fe());
782  const hp::MappingCollection<dim> mapping_collection (mapping);
783 
784  hp::FEValues<dim> x_fe_midpoint_value (mapping_collection, fe_collection,
785  q_collection,
786  DerivativeDescription::update_flags |
788 
789  // matrix Y=sum_i y_i y_i^T
790  Tensor<2,dim> Y;
791 
792 
793  // vector to hold iterators to all
794  // active neighbors of a cell
795  // reserve the maximal number of
796  // active neighbors
797  std::vector<TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
798  false> > > active_neighbors;
799 
800  active_neighbors.reserve (GeometryInfo<dim>::faces_per_cell *
802 
803  // vector
804  // g=sum_i y_i (f(x+y_i)-f(x))/|y_i|
805  // or related type for higher
806  // derivatives
807  typename DerivativeDescription::Derivative projected_derivative;
808 
809  // reinit fe values object...
810  x_fe_midpoint_value.reinit (cell);
811  const FEValues<dim> &fe_midpoint_value
812  = x_fe_midpoint_value.get_present_fe_values();
813 
814  // ...and get the value of the
815  // projected derivative...
816  const typename DerivativeDescription::ProjectedDerivative
817  this_midpoint_value
818  = DerivativeDescription::get_projected_derivative (fe_midpoint_value,
819  solution,
820  component);
821  // ...and the place where it lives
822  const Point<dim> this_center = fe_midpoint_value.quadrature_point(0);
823 
824  // loop over all neighbors and
825  // accumulate the difference
826  // quotients from them. note
827  // that things get a bit more
828  // complicated if the neighbor
829  // is more refined than the
830  // present one
831  //
832  // to make processing simpler,
833  // first collect all neighbor
834  // cells in a vector, and then
835  // collect the data from them
836  GridTools::get_active_neighbors<DoFHandlerType<dim,spacedim> >(cell, active_neighbors);
837 
838  // now loop over all active
839  // neighbors and collect the
840  // data we need
841  typename std::vector<TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
842  false> > >::const_iterator
843  neighbor_ptr = active_neighbors.begin();
844  for (; neighbor_ptr!=active_neighbors.end(); ++neighbor_ptr)
845  {
847  neighbor = *neighbor_ptr;
848 
849  // reinit fe values object...
850  x_fe_midpoint_value.reinit (neighbor);
851  const FEValues<dim> &neighbor_fe_midpoint_value
852  = x_fe_midpoint_value.get_present_fe_values();
853 
854  // ...and get the value of the
855  // solution...
856  const typename DerivativeDescription::ProjectedDerivative
857  neighbor_midpoint_value
858  = DerivativeDescription::get_projected_derivative (neighbor_fe_midpoint_value,
859  solution, component);
860 
861  // ...and the place where it lives
862  const Point<dim>
863  neighbor_center = neighbor_fe_midpoint_value.quadrature_point(0);
864 
865 
866  // vector for the
867  // normalized
868  // direction between
869  // the centers of two
870  // cells
871  Tensor<1,dim> y = neighbor_center - this_center;
872  const double distance = y.norm();
873  // normalize y
874  y /= distance;
875  // *** note that unlike in
876  // the docs, y denotes the
877  // normalized vector
878  // connecting the centers
879  // of the two cells, rather
880  // than the normal
881  // difference! ***
882 
883  // add up the
884  // contribution of
885  // this cell to Y
886  for (unsigned int i=0; i<dim; ++i)
887  for (unsigned int j=0; j<dim; ++j)
888  Y[i][j] += y[i] * y[j];
889 
890  // then update the sum
891  // of difference
892  // quotients
893  typename DerivativeDescription::ProjectedDerivative
894  projected_finite_difference
895  = (neighbor_midpoint_value -
896  this_midpoint_value);
897  projected_finite_difference /= distance;
898 
899  projected_derivative += outer_product(y, projected_finite_difference);
900  };
901 
902  // can we determine an
903  // approximation of the
904  // gradient for the present
905  // cell? if so, then we need to
906  // have passed over vectors y_i
907  // which span the whole space,
908  // otherwise we would not have
909  // all components of the
910  // gradient
911  AssertThrow (determinant(Y) != 0,
913 
914  // compute Y^-1 g
915  const Tensor<2,dim> Y_inverse = invert(Y);
916 
917  derivative = Y_inverse * projected_derivative;
918 
919  // finally symmetrize the derivative
920  DerivativeDescription::symmetrize (derivative);
921  }
922 
923 
924 
930  template <class DerivativeDescription, int dim,
931  template <int, int> class DoFHandlerType, class InputVector, int spacedim>
932  void
933  approximate
934  (SynchronousIterators<std_cxx11::tuple<TriaActiveIterator < ::DoFCellAccessor < DoFHandlerType < dim, spacedim >, false > >, Vector<float>::iterator> > const &cell,
935  const Mapping<dim,spacedim> &mapping,
936  const DoFHandlerType<dim,spacedim> &dof_handler,
937  const InputVector &solution,
938  const unsigned int component)
939  {
940  // if the cell is not locally owned, then there is nothing to do
941  if (std_cxx11::get<0>(*cell)->is_locally_owned() == false)
942  *std_cxx11::get<1>(*cell) = 0;
943  else
944  {
945  typename DerivativeDescription::Derivative derivative;
946  // call the function doing the actual
947  // work on this cell
948  approximate_cell<DerivativeDescription,dim,DoFHandlerType,InputVector, spacedim>
949  (mapping,dof_handler,solution,component,std_cxx11::get<0>(*cell),derivative);
950 
951  // evaluate the norm and fill the vector
952  //*derivative_norm_on_this_cell
953  *std_cxx11::get<1>(*cell) = DerivativeDescription::derivative_norm (derivative);
954  }
955  }
956 
957 
968  template <class DerivativeDescription, int dim,
969  template <int, int> class DoFHandlerType, class InputVector, int spacedim>
970  void
971  approximate_derivative (const Mapping<dim,spacedim> &mapping,
972  const DoFHandlerType<dim,spacedim> &dof_handler,
973  const InputVector &solution,
974  const unsigned int component,
976  {
977  Assert (derivative_norm.size() == dof_handler.get_triangulation().n_active_cells(),
979  dof_handler.get_triangulation().n_active_cells()));
980  Assert (component < dof_handler.get_fe().n_components(),
981  ExcIndexRange (component, 0, dof_handler.get_fe().n_components()));
982 
983  typedef std_cxx11::tuple<TriaActiveIterator<::DoFCellAccessor
984  <DoFHandlerType<dim, spacedim>, false> >,
985  Vector<float>::iterator> Iterators;
986  SynchronousIterators<Iterators> begin(Iterators(dof_handler.begin_active(),
987  derivative_norm.begin())),
988  end(Iterators(dof_handler.end(),
989  derivative_norm.end()));
990 
991  // There is no need for a copier because there is no conflict between threads
992  // to write in derivative_norm. Scratch and CopyData are also useless.
993  WorkStream::run(begin,
994  end,
995  static_cast<std_cxx11::function<void (SynchronousIterators<Iterators> const &,
996  Assembler::Scratch const &, Assembler::CopyData &)> >
997  (std_cxx11::bind(&approximate<DerivativeDescription,dim,DoFHandlerType,
998  InputVector,spacedim>,
999  std_cxx11::_1,
1000  std_cxx11::cref(mapping),
1001  std_cxx11::cref(dof_handler),
1002  std_cxx11::cref(solution),component)),
1003  std_cxx11::function<void (internal::Assembler::CopyData const &)> (),
1004  internal::Assembler::Scratch (),internal::Assembler::CopyData ());
1005  }
1006 
1007  } // namespace internal
1008 
1009 } // namespace DerivativeApproximation
1010 
1011 
1012 // ------------------------ finally for the public interface of this namespace
1013 
1014 namespace DerivativeApproximation
1015 {
1016  template <int dim, template <int, int> class DoFHandlerType, class InputVector, int spacedim>
1017  void
1019  const DoFHandlerType<dim,spacedim> &dof_handler,
1020  const InputVector &solution,
1022  const unsigned int component)
1023  {
1024  internal::approximate_derivative<internal::Gradient<dim>,dim> (mapping,
1025  dof_handler,
1026  solution,
1027  component,
1028  derivative_norm);
1029  }
1030 
1031 
1032  template <int dim, template <int, int> class DoFHandlerType, class InputVector, int spacedim>
1033  void
1034  approximate_gradient (const DoFHandlerType<dim,spacedim> &dof_handler,
1035  const InputVector &solution,
1037  const unsigned int component)
1038  {
1039  internal::approximate_derivative<internal::Gradient<dim>,dim> (StaticMappingQ1<dim>::mapping,
1040  dof_handler,
1041  solution,
1042  component,
1043  derivative_norm);
1044  }
1045 
1046 
1047  template <int dim, template <int, int> class DoFHandlerType, class InputVector, int spacedim>
1048  void
1050  const DoFHandlerType<dim,spacedim> &dof_handler,
1051  const InputVector &solution,
1053  const unsigned int component)
1054  {
1055  internal::approximate_derivative<internal::SecondDerivative<dim>,dim> (mapping,
1056  dof_handler,
1057  solution,
1058  component,
1059  derivative_norm);
1060  }
1061 
1062 
1063  template <int dim, template <int, int> class DoFHandlerType, class InputVector, int spacedim>
1064  void
1065  approximate_second_derivative (const DoFHandlerType<dim,spacedim> &dof_handler,
1066  const InputVector &solution,
1068  const unsigned int component)
1069  {
1070  internal::approximate_derivative<internal::SecondDerivative<dim>,dim> (StaticMappingQ1<dim>::mapping,
1071  dof_handler,
1072  solution,
1073  component,
1074  derivative_norm);
1075  }
1076 
1077 
1078  template <typename DoFHandlerType, class InputVector, int order>
1079  void
1082  const DoFHandlerType &dof,
1083  const InputVector &solution,
1084 #ifndef _MSC_VER
1085  const typename DoFHandlerType::active_cell_iterator &cell,
1086 #else
1088 #endif
1090  const unsigned int component)
1091  {
1092  internal::approximate_cell<typename internal::DerivativeSelector<order,DoFHandlerType::dimension>::DerivDescr>
1093  (mapping,
1094  dof,
1095  solution,
1096  component,
1097  cell,
1098  derivative);
1099  }
1100 
1101 
1102 
1103  template <typename DoFHandlerType, class InputVector, int order>
1104  void
1106  (const DoFHandlerType &dof,
1107  const InputVector &solution,
1108 #ifndef _MSC_VER
1109  const typename DoFHandlerType::active_cell_iterator &cell,
1110 #else
1112 #endif
1114  const unsigned int component)
1115  {
1116  // just call the respective function with Q1 mapping
1119  dof,
1120  solution,
1121  cell,
1122  derivative,
1123  component);
1124  }
1125 
1126 
1127 
1128 
1129 
1130  template <int dim, int order>
1131  double
1133  {
1134  return internal::DerivativeSelector<order,dim>::DerivDescr::derivative_norm(derivative);
1135  }
1136 
1137 }
1138 
1139 
1140 // --------------------------- explicit instantiations ---------------------
1141 #include "derivative_approximation.inst"
1142 
1143 
1144 DEAL_II_NAMESPACE_CLOSE
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
Definition: fe_values.cc:2850
Shape function values.
static ::ExceptionBase & ExcInsufficientDirections()
const FEValues< dim, spacedim > & get_present_fe_values() const
static void symmetrize(Derivative &derivative_tensor)
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Definition: fe_values.cc:2689
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:991
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
Definition: fe_values.cc:2973
static const double PI
Definition: numbers.h:94
const Point< spacedim > & quadrature_point(const unsigned int q) const
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
static void symmetrize(Derivative &derivative_tensor)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
Second derivatives of shape functions.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
double derivative_norm(const Tensor< order, dim > &derivative)
void approximate_derivative_tensor(const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof, const InputVector &solution, const typename DoFHandlerType::active_cell_iterator &cell, Tensor< order, DoFHandlerType::dimension > &derivative, const unsigned int component=0)
Shape function gradients.
Definition: fe.h:30
void run(const std::vector< std::vector< Iterator > > &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition: work_stream.h:1110
static ::ExceptionBase & ExcNotImplemented()
static double derivative_norm(const Derivative &d)
void approximate_second_derivative(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
void approximate_gradient(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcVectorLengthVsNActiveCells(int arg1, int arg2)