Reference documentation for deal.II version 8.5.1
tria_boundary.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 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/tensor.h>
17 #include <deal.II/grid/tria_boundary.h>
18 #include <deal.II/grid/tria.h>
19 #include <deal.II/grid/tria_iterator.h>
20 #include <deal.II/grid/tria_accessor.h>
21 #include <deal.II/fe/fe_q.h>
22 #include <cmath>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 
27 
28 /* -------------------------- Boundary --------------------- */
29 
30 
31 template <int dim, int spacedim>
33 {}
34 
35 
36 template <int dim, int spacedim>
37 void
39 get_intermediate_points_on_line (const typename Triangulation<dim, spacedim>::line_iterator &,
40  std::vector<Point<spacedim> > &) const
41 {
42  Assert (false, ExcPureFunctionCalled());
43 }
44 
45 
46 
47 template <int dim, int spacedim>
48 void
50 get_intermediate_points_on_quad (const typename Triangulation<dim, spacedim>::quad_iterator &,
51  std::vector<Point<spacedim> > &) const
52 {
53  Assert (false, ExcPureFunctionCalled());
54 }
55 
56 
57 template <int dim, int spacedim>
58 void
61  std::vector<Point<spacedim> > &points) const
62 {
63  Assert (dim>1, ExcImpossibleInDim(dim));
64 
65  switch (dim)
66  {
67  case 2:
68  get_intermediate_points_on_line (face, points);
69  break;
70  case 3:
71  get_intermediate_points_on_quad (face, points);
72  break;
73  default:
74  Assert (false, ExcNotImplemented());
75  }
76 }
77 
78 
79 template <>
80 void
83  std::vector<Point<1> > &) const
84 {
85  Assert (false, ExcImpossibleInDim(1));
86 }
87 
88 
89 template <>
90 void
93  std::vector<Point<2> > &) const
94 {
95  Assert (false, ExcImpossibleInDim(1));
96 }
97 
98 
99 template <>
100 void
103  std::vector<Point<3> > &) const
104 {
105  Assert (false, ExcImpossibleInDim(1));
106 }
107 
108 
109 template <int dim, int spacedim>
112 project_to_surface (const typename Triangulation<dim, spacedim>::line_iterator &,
113  const Point<spacedim> &trial_point) const
114 {
115  if (spacedim <= 1)
116  return trial_point;
117  else
118  {
119  Assert (false, ExcPureFunctionCalled());
120  return Point<spacedim>();
121  }
122 }
123 
124 
125 
126 template <int dim, int spacedim>
129 project_to_surface (const typename Triangulation<dim, spacedim>::quad_iterator &,
130  const Point<spacedim> &trial_point) const
131 {
132  if (spacedim <= 2)
133  return trial_point;
134  else
135  {
136  Assert (false, ExcPureFunctionCalled());
137  return Point<spacedim>();
138  }
139 }
140 
141 
142 
143 template <int dim, int spacedim>
146 project_to_surface (const typename Triangulation<dim, spacedim>::hex_iterator &,
147  const Point<spacedim> &trial_point) const
148 {
149  if (spacedim <= 3)
150  return trial_point;
151  else
152  {
153  Assert (false, ExcPureFunctionCalled());
154  return Point<spacedim>();
155  }
156 }
157 
158 
159 
160 template <int dim, int spacedim>
161 const std::vector<Point<1> > &
163 get_line_support_points (const unsigned int n_intermediate_points) const
164 {
165  if (points.size() <= n_intermediate_points ||
166  points[n_intermediate_points].get() == 0)
167  {
168  Threads::Mutex::ScopedLock lock(mutex);
169  if (points.size() <= n_intermediate_points)
170  points.resize(n_intermediate_points+1);
171 
172  // another thread might have created points in the meantime
173  if (points[n_intermediate_points].get() == 0)
174  {
175  std_cxx11::shared_ptr<QGaussLobatto<1> >
176  quadrature (new QGaussLobatto<1>(n_intermediate_points+2));
177  points[n_intermediate_points] = quadrature;
178  }
179  }
180  return points[n_intermediate_points]->get_points();
181 }
182 
183 
184 
185 
186 /* -------------------------- StraightBoundary --------------------- */
187 
188 
189 template <int dim, int spacedim>
191 {}
192 
193 
194 template <int dim, int spacedim>
197 get_new_point_on_line (const typename Triangulation<dim, spacedim>::line_iterator &line) const
198 {
199  return (line->vertex(0) + line->vertex(1)) / 2;
200 }
201 
202 
203 namespace
204 {
205  // compute the new midpoint of a quad --
206  // either of a 2d cell on a manifold in 3d
207  // or of a face of a 3d triangulation in 3d
208  template <int dim>
209  Point<3>
210  compute_new_point_on_quad (const typename Triangulation<dim, 3>::quad_iterator &quad)
211  {
212  // generate a new point in the middle of
213  // the face based on the points on the
214  // edges and the vertices.
215  //
216  // there is a pathological situation when
217  // this face is on a straight boundary, but
218  // one of its edges and the face behind it
219  // are not; if that face is refined first,
220  // the new point in the middle of that edge
221  // may not be at the same position as
222  // quad->line(.)->center() would have been,
223  // but would have been moved to the
224  // non-straight boundary. We cater to that
225  // situation by using existing edge
226  // midpoints if available, or center() if
227  // not
228  //
229  // note that this situation can not happen
230  // during mesh refinement, as there the
231  // edges are refined first and only then
232  // the face. thus, the check whether a line
233  // has children does not lead to the
234  // situation where the new face midpoints
235  // have different positions depending on
236  // which of the two cells is refined first.
237  //
238  // the situation where the edges aren't
239  // refined happens when a higher order
240  // MappingQ requests the midpoint of a
241  // face, though, and it is for these cases
242  // that we need to have the check available
243  //
244  // note that the factor of 1/8 for each
245  // of the 8 surrounding points isn't
246  // chosen arbitrarily. rather, we may ask
247  // where the harmonic map would place the
248  // point (0,0) if we map the square
249  // [-1,1]^2 onto the domain that is
250  // described using the 4 vertices and 4
251  // edge point points of this quad. we can
252  // then discretize the harmonic map using
253  // four cells and Q1 elements on each of
254  // the quadrants of the square [-1,1]^2
255  // and see where the midpoint would land
256  // (this is the procedure we choose, for
257  // example, in
258  // GridGenerator::laplace_solve) and it
259  // turns out that it will land at the
260  // mean of the 8 surrounding
261  // points. whether a discretization of
262  // the harmonic map with only 4 cells is
263  // adequate is a different question
264  // altogether, of course.
265  return (quad->vertex(0) + quad->vertex(1) +
266  quad->vertex(2) + quad->vertex(3) +
267  (quad->line(0)->has_children() ?
268  quad->line(0)->child(0)->vertex(1) :
269  quad->line(0)->center()) +
270  (quad->line(1)->has_children() ?
271  quad->line(1)->child(0)->vertex(1) :
272  quad->line(1)->center()) +
273  (quad->line(2)->has_children() ?
274  quad->line(2)->child(0)->vertex(1) :
275  quad->line(2)->center()) +
276  (quad->line(3)->has_children() ?
277  quad->line(3)->child(0)->vertex(1) :
278  quad->line(3)->center()) ) / 8;
279  }
280 }
281 
282 
283 
284 template <int dim, int spacedim>
287 get_new_point_on_quad (const typename Triangulation<dim, spacedim>::quad_iterator &quad) const
288 {
290 }
291 
292 
293 template <>
294 Point<3>
296 get_new_point_on_quad (const Triangulation<2,3>::quad_iterator &quad) const
297 {
298  return compute_new_point_on_quad<2> (quad);
299 }
300 
301 
302 
303 template <>
304 Point<3>
306 get_new_point_on_quad (const Triangulation<3>::quad_iterator &quad) const
307 {
308  return compute_new_point_on_quad<3> (quad);
309 }
310 
311 
312 
313 template <int dim, int spacedim>
314 void
316 get_intermediate_points_on_line (const typename Triangulation<dim, spacedim>::line_iterator &line,
317  std::vector<Point<spacedim> > &points) const
318 {
319  const unsigned int n=points.size();
320  Assert(n>0, ExcInternalError());
321 
322  // Use interior points of QGaussLobatto quadrature formula support points
323  // for consistency with MappingQ
324  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
325 
326  const Point<spacedim> vertices[2] = { line->vertex(0),
327  line->vertex(1)
328  };
329 
330  for (unsigned int i=0; i<n; ++i)
331  {
332  const double x = line_points[1+i][0];
333  points[i] = (1-x)*vertices[0] + x*vertices[1];
334  }
335 }
336 
337 
338 
339 
340 template <int dim, int spacedim>
341 void
343 get_intermediate_points_on_quad (const typename Triangulation<dim, spacedim>::quad_iterator &,
344  std::vector<Point<spacedim> > &) const
345 {
346  Assert(false, ExcImpossibleInDim(dim));
347 }
348 
349 
350 
351 template <>
352 void
354 get_intermediate_points_on_quad (const Triangulation<3>::quad_iterator &quad,
355  std::vector<Point<3> > &points) const
356 {
357  const unsigned int spacedim = 3;
358 
359  const unsigned int n=points.size(),
360  m=static_cast<unsigned int>(std::sqrt(static_cast<double>(n)));
361  // is n a square number
362  Assert(m*m==n, ExcInternalError());
363 
364  const std::vector<Point<1> > &line_points = this->get_line_support_points(m);
365 
366  const Point<spacedim> vertices[4] = { quad->vertex(0),
367  quad->vertex(1),
368  quad->vertex(2),
369  quad->vertex(3)
370  };
371 
372  for (unsigned int i=0; i<m; ++i)
373  {
374  const double y=line_points[1+i][0];
375  for (unsigned int j=0; j<m; ++j)
376  {
377  const double x=line_points[1+j][0];
378  points[i*m+j]=((1-x) * vertices[0] +
379  x * vertices[1]) * (1-y) +
380  ((1-x) * vertices[2] +
381  x * vertices[3]) * y;
382  }
383  }
384 }
385 
386 
387 
388 template <>
389 void
391 get_intermediate_points_on_quad (const Triangulation<2,3>::quad_iterator &quad,
392  std::vector<Point<3> > &points) const
393 {
394  const unsigned int spacedim = 3;
395 
396  const unsigned int n=points.size(),
397  m=static_cast<unsigned int>(std::sqrt(static_cast<double>(n)));
398  // is n a square number
399  Assert(m*m==n, ExcInternalError());
400 
401  const std::vector<Point<1> > &line_points = this->get_line_support_points(m);
402 
403  const Point<spacedim> vertices[4] = { quad->vertex(0),
404  quad->vertex(1),
405  quad->vertex(2),
406  quad->vertex(3)
407  };
408 
409  for (unsigned int i=0; i<m; ++i)
410  {
411  const double y=line_points[1+i][0];
412  for (unsigned int j=0; j<m; ++j)
413  {
414  const double x=line_points[1+j][0];
415  points[i*m+j]=((1-x) * vertices[0] +
416  x * vertices[1]) * (1-y) +
417  ((1-x) * vertices[2] +
418  x * vertices[3]) * y;
419  }
420  }
421 }
422 
423 
424 
425 template <>
429  const Point<1> &) const
430 {
431  Assert (false, ExcNotImplemented());
432  return Tensor<1,1>();
433 }
434 
435 
436 template <>
440  const Point<2> &) const
441 {
442  Assert (false, ExcNotImplemented());
443  return Tensor<1,2>();
444 }
445 
446 
447 template <>
451  const Point<3> &) const
452 {
453  Assert (false, ExcNotImplemented());
454  return Tensor<1,3>();
455 }
456 
457 
458 namespace internal
459 {
460  namespace
461  {
467  normalized_alternating_product (const Tensor<1,2> (&basis_vectors)[1])
468  {
469  Tensor<1,2> tmp = cross_product_2d (basis_vectors[0]);
470  return tmp/tmp.norm();
471  }
472 
473 
474 
476  normalized_alternating_product (const Tensor<1,3> ( &)[1])
477  {
478  // we get here from StraightBoundary<2,3>::normal_vector, but
479  // the implementation below is bogus for this case anyway
480  // (see the assert at the beginning of that function).
481  Assert (false, ExcNotImplemented());
482  return Tensor<1,3>();
483  }
484 
485 
486 
488  normalized_alternating_product (const Tensor<1,3> (&basis_vectors)[2])
489  {
490  Tensor<1,3> tmp = cross_product_3d (basis_vectors[0], basis_vectors[1]);
491  return tmp/tmp.norm();
492  }
493 
494  }
495 }
496 
497 
498 template <int dim, int spacedim>
502  const Point<spacedim> &p) const
503 {
504  // I don't think the implementation below will work when dim!=spacedim;
505  // in fact, I believe that we don't even have enough information here,
506  // because we would need to know not only about the tangent vectors
507  // of the face, but also of the cell, to compute the normal vector.
508  // Someone will have to think about this some more.
509  Assert (dim == spacedim, ExcNotImplemented());
510 
511  // in order to find out what the normal vector is, we first need to
512  // find the reference coordinates of the point p on the given face,
513  // or at least the reference coordinates of the closest point on the
514  // face
515  //
516  // in other words, we need to find a point xi so that f(xi)=||F(xi)-p||^2->min
517  // where F(xi) is the mapping. this algorithm is implemented in
518  // MappingQ1<dim,spacedim>::transform_real_to_unit_cell but only for cells,
519  // while we need it for faces here. it's also implemented in somewhat
520  // more generality there using the machinery of the MappingQ1 class
521  // while we really only need it for a specific case here
522  //
523  // in any case, the iteration we use here is a Gauss-Newton's iteration with
524  // xi^{n+1} = xi^n - H(xi^n)^{-1} J(xi^n)
525  // where
526  // J(xi) = (grad F(xi))^T (F(xi)-p)
527  // and
528  // H(xi) = [grad F(xi)]^T [grad F(xi)]
529  // In all this,
530  // F(xi) = sum_v vertex[v] phi_v(xi)
531  // We get the shape functions phi_v from an object of type FE_Q<dim-1>(1)
532 
533  // we start with the point xi=1/2, xi=(1/2,1/2), ...
534  const unsigned int facedim = dim-1;
535 
536  Point<facedim> xi;
537  for (unsigned int i=0; i<facedim; ++i)
538  xi[i] = 1./2;
539 
540  FE_Q<facedim> linear_fe(1);
541 
542  const double eps = 1e-12;
543  Tensor<1,spacedim> grad_F[facedim];
544  unsigned int iteration = 0;
545  while (true)
546  {
547  Point<spacedim> F;
548  for (unsigned int v=0; v<GeometryInfo<facedim>::vertices_per_cell; ++v)
549  F += face->vertex(v) * linear_fe.shape_value(v, xi);
550 
551  for (unsigned int i=0; i<facedim; ++i)
552  {
553  grad_F[i] = 0;
554  for (unsigned int v=0; v<GeometryInfo<facedim>::vertices_per_cell; ++v)
555  grad_F[i] += face->vertex(v) * linear_fe.shape_grad(v, xi)[i];
556  }
557 
559  for (unsigned int i=0; i<facedim; ++i)
560  for (unsigned int j=0; j<spacedim; ++j)
561  J[i] += grad_F[i][j] * (F-p)[j];
562 
564  for (unsigned int i=0; i<facedim; ++i)
565  for (unsigned int j=0; j<facedim; ++j)
566  for (unsigned int k=0; k<spacedim; ++k)
567  H[i][j] += grad_F[i][k] * grad_F[j][k];
568 
569  const Tensor<1,facedim> delta_xi = -invert(H) * J;
570  xi += delta_xi;
571  ++iteration;
572 
573  Assert (iteration<10,
574  ExcMessage("The Newton iteration to find the reference point "
575  "did not converge in 10 iterations. Do you have a "
576  "deformed cell? (See the glossary for a definition "
577  "of what a deformed cell is. You may want to output "
578  "the vertices of your cell."));
579 
580  if (delta_xi.norm() < eps)
581  break;
582  }
583 
584  // so now we have the reference coordinates xi of the point p.
585  // we then have to compute the normal vector, which we can do
586  // by taking the (normalize) alternating product of all the tangent
587  // vectors given by grad_F
588  return internal::normalized_alternating_product(grad_F);
589 }
590 
591 
592 
593 template <>
594 void
598 {
599  Assert (false, ExcImpossibleInDim(1));
600 }
601 
602 template <>
603 void
607 {
608  Assert (false, ExcNotImplemented());
609 }
610 
611 
612 template <>
613 void
617 {
618  Assert (false, ExcNotImplemented());
619 }
620 
621 
622 
623 template <>
624 void
627  Boundary<2,2>::FaceVertexNormals &face_vertex_normals) const
628 {
629  const Tensor<1,2> tangent = face->vertex(1) - face->vertex(0);
630  for (unsigned int vertex=0; vertex<GeometryInfo<2>::vertices_per_face; ++vertex)
631  // compute normals from tangent
632  face_vertex_normals[vertex] = Point<2>(tangent[1],
633  -tangent[0]);
634 }
635 
636 template <>
637 void
640  Boundary<2,3>::FaceVertexNormals &face_vertex_normals) const
641 {
642  const Tensor<1,3> tangent = face->vertex(1) - face->vertex(0);
643  for (unsigned int vertex=0; vertex<GeometryInfo<2>::vertices_per_face; ++vertex)
644  // compute normals from tangent
645  face_vertex_normals[vertex] = Point<3>(tangent[1],
646  -tangent[0],0);
647  Assert(false, ExcNotImplemented());
648 }
649 
650 
651 
652 
653 template <>
654 void
657  Boundary<3,3>::FaceVertexNormals &face_vertex_normals) const
658 {
659  const unsigned int vertices_per_face = GeometryInfo<3>::vertices_per_face;
660 
661  static const unsigned int neighboring_vertices[4][2]=
662  { {1,2},{3,0},{0,3},{2,1}};
663  for (unsigned int vertex=0; vertex<vertices_per_face; ++vertex)
664  {
665  // first define the two tangent vectors at the vertex by using the
666  // two lines radiating away from this vertex
667  const Tensor<1,3> tangents[2]
668  = { face->vertex(neighboring_vertices[vertex][0])
669  - face->vertex(vertex),
670  face->vertex(neighboring_vertices[vertex][1])
671  - face->vertex(vertex)
672  };
673 
674  // then compute the normal by taking the cross product. since the
675  // normal is not required to be normalized, no problem here
676  face_vertex_normals[vertex] = cross_product_3d(tangents[0], tangents[1]);
677  };
678 }
679 
680 
681 
682 template <int dim, int spacedim>
685 project_to_surface (const typename Triangulation<dim, spacedim>::line_iterator &line,
686  const Point<spacedim> &trial_point) const
687 {
688  if (spacedim <= 1)
689  return trial_point;
690  else
691  {
692  // find the point that lies on
693  // the line p1--p2. the
694  // formulas pan out to
695  // something rather simple
696  // because the mapping to the
697  // line is linear
698  const Point<spacedim> p1 = line->vertex(0),
699  p2 = line->vertex(1);
700  const double s = (trial_point-p1)*(p2-p1) / ((p2-p1)*(p2-p1));
701  return p1 + s*(p2-p1);
702  }
703 }
704 
705 
706 
707 namespace internal
708 {
709  template <typename Iterator, int spacedim, int dim>
711  compute_projection (const Iterator &object,
712  const Point<spacedim> &y,
714  {
715  // let's look at this for
716  // simplicity for a quad (dim==2)
717  // in a space with spacedim>2:
718 
719  // all points on the surface are given by
720  // x(\xi) = sum_i v_i phi_x(\xi)
721  // where v_i are the vertices of the quad,
722  // and \xi=(\xi_1,\xi_2) are the reference
723  // coordinates of the quad. so what we are
724  // trying to do is find a point x on
725  // the surface that is closest to the point
726  // y. there are different ways
727  // to solve this problem, but in the end
728  // it's a nonlinear problem and we have to
729  // find reference coordinates \xi so that
730  // J(\xi) = 1/2 || x(\xi)-y ||^2
731  // is minimal. x(\xi) is a function that
732  // is dim-linear in \xi, so J(\xi) is
733  // a polynomial of degree 2*dim that
734  // we'd like to minimize. unless dim==1,
735  // we'll have to use a Newton
736  // method to find the
737  // answer. This leads to the
738  // following formulation of
739  // Newton steps:
740  //
741  // Given \xi_k, find \delta\xi_k so that
742  // H_k \delta\xi_k = - F_k
743  // where H_k is an approximation to the
744  // second derivatives of J at \xi_k, and
745  // F_k is the first derivative of J.
746  // We'll iterate this a number of times
747  // until the right hand side is small
748  // enough. As a stopping criterion, we
749  // terminate if ||\delta\xi||<eps.
750  //
751  // As for the Hessian, the best choice
752  // would be
753  // H_k = J''(\xi_k)
754  // but we'll opt for the simpler
755  // Gauss-Newton form
756  // H_k = A^T A
757  // i.e.
758  // (H_k)_{nm} = \sum_{i,j} v_i*v_j *
759  // \partial_n phi_i *
760  // \partial_m phi_j
761  // we start at xi=(0.5,0.5).
762  Point<dim> xi;
763  for (unsigned int d=0; d<dim; ++d)
764  xi[d] = 0.5;
765 
766  Point<spacedim> x_k;
767  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
768  x_k += object->vertex(i) *
770 
771  do
772  {
773  Tensor<1,dim> F_k;
774  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
775  F_k += (x_k-y)*object->vertex(i) *
777 
778  Tensor<2,dim> H_k;
779  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
780  for (unsigned int j=0; j<GeometryInfo<dim>::vertices_per_cell; ++j)
781  {
782  Tensor<2, dim> tmp = outer_product(
785  H_k += (object->vertex(i) * object->vertex(j)) * tmp;
786  }
787 
788  const Tensor<1,dim> delta_xi = - invert(H_k) * F_k;
789  xi += delta_xi;
790 
791  x_k = Point<spacedim>();
792  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
793  x_k += object->vertex(i) *
795 
796  if (delta_xi.norm() < 1e-5)
797  break;
798  }
799  while (true);
800 
801  return x_k;
802  }
803 
804 
805  // specialization for a quad in 1d
806  template <typename Iterator>
807  Point<1>
808  compute_projection (const Iterator &,
809  const Point<1> &y,
810  /* it's a quad: */internal::int2type<2>)
811  {
812  return y;
813  }
814 
815  // specialization for a quad in 2d
816  template <typename Iterator>
817  Point<2>
818  compute_projection (const Iterator &,
819  const Point<2> &y,
820  /* it's a quad: */internal::int2type<2>)
821  {
822  return y;
823  }
824 }
825 
826 
827 
828 
829 
830 template <>
831 Point<3>
833 project_to_surface (const Triangulation<1, 3>::quad_iterator &,
834  const Point<3> &y) const
835 {
836  return y;
837 }
838 
839 //TODO[SP]: This is just a horrible way out to make it compile in codim 2.
840 template <int dim, int spacedim>
843 project_to_surface (const typename Triangulation<dim, spacedim>::quad_iterator &quad,
844  const Point<spacedim> &y) const
845 {
846  if (spacedim <= 2)
847  return y;
848  else
849  return internal::compute_projection (quad, y,
850  /* it's a quad */internal::int2type<2>());
851 }
852 
853 
854 
855 template <int dim, int spacedim>
858 project_to_surface (const typename Triangulation<dim, spacedim>::hex_iterator &,
859  const Point<spacedim> &trial_point) const
860 {
861  if (spacedim <= 3)
862  return trial_point;
863  else
864  {
865  // we can presumably call the
866  // same function as above (it's
867  // written in a generic way)
868  // but someone needs to check
869  // whether that actually yields
870  // the correct result
871  Assert (false, ExcNotImplemented());
872  return Point<spacedim>();
873  }
874 }
875 
876 
877 
878 // explicit instantiations
879 #include "tria_boundary.inst"
880 
881 DEAL_II_NAMESPACE_CLOSE
882 
const std::vector< Point< 1 > > & get_line_support_points(const unsigned int n_intermediate_points) const
static ::ExceptionBase & ExcPureFunctionCalled()
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:991
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
virtual Point< spacedim > project_to_surface(const typename Triangulation< dim, spacedim >::line_iterator &line, const Point< spacedim > &candidate) const
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const
Definition: tria.h:52
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
Definition: fe_q.h:552
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
virtual ~Boundary()
#define Assert(cond, exc)
Definition: exceptions.h:313
void get_intermediate_points_on_face(const typename Triangulation< dim, spacedim >::face_iterator &face, std::vector< Point< spacedim > > &points) const
Tensor< 1, dim, Number > cross_product_3d(const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2)
Definition: tensor.h:1698
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: mpi.h:41
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
Definition: manifold.cc:386
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const
static ::ExceptionBase & ExcNotImplemented()
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
virtual Point< spacedim > project_to_surface(const typename Triangulation< dim, spacedim >::line_iterator &line, const Point< spacedim > &candidate) const
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static ::ExceptionBase & ExcInternalError()