43struct MeshMovementData
46 : temporal(MeshMovementAdvanceInTime::Undefined),
47 shape(MeshMovementShape::Undefined),
52 spatial_number_of_oscillations(1)
56 MeshMovementAdvanceInTime temporal;
57 MeshMovementShape shape;
58 dealii::Tensor<1, dim> dimensions;
63 double spatial_number_of_oscillations;
67class CubeMeshMovementFunctions :
public dealii::Function<dim>
71 : dealii::Function<dim>(dim),
73 width(data_in.dimensions[0]),
74 left(-1.0 / 2.0 * width),
76 runtime(data_in.t_end - data_in.t_start),
77 time_period(data_in.period)
82 value(dealii::Point<dim>
const & x,
unsigned int const coordinate_direction = 0)
const override
84 double displacement = 0.0;
86 displacement = compute_displacement_share(x, coordinate_direction) * compute_time_share();
93 compute_displacement_share(dealii::Point<dim>
const & x,
94 unsigned int const coordinate_direction = 0)
const
96 double solution = 0.0;
100 case MeshMovementShape::Undefined:
101 AssertThrow(
false, dealii::ExcMessage(
"Undefined parameter MeshMovementShape."));
104 case MeshMovementShape::Sin:
105 if(coordinate_direction == 0)
107 std::sin(2.0 * pi * (x(1) - left) * data.spatial_number_of_oscillations / width) *
110 std::sin(2.0 * pi * (x(2) - left) * data.spatial_number_of_oscillations / width) :
112 else if(coordinate_direction == 1)
114 std::sin(2.0 * pi * (x(0) - left) * data.spatial_number_of_oscillations / width) *
117 std::sin(2.0 * pi * (x(2) - left) * data.spatial_number_of_oscillations / width) :
119 else if(coordinate_direction == 2)
121 std::sin(2.0 * pi * (x(0) - left) * data.spatial_number_of_oscillations / width) *
123 std::sin(2.0 * pi * (x(1) - left) * data.spatial_number_of_oscillations / width);
126 case MeshMovementShape::SineZeroAtBoundary:
127 if(coordinate_direction == 0)
129 std::sin(pi * (x(0) - left) * data.spatial_number_of_oscillations / width) *
130 std::sin(2.0 * pi * (x(1) - left) * data.spatial_number_of_oscillations / width) *
133 std::sin(2.0 * pi * (x(2) - left) * data.spatial_number_of_oscillations / width) :
135 else if(coordinate_direction == 1)
137 std::sin(pi * (x(1) - left) * data.spatial_number_of_oscillations / width) *
138 std::sin(2.0 * pi * (x(0) - left) * data.spatial_number_of_oscillations / width) *
141 std::sin(2.0 * pi * (x(2) - left) * data.spatial_number_of_oscillations / width) :
143 else if(coordinate_direction == 2)
145 std::sin(pi * (x(2) - left) * data.spatial_number_of_oscillations / width) *
146 std::sin(2.0 * pi * (x(0) - left) * data.spatial_number_of_oscillations / width) *
148 std::sin(2.0 * pi * (x(1) - left) * data.spatial_number_of_oscillations / width);
151 case MeshMovementShape::SineAligned:
152 solution = std::sin(2.0 * pi * (x(coordinate_direction) - left) *
153 data.spatial_number_of_oscillations / width) *
158 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
166 compute_time_share()
const
168 double solution = 0.0;
170 switch(data.temporal)
172 case MeshMovementAdvanceInTime::Undefined:
173 AssertThrow(
false, dealii::ExcMessage(
"Undefined parameter MeshMovementAdvanceInTime."));
176 case MeshMovementAdvanceInTime::SinSquared:
177 solution = std::pow(std::sin(2.0 * pi * this->get_time() / time_period), 2);
180 case MeshMovementAdvanceInTime::Sin:
181 solution = std::sin(2.0 * pi * this->get_time() / time_period);
185 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
192 double const pi = dealii::numbers::PI;
197 double const runtime;
198 double const time_period;
202class RectangleMeshMovementFunctions :
public dealii::Function<dim>
206 : dealii::Function<dim>(dim),
208 length(data_in.dimensions[0]),
209 height(data_in.dimensions[1]),
210 depth(dim == 3 ? data_in.dimensions[2] : 1.0),
211 runtime(data_in.t_end - data_in.t_start),
212 time_period(data_in.period)
217 value(dealii::Point<dim>
const & x_in,
unsigned int const coordinate_direction = 0)
const override
220 dealii::Point<dim> x = x_in;
222 x[0] -= length / 2.0;
224 double displacement = 0.0;
226 displacement = compute_displacement_share(x, coordinate_direction) * compute_time_share();
233 compute_displacement_share(dealii::Point<dim>
const & x,
234 unsigned int const coordinate_direction = 0)
const
236 double solution = 0.0;
240 case MeshMovementShape::Undefined:
241 AssertThrow(
false, dealii::ExcMessage(
"Undefined parameter MeshMovementShape."));
244 case MeshMovementShape::Sin:
245 if(coordinate_direction == 0)
246 solution = std::sin(2.0 * pi * (x(1) - (height / 2.0)) / height) * data.amplitude *
247 (dim == 3 ? std::sin(2 * pi * (x(2.0) - (depth / 2)) / depth) : 1.0);
248 else if(coordinate_direction == 1)
249 solution = std::sin(2.0 * pi * (x(0) - (length / 2.0)) / length) * data.amplitude *
250 (dim == 3 ? std::sin(2 * pi * (x(2) - (depth / 2.0)) / depth) : 1.0);
251 else if(coordinate_direction == 2)
252 solution = std::sin(2.0 * pi * (x(1) - (height / 2.0)) / height) * data.amplitude *
253 std::sin(2.0 * pi * (x(0) - (length / 2.0)) / length);
257 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
265 compute_time_share()
const
267 double solution = 0.0;
269 switch(data.temporal)
271 case MeshMovementAdvanceInTime::Undefined:
272 AssertThrow(
false, dealii::ExcMessage(
"Undefined parameter MeshMovementAdvanceInTime."));
275 case MeshMovementAdvanceInTime::SinSquared:
276 solution = std::pow(std::sin(2.0 * pi * this->get_time() / time_period), 2);
279 case MeshMovementAdvanceInTime::Sin:
280 solution = std::sin(2.0 * pi * this->get_time() / time_period);
284 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
291 double const pi = dealii::numbers::PI;
296 double const runtime;
297 double const time_period;