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> 24 DEAL_II_NAMESPACE_OPEN
31 template <
int dim,
int spacedim>
36 template <
int dim,
int spacedim>
47 template <
int dim,
int spacedim>
57 template <
int dim,
int spacedim>
68 get_intermediate_points_on_line (face, points);
71 get_intermediate_points_on_quad (face, points);
109 template <
int dim,
int spacedim>
126 template <
int dim,
int spacedim>
143 template <
int dim,
int spacedim>
160 template <
int dim,
int spacedim>
161 const std::vector<Point<1> > &
165 if (points.size() <= n_intermediate_points ||
166 points[n_intermediate_points].get() == 0)
169 if (points.size() <= n_intermediate_points)
170 points.resize(n_intermediate_points+1);
173 if (points[n_intermediate_points].
get() == 0)
175 std_cxx11::shared_ptr<QGaussLobatto<1> >
177 points[n_intermediate_points] = quadrature;
180 return points[n_intermediate_points]->get_points();
189 template <
int dim,
int spacedim>
194 template <
int dim,
int spacedim>
199 return (line->vertex(0) + line->vertex(1)) / 2;
210 compute_new_point_on_quad (
const typename Triangulation<dim, 3>::quad_iterator &quad)
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;
284 template <
int dim,
int spacedim>
298 return compute_new_point_on_quad<2> (quad);
308 return compute_new_point_on_quad<3> (quad);
313 template <
int dim,
int spacedim>
319 const unsigned int n=points.size();
324 const std::vector<Point<1> > &line_points = this->get_line_support_points(n);
330 for (
unsigned int i=0; i<n; ++i)
332 const double x = line_points[1+i][0];
333 points[i] = (1-x)*vertices[0] + x*vertices[1];
340 template <
int dim,
int spacedim>
355 std::vector<
Point<3> > &points)
const 357 const unsigned int spacedim = 3;
359 const unsigned int n=points.size(),
360 m=
static_cast<unsigned int>(std::sqrt(static_cast<double>(n)));
364 const std::vector<Point<1> > &line_points = this->get_line_support_points(m);
372 for (
unsigned int i=0; i<m; ++i)
374 const double y=line_points[1+i][0];
375 for (
unsigned int j=0; j<m; ++j)
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;
392 std::vector<
Point<3> > &points)
const 394 const unsigned int spacedim = 3;
396 const unsigned int n=points.size(),
397 m=
static_cast<unsigned int>(std::sqrt(static_cast<double>(n)));
401 const std::vector<Point<1> > &line_points = this->get_line_support_points(m);
409 for (
unsigned int i=0; i<m; ++i)
411 const double y=line_points[1+i][0];
412 for (
unsigned int j=0; j<m; ++j)
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;
467 normalized_alternating_product (
const Tensor<1,2> (&basis_vectors)[1])
469 Tensor<1,2> tmp = cross_product_2d (basis_vectors[0]);
470 return tmp/tmp.
norm();
476 normalized_alternating_product (
const Tensor<1,3> ( &)[1])
488 normalized_alternating_product (
const Tensor<1,3> (&basis_vectors)[2])
490 Tensor<1,3> tmp = cross_product_3d (basis_vectors[0], basis_vectors[1]);
491 return tmp/tmp.
norm();
498 template <
int dim,
int spacedim>
534 const unsigned int facedim = dim-1;
537 for (
unsigned int i=0; i<facedim; ++i)
542 const double eps = 1e-12;
544 unsigned int iteration = 0;
548 for (
unsigned int v=0; v<GeometryInfo<facedim>::vertices_per_cell; ++v)
549 F += face->vertex(v) * linear_fe.
shape_value(v, xi);
551 for (
unsigned int i=0; i<facedim; ++i)
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];
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];
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];
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."));
580 if (delta_xi.
norm() < eps)
588 return internal::normalized_alternating_product(grad_F);
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)
632 face_vertex_normals[vertex] =
Point<2>(tangent[1],
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)
645 face_vertex_normals[vertex] =
Point<3>(tangent[1],
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)
668 = { face->vertex(neighboring_vertices[vertex][0])
669 - face->vertex(vertex),
670 face->vertex(neighboring_vertices[vertex][1])
671 - face->vertex(vertex)
682 template <
int dim,
int spacedim>
699 p2 = line->vertex(1);
700 const double s = (trial_point-p1)*(p2-p1) / ((p2-p1)*(p2-p1));
701 return p1 + s*(p2-p1);
709 template <
typename Iterator,
int spacedim,
int dim>
711 compute_projection (
const Iterator &
object,
763 for (
unsigned int d=0; d<dim; ++d)
767 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
768 x_k += object->vertex(i) *
774 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
775 F_k += (x_k-y)*
object->vertex(i) *
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)
785 H_k += (
object->vertex(i) *
object->vertex(j)) * tmp;
792 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
793 x_k += object->vertex(i) *
796 if (delta_xi.
norm() < 1
e-5)
806 template <
typename Iterator>
808 compute_projection (
const Iterator &,
816 template <
typename Iterator>
818 compute_projection (
const Iterator &,
840 template <
int dim,
int spacedim>
849 return internal::compute_projection (quad, y,
855 template <
int dim,
int spacedim>
879 #include "tria_boundary.inst" 881 DEAL_II_NAMESPACE_CLOSE
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
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
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define Assert(cond, exc)
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)
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
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()