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