16#ifndef dealii_tensor_product_manifold_h
17#define dealii_tensor_product_manifold_h
71 :
public ChartManifold<dim, spacedim_A + spacedim_B, chartdim_A + chartdim_B>
78 static const unsigned int chartdim = chartdim_A + chartdim_B;
83 static const unsigned int spacedim = spacedim_A + spacedim_B;
95 virtual std::unique_ptr<Manifold<dim, spacedim_A + spacedim_B>>
96 clone()
const override;
117 std::unique_ptr<const ChartManifold<dim_A, spacedim_A, chartdim_A>>
120 std::unique_ptr<const ChartManifold<dim_B, spacedim_B, chartdim_B>>
132 namespace TensorProductManifoldImplementation
134 template <
int dim1,
int dim2>
139 for (
unsigned int d = 0; d < dim1; ++d)
141 for (
unsigned int d = 0; d < dim2; ++d)
146 template <
int dim1,
int dim2>
151 for (
unsigned int d = 0; d < dim1; ++d)
153 for (
unsigned int d = 0; d < dim2; ++d)
158 template <
int dim1,
int dim2>
164 for (
unsigned int d = 0; d < dim1; ++d)
166 for (
unsigned int d = 0; d < dim2; ++d)
167 p2[d] = source[dim1 + d];
187 TensorProductManifold(
190 :
ChartManifold<dim, spacedim_A + spacedim_B, chartdim_A + chartdim_B>(
191 internal::TensorProductManifoldImplementation::concat(
192 manifold_A.get_periodicity(),
193 manifold_B.get_periodicity()))
194 , manifold_A(
Utilities::dynamic_unique_cast<
196 Manifold<dim_A, spacedim_A>>(manifold_A.clone()))
197 , manifold_B(
Utilities::dynamic_unique_cast<
199 Manifold<dim_B, spacedim_B>>(manifold_B.clone()))
209std::unique_ptr<Manifold<dim, spacedim_A + spacedim_B>>
216 chartdim_B>::clone()
const
224 chartdim_B>>(*manifold_A,
241 chartdim_B>::chartdim>
256 chartdim_B>::spacedim> &space_point)
const
284 chartdim_B>::spacedim>
299 chartdim_B>::chartdim> &chart_point)
const
328 chartdim_B>::chartdim,
335 chartdim_B>::spacedim>
344 push_forward_gradient(
351 chartdim_B>::chartdim> &chart_point)
const
360 manifold_A->push_forward_gradient(chart_point_A);
362 manifold_B->push_forward_gradient(chart_point_B);
366 for (
unsigned int i = 0; i < chartdim_A; ++i)
367 for (
unsigned int j = 0; j < spacedim_A; ++j)
368 result[j][i] = result_A[j][i];
369 for (
unsigned int i = 0; i < chartdim_B; ++i)
370 for (
unsigned int j = 0; j < spacedim_B; ++j)
371 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)