17 #include <deal.II/base/tensor.h> 18 #include <deal.II/base/quadrature.h> 19 #include <deal.II/base/qprojector.h> 20 #include <deal.II/base/memory_consumption.h> 21 #include <deal.II/lac/full_matrix.h> 22 #include <deal.II/grid/tria.h> 23 #include <deal.II/grid/tria_iterator.h> 24 #include <deal.II/dofs/dof_accessor.h> 25 #include <deal.II/fe/mapping_cartesian.h> 26 #include <deal.II/fe/fe_values.h> 32 DEAL_II_NAMESPACE_OPEN
35 template<
int dim,
int spacedim>
40 template<
int dim,
int spacedim>
43 quadrature_points (q.get_points ())
48 template<
int dim,
int spacedim>
60 template <
int dim,
int spacedim>
69 template<
int dim,
int spacedim>
87 template<
int dim,
int spacedim>
104 template<
int dim,
int spacedim>
126 template<
int dim,
int spacedim>
149 template<
int dim,
int spacedim>
152 const unsigned int face_no,
153 const unsigned int sub_no,
154 const CellSimilarity::Similarity cell_similarity,
180 cell->face(face_no)->has_children() ||
182 ExcIndexRange (sub_no, 0,
185 !cell->face(face_no)->has_children() ||
186 (sub_no < cell->face(face_no)->n_children()),
187 ExcIndexRange (sub_no, 0, cell->face(face_no)->n_children()));
205 if (cell_similarity != CellSimilarity::translation)
222 Assert(
false, ExcNotImplemented());
242 cell->face_orientation(face_no),
243 cell->face_flip(face_no),
244 cell->face_rotation(face_no),
249 cell->face_orientation(face_no),
250 cell->face_flip(face_no),
251 cell->face_rotation(face_no),
253 cell->subface_case(face_no))
259 for (
unsigned int d=0; d<dim; ++d)
287 std::fill (normal_vectors.begin(),
288 normal_vectors.end(),
302 std::fill (normal_vectors.begin(),
303 normal_vectors.end(),
319 std::fill (normal_vectors.begin(),
320 normal_vectors.end(),
326 Assert (
false, ExcNotImplemented());
333 template<
int dim,
int spacedim>
334 CellSimilarity::Similarity
337 const CellSimilarity::Similarity cell_similarity,
344 Assert (dynamic_cast<const InternalData *> (&internal_data) != 0, ExcInternalError());
347 std::vector<Tensor<1,dim> > dummy;
357 if (cell_similarity != CellSimilarity::translation)
360 for (
unsigned int d=1; d<dim; ++d)
364 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
370 if (cell_similarity != CellSimilarity::translation)
371 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
374 for (
unsigned int j=0; j<dim; ++j)
380 if (cell_similarity != CellSimilarity::translation)
385 if (cell_similarity != CellSimilarity::translation)
392 if (cell_similarity != CellSimilarity::translation)
397 if (cell_similarity != CellSimilarity::translation)
402 if (cell_similarity != CellSimilarity::translation)
407 if (cell_similarity != CellSimilarity::translation)
414 if (cell_similarity != CellSimilarity::translation)
418 for (
unsigned int j=0; j<dim; ++j)
422 return cell_similarity;
427 template<
int dim,
int spacedim>
431 const unsigned int face_no,
440 Assert (dynamic_cast<const InternalData *> (&internal_data) != 0,
445 CellSimilarity::none,
453 for (
unsigned int d=0; d<dim; ++d)
458 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
468 for (
unsigned int d=1; d<dim; ++d)
474 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
477 for (
unsigned int d=0; d<dim; ++d)
509 for (
unsigned int d=0; d<dim; ++d)
516 template<
int dim,
int spacedim>
520 const unsigned int face_no,
521 const unsigned int subface_no,
528 Assert (dynamic_cast<const InternalData *> (&internal_data) != 0, ExcInternalError());
531 compute_fill (cell, face_no, subface_no, CellSimilarity::none,
539 for (
unsigned int d=0; d<dim; ++d)
549 const unsigned int n_subfaces=
550 cell->face(face_no)->has_children() ?
551 cell->face(face_no)->n_children() :
553 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
564 for (
unsigned int d=1; d<dim; ++d)
570 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
573 for (
unsigned int d=0; d<dim; ++d)
605 for (
unsigned int d=0; d<dim; ++d)
613 template<
int dim,
int spacedim>
622 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
626 switch (mapping_type)
633 for (
unsigned int i=0; i<output.size(); ++i)
634 for (
unsigned int d=0; d<dim; ++d)
644 for (
unsigned int i=0; i<output.size(); ++i)
645 for (
unsigned int d=0; d<dim; ++d)
656 for (
unsigned int i=0; i<output.size(); ++i)
657 for (
unsigned int d=0; d<dim; ++d)
662 Assert(
false, ExcNotImplemented());
668 template<
int dim,
int spacedim>
677 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
681 switch (mapping_type)
688 for (
unsigned int i=0; i<output.size(); ++i)
689 for (
unsigned int d1=0; d1<dim; ++d1)
690 for (
unsigned int d2=0; d2<dim; ++d2)
691 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
700 for (
unsigned int i=0; i<output.size(); ++i)
701 for (
unsigned int d1=0; d1<dim; ++d1)
702 for (
unsigned int d2=0; d2<dim; ++d2)
703 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
712 for (
unsigned int i=0; i<output.size(); ++i)
713 for (
unsigned int d1=0; d1<dim; ++d1)
714 for (
unsigned int d2=0; d2<dim; ++d2)
724 for (
unsigned int i=0; i<output.size(); ++i)
725 for (
unsigned int d1=0; d1<dim; ++d1)
726 for (
unsigned int d2=0; d2<dim; ++d2)
738 for (
unsigned int i=0; i<output.size(); ++i)
739 for (
unsigned int d1=0; d1<dim; ++d1)
740 for (
unsigned int d2=0; d2<dim; ++d2)
741 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
753 for (
unsigned int i=0; i<output.size(); ++i)
754 for (
unsigned int d1=0; d1<dim; ++d1)
755 for (
unsigned int d2=0; d2<dim; ++d2)
756 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
762 Assert(
false, ExcNotImplemented());
769 template<
int dim,
int spacedim>
779 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
783 switch (mapping_type)
790 for (
unsigned int i=0; i<output.size(); ++i)
791 for (
unsigned int d1=0; d1<dim; ++d1)
792 for (
unsigned int d2=0; d2<dim; ++d2)
793 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
802 for (
unsigned int i=0; i<output.size(); ++i)
803 for (
unsigned int d1=0; d1<dim; ++d1)
804 for (
unsigned int d2=0; d2<dim; ++d2)
805 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
814 for (
unsigned int i=0; i<output.size(); ++i)
815 for (
unsigned int d1=0; d1<dim; ++d1)
816 for (
unsigned int d2=0; d2<dim; ++d2)
826 for (
unsigned int i=0; i<output.size(); ++i)
827 for (
unsigned int d1=0; d1<dim; ++d1)
828 for (
unsigned int d2=0; d2<dim; ++d2)
840 for (
unsigned int i=0; i<output.size(); ++i)
841 for (
unsigned int d1=0; d1<dim; ++d1)
842 for (
unsigned int d2=0; d2<dim; ++d2)
843 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
855 for (
unsigned int i=0; i<output.size(); ++i)
856 for (
unsigned int d1=0; d1<dim; ++d1)
857 for (
unsigned int d2=0; d2<dim; ++d2)
858 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
864 Assert(
false, ExcNotImplemented());
870 template<
int dim,
int spacedim>
880 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
884 switch (mapping_type)
891 for (
unsigned int q=0; q<output.size(); ++q)
892 for (
unsigned int i=0; i<spacedim; ++i)
893 for (
unsigned int j=0; j<spacedim; ++j)
894 for (
unsigned int k=0; k<spacedim; ++k)
903 Assert(
false, ExcNotImplemented());
908 template<
int dim,
int spacedim>
918 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
922 switch (mapping_type)
931 for (
unsigned int q=0; q<output.size(); ++q)
932 for (
unsigned int i=0; i<spacedim; ++i)
933 for (
unsigned int j=0; j<spacedim; ++j)
934 for (
unsigned int k=0; k<spacedim; ++k)
936 output[q][i][j][k] = input[q][i][j][k]
949 for (
unsigned int q=0; q<output.size(); ++q)
950 for (
unsigned int i=0; i<spacedim; ++i)
951 for (
unsigned int j=0; j<spacedim; ++j)
952 for (
unsigned int k=0; k<spacedim; ++k)
954 output[q][i][j][k] = input[q][i][j][k]
972 for (
unsigned int q=0; q<output.size(); ++q)
973 for (
unsigned int i=0; i<spacedim; ++i)
974 for (
unsigned int j=0; j<spacedim; ++j)
975 for (
unsigned int k=0; k<spacedim; ++k)
977 output[q][i][j][k] = input[q][i][j][k]
988 Assert(
false, ExcNotImplemented());
993 template <
int dim,
int spacedim>
1004 length[0] = cell->vertex(1)(0) - start(0);
1007 length[0] = cell->vertex(1)(0) - start(0);
1008 length[1] = cell->vertex(2)(1) - start(1);
1011 length[0] = cell->vertex(1)(0) - start(0);
1012 length[1] = cell->vertex(2)(1) - start(1);
1013 length[2] = cell->vertex(4)(2) - start(2);
1016 Assert(
false, ExcNotImplemented());
1020 for (
unsigned int d=0; d<dim; ++d)
1021 p_real(d) += length[d]*p(d);
1028 template<
int dim,
int spacedim>
1035 if (dim != spacedim)
1036 Assert(
false, ExcNotImplemented());
1044 real(0) /= cell->vertex(1)(0) - start(0);
1047 real(0) /= cell->vertex(1)(0) - start(0);
1048 real(1) /= cell->vertex(2)(1) - start(1);
1051 real(0) /= cell->vertex(1)(0) - start(0);
1052 real(1) /= cell->vertex(2)(1) - start(1);
1053 real(2) /= cell->vertex(4)(2) - start(2);
1056 Assert(
false, ExcNotImplemented());
1062 template<
int dim,
int spacedim>
1072 #include "mapping_cartesian.inst" 1075 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
void compute_fill(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const CellSimilarity::Similarity cell_similarity, const InternalData &data, std::vector< Point< dim > > &quadrature_points, std::vector< Tensor< 1, dim > > &normal_vectors) const
#define AssertDimension(dim1, dim2)
Contravariant transformation.
virtual Mapping< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const
Outer normal vector, not normalized.
bool preserves_vertex_locations() const
std::vector< Point< dim > > quadrature_points
virtual Mapping< dim, spacedim >::InternalDataBase * get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
Determinant of the Jacobian.
Transformed quadrature points.
static DataSetDescriptor cell()
virtual std::size_t memory_consumption() const
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
#define Assert(cond, exc)
Abstract base class for mapping classes.
Gradient of volume element.
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const
virtual Mapping< dim, spacedim >::InternalDataBase * get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
Tensor< 1, dim > cell_extents
static const unsigned int invalid_face_number
virtual Mapping< dim, spacedim > * clone() const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
InternalData(const Quadrature< dim > &quadrature)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Covariant transformation.
double weight(const unsigned int i) const
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)