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/signaling_nan.h> 21 #include <deal.II/base/memory_consumption.h> 22 #include <deal.II/lac/full_matrix.h> 23 #include <deal.II/grid/tria.h> 24 #include <deal.II/grid/tria_iterator.h> 25 #include <deal.II/dofs/dof_accessor.h> 26 #include <deal.II/fe/mapping_cartesian.h> 27 #include <deal.II/fe/fe_values.h> 33 DEAL_II_NAMESPACE_OPEN
36 template<
int dim,
int spacedim>
41 template<
int dim,
int spacedim>
45 volume_element (
numbers::signaling_nan<double>()),
46 quadrature_points (q.get_points ())
51 template<
int dim,
int spacedim>
63 template <
int dim,
int spacedim>
72 template<
int dim,
int spacedim>
90 template<
int dim,
int spacedim>
107 template<
int dim,
int spacedim>
129 template<
int dim,
int spacedim>
152 template<
int dim,
int spacedim>
155 const unsigned int face_no,
156 const unsigned int sub_no,
183 cell->face(face_no)->has_children() ||
188 !cell->face(face_no)->has_children() ||
189 (sub_no < cell->face(face_no)->n_children()),
190 ExcIndexRange (sub_no, 0, cell->face(face_no)->n_children()));
245 cell->face_orientation(face_no),
246 cell->face_flip(face_no),
247 cell->face_rotation(face_no),
248 quadrature_points.size())
252 cell->face_orientation(face_no),
253 cell->face_flip(face_no),
254 cell->face_rotation(face_no),
255 quadrature_points.size(),
256 cell->subface_case(face_no))
259 for (
unsigned int i=0; i<quadrature_points.size(); ++i)
261 quadrature_points[i] = start;
262 for (
unsigned int d=0; d<dim; ++d)
290 std::fill (normal_vectors.begin(),
291 normal_vectors.end(),
305 std::fill (normal_vectors.begin(),
306 normal_vectors.end(),
322 std::fill (normal_vectors.begin(),
323 normal_vectors.end(),
336 template<
int dim,
int spacedim>
348 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
350 std::vector<Tensor<1,dim> > dummy;
362 double J = data.cell_extents[0];
363 for (
unsigned int d=1; d<dim; ++d)
364 J *= data.cell_extents[d];
365 data.volume_element = J;
367 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
374 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
377 for (
unsigned int j=0; j<dim; ++j)
378 output_data.
jacobians[i][j][j] = data.cell_extents[j];
421 for (
unsigned int j=0; j<dim; ++j)
425 return cell_similarity;
430 template<
int dim,
int spacedim>
434 const unsigned int face_no,
443 Assert (dynamic_cast<const InternalData *> (&internal_data) != 0,
445 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
456 for (
unsigned int d=0;
d<dim; ++
d)
458 J *= data.cell_extents[
d];
461 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
470 J = data.cell_extents[0];
471 for (
unsigned int d=1;
d<dim; ++
d)
472 J *= data.cell_extents[d];
473 data.volume_element = J;
477 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
480 for (
unsigned int d=0;
d<dim; ++
d)
481 output_data.
jacobians[i][d][d] = data.cell_extents[d];
512 for (
unsigned int d=0;
d<dim; ++
d)
519 template<
int dim,
int spacedim>
523 const unsigned int face_no,
524 const unsigned int subface_no,
532 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
542 for (
unsigned int d=0;
d<dim; ++
d)
544 J *= data.cell_extents[
d];
552 const unsigned int n_subfaces=
553 cell->face(face_no)->has_children() ?
554 cell->face(face_no)->n_children() :
556 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
566 J = data.cell_extents[0];
567 for (
unsigned int d=1;
d<dim; ++
d)
568 J *= data.cell_extents[d];
569 data.volume_element = J;
573 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
576 for (
unsigned int d=0;
d<dim; ++
d)
577 output_data.
jacobians[i][d][d] = data.cell_extents[d];
608 for (
unsigned int d=0;
d<dim; ++
d)
616 template<
int dim,
int spacedim>
625 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
629 switch (mapping_type)
636 for (
unsigned int i=0; i<output.size(); ++i)
637 for (
unsigned int d=0; d<dim; ++d)
647 for (
unsigned int i=0; i<output.size(); ++i)
648 for (
unsigned int d=0; d<dim; ++d)
659 for (
unsigned int i=0; i<output.size(); ++i)
660 for (
unsigned int d=0; d<dim; ++d)
671 template<
int dim,
int spacedim>
680 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
684 switch (mapping_type)
691 for (
unsigned int i=0; i<output.size(); ++i)
692 for (
unsigned int d1=0; d1<dim; ++d1)
693 for (
unsigned int d2=0; d2<dim; ++d2)
694 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
703 for (
unsigned int i=0; i<output.size(); ++i)
704 for (
unsigned int d1=0; d1<dim; ++d1)
705 for (
unsigned int d2=0; d2<dim; ++d2)
706 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
715 for (
unsigned int i=0; i<output.size(); ++i)
716 for (
unsigned int d1=0; d1<dim; ++d1)
717 for (
unsigned int d2=0; d2<dim; ++d2)
727 for (
unsigned int i=0; i<output.size(); ++i)
728 for (
unsigned int d1=0; d1<dim; ++d1)
729 for (
unsigned int d2=0; d2<dim; ++d2)
741 for (
unsigned int i=0; i<output.size(); ++i)
742 for (
unsigned int d1=0; d1<dim; ++d1)
743 for (
unsigned int d2=0; d2<dim; ++d2)
744 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
756 for (
unsigned int i=0; i<output.size(); ++i)
757 for (
unsigned int d1=0; d1<dim; ++d1)
758 for (
unsigned int d2=0; d2<dim; ++d2)
759 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
772 template<
int dim,
int spacedim>
782 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
786 switch (mapping_type)
793 for (
unsigned int i=0; i<output.size(); ++i)
794 for (
unsigned int d1=0; d1<dim; ++d1)
795 for (
unsigned int d2=0; d2<dim; ++d2)
796 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
805 for (
unsigned int i=0; i<output.size(); ++i)
806 for (
unsigned int d1=0; d1<dim; ++d1)
807 for (
unsigned int d2=0; d2<dim; ++d2)
808 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
817 for (
unsigned int i=0; i<output.size(); ++i)
818 for (
unsigned int d1=0; d1<dim; ++d1)
819 for (
unsigned int d2=0; d2<dim; ++d2)
829 for (
unsigned int i=0; i<output.size(); ++i)
830 for (
unsigned int d1=0; d1<dim; ++d1)
831 for (
unsigned int d2=0; d2<dim; ++d2)
843 for (
unsigned int i=0; i<output.size(); ++i)
844 for (
unsigned int d1=0; d1<dim; ++d1)
845 for (
unsigned int d2=0; d2<dim; ++d2)
846 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
858 for (
unsigned int i=0; i<output.size(); ++i)
859 for (
unsigned int d1=0; d1<dim; ++d1)
860 for (
unsigned int d2=0; d2<dim; ++d2)
861 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
873 template<
int dim,
int spacedim>
883 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
887 switch (mapping_type)
894 for (
unsigned int q=0; q<output.size(); ++q)
895 for (
unsigned int i=0; i<spacedim; ++i)
896 for (
unsigned int j=0; j<spacedim; ++j)
897 for (
unsigned int k=0; k<spacedim; ++k)
911 template<
int dim,
int spacedim>
921 Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
925 switch (mapping_type)
934 for (
unsigned int q=0; q<output.size(); ++q)
935 for (
unsigned int i=0; i<spacedim; ++i)
936 for (
unsigned int j=0; j<spacedim; ++j)
937 for (
unsigned int k=0; k<spacedim; ++k)
939 output[q][i][j][k] = input[q][i][j][k]
952 for (
unsigned int q=0; q<output.size(); ++q)
953 for (
unsigned int i=0; i<spacedim; ++i)
954 for (
unsigned int j=0; j<spacedim; ++j)
955 for (
unsigned int k=0; k<spacedim; ++k)
957 output[q][i][j][k] = input[q][i][j][k]
975 for (
unsigned int q=0; q<output.size(); ++q)
976 for (
unsigned int i=0; i<spacedim; ++i)
977 for (
unsigned int j=0; j<spacedim; ++j)
978 for (
unsigned int k=0; k<spacedim; ++k)
980 output[q][i][j][k] = input[q][i][j][k]
996 template <
int dim,
int spacedim>
1007 length[0] = cell->vertex(1)(0) - start(0);
1010 length[0] = cell->vertex(1)(0) - start(0);
1011 length[1] = cell->vertex(2)(1) - start(1);
1014 length[0] = cell->vertex(1)(0) - start(0);
1015 length[1] = cell->vertex(2)(1) - start(1);
1016 length[2] = cell->vertex(4)(2) - start(2);
1023 for (
unsigned int d=0; d<dim; ++d)
1024 p_real(d) += length[d]*p(d);
1031 template<
int dim,
int spacedim>
1038 if (dim != spacedim)
1047 real(0) /= cell->vertex(1)(0) - start(0);
1050 real(0) /= cell->vertex(1)(0) - start(0);
1051 real(1) /= cell->vertex(2)(1) - start(1);
1054 real(0) /= cell->vertex(1)(0) - start(0);
1055 real(1) /= cell->vertex(2)(1) - start(1);
1056 real(2) /= cell->vertex(4)(2) - start(2);
1065 template<
int dim,
int spacedim>
1075 #include "mapping_cartesian.inst" 1078 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 ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
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
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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)
static ::ExceptionBase & ExcNotImplemented()
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Covariant transformation.
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()
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)