15#ifndef dealii_tensor_product_manifold_h
16#define dealii_tensor_product_manifold_h
70 :
public ChartManifold<dim, spacedim_A + spacedim_B, chartdim_A + chartdim_B>
77 static const unsigned int chartdim = chartdim_A + chartdim_B;
82 static const unsigned int spacedim = spacedim_A + spacedim_B;
94 virtual std::unique_ptr<Manifold<dim, spacedim_A + spacedim_B>>
95 clone()
const override;
116 std::unique_ptr<const ChartManifold<dim_A, spacedim_A, chartdim_A>>
119 std::unique_ptr<const ChartManifold<dim_B, spacedim_B, chartdim_B>>
131 namespace TensorProductManifoldImplementation
133 template <
int dim1,
int dim2>
138 for (
unsigned int d = 0; d < dim1; ++d)
140 for (
unsigned int d = 0; d < dim2; ++d)
145 template <
int dim1,
int dim2>
150 for (
unsigned int d = 0; d < dim1; ++d)
152 for (
unsigned int d = 0; d < dim2; ++d)
157 template <
int dim1,
int dim2>
163 for (
unsigned int d = 0; d < dim1; ++d)
165 for (
unsigned int d = 0; d < dim2; ++d)
166 p2[d] = source[dim1 + d];
186 TensorProductManifold(
189 :
ChartManifold<dim, spacedim_A + spacedim_B, chartdim_A + chartdim_B>(
190 internal::TensorProductManifoldImplementation::concat(
191 manifold_A.get_periodicity(),
192 manifold_B.get_periodicity()))
193 , manifold_A(
Utilities::dynamic_unique_cast<
195 Manifold<dim_A, spacedim_A>>(manifold_A.clone()))
196 , manifold_B(
Utilities::dynamic_unique_cast<
198 Manifold<dim_B, spacedim_B>>(manifold_B.clone()))
208std::unique_ptr<Manifold<dim, spacedim_A + spacedim_B>>
215 chartdim_B>::clone()
const
223 chartdim_B>>(*manifold_A,
240 chartdim_B>::chartdim>
255 chartdim_B>::spacedim> &space_point)
const
283 chartdim_B>::spacedim>
298 chartdim_B>::chartdim> &chart_point)
const
327 chartdim_B>::chartdim,
334 chartdim_B>::spacedim>
343 push_forward_gradient(
350 chartdim_B>::chartdim> &chart_point)
const
359 manifold_A->push_forward_gradient(chart_point_A);
361 manifold_B->push_forward_gradient(chart_point_B);
365 for (
unsigned int i = 0; i < chartdim_A; ++i)
366 for (
unsigned int j = 0; j < spacedim_A; ++j)
367 result[j][i] = result_A[j][i];
368 for (
unsigned int i = 0; i < chartdim_B; ++i)
369 for (
unsigned int j = 0; j < spacedim_B; ++j)
370 result[j + spacedim_A][i + chartdim_A] = result_B[j][i];
Tensor product manifold of two ChartManifolds.
static const unsigned int chartdim
virtual Point< spacedim > push_forward(const Point< chartdim > &chart_point) const override
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const override
static const unsigned int spacedim
virtual Point< chartdim > pull_back(const Point< spacedim > &space_point) const override
virtual std::unique_ptr< Manifold< dim, spacedim_A+spacedim_B > > clone() const override
std::unique_ptr< const ChartManifold< dim_B, spacedim_B, chartdim_B > > manifold_B
std::unique_ptr< const ChartManifold< dim_A, spacedim_A, chartdim_A > > manifold_A
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
Tensor< 1, dim1+dim2 > concat(const Tensor< 1, dim1 > &p1, const Tensor< 1, dim2 > &p2)
void split_point(const Point< dim1+dim2 > &source, Point< dim1 > &p1, Point< dim2 > &p2)