16 #include <deal.II/base/std_cxx14/memory.h> 17 #include <deal.II/base/array_view.h> 18 #include <deal.II/base/tensor.h> 19 #include <deal.II/base/quadrature.h> 20 #include <deal.II/base/qprojector.h> 21 #include <deal.II/base/signaling_nan.h> 22 #include <deal.II/base/memory_consumption.h> 23 #include <deal.II/lac/full_matrix.h> 24 #include <deal.II/grid/tria.h> 25 #include <deal.II/grid/tria_iterator.h> 26 #include <deal.II/dofs/dof_accessor.h> 27 #include <deal.II/fe/mapping_cartesian.h> 28 #include <deal.II/fe/fe_values.h> 34 DEAL_II_NAMESPACE_OPEN
37 template <
int dim,
int spacedim>
42 template <
int dim,
int spacedim>
46 volume_element (
numbers::signaling_nan<double>()),
47 quadrature_points (q.get_points ())
52 template <
int dim,
int spacedim>
64 template <
int dim,
int spacedim>
73 template <
int dim,
int spacedim>
91 template <
int dim,
int spacedim>
92 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
96 auto data = std_cxx14::make_unique<InternalData> (q);
103 return std::move(data);
108 template <
int dim,
int spacedim>
109 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
122 data->update_each = update_flags;
124 return std::move(data);
129 template <
int dim,
int spacedim>
130 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
143 data->update_each = update_flags;
145 return std::move(data);
151 template <
int dim,
int spacedim>
154 const unsigned int face_no,
155 const unsigned int sub_no,
182 cell->face(face_no)->has_children() ||
187 !cell->face(face_no)->has_children() ||
188 (sub_no < cell->face(face_no)->n_children()),
189 ExcIndexRange (sub_no, 0, cell->face(face_no)->n_children()));
244 cell->face_orientation(face_no),
245 cell->face_flip(face_no),
246 cell->face_rotation(face_no),
247 quadrature_points.size())
251 cell->face_orientation(face_no),
252 cell->face_flip(face_no),
253 cell->face_rotation(face_no),
254 quadrature_points.size(),
255 cell->subface_case(face_no))
258 for (
unsigned int i=0; i<quadrature_points.size(); ++i)
260 quadrature_points[i] = start;
261 for (
unsigned int d=0; d<dim; ++d)
289 std::fill (normal_vectors.begin(),
290 normal_vectors.end(),
304 std::fill (normal_vectors.begin(),
305 normal_vectors.end(),
321 std::fill (normal_vectors.begin(),
322 normal_vectors.end(),
335 template <
int dim,
int spacedim>
347 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
349 std::vector<Tensor<1,dim> > dummy;
361 double J = data.cell_extents[0];
362 for (
unsigned int d=1; d<dim; ++d)
363 J *= data.cell_extents[d];
364 data.volume_element = J;
366 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
373 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
376 for (
unsigned int j=0; j<dim; ++j)
377 output_data.
jacobians[i][j][j] = data.cell_extents[j];
420 for (
unsigned int j=0; j<dim; ++j)
424 return cell_similarity;
429 template <
int dim,
int spacedim>
433 const unsigned int face_no,
442 Assert (dynamic_cast<const InternalData *> (&internal_data) !=
nullptr,
444 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
455 for (
unsigned int d=0;
d<dim; ++
d)
457 J *= data.cell_extents[
d];
460 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
469 J = data.cell_extents[0];
470 for (
unsigned int d=1;
d<dim; ++
d)
471 J *= data.cell_extents[d];
472 data.volume_element = J;
476 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
479 for (
unsigned int d=0;
d<dim; ++
d)
480 output_data.
jacobians[i][d][d] = data.cell_extents[d];
511 for (
unsigned int d=0;
d<dim; ++
d)
518 template <
int dim,
int spacedim>
522 const unsigned int face_no,
523 const unsigned int subface_no,
531 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
541 for (
unsigned int d=0;
d<dim; ++
d)
543 J *= data.cell_extents[
d];
551 const unsigned int n_subfaces=
552 cell->face(face_no)->has_children() ?
553 cell->face(face_no)->n_children() :
555 for (
unsigned int i=0; i<output_data.
JxW_values.size(); ++i)
565 J = data.cell_extents[0];
566 for (
unsigned int d=1;
d<dim; ++
d)
567 J *= data.cell_extents[d];
568 data.volume_element = J;
572 for (
unsigned int i=0; i<output_data.
jacobians.size(); ++i)
575 for (
unsigned int d=0;
d<dim; ++
d)
576 output_data.
jacobians[i][d][d] = data.cell_extents[d];
607 for (
unsigned int d=0;
d<dim; ++
d)
615 template <
int dim,
int spacedim>
624 Assert (dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
628 switch (mapping_type)
635 for (
unsigned int i=0; i<output.size(); ++i)
636 for (
unsigned int d=0; d<dim; ++d)
646 for (
unsigned int i=0; i<output.size(); ++i)
647 for (
unsigned int d=0; d<dim; ++d)
658 for (
unsigned int i=0; i<output.size(); ++i)
659 for (
unsigned int d=0; d<dim; ++d)
670 template <
int dim,
int spacedim>
679 Assert (dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
683 switch (mapping_type)
690 for (
unsigned int i=0; i<output.size(); ++i)
691 for (
unsigned int d1=0; d1<dim; ++d1)
692 for (
unsigned int d2=0; d2<dim; ++d2)
693 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
702 for (
unsigned int i=0; i<output.size(); ++i)
703 for (
unsigned int d1=0; d1<dim; ++d1)
704 for (
unsigned int d2=0; d2<dim; ++d2)
705 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
714 for (
unsigned int i=0; i<output.size(); ++i)
715 for (
unsigned int d1=0; d1<dim; ++d1)
716 for (
unsigned int d2=0; d2<dim; ++d2)
726 for (
unsigned int i=0; i<output.size(); ++i)
727 for (
unsigned int d1=0; d1<dim; ++d1)
728 for (
unsigned int d2=0; d2<dim; ++d2)
740 for (
unsigned int i=0; i<output.size(); ++i)
741 for (
unsigned int d1=0; d1<dim; ++d1)
742 for (
unsigned int d2=0; d2<dim; ++d2)
743 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
755 for (
unsigned int i=0; i<output.size(); ++i)
756 for (
unsigned int d1=0; d1<dim; ++d1)
757 for (
unsigned int d2=0; d2<dim; ++d2)
758 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
771 template <
int dim,
int spacedim>
781 Assert (dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
785 switch (mapping_type)
792 for (
unsigned int i=0; i<output.size(); ++i)
793 for (
unsigned int d1=0; d1<dim; ++d1)
794 for (
unsigned int d2=0; d2<dim; ++d2)
795 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
804 for (
unsigned int i=0; i<output.size(); ++i)
805 for (
unsigned int d1=0; d1<dim; ++d1)
806 for (
unsigned int d2=0; d2<dim; ++d2)
807 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
816 for (
unsigned int i=0; i<output.size(); ++i)
817 for (
unsigned int d1=0; d1<dim; ++d1)
818 for (
unsigned int d2=0; d2<dim; ++d2)
828 for (
unsigned int i=0; i<output.size(); ++i)
829 for (
unsigned int d1=0; d1<dim; ++d1)
830 for (
unsigned int d2=0; d2<dim; ++d2)
842 for (
unsigned int i=0; i<output.size(); ++i)
843 for (
unsigned int d1=0; d1<dim; ++d1)
844 for (
unsigned int d2=0; d2<dim; ++d2)
845 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
857 for (
unsigned int i=0; i<output.size(); ++i)
858 for (
unsigned int d1=0; d1<dim; ++d1)
859 for (
unsigned int d2=0; d2<dim; ++d2)
860 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2]
872 template <
int dim,
int spacedim>
882 Assert (dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
886 switch (mapping_type)
893 for (
unsigned int q=0; q<output.size(); ++q)
894 for (
unsigned int i=0; i<spacedim; ++i)
895 for (
unsigned int j=0; j<spacedim; ++j)
896 for (
unsigned int k=0; k<spacedim; ++k)
910 template <
int dim,
int spacedim>
920 Assert (dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
924 switch (mapping_type)
933 for (
unsigned int q=0; q<output.size(); ++q)
934 for (
unsigned int i=0; i<spacedim; ++i)
935 for (
unsigned int j=0; j<spacedim; ++j)
936 for (
unsigned int k=0; k<spacedim; ++k)
938 output[q][i][j][k] = input[q][i][j][k]
951 for (
unsigned int q=0; q<output.size(); ++q)
952 for (
unsigned int i=0; i<spacedim; ++i)
953 for (
unsigned int j=0; j<spacedim; ++j)
954 for (
unsigned int k=0; k<spacedim; ++k)
956 output[q][i][j][k] = input[q][i][j][k]
974 for (
unsigned int q=0; q<output.size(); ++q)
975 for (
unsigned int i=0; i<spacedim; ++i)
976 for (
unsigned int j=0; j<spacedim; ++j)
977 for (
unsigned int k=0; k<spacedim; ++k)
979 output[q][i][j][k] = input[q][i][j][k]
995 template <
int dim,
int spacedim>
1006 length[0] = cell->vertex(1)(0) - start(0);
1009 length[0] = cell->vertex(1)(0) - start(0);
1010 length[1] = cell->vertex(2)(1) - start(1);
1013 length[0] = cell->vertex(1)(0) - start(0);
1014 length[1] = cell->vertex(2)(1) - start(1);
1015 length[2] = cell->vertex(4)(2) - start(2);
1022 for (
unsigned int d=0; d<dim; ++d)
1023 p_real(d) += length[d]*p(d);
1030 template <
int dim,
int spacedim>
1037 if (dim != spacedim)
1046 real(0) /= cell->vertex(1)(0) - start(0);
1049 real(0) /= cell->vertex(1)(0) - start(0);
1050 real(1) /= cell->vertex(2)(1) - start(1);
1053 real(0) /= cell->vertex(1)(0) - start(0);
1054 real(1) /= cell->vertex(2)(1) - start(1);
1055 real(2) /= cell->vertex(4)(2) - start(2);
1064 template <
int dim,
int spacedim>
1065 std::unique_ptr<Mapping<dim, spacedim> >
1068 return std_cxx14::make_unique<MappingCartesian<dim,spacedim>>(*this);
1074 #include "mapping_cartesian.inst" 1077 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
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.
Outer normal vector, not normalized.
std::vector< Point< dim > > quadrature_points
Determinant of the Jacobian.
Transformed quadrature points.
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static DataSetDescriptor cell()
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
virtual std::size_t memory_consumption() const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
#define Assert(cond, exc)
Abstract base class for mapping classes.
Gradient of volume element.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 1, dim > cell_extents
static const unsigned int invalid_face_number
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
InternalData(const Quadrature< dim > &quadrature)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual bool preserves_vertex_locations() const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
static ::ExceptionBase & ExcNotImplemented()
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 override
Covariant transformation.
double weight(const unsigned int i) const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
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)