Reference documentation for deal.II version 9.0.0
tria_boundary_lib.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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/grid/tria_boundary_lib.h>
17 #include <deal.II/grid/tria.h>
18 #include <deal.II/grid/tria_iterator.h>
19 #include <deal.II/grid/tria_accessor.h>
20 #include <deal.II/base/tensor.h>
21 #include <deal.II/base/std_cxx14/memory.h>
22 #include <cmath>
23 
24 
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 
29 
30 template <int dim, int spacedim>
32  const unsigned int axis)
33  :
34  radius(radius),
35  direction (get_axis_vector (axis)),
36  point_on_axis (Point<spacedim>())
37 {}
38 
39 
40 template <int dim, int spacedim>
42  const Point<spacedim> &direction,
43  const Point<spacedim> &point_on_axis)
44  :
45  radius(radius),
46  direction (direction / direction.norm()),
47  point_on_axis (point_on_axis)
48 {}
49 
50 
51 
52 template<int dim, int spacedim>
53 std::unique_ptr<Manifold<dim, spacedim> >
55 {
56  return std_cxx14::make_unique<CylinderBoundary<dim, spacedim> >(radius, direction, point_on_axis);
57 }
58 
59 
60 
61 template <int dim, int spacedim>
64 {
65  Assert (axis < spacedim, ExcIndexRange (axis, 0, spacedim));
66 
67  Point<spacedim> axis_vector;
68  axis_vector[axis] = 1;
69  return axis_vector;
70 }
71 
72 
73 
74 template <int dim, int spacedim>
77 get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const
78 {
79  // compute a proposed new point
81 
82  // we then have to project this
83  // point out to the given radius
84  // from the axis. to this end, we
85  // have to take into account the
86  // offset point_on_axis and the
87  // direction of the axis
88  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
89  ((middle-point_on_axis) * direction) * direction;
90  // scale it to the desired length
91  // and put everything back
92  // together, unless we have a point
93  // on the axis
94  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
95  return middle;
96  else
97  return Point<spacedim>(vector_from_axis / vector_from_axis.norm() * radius +
98  ((middle-point_on_axis) * direction) * direction +
99  point_on_axis);
100 }
101 
102 
103 
104 template <>
105 Point<3>
107 get_new_point_on_quad (const Triangulation<3>::quad_iterator &quad) const
108 {
110 
111  // same algorithm as above
112  const unsigned int spacedim = 3;
113 
114  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
115  ((middle-point_on_axis) * direction) * direction;
116  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
117  return middle;
118  else
119  return Point<3>(vector_from_axis / vector_from_axis.norm() * radius +
120  ((middle-point_on_axis) * direction) * direction +
121  point_on_axis);
122 }
123 
124 template <>
125 Point<3>
127 get_new_point_on_quad (const Triangulation<2,3>::quad_iterator &quad) const
128 {
130 
131  // same algorithm as above
132  const unsigned int spacedim = 3;
133  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
134  ((middle-point_on_axis) * direction) * direction;
135  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
136  return middle;
137  else
138  return Point<3>(vector_from_axis / vector_from_axis.norm() * radius +
139  ((middle-point_on_axis) * direction) * direction +
140  point_on_axis);
141 }
142 
143 
144 template <int dim, int spacedim>
147 get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &) const
148 {
149  Assert (false, ExcImpossibleInDim(dim));
150  return Point<spacedim>();
151 }
152 
153 
154 
155 template <int dim, int spacedim>
156 void
158  const typename Triangulation<dim,spacedim>::line_iterator &line,
159  std::vector<Point<spacedim> > &points) const
160 {
161  if (points.size()==1)
162  points[0]=get_new_point_on_line(line);
163  else
164  get_intermediate_points_between_points(line->vertex(0), line->vertex(1), points);
165 }
166 
167 
168 template <int dim, int spacedim>
169 void
171  const Point<spacedim> &v0,
172  const Point<spacedim> &v1,
173  std::vector<Point<spacedim> > &points) const
174 {
175  const unsigned int n=points.size();
176  Assert(n>0, ExcInternalError());
177 
178  // Do a simple linear interpolation followed by projection, using the same
179  // algorithm as above
180  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
181 
182  for (unsigned int i=0; i<n; ++i)
183  {
184  const double x = line_points[i+1][0];
185  const Point<spacedim> middle = (1-x)*v0 + x*v1;
186 
187  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
188  ((middle-point_on_axis) * direction) * direction;
189  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
190  points[i] = middle;
191  else
192  points[i] = Point<spacedim>(vector_from_axis / vector_from_axis.norm() * radius +
193  ((middle-point_on_axis) * direction) * direction +
194  point_on_axis);
195  }
196 }
197 
198 
199 
200 template <>
201 void
203  const Triangulation<3>::quad_iterator &quad,
204  std::vector<Point<3> > &points) const
205 {
206  if (points.size()==1)
207  points[0]=get_new_point_on_quad(quad);
208  else
209  {
210  unsigned int m=static_cast<unsigned int> (std::sqrt(static_cast<double>(points.size())));
211  Assert(points.size()==m*m, ExcInternalError());
212 
213  std::vector<Point<3> > lp0(m);
214  std::vector<Point<3> > lp1(m);
215 
216  get_intermediate_points_on_line(quad->line(0), lp0);
217  get_intermediate_points_on_line(quad->line(1), lp1);
218 
219  std::vector<Point<3> > lps(m);
220  for (unsigned int i=0; i<m; ++i)
221  {
222  get_intermediate_points_between_points(lp0[i], lp1[i], lps);
223 
224  for (unsigned int j=0; j<m; ++j)
225  points[i*m+j]=lps[j];
226  }
227  }
228 }
229 
230 
231 
232 template <int dim, int spacedim>
233 void
235  const typename Triangulation<dim,spacedim>::quad_iterator &,
236  std::vector<Point<spacedim> > &) const
237 {
238  Assert (false, ExcImpossibleInDim(dim));
239 }
240 
241 
242 
243 
244 template <>
245 void
249 {
250  Assert (false, ExcImpossibleInDim(1));
251 }
252 
253 
254 
255 
256 template <int dim, int spacedim>
257 void
260  typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const
261 {
262  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
263  {
264  const Point<spacedim> vertex = face->vertex(v);
265 
266  const Tensor<1,spacedim> vector_from_axis = (vertex-point_on_axis) -
267  ((vertex-point_on_axis) * direction) * direction;
268 
269  face_vertex_normals[v] = (vector_from_axis / vector_from_axis.norm());
270  }
271 }
272 
273 
274 
275 template <int dim, int spacedim>
276 double
278 {
279  return radius;
280 }
281 
282 
283 //======================================================================//
284 
285 template <int dim>
286 ConeBoundary<dim>::ConeBoundary (const double radius_0,
287  const double radius_1,
288  const Point<dim> x_0,
289  const Point<dim> x_1)
290  :
291  radius_0 (radius_0),
292  radius_1 (radius_1),
293  x_0 (x_0),
294  x_1 (x_1)
295 {}
296 
297 
298 
299 template<int dim>
300 std::unique_ptr<Manifold<dim,dim> > ConeBoundary<dim>::clone() const
301 {
302  return std_cxx14::make_unique<ConeBoundary<dim> >(radius_0,radius_1,x_0,x_1);
303 }
304 
305 
306 
307 template <int dim>
309 {
310  for (unsigned int i = 0; i < dim; ++i)
311  if ((x_1 (i) - x_0 (i)) != 0)
312  return (radius_1 - radius_0) * (x (i) - x_0 (i)) / (x_1 (i) - x_0 (i)) + radius_0;
313 
314  return 0;
315 }
316 
317 
318 
319 template <int dim>
320 void
323  const Point<dim> &p1,
324  std::vector<Point<dim> > &points) const
325 {
326  const unsigned int n = points.size ();
327  const Tensor<1,dim> axis = x_1 - x_0;
328 
329  Assert (n > 0, ExcInternalError ());
330 
331  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
332 
333  for (unsigned int i=0; i<n; ++i)
334  {
335  const double x = line_points[i+1][0];
336 
337  // Compute the current point.
338  const Point<dim> x_i = (1-x)*p0 + x*p1;
339  // To project this point on the boundary of the cone we first compute
340  // the orthogonal projection of this point onto the axis of the cone.
341  const double c = (x_i - x_0) * axis / (axis*axis);
342  const Point<dim> x_ip = x_0 + c * axis;
343  // Compute the projection of the middle point on the boundary of the
344  // cone.
345  points[i] = x_ip + get_radius (x_ip) * (x_i - x_ip) / (x_i - x_ip).norm ();
346  }
347 }
348 
349 template <int dim>
353 {
354  const Tensor<1,dim> axis = x_1 - x_0;
355  // Compute the middle point of the line.
357  // To project it on the boundary of the cone we first compute the orthogonal
358  // projection of the middle point onto the axis of the cone.
359  const double c = (middle - x_0) * axis / (axis*axis);
360  const Point<dim> middle_p = x_0 + c * axis;
361  // Compute the projection of the middle point on the boundary of the cone.
362  return middle_p + get_radius (middle_p) * (middle - middle_p) / (middle - middle_p).norm ();
363 }
364 
365 
366 
367 template <>
368 Point<3>
370 get_new_point_on_quad (const Triangulation<3>::quad_iterator &quad) const
371 {
372  const int dim = 3;
373 
374  const Tensor<1,dim> axis = x_1 - x_0;
375  // Compute the middle point of the quad.
377  // Same algorithm as above: To project it on the boundary of the cone we
378  // first compute the orthogonal projection of the middle point onto the axis
379  // of the cone.
380  const double c = (middle - x_0) * axis / (axis*axis);
381  const Point<dim> middle_p = x_0 + c * axis;
382  // Compute the projection of the middle point on the boundary of the cone.
383  return middle_p + get_radius (middle_p) * (middle - middle_p) / (middle - middle_p).norm ();
384 }
385 
386 
387 
388 template <int dim>
392 {
393  Assert (false, ExcImpossibleInDim (dim));
394 
395  return Point<dim>();
396 }
397 
398 
399 
400 template <int dim>
401 void
404  std::vector<Point<dim> > &points) const
405 {
406  if (points.size () == 1)
407  points[0] = get_new_point_on_line (line);
408  else
409  get_intermediate_points_between_points (line->vertex (0), line->vertex (1), points);
410 }
411 
412 
413 
414 
415 template <>
416 void
418 get_intermediate_points_on_quad (const Triangulation<3>::quad_iterator &quad,
419  std::vector<Point<3> > &points) const
420 {
421  if (points.size () == 1)
422  points[0] = get_new_point_on_quad (quad);
423  else
424  {
425  unsigned int n = static_cast<unsigned int> (std::sqrt (static_cast<double> (points.size ())));
426 
427  Assert (points.size () == n * n, ExcInternalError ());
428 
429  std::vector<Point<3> > points_line_0 (n);
430  std::vector<Point<3> > points_line_1 (n);
431 
432  get_intermediate_points_on_line (quad->line (0), points_line_0);
433  get_intermediate_points_on_line (quad->line (1), points_line_1);
434 
435  std::vector<Point<3> > points_line_segment (n);
436 
437  for (unsigned int i = 0; i < n; ++i)
438  {
439  get_intermediate_points_between_points (points_line_0[i],
440  points_line_1[i],
441  points_line_segment);
442 
443  for (unsigned int j = 0; j < n; ++j)
444  points[i * n + j] = points_line_segment[j];
445  }
446  }
447 }
448 
449 
450 
451 template <int dim>
452 void
455  std::vector<Point<dim> > &) const
456 {
457  Assert (false, ExcImpossibleInDim (dim));
458 }
459 
460 
461 
462 
463 template <>
464 void
468 {
469  Assert (false, ExcImpossibleInDim (1));
470 }
471 
472 
473 
474 template <int dim>
475 void
478  typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const
479 {
480  const Tensor<1,dim> axis = x_1 - x_0;
481 
482  for (unsigned int vertex = 0; vertex < GeometryInfo<dim>::vertices_per_face; ++vertex)
483  {
484  // Compute the orthogonal projection of the vertex onto the axis of the
485  // cone.
486  const double c = (face->vertex (vertex) - x_0) * axis / (axis*axis);
487  const Point<dim> vertex_p = x_0 + c * axis;
488  // Then compute the vector pointing from the point <tt>vertex_p</tt> on
489  // the axis to the vertex.
490  const Tensor<1,dim> axis_to_vertex = face->vertex (vertex) - vertex_p;
491 
492  face_vertex_normals[vertex] = axis_to_vertex / axis_to_vertex.norm ();
493  }
494 }
495 
496 
497 template <int dim>
498 inline
502  const Point<dim> &p) const
503 {
504 // TODO only for cone opening along z-axis
505  Assert (dim == 3, ExcNotImplemented());
506  Assert (this->radius_0 == 0., ExcNotImplemented());
507  Assert (this->x_0 == Point<dim>(), ExcNotImplemented());
508  for (unsigned int d=0; d<dim; ++d)
509  if (d != dim-1) // don't test the last component of the vector
510  Assert (this->x_1[d] == 0., ExcNotImplemented());
511  Assert (this->x_1[dim-1] > 0, ExcNotImplemented());
512 
513  const double c_squared = (this->radius_1 / this->x_1[dim-1])*(this->radius_1 / this->x_1[dim-1]);
514  Tensor<1,dim> normal = p;
515  normal[0] *= -2.0/c_squared;
516  normal[1] *= -2.0/c_squared;
517  normal[dim-1] *= 2.0;
518 
519  return normal/normal.norm();
520 }
521 
522 
523 
524 //======================================================================//
525 
526 template <int dim, int spacedim>
528  const double radius)
529  :
530  center(p),
531  radius(radius),
532  compute_radius_automatically(false)
533 {}
534 
535 
536 
537 template<int dim, int spacedim>
538 std::unique_ptr<Manifold<dim, spacedim> >
540 {
541  return std_cxx14::make_unique<HyperBallBoundary<dim,spacedim> >(center, radius);
542 }
543 
544 
545 
546 template <int dim, int spacedim>
548 HyperBallBoundary<dim,spacedim>::get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const
549 {
551 
552  middle -= center;
553 
554  double r=0;
555  if (compute_radius_automatically)
556  r = (line->vertex(0) - center).norm();
557  else
558  r = radius;
559 
560  // project to boundary
561  middle *= r / std::sqrt(middle.square());
562  middle += center;
563  return middle;
564 }
565 
566 
567 
568 template <>
569 Point<1>
571 get_new_point_on_quad (const Triangulation<1,1>::quad_iterator &) const
572 {
573  Assert (false, ExcInternalError());
574  return Point<1>();
575 }
576 
577 
578 template <>
579 Point<2>
581 get_new_point_on_quad (const Triangulation<1,2>::quad_iterator &) const
582 {
583  Assert (false, ExcInternalError());
584  return Point<2>();
585 }
586 
587 
588 
589 template <int dim, int spacedim>
592 get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const
593 {
595 
596  middle -= center;
597 
598  double r=0;
599  if (compute_radius_automatically)
600  r = (quad->vertex(0) - center).norm();
601  else
602  r = radius;
603 
604  // project to boundary
605  middle *= r / std::sqrt(middle.square());
606 
607  middle += center;
608  return middle;
609 }
610 
611 
612 
613 template <>
614 void
616  const Triangulation<1>::line_iterator &,
617  std::vector<Point<1> > &) const
618 {
619  Assert (false, ExcImpossibleInDim(1));
620 }
621 
622 
623 
624 template <int dim, int spacedim>
625 void
627  const typename Triangulation<dim,spacedim>::line_iterator &line,
628  std::vector<Point<spacedim> > &points) const
629 {
630  if (points.size()==1)
631  points[0]=get_new_point_on_line(line);
632  else
633  get_intermediate_points_between_points(line->vertex(0), line->vertex(1), points);
634 }
635 
636 
637 
638 template <int dim, int spacedim>
639 void
641  const Point<spacedim> &p0, const Point<spacedim> &p1,
642  std::vector<Point<spacedim> > &points) const
643 {
644  const unsigned int n=points.size();
645  Assert(n>0, ExcInternalError());
646 
647  const Tensor<1,spacedim> v0=p0-center,
648  v1=p1-center;
649  const double length=(v1-v0).norm();
650 
651  double eps=1e-12;
652  (void)eps;
653  double r=0;
654  if (compute_radius_automatically)
655  r = (p0 - center).norm();
656  else
657  r = radius;
658 
659  Assert((std::fabs(v0*v0-r*r)<eps*r*r)
660  &&
661  (std::fabs(v1*v1-r*r)<eps*r*r),
662  ExcMessage("In computing the location of midpoints of an edge of a "
663  "HyperBallBoundary, one of the vertices of the edge "
664  "does not have the expected distance from the center "
665  "of the sphere. This may happen if the geometry has "
666  "been deformed, or if the HyperBallBoundary object has "
667  "been assigned to a part of the boundary that is not "
668  "in fact a sphere (e.g., the sides of a quarter shell "
669  "or one octant of a ball)."));
670 
671  const double alpha=std::acos((v0*v1)/std::sqrt((v0*v0)*(v1*v1)));
672  const Tensor<1,spacedim> pm=0.5*(v0+v1);
673 
674  const double h=pm.norm();
675 
676  // n even: m=n/2,
677  // n odd: m=(n-1)/2
678  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
679  const unsigned int m=n/2;
680  for (unsigned int i=0; i<m ; ++i)
681  {
682  const double beta = alpha * (line_points[i+1][0]-0.5);
683  const double d = h*std::tan(beta);
684  points[i] = Point<spacedim>(pm+d/length*(v1-v0));
685  points[n-1-i] = Point<spacedim>(pm-d/length*(v1-v0));
686  }
687 
688  if ((n+1)%2==0)
689  // if the number of parts is even insert the midpoint
690  points[(n-1)/2] = Point<spacedim>(pm);
691 
692 
693  // project the points from the straight line to the HyperBallBoundary
694  for (unsigned int i=0; i<n; ++i)
695  {
696  points[i] *= r / std::sqrt(points[i].square());
697  points[i] += center;
698  }
699 }
700 
701 
702 
703 template <>
704 void
706  const Triangulation<3>::quad_iterator &quad,
707  std::vector<Point<3> > &points) const
708 {
709  if (points.size()==1)
710  points[0]=get_new_point_on_quad(quad);
711  else
712  {
713  unsigned int m=static_cast<unsigned int> (std::sqrt(static_cast<double>(points.size())));
714  Assert(points.size()==m*m, ExcInternalError());
715 
716  std::vector<Point<3> > lp0(m);
717  std::vector<Point<3> > lp1(m);
718 
719  get_intermediate_points_on_line(quad->line(0), lp0);
720  get_intermediate_points_on_line(quad->line(1), lp1);
721 
722  std::vector<Point<3> > lps(m);
723  for (unsigned int i=0; i<m; ++i)
724  {
725  get_intermediate_points_between_points(lp0[i], lp1[i], lps);
726 
727  for (unsigned int j=0; j<m; ++j)
728  points[i*m+j]=lps[j];
729  }
730  }
731 }
732 
733 
734 
735 template <>
736 void
738  const Triangulation<2,3>::quad_iterator &quad,
739  std::vector<Point<3> > &points) const
740 {
741  if (points.size()==1)
742  points[0]=get_new_point_on_quad(quad);
743  else
744  {
745  unsigned int m=static_cast<unsigned int> (std::sqrt(static_cast<double>(points.size())));
746  Assert(points.size()==m*m, ExcInternalError());
747 
748  std::vector<Point<3> > lp0(m);
749  std::vector<Point<3> > lp1(m);
750 
751  get_intermediate_points_on_line(quad->line(0), lp0);
752  get_intermediate_points_on_line(quad->line(1), lp1);
753 
754  std::vector<Point<3> > lps(m);
755  for (unsigned int i=0; i<m; ++i)
756  {
757  get_intermediate_points_between_points(lp0[i], lp1[i], lps);
758 
759  for (unsigned int j=0; j<m; ++j)
760  points[i*m+j]=lps[j];
761  }
762  }
763 }
764 
765 
766 
767 template <int dim, int spacedim>
768 void
770  const typename Triangulation<dim,spacedim>::quad_iterator &,
771  std::vector<Point<spacedim> > &) const
772 {
773  Assert(false, ExcImpossibleInDim(dim));
774 }
775 
776 
777 
778 template <int dim, int spacedim>
782  const Point<spacedim> &p) const
783 {
784  const Tensor<1,spacedim> unnormalized_normal = p-center;
785  return unnormalized_normal/unnormalized_normal.norm();
786 }
787 
788 
789 
790 template <>
791 void
795 {
796  Assert (false, ExcImpossibleInDim(1));
797 }
798 
799 template <>
800 void
804 {
805  Assert (false, ExcImpossibleInDim(1));
806 }
807 
808 
809 
810 template <int dim, int spacedim>
811 void
814  typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const
815 {
816  for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
817  face_vertex_normals[vertex] = face->vertex(vertex)-center;
818 }
819 
820 
821 
822 template <int dim, int spacedim>
825 {
826  return center;
827 }
828 
829 
830 
831 template <int dim, int spacedim>
832 double
834 {
835  Assert(!compute_radius_automatically, ExcRadiusNotSet());
836  return radius;
837 }
838 
839 
840 /* ---------------------------------------------------------------------- */
841 
842 
843 template <int dim>
845  const double radius) :
846  HyperBallBoundary<dim> (center, radius)
847 {}
848 
849 
850 
851 template<int dim>
852 std::unique_ptr<Manifold<dim,dim> > HalfHyperBallBoundary<dim>::clone() const
853 {
854  return std_cxx14::make_unique<HalfHyperBallBoundary<dim> >(this->get_center(),
855  this->get_radius());
856 }
857 
858 
859 
860 template <int dim>
864 {
865  // check whether center of object is at x==x_center, since then it belongs
866  // to the plane part of the boundary. however, this is not the case if it is
867  // at the outer perimeter
868  const Point<dim> line_center = line->center();
869  const Point<dim> vertices[2] = { line->vertex(0), line->vertex(1) };
870 
871  if ((line_center(0) == this->center(0))
872  &&
873  ((std::fabs(vertices[0].distance(this->center)-this->radius) >
874  1e-5*this->radius)
875  ||
876  (std::fabs(vertices[1].distance(this->center)-this->radius) >
877  1e-5*this->radius)))
878  return line_center;
879  else
881 }
882 
883 
884 
885 template <>
886 Point<1>
888 get_new_point_on_quad (const Triangulation<1>::quad_iterator &) const
889 {
890  Assert (false, ExcInternalError());
891  return Point<1>();
892 }
893 
894 
895 
896 template <int dim>
900 {
901  const Point<dim> quad_center = quad->center();
902  if (quad_center(0) == this->center(0))
903  return quad_center;
904  else
906 }
907 
908 
909 
910 template <int dim>
911 void
914  std::vector<Point<dim> > &points) const
915 {
916  // check whether center of object is at x==0, since then it belongs to the
917  // plane part of the boundary
918  const Point<dim> line_center = line->center();
919  if (line_center(0) == this->center(0))
921  else
923 }
924 
925 
926 
927 template <int dim>
928 void
931  std::vector<Point<dim> > &points) const
932 {
933  if (points.size()==1)
934  points[0]=get_new_point_on_quad(quad);
935  else
936  {
937  // check whether center of object is at x==0, since then it belongs to
938  // the plane part of the boundary
939  const Point<dim> quad_center = quad->center();
940  if (quad_center(0) == this->center(0))
942  else
944  }
945 }
946 
947 
948 
949 template <>
950 void
952 get_intermediate_points_on_quad (const Triangulation<1>::quad_iterator &,
953  std::vector<Point<1> > &) const
954 {
955  Assert (false, ExcInternalError());
956 }
957 
958 
959 
960 template <>
961 void
965 {
966  Assert (false, ExcImpossibleInDim(1));
967 }
968 
969 
970 
971 template <int dim>
972 void
975  typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const
976 {
977  // check whether center of object is at x==0, since then it belongs to the
978  // plane part of the boundary
979  const Point<dim> quad_center = face->center();
980  if (quad_center(0) == this->center(0))
981  StraightBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
982  else
983  HyperBallBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
984 }
985 
986 
987 /* ---------------------------------------------------------------------- */
988 
989 
990 
991 template <int dim>
993  :
994  HyperBallBoundary<dim>(center, 0.)
995 {
996  this->compute_radius_automatically=true;
997 }
998 
999 
1000 
1001 template<int dim>
1002 std::unique_ptr<Manifold<dim,dim> >
1004 {
1005  return std_cxx14::make_unique<HyperShellBoundary<dim> >(this->get_center());
1006 }
1007 
1008 
1009 
1010 /* ---------------------------------------------------------------------- */
1011 
1012 
1013 
1014 
1015 template <int dim>
1017  const double inner_radius,
1018  const double outer_radius)
1019  :
1020  HyperShellBoundary<dim> (center),
1021  inner_radius (inner_radius),
1022  outer_radius (outer_radius)
1023 {
1024  if (dim > 2)
1025  Assert ((inner_radius >= 0) &&
1026  (outer_radius > 0) &&
1027  (outer_radius > inner_radius),
1028  ExcMessage ("Inner and outer radii must be specified explicitly in 3d."));
1029 }
1030 
1031 
1032 
1033 template<int dim>
1034 std::unique_ptr<Manifold<dim,dim> >
1036 {
1037  return std_cxx14::make_unique<HalfHyperShellBoundary<dim> >(this->get_center(),
1038  inner_radius,
1039  outer_radius);
1040 }
1041 
1042 
1043 
1044 template <int dim>
1045 Point<dim>
1048 {
1049  switch (dim)
1050  {
1051  // in 2d, first check whether the two end points of the line are on the
1052  // axis of symmetry. if so, then return the mid point
1053  case 2:
1054  {
1055  if ((line->vertex(0)(0) == this->center(0))
1056  &&
1057  (line->vertex(1)(0) == this->center(0)))
1058  return (line->vertex(0) + line->vertex(1))/2;
1059  else
1060  // otherwise we are on the outer or inner part of the shell. proceed
1061  // as in the base class
1063  }
1064 
1065  // in 3d, a line is a straight line if it is on the symmetry plane and if
1066  // not both of its end points are on either the inner or outer sphere
1067  case 3:
1068  {
1069 
1070  if (((line->vertex(0)(0) == this->center(0))
1071  &&
1072  (line->vertex(1)(0) == this->center(0)))
1073  &&
1074  !(((std::fabs (line->vertex(0).distance (this->center)
1075  - inner_radius) < 1e-12 * outer_radius)
1076  &&
1077  (std::fabs (line->vertex(1).distance (this->center)
1078  - inner_radius) < 1e-12 * outer_radius))
1079  ||
1080  ((std::fabs (line->vertex(0).distance (this->center)
1081  - outer_radius) < 1e-12 * outer_radius)
1082  &&
1083  (std::fabs (line->vertex(1).distance (this->center)
1084  - outer_radius) < 1e-12 * outer_radius))))
1085  return (line->vertex(0) + line->vertex(1))/2;
1086  else
1087  // otherwise we are on the outer or inner part of the shell. proceed
1088  // as in the base class
1090  }
1091 
1092  default:
1093  Assert (false, ExcNotImplemented());
1094  }
1095 
1096  return Point<dim>();
1097 }
1098 
1099 
1100 
1101 template <>
1102 Point<1>
1104 get_new_point_on_quad (const Triangulation<1>::quad_iterator &) const
1105 {
1106  Assert (false, ExcInternalError());
1107  return Point<1>();
1108 }
1109 
1110 
1111 
1112 
1113 template <int dim>
1114 Point<dim>
1117 {
1118  // if this quad is on the symmetry plane, take the center point and project
1119  // it outward to the same radius as the centers of the two radial lines
1120  if ((quad->vertex(0)(0) == this->center(0)) &&
1121  (quad->vertex(1)(0) == this->center(0)) &&
1122  (quad->vertex(2)(0) == this->center(0)) &&
1123  (quad->vertex(3)(0) == this->center(0)))
1124  {
1125  const Point<dim> quad_center = (quad->vertex(0) + quad->vertex(1) +
1126  quad->vertex(2) + quad->vertex(3) )/4;
1127  const Tensor<1,dim> quad_center_offset = quad_center - this->center;
1128 
1129 
1130  if (std::fabs (quad->line(0)->center().distance(this->center) -
1131  quad->line(1)->center().distance(this->center))
1132  < 1e-12 * outer_radius)
1133  {
1134  // lines 0 and 1 are radial
1135  const double needed_radius
1136  = quad->line(0)->center().distance(this->center);
1137 
1138  return (this->center +
1139  quad_center_offset/quad_center_offset.norm() * needed_radius);
1140  }
1141  else if (std::fabs (quad->line(2)->center().distance(this->center) -
1142  quad->line(3)->center().distance(this->center))
1143  < 1e-12 * outer_radius)
1144  {
1145  // lines 2 and 3 are radial
1146  const double needed_radius
1147  = quad->line(2)->center().distance(this->center);
1148 
1149  return (this->center +
1150  quad_center_offset/quad_center_offset.norm() * needed_radius);
1151  }
1152  else
1153  Assert (false, ExcInternalError());
1154  }
1155 
1156  // otherwise we are on the outer or inner part of the shell. proceed as in
1157  // the base class
1159 }
1160 
1161 
1162 
1163 template <int dim>
1164 void
1167  std::vector<Point<dim> > &points) const
1168 {
1169  switch (dim)
1170  {
1171  // in 2d, first check whether the two end points of the line are on the
1172  // axis of symmetry. if so, then return the mid point
1173  case 2:
1174  {
1175  if ((line->vertex(0)(0) == this->center(0))
1176  &&
1177  (line->vertex(1)(0) == this->center(0)))
1179  else
1180  // otherwise we are on the outer or inner part of the shell. proceed
1181  // as in the base class
1183  break;
1184  }
1185 
1186  // in 3d, a line is a straight line if it is on the symmetry plane and if
1187  // not both of its end points are on either the inner or outer sphere
1188  case 3:
1189  {
1190  if (((line->vertex(0)(0) == this->center(0))
1191  &&
1192  (line->vertex(1)(0) == this->center(0)))
1193  &&
1194  !(((std::fabs (line->vertex(0).distance (this->center)
1195  - inner_radius) < 1e-12 * outer_radius)
1196  &&
1197  (std::fabs (line->vertex(1).distance (this->center)
1198  - inner_radius) < 1e-12 * outer_radius))
1199  ||
1200  ((std::fabs (line->vertex(0).distance (this->center)
1201  - outer_radius) < 1e-12 * outer_radius)
1202  &&
1203  (std::fabs (line->vertex(1).distance (this->center)
1204  - outer_radius) < 1e-12 * outer_radius))))
1206  else
1207  // otherwise we are on the outer or inner part of the shell. proceed
1208  // as in the base class
1210 
1211  break;
1212  }
1213 
1214  default:
1215  Assert (false, ExcNotImplemented());
1216  }
1217 }
1218 
1219 
1220 
1221 template <int dim>
1222 void
1225  std::vector<Point<dim> > &points) const
1226 {
1227  Assert (dim < 3, ExcNotImplemented());
1228 
1229  // check whether center of object is at x==0, since then it belongs to the
1230  // plane part of the boundary
1231  const Point<dim> quad_center = quad->center();
1232  if (quad_center(0) == this->center(0))
1234  else
1236 }
1237 
1238 
1239 
1240 template <>
1241 void
1243 get_intermediate_points_on_quad (const Triangulation<1>::quad_iterator &,
1244  std::vector<Point<1> > &) const
1245 {
1246  Assert (false, ExcInternalError());
1247 }
1248 
1249 
1250 
1251 template <>
1252 void
1256 {
1257  Assert (false, ExcImpossibleInDim(1));
1258 }
1259 
1260 
1261 
1262 
1263 
1264 template <int dim>
1265 void
1268  typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const
1269 {
1270  if (face->center()(0) == this->center(0))
1271  StraightBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
1272  else
1273  HyperShellBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
1274 }
1275 
1276 
1277 
1278 
1279 template <int dim, int spacedim>
1281  const double r_)
1282  :
1283  R(R_),
1284  r(r_)
1285 {
1286  Assert (false, ExcNotImplemented());
1287 }
1288 
1289 
1290 
1291 template <>
1292 TorusBoundary<2,3>::TorusBoundary (const double R_,
1293  const double r_)
1294  :
1295  R(R_),
1296  r(r_)
1297 {
1298  Assert (R>r, ExcMessage("Outer radius must be greater than inner radius."));
1299 }
1300 
1301 template<int dim, int spacedim>
1302 std::unique_ptr<Manifold<dim, spacedim> >
1304 {
1305  return std_cxx14::make_unique<TorusBoundary<dim,spacedim> >(R,r);
1306 }
1307 
1308 
1309 
1310 template <int dim, int spacedim>
1311 double
1313  const double x,
1314  const double y) const
1315 {
1316  if (y>=0)
1317  {
1318  if (x >=0)
1319  return angle;
1320 
1321  return numbers::PI-angle;
1322  }
1323 
1324  if (x <=0)
1325  return numbers::PI+angle;
1326 
1327  return 2.0*numbers::PI-angle;
1328 }
1329 
1330 
1331 
1332 template <>
1333 Point<3>
1334 TorusBoundary<2,3>::get_real_coord (const Point<2> &surfP) const
1335 {
1336  const double theta=surfP(0);
1337  const double phi=surfP(1);
1338 
1339  return Point<3> ((R+r*std::cos(phi))*std::cos(theta),
1340  r*std::sin(phi),
1341  (R+r*std::cos(phi))*std::sin(theta));
1342 }
1343 
1344 
1345 
1346 template <>
1347 Point<2>
1349 {
1350  const double phi=std::asin(std::abs(p(1))/r);
1351  const double Rr_2=p(0)*p(0)+p(2)*p(2);
1352 
1353  Point<2> surfP;
1354  surfP(1)=get_correct_angle(phi,Rr_2-R*R,p(1));//phi
1355 
1356  if (std::abs(p(0))<1.E-5)
1357  {
1358  if (p(2)>=0)
1359  surfP(0) = numbers::PI*0.5;
1360  else
1361  surfP(0) = -numbers::PI*0.5;
1362  }
1363  else
1364  {
1365  const double theta = std::atan(std::abs(p(2)/p(0)));
1366  surfP(0)= get_correct_angle(theta,p(0),p(2));
1367  }
1368 
1369  return surfP;
1370 }
1371 
1372 
1373 
1374 template <>
1375 Point<3>
1376 TorusBoundary<2,3>::get_new_point_on_line (const Triangulation<2,3>::line_iterator &line) const
1377 {
1378  //Just get the average
1379  Point<2> p0=get_surf_coord(line->vertex(0));
1380  Point<2> p1=get_surf_coord(line->vertex(1));
1381 
1382  Point<2> middle(0,0);
1383 
1384  //Take care for periodic conditions, For instance phi0= 0, phi1= 3/2*Pi
1385  //middle has to be 7/4*Pi not 3/4*Pi. This also works for -Pi/2 + Pi, middle
1386  //is 5/4*Pi
1387  for (unsigned int i=0; i<2; i++)
1388  if (std::abs(p0(i)-p1(i))> numbers::PI)
1389  middle(i)=2*numbers::PI;
1390 
1391  middle+= p0 + p1;
1392  middle*=0.5;
1393 
1394  Point<3> midReal=get_real_coord(middle);
1395  return midReal;
1396 }
1397 
1398 
1399 
1400 template <>
1401 Point<3>
1402 TorusBoundary<2,3>::get_new_point_on_quad (const Triangulation<2,3>::quad_iterator &quad) const
1403 {
1404  //Just get the average
1405  Point<2> p[4];
1406 
1407  for (unsigned int i=0; i<4; i++)
1408  p[i]=get_surf_coord(quad->vertex(i));
1409 
1410  Point<2> middle(0,0);
1411 
1412  //Take care for periodic conditions, see get_new_point_on_line() above
1413  //For instance phi0= 0, phi1= 3/2*Pi middle has to be 7/4*Pi not 3/4*Pi
1414  //This also works for -Pi/2 + Pi + Pi- Pi/2, middle is 5/4*Pi
1415  for (unsigned int i=0; i<2; i++)
1416  for (unsigned int j=1; j<4; j++)
1417  {
1418  if (std::abs(p[0](i)-p[j](i))> numbers::PI)
1419  {
1420  middle(i)+=2*numbers::PI;
1421  }
1422  }
1423 
1424  for (unsigned int i=0; i<4; i++)
1425  middle+=p[i];
1426 
1427  middle*= 0.25;
1428 
1429  return get_real_coord(middle);
1430 }
1431 
1432 
1433 
1434 //Normal field without unit length
1435 template <>
1436 Point<3>
1438 {
1439 
1440  Point<3> n;
1441  double theta=surfP[0];
1442  double phi=surfP[1];
1443 
1444  double f=R+r*std::cos(phi);
1445 
1446  n[0]=r*std::cos(phi)*std::cos(theta)*f;
1447  n[1]=r*std::sin(phi)*f;
1448  n[2]=r*std::sin(theta)*std::cos(phi)*f;
1449 
1450  return n;
1451 }
1452 
1453 
1454 
1455 //Normal field without unit length
1456 template <>
1457 Point<3>
1459 {
1460 
1461  Point<2> surfP=get_surf_coord(p);
1462  return get_surf_norm_from_sp(surfP);
1463 
1464 }
1465 
1466 
1467 
1468 template <>
1469 void
1471 get_intermediate_points_on_line (const Triangulation<2, 3>::line_iterator &line,
1472  std::vector< Point< 3 > > &points) const
1473 {
1474  //Almost the same implementation as StraightBoundary<2,3>
1475  unsigned int npoints=points.size();
1476  if (npoints==0) return;
1477 
1478  Point<2> p[2];
1479 
1480  for (unsigned int i=0; i<2; i++)
1481  p[i]=get_surf_coord(line->vertex(i));
1482 
1483  unsigned int offset[2];
1484  offset[0]=0;
1485  offset[1]=0;
1486 
1487  //Take care for periodic conditions & negative angles, see
1488  //get_new_point_on_line() above. Because we don't have a symmetric
1489  //interpolation (just the middle) we need to add 2*Pi to each almost zero
1490  //and negative angles.
1491  for (unsigned int i=0; i<2; i++)
1492  for (unsigned int j=1; j<2; j++)
1493  {
1494  if (std::abs(p[0](i)-p[j](i))> numbers::PI)
1495  {
1496  offset[i]++;
1497  break;
1498  }
1499  }
1500 
1501  for (unsigned int i=0; i<2; i++)
1502  for (unsigned int j=0; j<2; j++)
1503  if (p[j](i)<1.E-12 ) //Take care for periodic conditions & negative angles
1504  p[j](i)+=2*numbers::PI*offset[i];
1505 
1506 
1507  Point<2> target;
1508  const std::vector<Point<1> > &line_points = this->get_line_support_points(npoints);
1509  for (unsigned int i=0; i<npoints; i++)
1510  {
1511  const double x = line_points[i+1][0];
1512  target= (1-x)*p[0] + x*p[1];
1513  points[i]=get_real_coord(target);
1514  }
1515 }
1516 
1517 
1518 
1519 template <>
1520 void
1522 get_intermediate_points_on_quad (const Triangulation< 2, 3 >::quad_iterator &quad,
1523  std::vector< Point< 3 > > &points )const
1524 {
1525  //Almost the same implementation as StraightBoundary<2,3>
1526  const unsigned int n=points.size(),
1527  m=static_cast<unsigned int>(std::sqrt(static_cast<double>(n)));
1528  // is n a square number
1529  Assert(m*m==n, ExcInternalError());
1530 
1531  Point<2> p[4];
1532 
1533  for (unsigned int i=0; i<4; i++)
1534  p[i]=get_surf_coord(quad->vertex(i));
1535 
1536  Point<2> target;
1537  unsigned int offset[2];
1538  offset[0]=0;
1539  offset[1]=0;
1540 
1541  //Take care for periodic conditions & negative angles, see
1542  //get_new_point_on_line() above. Because we don't have a symmetric
1543  //interpolation (just the middle) we need to add 2*Pi to each almost zero
1544  //and negative angles.
1545  for (unsigned int i=0; i<2; i++)
1546  for (unsigned int j=1; j<4; j++)
1547  {
1548  if (std::abs(p[0](i)-p[j](i))> numbers::PI)
1549  {
1550  offset[i]++;
1551  break;
1552  }
1553  }
1554 
1555  for (unsigned int i=0; i<2; i++)
1556  for (unsigned int j=0; j<4; j++)
1557  if (p[j](i)<1.E-12 ) //Take care for periodic conditions & negative angles
1558  p[j](i)+=2*numbers::PI*offset[i];
1559 
1560  const std::vector<Point<1> > &line_points = this->get_line_support_points(m);
1561  for (unsigned int i=0; i<m; ++i)
1562  {
1563  const double y=line_points[i+1][0];
1564  for (unsigned int j=0; j<m; ++j)
1565  {
1566  const double x=line_points[j+1][0];
1567  target=((1-x) * p[0] +
1568  x * p[1]) * (1-y) +
1569  ((1-x) * p[2] +
1570  x * p[3]) * y;
1571 
1572  points[i*m+j]=get_real_coord(target);
1573  }
1574  }
1575 }
1576 
1577 
1578 
1579 template <>
1580 void
1583  Boundary<2,3>::FaceVertexNormals &face_vertex_normals) const
1584 {
1585  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_face; i++)
1586  face_vertex_normals[i]=get_surf_norm(face->vertex(i));
1587 }
1588 
1589 
1590 
1591 // explicit instantiations
1592 #include "tria_boundary_lib.inst"
1593 
1594 DEAL_II_NAMESPACE_CLOSE
virtual void get_normals_at_vertices(const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const override
void get_intermediate_points_between_points(const Point< spacedim > &p0, const Point< spacedim > &p1, std::vector< Point< spacedim > > &points) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const override
void get_intermediate_points_between_points(const Point< dim > &p0, const Point< dim > &p1, std::vector< Point< dim > > &points) const
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const override
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const override
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const override
virtual void get_intermediate_points_on_line(const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const override
virtual std::unique_ptr< Manifold< dim, dim > > clone() const override
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const override
const double R
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:1273
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const override
Definition: point.h:104
double get_radius(const Point< dim > x) const
Point< dim > get_surf_coord(const Point< spacedim > &p) const
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const override
double get_correct_angle(const double angle, const double x, const double y) const
virtual void get_normals_at_vertices(const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const override
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const override
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const override
static Point< spacedim > get_axis_vector(const unsigned int axis)
static const double PI
Definition: numbers.h:127
HalfHyperShellBoundary(const Point< dim > &center=Point< dim >(), const double inner_radius=-1, const double outer_radius=-1)
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
Definition: tria.h:45
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const override
virtual Point< dim > get_new_point_on_line(const typename Triangulation< dim >::line_iterator &line) const override
virtual void get_intermediate_points_on_line(const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const override
virtual std::unique_ptr< Manifold< dim, dim > > clone() const override
HalfHyperBallBoundary(const Point< dim > p=Point< dim >(), const double radius=1.0)
#define Assert(cond, exc)
Definition: exceptions.h:1142
Point< spacedim > get_surf_norm(const Point< spacedim > &p) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const override
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const override
Point< spacedim > get_real_coord(const Point< dim > &surfP) const
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const override
virtual Point< dim > get_new_point_on_line(const typename Triangulation< dim >::line_iterator &line) const override
virtual Point< dim > get_new_point_on_quad(const typename Triangulation< dim >::quad_iterator &quad) const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual Point< dim > get_new_point_on_quad(const typename Triangulation< dim >::quad_iterator &quad) const override
ConeBoundary(const double radius_0, const double radius_1, const Point< dim > x_0, const Point< dim > x_1)
virtual Tensor< 1, dim > normal_vector(const typename Triangulation< dim >::face_iterator &face, const Point< dim > &p) const override
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const override
virtual Point< dim > get_new_point_on_line(const typename Triangulation< dim >::line_iterator &line) const override
virtual void get_normals_at_vertices(const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const override
double get_radius() const
HyperShellBoundary(const Point< dim > &center=Point< dim >())
Point< spacedim > get_surf_norm_from_sp(const Point< dim > &surfP) const
numbers::NumberTraits< Number >::real_type square() const
double get_radius() const
void get_intermediate_points_between_points(const Point< spacedim > &p0, const Point< spacedim > &p1, std::vector< Point< spacedim > > &points) const
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
virtual Point< dim > get_new_point_on_quad(const typename Triangulation< dim >::quad_iterator &quad) const override
TorusBoundary(const double R, const double r)
static ::ExceptionBase & ExcNotImplemented()
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const override
virtual std::unique_ptr< Manifold< dim, dim > > clone() const override
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const override
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const override
HyperBallBoundary(const Point< spacedim > p=Point< spacedim >(), const double radius=1.0)
CylinderBoundary(const double radius=1.0, const unsigned int axis=0)
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Point< spacedim > get_center() const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const override
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const override
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const override
virtual std::unique_ptr< Manifold< dim, dim > > clone() const override
static ::ExceptionBase & ExcInternalError()