Reference documentation for deal.II version 8.4.2
tria_boundary_lib.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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/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 <cmath>
22 
23 
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 
29 template <int dim, int spacedim>
31  const unsigned int axis)
32  :
33  radius(radius),
34  direction (get_axis_vector (axis)),
35  point_on_axis (Point<spacedim>())
36 {}
37 
38 
39 template <int dim, int spacedim>
43  :
44  radius(radius),
45  direction (direction / direction.norm()),
46  point_on_axis (point_on_axis)
47 {}
48 
49 
50 template <int dim, int spacedim>
53 {
54  Assert (axis < spacedim, ExcIndexRange (axis, 0, spacedim));
55 
56  Point<spacedim> axis_vector;
57  axis_vector[axis] = 1;
58  return axis_vector;
59 }
60 
61 
62 
63 template <int dim, int spacedim>
66 get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const
67 {
68  // compute a proposed new point
70 
71  // we then have to project this
72  // point out to the given radius
73  // from the axis. to this end, we
74  // have to take into account the
75  // offset point_on_axis and the
76  // direction of the axis
77  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
78  ((middle-point_on_axis) * direction) * direction;
79  // scale it to the desired length
80  // and put everything back
81  // together, unless we have a point
82  // on the axis
83  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
84  return middle;
85  else
86  return Point<spacedim>(vector_from_axis / vector_from_axis.norm() * radius +
87  ((middle-point_on_axis) * direction) * direction +
89 }
90 
91 
92 
93 template<>
96 get_new_point_on_quad (const Triangulation<3>::quad_iterator &quad) const
97 {
99 
100  // same algorithm as above
101  const unsigned int spacedim = 3;
102 
103  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
104  ((middle-point_on_axis) * direction) * direction;
105  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
106  return middle;
107  else
108  return Point<3>(vector_from_axis / vector_from_axis.norm() * radius +
109  ((middle-point_on_axis) * direction) * direction +
110  point_on_axis);
111 }
112 
113 template<>
114 Point<3>
116 get_new_point_on_quad (const Triangulation<2,3>::quad_iterator &quad) const
117 {
119 
120  // same algorithm as above
121  const unsigned int spacedim = 3;
122  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
123  ((middle-point_on_axis) * direction) * direction;
124  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
125  return middle;
126  else
127  return Point<3>(vector_from_axis / vector_from_axis.norm() * radius +
128  ((middle-point_on_axis) * direction) * direction +
129  point_on_axis);
130 }
131 
132 
133 template <int dim, int spacedim>
136 get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &) const
137 {
138  Assert (false, ExcImpossibleInDim(dim));
139  return Point<spacedim>();
140 }
141 
142 
143 
144 template <int dim, int spacedim>
145 void
147  const typename Triangulation<dim,spacedim>::line_iterator &line,
148  std::vector<Point<spacedim> > &points) const
149 {
150  if (points.size()==1)
151  points[0]=get_new_point_on_line(line);
152  else
153  get_intermediate_points_between_points(line->vertex(0), line->vertex(1), points);
154 }
155 
156 
157 template <int dim, int spacedim>
158 void
160  const Point<spacedim> &v0,
161  const Point<spacedim> &v1,
162  std::vector<Point<spacedim> > &points) const
163 {
164  const unsigned int n=points.size();
165  Assert(n>0, ExcInternalError());
166 
167  // Do a simple linear interpolation followed by projection, using the same
168  // algorithm as above
169  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
170 
171  for (unsigned int i=0; i<n; ++i)
172  {
173  const double x = line_points[i+1][0];
174  const Point<spacedim> middle = (1-x)*v0 + x*v1;
175 
176  const Tensor<1,spacedim> vector_from_axis = (middle-point_on_axis) -
177  ((middle-point_on_axis) * direction) * direction;
178  if (vector_from_axis.norm() <= 1e-10 * middle.norm())
179  points[i] = middle;
180  else
181  points[i] = Point<spacedim>(vector_from_axis / vector_from_axis.norm() * radius +
182  ((middle-point_on_axis) * direction) * direction +
183  point_on_axis);
184  }
185 }
186 
187 
188 
189 template <>
190 void
192  const Triangulation<3>::quad_iterator &quad,
193  std::vector<Point<3> > &points) const
194 {
195  if (points.size()==1)
196  points[0]=get_new_point_on_quad(quad);
197  else
198  {
199  unsigned int m=static_cast<unsigned int> (std::sqrt(static_cast<double>(points.size())));
200  Assert(points.size()==m*m, ExcInternalError());
201 
202  std::vector<Point<3> > lp0(m);
203  std::vector<Point<3> > lp1(m);
204 
205  get_intermediate_points_on_line(quad->line(0), lp0);
206  get_intermediate_points_on_line(quad->line(1), lp1);
207 
208  std::vector<Point<3> > lps(m);
209  for (unsigned int i=0; i<m; ++i)
210  {
211  get_intermediate_points_between_points(lp0[i], lp1[i], lps);
212 
213  for (unsigned int j=0; j<m; ++j)
214  points[i*m+j]=lps[j];
215  }
216  }
217 }
218 
219 
220 
221 template <int dim, int spacedim>
222 void
224  const typename Triangulation<dim,spacedim>::quad_iterator &,
225  std::vector<Point<spacedim> > &) const
226 {
227  Assert (false, ExcImpossibleInDim(dim));
228 }
229 
230 
231 
232 
233 template <>
234 void
238 {
239  Assert (false, ExcImpossibleInDim(1));
240 }
241 
242 
243 
244 
245 template <int dim, int spacedim>
246 void
249  typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const
250 {
251  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
252  {
253  const Point<spacedim> vertex = face->vertex(v);
254 
255  const Tensor<1,spacedim> vector_from_axis = (vertex-point_on_axis) -
256  ((vertex-point_on_axis) * direction) * direction;
257 
258  face_vertex_normals[v] = (vector_from_axis / vector_from_axis.norm());
259  }
260 }
261 
262 
263 
264 template <int dim, int spacedim>
265 double
267 {
268  return radius;
269 }
270 
271 
272 //======================================================================//
273 
274 template<int dim>
275 ConeBoundary<dim>::ConeBoundary (const double radius_0,
276  const double radius_1,
277  const Point<dim> x_0,
278  const Point<dim> x_1)
279  :
280  radius_0 (radius_0),
281  radius_1 (radius_1),
282  x_0 (x_0),
283  x_1 (x_1)
284 {}
285 
286 
287 
288 template<int dim>
290 {
291  for (unsigned int i = 0; i < dim; ++i)
292  if ((x_1 (i) - x_0 (i)) != 0)
293  return (radius_1 - radius_0) * (x (i) - x_0 (i)) / (x_1 (i) - x_0 (i)) + radius_0;
294 
295  return 0;
296 }
297 
298 
299 
300 template<int dim>
301 void
304  const Point<dim> &p1,
305  std::vector<Point<dim> > &points) const
306 {
307  const unsigned int n = points.size ();
308  const Tensor<1,dim> axis = x_1 - x_0;
309 
310  Assert (n > 0, ExcInternalError ());
311 
312  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
313 
314  for (unsigned int i=0; i<n; ++i)
315  {
316  const double x = line_points[i+1][0];
317 
318  // Compute the current point.
319  const Point<dim> x_i = (1-x)*p0 + x*p1;
320  // To project this point on the boundary of the cone we first compute
321  // the orthogonal projection of this point onto the axis of the cone.
322  const double c = (x_i - x_0) * axis / (axis*axis);
323  const Point<dim> x_ip = x_0 + c * axis;
324  // Compute the projection of the middle point on the boundary of the
325  // cone.
326  points[i] = x_ip + get_radius (x_ip) * (x_i - x_ip) / (x_i - x_ip).norm ();
327  }
328 }
329 
330 template<int dim>
334 {
335  const Tensor<1,dim> axis = x_1 - x_0;
336  // Compute the middle point of the line.
338  // To project it on the boundary of the cone we first compute the orthogonal
339  // projection of the middle point onto the axis of the cone.
340  const double c = (middle - x_0) * axis / (axis*axis);
341  const Point<dim> middle_p = x_0 + c * axis;
342  // Compute the projection of the middle point on the boundary of the cone.
343  return middle_p + get_radius (middle_p) * (middle - middle_p) / (middle - middle_p).norm ();
344 }
345 
346 
347 
348 template <>
349 Point<3>
351 get_new_point_on_quad (const Triangulation<3>::quad_iterator &quad) const
352 {
353  const int dim = 3;
354 
355  const Tensor<1,dim> axis = x_1 - x_0;
356  // Compute the middle point of the quad.
358  // Same algorithm as above: To project it on the boundary of the cone we
359  // first compute the orthogonal projection of the middle point onto the axis
360  // of the cone.
361  const double c = (middle - x_0) * axis / (axis*axis);
362  const Point<dim> middle_p = x_0 + c * axis;
363  // Compute the projection of the middle point on the boundary of the cone.
364  return middle_p + get_radius (middle_p) * (middle - middle_p) / (middle - middle_p).norm ();
365 }
366 
367 
368 
369 template<int dim>
373 {
374  Assert (false, ExcImpossibleInDim (dim));
375 
376  return Point<dim>();
377 }
378 
379 
380 
381 template<int dim>
382 void
385  std::vector<Point<dim> > &points) const
386 {
387  if (points.size () == 1)
388  points[0] = get_new_point_on_line (line);
389  else
390  get_intermediate_points_between_points (line->vertex (0), line->vertex (1), points);
391 }
392 
393 
394 
395 
396 template<>
397 void
399 get_intermediate_points_on_quad (const Triangulation<3>::quad_iterator &quad,
400  std::vector<Point<3> > &points) const
401 {
402  if (points.size () == 1)
403  points[0] = get_new_point_on_quad (quad);
404  else
405  {
406  unsigned int n = static_cast<unsigned int> (std::sqrt (static_cast<double> (points.size ())));
407 
408  Assert (points.size () == n * n, ExcInternalError ());
409 
410  std::vector<Point<3> > points_line_0 (n);
411  std::vector<Point<3> > points_line_1 (n);
412 
413  get_intermediate_points_on_line (quad->line (0), points_line_0);
414  get_intermediate_points_on_line (quad->line (1), points_line_1);
415 
416  std::vector<Point<3> > points_line_segment (n);
417 
418  for (unsigned int i = 0; i < n; ++i)
419  {
420  get_intermediate_points_between_points (points_line_0[i],
421  points_line_1[i],
422  points_line_segment);
423 
424  for (unsigned int j = 0; j < n; ++j)
425  points[i * n + j] = points_line_segment[j];
426  }
427  }
428 }
429 
430 
431 
432 template <int dim>
433 void
436  std::vector<Point<dim> > &) const
437 {
438  Assert (false, ExcImpossibleInDim (dim));
439 }
440 
441 
442 
443 
444 template<>
445 void
449 {
450  Assert (false, ExcImpossibleInDim (1));
451 }
452 
453 
454 
455 template<int dim>
456 void
459  typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const
460 {
461  const Tensor<1,dim> axis = x_1 - x_0;
462 
463  for (unsigned int vertex = 0; vertex < GeometryInfo<dim>::vertices_per_face; ++vertex)
464  {
465  // Compute the orthogonal projection of the vertex onto the axis of the
466  // cone.
467  const double c = (face->vertex (vertex) - x_0) * axis / (axis*axis);
468  const Point<dim> vertex_p = x_0 + c * axis;
469  // Then compute the vector pointing from the point <tt>vertex_p</tt> on
470  // the axis to the vertex.
471  const Tensor<1,dim> axis_to_vertex = face->vertex (vertex) - vertex_p;
472 
473  face_vertex_normals[vertex] = axis_to_vertex / axis_to_vertex.norm ();
474  }
475 }
476 
477 
478 //======================================================================//
479 
480 template <int dim, int spacedim>
482  const double radius)
483  :
484  center(p),
485  radius(radius),
486  compute_radius_automatically(false)
487 {}
488 
489 
490 
491 template <int dim, int spacedim>
493 HyperBallBoundary<dim,spacedim>::get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const
494 {
496 
497  middle -= center;
498 
499  double r=0;
501  r = (line->vertex(0) - center).norm();
502  else
503  r = radius;
504 
505  // project to boundary
506  middle *= r / std::sqrt(middle.square());
507  middle += center;
508  return middle;
509 }
510 
511 
512 
513 template <>
514 Point<1>
516 get_new_point_on_quad (const Triangulation<1,1>::quad_iterator &) const
517 {
518  Assert (false, ExcInternalError());
519  return Point<1>();
520 }
521 
522 
523 template <>
524 Point<2>
526 get_new_point_on_quad (const Triangulation<1,2>::quad_iterator &) const
527 {
528  Assert (false, ExcInternalError());
529  return Point<2>();
530 }
531 
532 
533 
534 template <int dim, int spacedim>
537 get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const
538 {
540 
541  middle -= center;
542 
543  double r=0;
545  r = (quad->vertex(0) - center).norm();
546  else
547  r = radius;
548 
549  // project to boundary
550  middle *= r / std::sqrt(middle.square());
551 
552  middle += center;
553  return middle;
554 }
555 
556 
557 
558 template <>
559 void
561  const Triangulation<1>::line_iterator &,
562  std::vector<Point<1> > &) const
563 {
564  Assert (false, ExcImpossibleInDim(1));
565 }
566 
567 
568 
569 template <int dim, int spacedim>
570 void
572  const typename Triangulation<dim,spacedim>::line_iterator &line,
573  std::vector<Point<spacedim> > &points) const
574 {
575  if (points.size()==1)
576  points[0]=get_new_point_on_line(line);
577  else
578  get_intermediate_points_between_points(line->vertex(0), line->vertex(1), points);
579 }
580 
581 
582 
583 template <int dim, int spacedim>
584 void
586  const Point<spacedim> &p0, const Point<spacedim> &p1,
587  std::vector<Point<spacedim> > &points) const
588 {
589  const unsigned int n=points.size();
590  Assert(n>0, ExcInternalError());
591 
592  const Tensor<1,spacedim> v0=p0-center,
593  v1=p1-center;
594  const double length=(v1-v0).norm();
595 
596  double eps=1e-12;
597  (void)eps;
598  double r=0;
600  r = (p0 - center).norm();
601  else
602  r = radius;
603 
604  Assert(std::fabs(v0*v0-r*r)<eps*r*r, ExcInternalError());
605  Assert(std::fabs(v1*v1-r*r)<eps*r*r, ExcInternalError());
606 
607  const double alpha=std::acos((v0*v1)/std::sqrt((v0*v0)*(v1*v1)));
608  const Tensor<1,spacedim> pm=0.5*(v0+v1);
609 
610  const double h=pm.norm();
611 
612  // n even: m=n/2,
613  // n odd: m=(n-1)/2
614  const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
615  const unsigned int m=n/2;
616  for (unsigned int i=0; i<m ; ++i)
617  {
618  const double beta = alpha * (line_points[i+1][0]-0.5);
619  const double d = h*std::tan(beta);
620  points[i] = Point<spacedim>(pm+d/length*(v1-v0));
621  points[n-1-i] = Point<spacedim>(pm-d/length*(v1-v0));
622  }
623 
624  if ((n+1)%2==0)
625  // if the number of parts is even insert the midpoint
626  points[(n-1)/2] = Point<spacedim>(pm);
627 
628 
629  // project the points from the straight line to the HyperBallBoundary
630  for (unsigned int i=0; i<n; ++i)
631  {
632  points[i] *= r / std::sqrt(points[i].square());
633  points[i] += center;
634  }
635 }
636 
637 
638 
639 template <>
640 void
642  const Triangulation<3>::quad_iterator &quad,
643  std::vector<Point<3> > &points) const
644 {
645  if (points.size()==1)
646  points[0]=get_new_point_on_quad(quad);
647  else
648  {
649  unsigned int m=static_cast<unsigned int> (std::sqrt(static_cast<double>(points.size())));
650  Assert(points.size()==m*m, ExcInternalError());
651 
652  std::vector<Point<3> > lp0(m);
653  std::vector<Point<3> > lp1(m);
654 
655  get_intermediate_points_on_line(quad->line(0), lp0);
656  get_intermediate_points_on_line(quad->line(1), lp1);
657 
658  std::vector<Point<3> > lps(m);
659  for (unsigned int i=0; i<m; ++i)
660  {
661  get_intermediate_points_between_points(lp0[i], lp1[i], lps);
662 
663  for (unsigned int j=0; j<m; ++j)
664  points[i*m+j]=lps[j];
665  }
666  }
667 }
668 
669 
670 
671 template <>
672 void
674  const Triangulation<2,3>::quad_iterator &quad,
675  std::vector<Point<3> > &points) const
676 {
677  if (points.size()==1)
678  points[0]=get_new_point_on_quad(quad);
679  else
680  {
681  unsigned int m=static_cast<unsigned int> (std::sqrt(static_cast<double>(points.size())));
682  Assert(points.size()==m*m, ExcInternalError());
683 
684  std::vector<Point<3> > lp0(m);
685  std::vector<Point<3> > lp1(m);
686 
687  get_intermediate_points_on_line(quad->line(0), lp0);
688  get_intermediate_points_on_line(quad->line(1), lp1);
689 
690  std::vector<Point<3> > lps(m);
691  for (unsigned int i=0; i<m; ++i)
692  {
693  get_intermediate_points_between_points(lp0[i], lp1[i], lps);
694 
695  for (unsigned int j=0; j<m; ++j)
696  points[i*m+j]=lps[j];
697  }
698  }
699 }
700 
701 
702 
703 template <int dim, int spacedim>
704 void
706  const typename Triangulation<dim,spacedim>::quad_iterator &,
707  std::vector<Point<spacedim> > &) const
708 {
709  Assert(false, ExcImpossibleInDim(dim));
710 }
711 
712 
713 
714 template <int dim, int spacedim>
718  const Point<spacedim> &p) const
719 {
720  const Tensor<1,spacedim> unnormalized_normal = p-center;
721  return unnormalized_normal/unnormalized_normal.norm();
722 }
723 
724 
725 
726 template <>
727 void
731 {
732  Assert (false, ExcImpossibleInDim(1));
733 }
734 
735 template <>
736 void
740 {
741  Assert (false, ExcImpossibleInDim(1));
742 }
743 
744 
745 
746 template <int dim, int spacedim>
747 void
750  typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const
751 {
752  for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
753  face_vertex_normals[vertex] = face->vertex(vertex)-center;
754 }
755 
756 
757 
758 template <int dim, int spacedim>
761 {
762  return center;
763 }
764 
765 
766 
767 template <int dim, int spacedim>
768 double
770 {
771  Assert(!compute_radius_automatically, ExcRadiusNotSet());
772  return radius;
773 }
774 
775 
776 /* ---------------------------------------------------------------------- */
777 
778 
779 template <int dim>
781  const double radius) :
782  HyperBallBoundary<dim> (center, radius)
783 {}
784 
785 
786 
787 template <int dim>
791 {
792  // check whether center of object is at x==x_center, since then it belongs
793  // to the plane part of the boundary. however, this is not the case if it is
794  // at the outer perimeter
795  const Point<dim> line_center = line->center();
796  const Point<dim> vertices[2] = { line->vertex(0), line->vertex(1) };
797 
798  if ((line_center(0) == this->center(0))
799  &&
800  ((std::fabs(vertices[0].distance(this->center)-this->radius) >
801  1e-5*this->radius)
802  ||
803  (std::fabs(vertices[1].distance(this->center)-this->radius) >
804  1e-5*this->radius)))
805  return line_center;
806  else
808 }
809 
810 
811 
812 template <>
813 Point<1>
815 get_new_point_on_quad (const Triangulation<1>::quad_iterator &) const
816 {
817  Assert (false, ExcInternalError());
818  return Point<1>();
819 }
820 
821 
822 
823 template <int dim>
827 {
828  const Point<dim> quad_center = quad->center();
829  if (quad_center(0) == this->center(0))
830  return quad_center;
831  else
833 }
834 
835 
836 
837 template <int dim>
838 void
841  std::vector<Point<dim> > &points) const
842 {
843  // check whether center of object is at x==0, since then it belongs to the
844  // plane part of the boundary
845  const Point<dim> line_center = line->center();
846  if (line_center(0) == this->center(0))
848  else
850 }
851 
852 
853 
854 template <int dim>
855 void
858  std::vector<Point<dim> > &points) const
859 {
860  if (points.size()==1)
861  points[0]=get_new_point_on_quad(quad);
862  else
863  {
864  // check whether center of object is at x==0, since then it belongs to
865  // the plane part of the boundary
866  const Point<dim> quad_center = quad->center();
867  if (quad_center(0) == this->center(0))
869  else
871  }
872 }
873 
874 
875 
876 template <>
877 void
879 get_intermediate_points_on_quad (const Triangulation<1>::quad_iterator &,
880  std::vector<Point<1> > &) const
881 {
882  Assert (false, ExcInternalError());
883 }
884 
885 
886 
887 template <>
888 void
892 {
893  Assert (false, ExcImpossibleInDim(1));
894 }
895 
896 
897 
898 template <int dim>
899 void
902  typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const
903 {
904  // check whether center of object is at x==0, since then it belongs to the
905  // plane part of the boundary
906  const Point<dim> quad_center = face->center();
907  if (quad_center(0) == this->center(0))
908  StraightBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
909  else
910  HyperBallBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
911 }
912 
913 
914 /* ---------------------------------------------------------------------- */
915 
916 
917 
918 template <int dim>
920  :
921  HyperBallBoundary<dim>(center, 0.)
922 {
923  this->compute_radius_automatically=true;
924 }
925 
926 
927 /* ---------------------------------------------------------------------- */
928 
929 
930 
931 
932 template <int dim>
934  const double inner_radius,
935  const double outer_radius)
936  :
937  HyperShellBoundary<dim> (center),
938  inner_radius (inner_radius),
939  outer_radius (outer_radius)
940 {
941  if (dim > 2)
942  Assert ((inner_radius >= 0) &&
943  (outer_radius > 0) &&
944  (outer_radius > inner_radius),
945  ExcMessage ("Inner and outer radii must be specified explicitly in 3d."));
946 }
947 
948 
949 
950 template <int dim>
954 {
955  switch (dim)
956  {
957  // in 2d, first check whether the two end points of the line are on the
958  // axis of symmetry. if so, then return the mid point
959  case 2:
960  {
961  if ((line->vertex(0)(0) == this->center(0))
962  &&
963  (line->vertex(1)(0) == this->center(0)))
964  return (line->vertex(0) + line->vertex(1))/2;
965  else
966  // otherwise we are on the outer or inner part of the shell. proceed
967  // as in the base class
969  }
970 
971  // in 3d, a line is a straight line if it is on the symmetry plane and if
972  // not both of its end points are on either the inner or outer sphere
973  case 3:
974  {
975 
976  if (((line->vertex(0)(0) == this->center(0))
977  &&
978  (line->vertex(1)(0) == this->center(0)))
979  &&
980  !(((std::fabs (line->vertex(0).distance (this->center)
981  - inner_radius) < 1e-12 * outer_radius)
982  &&
983  (std::fabs (line->vertex(1).distance (this->center)
984  - inner_radius) < 1e-12 * outer_radius))
985  ||
986  ((std::fabs (line->vertex(0).distance (this->center)
987  - outer_radius) < 1e-12 * outer_radius)
988  &&
989  (std::fabs (line->vertex(1).distance (this->center)
990  - outer_radius) < 1e-12 * outer_radius))))
991  return (line->vertex(0) + line->vertex(1))/2;
992  else
993  // otherwise we are on the outer or inner part of the shell. proceed
994  // as in the base class
996  }
997 
998  default:
999  Assert (false, ExcNotImplemented());
1000  }
1001 
1002  return Point<dim>();
1003 }
1004 
1005 
1006 
1007 template <>
1008 Point<1>
1010 get_new_point_on_quad (const Triangulation<1>::quad_iterator &) const
1011 {
1012  Assert (false, ExcInternalError());
1013  return Point<1>();
1014 }
1015 
1016 
1017 
1018 
1019 template <int dim>
1020 Point<dim>
1023 {
1024  // if this quad is on the symmetry plane, take the center point and project
1025  // it outward to the same radius as the centers of the two radial lines
1026  if ((quad->vertex(0)(0) == this->center(0)) &&
1027  (quad->vertex(1)(0) == this->center(0)) &&
1028  (quad->vertex(2)(0) == this->center(0)) &&
1029  (quad->vertex(3)(0) == this->center(0)))
1030  {
1031  const Point<dim> quad_center = (quad->vertex(0) + quad->vertex(1) +
1032  quad->vertex(2) + quad->vertex(3) )/4;
1033  const Tensor<1,dim> quad_center_offset = quad_center - this->center;
1034 
1035 
1036  if (std::fabs (quad->line(0)->center().distance(this->center) -
1037  quad->line(1)->center().distance(this->center))
1038  < 1e-12 * outer_radius)
1039  {
1040  // lines 0 and 1 are radial
1041  const double needed_radius
1042  = quad->line(0)->center().distance(this->center);
1043 
1044  return (this->center +
1045  quad_center_offset/quad_center_offset.norm() * needed_radius);
1046  }
1047  else if (std::fabs (quad->line(2)->center().distance(this->center) -
1048  quad->line(3)->center().distance(this->center))
1049  < 1e-12 * outer_radius)
1050  {
1051  // lines 2 and 3 are radial
1052  const double needed_radius
1053  = quad->line(2)->center().distance(this->center);
1054 
1055  return (this->center +
1056  quad_center_offset/quad_center_offset.norm() * needed_radius);
1057  }
1058  else
1059  Assert (false, ExcInternalError());
1060  }
1061 
1062  // otherwise we are on the outer or inner part of the shell. proceed as in
1063  // the base class
1065 }
1066 
1067 
1068 
1069 template <int dim>
1070 void
1073  std::vector<Point<dim> > &points) const
1074 {
1075  switch (dim)
1076  {
1077  // in 2d, first check whether the two end points of the line are on the
1078  // axis of symmetry. if so, then return the mid point
1079  case 2:
1080  {
1081  if ((line->vertex(0)(0) == this->center(0))
1082  &&
1083  (line->vertex(1)(0) == this->center(0)))
1085  else
1086  // otherwise we are on the outer or inner part of the shell. proceed
1087  // as in the base class
1089  break;
1090  }
1091 
1092  // in 3d, a line is a straight line if it is on the symmetry plane and if
1093  // not both of its end points are on either the inner or outer sphere
1094  case 3:
1095  {
1096  if (((line->vertex(0)(0) == this->center(0))
1097  &&
1098  (line->vertex(1)(0) == this->center(0)))
1099  &&
1100  !(((std::fabs (line->vertex(0).distance (this->center)
1101  - inner_radius) < 1e-12 * outer_radius)
1102  &&
1103  (std::fabs (line->vertex(1).distance (this->center)
1104  - inner_radius) < 1e-12 * outer_radius))
1105  ||
1106  ((std::fabs (line->vertex(0).distance (this->center)
1107  - outer_radius) < 1e-12 * outer_radius)
1108  &&
1109  (std::fabs (line->vertex(1).distance (this->center)
1110  - outer_radius) < 1e-12 * outer_radius))))
1112  else
1113  // otherwise we are on the outer or inner part of the shell. proceed
1114  // as in the base class
1116 
1117  break;
1118  }
1119 
1120  default:
1121  Assert (false, ExcNotImplemented());
1122  }
1123 }
1124 
1125 
1126 
1127 template <int dim>
1128 void
1131  std::vector<Point<dim> > &points) const
1132 {
1133  Assert (dim < 3, ExcNotImplemented());
1134 
1135  // check whether center of object is at x==0, since then it belongs to the
1136  // plane part of the boundary
1137  const Point<dim> quad_center = quad->center();
1138  if (quad_center(0) == this->center(0))
1140  else
1142 }
1143 
1144 
1145 
1146 template <>
1147 void
1149 get_intermediate_points_on_quad (const Triangulation<1>::quad_iterator &,
1150  std::vector<Point<1> > &) const
1151 {
1152  Assert (false, ExcInternalError());
1153 }
1154 
1155 
1156 
1157 template <>
1158 void
1162 {
1163  Assert (false, ExcImpossibleInDim(1));
1164 }
1165 
1166 
1167 
1168 
1169 
1170 template <int dim>
1171 void
1174  typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const
1175 {
1176  if (face->center()(0) == this->center(0))
1177  StraightBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
1178  else
1179  HyperShellBoundary<dim>::get_normals_at_vertices (face, face_vertex_normals);
1180 }
1181 
1182 
1183 
1184 
1185 template <int dim, int spacedim>
1187  const double r__)
1188  :
1189  R(R__),
1190  r(r__)
1191 {
1192  Assert (false, ExcNotImplemented());
1193 }
1194 
1195 
1196 
1197 template <>
1198 TorusBoundary<2,3>::TorusBoundary (const double R__,
1199  const double r__)
1200  :
1201  R(R__),
1202  r(r__)
1203 {
1204  Assert (R>r, ExcMessage("Outer radius must be greater than inner radius."));
1205 }
1206 
1207 
1208 
1209 template <int dim, int spacedim>
1210 double
1212  const double x,
1213  const double y) const
1214 {
1215  if (y>=0)
1216  {
1217  if (x >=0)
1218  return angle;
1219 
1220  return numbers::PI-angle;
1221  }
1222 
1223  if (x <=0)
1224  return numbers::PI+angle;
1225 
1226  return 2.0*numbers::PI-angle;
1227 }
1228 
1229 
1230 
1231 template <>
1232 Point<3>
1233 TorusBoundary<2,3>::get_real_coord (const Point<2> &surfP) const
1234 {
1235  const double theta=surfP(0);
1236  const double phi=surfP(1);
1237 
1238  return Point<3> ((R+r*std::cos(phi))*std::cos(theta),
1239  r*std::sin(phi),
1240  (R+r*std::cos(phi))*std::sin(theta));
1241 }
1242 
1243 
1244 
1245 template <>
1246 Point<2>
1248 {
1249  const double phi=std::asin(std::abs(p(1))/r);
1250  const double Rr_2=p(0)*p(0)+p(2)*p(2);
1251 
1252  Point<2> surfP;
1253  surfP(1)=get_correct_angle(phi,Rr_2-R*R,p(1));//phi
1254 
1255  if (std::abs(p(0))<1.E-5)
1256  {
1257  if (p(2)>=0)
1258  surfP(0) = numbers::PI*0.5;
1259  else
1260  surfP(0) = -numbers::PI*0.5;
1261  }
1262  else
1263  {
1264  const double theta = std::atan(std::abs(p(2)/p(0)));
1265  surfP(0)= get_correct_angle(theta,p(0),p(2));
1266  }
1267 
1268  return surfP;
1269 }
1270 
1271 
1272 
1273 template <>
1274 Point<3>
1275 TorusBoundary<2,3>::get_new_point_on_line (const Triangulation<2,3>::line_iterator &line) const
1276 {
1277  //Just get the average
1278  Point<2> p0=get_surf_coord(line->vertex(0));
1279  Point<2> p1=get_surf_coord(line->vertex(1));
1280 
1281  Point<2> middle(0,0);
1282 
1283  //Take care for periodic conditions, For instance phi0= 0, phi1= 3/2*Pi
1284  //middle has to be 7/4*Pi not 3/4*Pi. This also works for -Pi/2 + Pi, middle
1285  //is 5/4*Pi
1286  for (unsigned int i=0; i<2; i++)
1287  if (std::abs(p0(i)-p1(i))> numbers::PI)
1288  middle(i)=2*numbers::PI;
1289 
1290  middle+= p0 + p1;
1291  middle*=0.5;
1292 
1293  Point<3> midReal=get_real_coord(middle);
1294  return midReal;
1295 }
1296 
1297 
1298 
1299 template <>
1300 Point<3>
1301 TorusBoundary<2,3>::get_new_point_on_quad (const Triangulation<2,3>::quad_iterator &quad) const
1302 {
1303  //Just get the average
1304  Point<2> p[4];
1305 
1306  for (unsigned int i=0; i<4; i++)
1307  p[i]=get_surf_coord(quad->vertex(i));
1308 
1309  Point<2> middle(0,0);
1310 
1311  //Take care for periodic conditions, see get_new_point_on_line() above
1312  //For instance phi0= 0, phi1= 3/2*Pi middle has to be 7/4*Pi not 3/4*Pi
1313  //This also works for -Pi/2 + Pi + Pi- Pi/2, middle is 5/4*Pi
1314  for (unsigned int i=0; i<2; i++)
1315  for (unsigned int j=1; j<4; j++)
1316  {
1317  if (std::abs(p[0](i)-p[j](i))> numbers::PI)
1318  {
1319  middle(i)+=2*numbers::PI;
1320  }
1321  }
1322 
1323  for (unsigned int i=0; i<4; i++)
1324  middle+=p[i];
1325 
1326  middle*= 0.25;
1327 
1328  return get_real_coord(middle);
1329 }
1330 
1331 
1332 
1333 //Normal field without unit length
1334 template <>
1335 Point<3>
1337 {
1338 
1339  Point<3> n;
1340  double theta=surfP[0];
1341  double phi=surfP[1];
1342 
1343  double f=R+r*std::cos(phi);
1344 
1345  n[0]=r*std::cos(phi)*std::cos(theta)*f;
1346  n[1]=r*std::sin(phi)*f;
1347  n[2]=r*std::sin(theta)*std::cos(phi)*f;
1348 
1349  return n;
1350 }
1351 
1352 
1353 
1354 //Normal field without unit length
1355 template <>
1356 Point<3>
1358 {
1359 
1360  Point<2> surfP=get_surf_coord(p);
1361  return get_surf_norm_from_sp(surfP);
1362 
1363 }
1364 
1365 
1366 
1367 template<>
1368 void
1370 get_intermediate_points_on_line (const Triangulation<2, 3>::line_iterator &line,
1371  std::vector< Point< 3 > > &points) const
1372 {
1373  //Almost the same implementation as StraightBoundary<2,3>
1374  unsigned int npoints=points.size();
1375  if (npoints==0) return;
1376 
1377  Point<2> p[2];
1378 
1379  for (unsigned int i=0; i<2; i++)
1380  p[i]=get_surf_coord(line->vertex(i));
1381 
1382  unsigned int offset[2];
1383  offset[0]=0;
1384  offset[1]=0;
1385 
1386  //Take care for periodic conditions & negative angles, see
1387  //get_new_point_on_line() above. Because we dont have a symmetric
1388  //interpolation (just the middle) we need to add 2*Pi to each almost zero
1389  //and negative angles.
1390  for (unsigned int i=0; i<2; i++)
1391  for (unsigned int j=1; j<2; j++)
1392  {
1393  if (std::abs(p[0](i)-p[j](i))> numbers::PI)
1394  {
1395  offset[i]++;
1396  break;
1397  }
1398  }
1399 
1400  for (unsigned int i=0; i<2; i++)
1401  for (unsigned int j=0; j<2; j++)
1402  if (p[j](i)<1.E-12 ) //Take care for periodic conditions & negative angles
1403  p[j](i)+=2*numbers::PI*offset[i];
1404 
1405 
1406  Point<2> target;
1407  const std::vector<Point<1> > &line_points = this->get_line_support_points(npoints);
1408  for (unsigned int i=0; i<npoints; i++)
1409  {
1410  const double x = line_points[i+1][0];
1411  target= (1-x)*p[0] + x*p[1];
1412  points[i]=get_real_coord(target);
1413  }
1414 }
1415 
1416 
1417 
1418 template<>
1419 void
1421 get_intermediate_points_on_quad (const Triangulation< 2, 3 >::quad_iterator &quad,
1422  std::vector< Point< 3 > > &points )const
1423 {
1424  //Almost the same implementation as StraightBoundary<2,3>
1425  const unsigned int n=points.size(),
1426  m=static_cast<unsigned int>(std::sqrt(static_cast<double>(n)));
1427  // is n a square number
1428  Assert(m*m==n, ExcInternalError());
1429 
1430  Point<2> p[4];
1431 
1432  for (unsigned int i=0; i<4; i++)
1433  p[i]=get_surf_coord(quad->vertex(i));
1434 
1435  Point<2> target;
1436  unsigned int offset[2];
1437  offset[0]=0;
1438  offset[1]=0;
1439 
1440  //Take care for periodic conditions & negative angles, see
1441  //get_new_point_on_line() above. Because we dont have a symmetric
1442  //interpolation (just the middle) we need to add 2*Pi to each almost zero
1443  //and negative angles.
1444  for (unsigned int i=0; i<2; i++)
1445  for (unsigned int j=1; j<4; j++)
1446  {
1447  if (std::abs(p[0](i)-p[j](i))> numbers::PI)
1448  {
1449  offset[i]++;
1450  break;
1451  }
1452  }
1453 
1454  for (unsigned int i=0; i<2; i++)
1455  for (unsigned int j=0; j<4; j++)
1456  if (p[j](i)<1.E-12 ) //Take care for periodic conditions & negative angles
1457  p[j](i)+=2*numbers::PI*offset[i];
1458 
1459  const std::vector<Point<1> > &line_points = this->get_line_support_points(m);
1460  for (unsigned int i=0; i<m; ++i)
1461  {
1462  const double y=line_points[i+1][0];
1463  for (unsigned int j=0; j<m; ++j)
1464  {
1465  const double x=line_points[j+1][0];
1466  target=((1-x) * p[0] +
1467  x * p[1]) * (1-y) +
1468  ((1-x) * p[2] +
1469  x * p[3]) * y;
1470 
1471  points[i*m+j]=get_real_coord(target);
1472  }
1473  }
1474 }
1475 
1476 
1477 
1478 template<>
1479 void
1482  Boundary<2,3>::FaceVertexNormals &face_vertex_normals) const
1483 {
1484  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_face; i++)
1485  face_vertex_normals[i]=get_surf_norm(face->vertex(i));
1486 }
1487 
1488 
1489 
1490 // explicit instantiations
1491 #include "tria_boundary_lib.inst"
1492 
1493 DEAL_II_NAMESPACE_CLOSE
const std::vector< Point< 1 > > & get_line_support_points(const unsigned int n_intermediate_points) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
void get_intermediate_points_between_points(const Point< spacedim > &p0, const Point< spacedim > &p1, std::vector< Point< spacedim > > &points) const
void get_intermediate_points_between_points(const Point< dim > &p0, const Point< dim > &p1, std::vector< Point< dim > > &points) const
virtual Point< dim > get_new_point_on_line(const typename Triangulation< dim >::line_iterator &line) const
virtual void get_normals_at_vertices(const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const
::ExceptionBase & ExcMessage(std::string arg1)
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
const Point< spacedim > point_on_axis
const Point< spacedim > direction
const double R
const Point< dim > x_1
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
virtual Point< dim > get_new_point_on_quad(const typename Triangulation< dim >::quad_iterator &quad) const
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:1047
virtual void get_intermediate_points_on_line(const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const
Definition: point.h:89
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const
double get_radius(const Point< dim > x) const
Point< dim > get_surf_coord(const Point< spacedim > &p) const
double get_correct_angle(const double angle, const double x, const double y) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
virtual Point< dim > get_new_point_on_line(const typename Triangulation< dim >::line_iterator &line) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
virtual Point< dim > get_new_point_on_line(const typename Triangulation< dim >::line_iterator &line) const
static Point< spacedim > get_axis_vector(const unsigned int axis)
static const double PI
Definition: numbers.h:74
HalfHyperShellBoundary(const Point< dim > &center=Point< dim >(), const double inner_radius=-1, const double outer_radius=-1)
Definition: tria.h:46
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
std::vector< std_cxx11::shared_ptr< QGaussLobatto< 1 > > > points
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const
virtual Point< dim > get_new_point_on_quad(const typename Triangulation< dim >::quad_iterator &quad) const
HalfHyperBallBoundary(const Point< dim > p=Point< dim >(), const double radius=1.0)
#define Assert(cond, exc)
Definition: exceptions.h:294
Point< spacedim > get_surf_norm(const Point< spacedim > &p) const
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
ConeBoundary(const double radius_0, const double radius_1, const Point< dim > x_0, const Point< dim > x_1)
const Point< spacedim > center
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const
const double radius_1
virtual Point< dim > get_new_point_on_quad(const typename Triangulation< dim >::quad_iterator &quad) const
virtual void get_normals_at_vertices(const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const
double get_radius() const
HyperShellBoundary(const Point< dim > &center=Point< dim >())
Point< spacedim > get_surf_norm_from_sp(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
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 void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const
TorusBoundary(const double R, const double r)
const Point< dim > x_0
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const
HyperBallBoundary(const Point< spacedim > p=Point< spacedim >(), const double radius=1.0)
CylinderBoundary(const double radius=1.0, const unsigned int axis=0)
const double radius_0
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const
Point< spacedim > get_center() const
virtual void get_intermediate_points_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
virtual void get_normals_at_vertices(const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const