36class DeformedCubeManifold :
public dealii::ChartManifold<dim, dim, dim>
39 DeformedCubeManifold(
double const left,
41 double const deformation,
42 unsigned int const frequency = 1)
43 : left(left), right(right), deformation(deformation), frequency(frequency)
48 push_forward(dealii::Point<dim>
const & chart_point)
const override
50 double sinval = deformation;
51 for(
unsigned int d = 0; d < dim; ++d)
53 std::sin(frequency * dealii::numbers::PI * (chart_point(d) - left) / (right - left));
54 dealii::Point<dim> space_point;
55 for(
unsigned int d = 0; d < dim; ++d)
56 space_point(d) = chart_point(d) + sinval;
62 pull_back(dealii::Point<dim>
const & space_point)
const override
64 dealii::Point<dim> x = space_point;
65 dealii::Point<dim> one;
66 for(
unsigned int d = 0; d < dim; ++d)
70 dealii::Tensor<1, dim> sinvals, residual;
72 for(
unsigned int d = 0; d < dim; ++d)
73 sinvals[d] = std::sin(frequency * dealii::numbers::PI * (x(d) - left) / (right - left));
74 double sinval = deformation;
75 for(
unsigned int d = 0; d < dim; ++d)
77 residual = space_point - (x + sinval * one);
80 while(residual.norm() > 1e-12 and its < 100)
82 dealii::Tensor<2, dim> jacobian;
83 for(
unsigned int d = 0; d < dim; ++d)
85 for(
unsigned int d = 0; d < dim; ++d)
88 deformation * frequency / (right - left) * dealii::numbers::PI *
89 std::cos(frequency * dealii::numbers::PI * (x(d) - left) / (right - left));
90 for(
unsigned int e = 0; e < dim; ++e)
92 sinval_der *= sinvals[e];
93 for(
unsigned int e = 0; e < dim; ++e)
94 jacobian[e][d] += sinval_der;
97 x += invert(jacobian) * residual;
99 for(
unsigned int d = 0; d < dim; ++d)
100 sinvals[d] = std::sin(frequency * dealii::numbers::PI * (x(d) - left) / (right - left));
101 sinval = deformation;
102 for(
unsigned int d = 0; d < dim; ++d)
103 sinval *= sinvals[d];
104 residual = space_point - (x + sinval * one);
109 AssertThrow(residual.norm() < 1e-12, dealii::ExcMessage(
"Newton for point did not converge."));
114 std::unique_ptr<dealii::Manifold<dim>>
115 clone()
const override
117 return std::make_unique<DeformedCubeManifold<dim>>(left, right, deformation, frequency);
123 double const deformation;
124 unsigned int const frequency;