176 using This = Operator<dim, Number>;
178 using BlockVectorType = dealii::LinearAlgebra::distributed::BlockVector<Number>;
180 using scalar = dealii::VectorizedArray<Number>;
181 using vector = dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>;
183 using Range = std::pair<unsigned int, unsigned int>;
185 using CellIntegratorU = CellIntegrator<dim, dim, Number>;
186 using CellIntegratorP = CellIntegrator<dim, 1, Number>;
188 using FaceIntegratorU = FaceIntegrator<dim, dim, Number>;
189 using FaceIntegratorP = FaceIntegrator<dim, 1, Number>;
192 Operator() : evaluation_time(Number{0.0}), tau(Number{0.0}), gamma(Number{0.0})
197 initialize(dealii::MatrixFree<dim, Number>
const & matrix_free_in,
201 matrix_free = &matrix_free_in;
205 tau = (Number)(0.5 * data_in.speed_of_sound);
206 gamma = (Number)(0.5 / data_in.speed_of_sound);
209 if(data_in.formulation == Formulation::Weak)
211 integrator_flags_p.cell_evaluate = dealii::EvaluationFlags::values;
212 integrator_flags_p.cell_integrate = dealii::EvaluationFlags::gradients;
213 integrator_flags_p.face_evaluate = dealii::EvaluationFlags::values;
214 integrator_flags_p.face_integrate = dealii::EvaluationFlags::values;
216 integrator_flags_u.cell_evaluate = dealii::EvaluationFlags::values;
217 integrator_flags_u.cell_integrate = dealii::EvaluationFlags::gradients;
218 integrator_flags_u.face_evaluate = dealii::EvaluationFlags::values;
219 integrator_flags_u.face_integrate = dealii::EvaluationFlags::values;
221 else if(data_in.formulation == Formulation::Strong)
223 integrator_flags_p.cell_evaluate = dealii::EvaluationFlags::gradients;
224 integrator_flags_p.cell_integrate = dealii::EvaluationFlags::values;
225 integrator_flags_p.face_evaluate = dealii::EvaluationFlags::values;
226 integrator_flags_p.face_integrate = dealii::EvaluationFlags::values;
228 integrator_flags_u.cell_evaluate = dealii::EvaluationFlags::gradients;
229 integrator_flags_u.cell_integrate = dealii::EvaluationFlags::values;
230 integrator_flags_u.face_evaluate = dealii::EvaluationFlags::values;
231 integrator_flags_u.face_integrate = dealii::EvaluationFlags::values;
233 else if(data_in.formulation == Formulation::SkewSymmetric)
235 integrator_flags_p.cell_evaluate = dealii::EvaluationFlags::gradients;
236 integrator_flags_p.cell_integrate = dealii::EvaluationFlags::gradients;
237 integrator_flags_p.face_evaluate = dealii::EvaluationFlags::values;
238 integrator_flags_p.face_integrate = dealii::EvaluationFlags::values;
240 integrator_flags_u.cell_evaluate = dealii::EvaluationFlags::values;
241 integrator_flags_u.cell_integrate = dealii::EvaluationFlags::values;
242 integrator_flags_u.face_evaluate = dealii::EvaluationFlags::values;
243 integrator_flags_u.face_integrate = dealii::EvaluationFlags::values;
247 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
252 evaluate(BlockVectorType & dst, BlockVectorType
const & src,
double const time)
const
254 do_evaluate(dst, src, time,
true);
258 evaluate_add(BlockVectorType & dst, BlockVectorType
const & src,
double const time)
const
260 do_evaluate(dst, src, time,
true);
265 do_evaluate(BlockVectorType & dst,
266 BlockVectorType
const & src,
268 bool const zero_dst_vector)
const
270 evaluation_time = (Number)time;
272 matrix_free->loop(&This::cell_loop,
274 &This::boundary_face_loop,
279 dealii::MatrixFree<dim, Number>::DataAccessOnFaces::values,
280 dealii::MatrixFree<dim, Number>::DataAccessOnFaces::values);
284 do_cell_integral(CellIntegratorP & pressure, CellIntegratorU & velocity)
const
286 if(data.formulation == Formulation::Weak)
288 for(
unsigned int q : pressure.quadrature_point_indices())
290 vector
const flux_momentum = kernel.get_volume_flux_weak_momentum(velocity, q);
291 scalar
const flux_mass = kernel.get_volume_flux_weak_mass(pressure, q);
293 pressure.submit_gradient(flux_momentum, q);
294 velocity.submit_divergence(flux_mass, q);
297 else if(data.formulation == Formulation::Strong)
299 for(
unsigned int q : pressure.quadrature_point_indices())
301 scalar
const flux_momentum = kernel.get_volume_flux_strong_momentum(velocity, q);
302 vector
const flux_mass = kernel.get_volume_flux_strong_mass(pressure, q);
304 pressure.submit_value(flux_momentum, q);
305 velocity.submit_value(flux_mass, q);
308 else if(data.formulation == Formulation::SkewSymmetric)
310 for(
unsigned int q : pressure.quadrature_point_indices())
312 vector
const flux_momentum = kernel.get_volume_flux_weak_momentum(velocity, q);
313 vector
const flux_mass = kernel.get_volume_flux_strong_mass(pressure, q);
315 pressure.submit_gradient(flux_momentum, q);
316 velocity.submit_value(flux_mass, q);
321 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
326 template<
bool weight_neighbor,
327 typename ExteriorFaceIntegratorP,
328 typename ExteriorFaceIntegratorU>
330 do_face_integral(FaceIntegratorP & pressure_m,
331 ExteriorFaceIntegratorP & pressure_p,
332 FaceIntegratorU & velocity_m,
333 ExteriorFaceIntegratorU & velocity_p)
const
335 for(
unsigned int q : pressure_m.quadrature_point_indices())
337 scalar
const pm = pressure_m.get_value(q);
338 scalar
const pp = pressure_p.get_value(q);
339 vector
const rho_um = velocity_m.get_value(q);
340 vector
const rho_up = velocity_p.get_value(q);
341 vector
const n = pressure_m.normal_vector(q);
343 vector
const flux_momentum =
344 kernel.calculate_lax_friedrichs_flux_momentum(rho_um, rho_up, gamma, pm, pp, n);
346 scalar
const flux_mass =
347 kernel.calculate_lax_friedrichs_flux_mass(pm, pp, tau, rho_um, rho_up, n);
349 if(data.formulation == Formulation::Weak)
351 scalar
const flux_momentum_weak = flux_momentum * n;
352 vector
const flux_mass_weak = flux_mass * n;
354 pressure_m.submit_value(flux_momentum_weak, q);
355 velocity_m.submit_value(flux_mass_weak, q);
357 if constexpr(weight_neighbor)
360 pressure_p.submit_value(-flux_momentum_weak, q);
361 velocity_p.submit_value(-flux_mass_weak, q);
364 else if(data.formulation == Formulation::Strong)
366 pressure_m.submit_value((flux_momentum - rho_um) * n, q);
367 velocity_m.submit_value((flux_mass - pm) * n, q);
369 if constexpr(weight_neighbor)
372 pressure_p.submit_value((rho_up - flux_momentum) * n, q);
373 velocity_p.submit_value((pp - flux_mass) * n, q);
376 else if(data.formulation == Formulation::SkewSymmetric)
378 scalar
const flux_momentum_weak = flux_momentum * n;
380 pressure_m.submit_value(flux_momentum_weak, q);
381 velocity_m.submit_value((flux_mass - pm) * n, q);
383 if constexpr(weight_neighbor)
386 pressure_p.submit_value(-flux_momentum_weak, q);
387 velocity_p.submit_value((pp - flux_mass) * n, q);
392 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
398 cell_loop(dealii::MatrixFree<dim, Number>
const & matrix_free_in,
399 BlockVectorType & dst,
400 BlockVectorType
const & src,
401 Range
const & cell_range)
const
403 CellIntegratorP pressure(matrix_free_in, data.dof_index_pressure, data.quad_index);
404 CellIntegratorU velocity(matrix_free_in, data.dof_index_velocity, data.quad_index);
406 for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
408 pressure.reinit(cell);
409 pressure.gather_evaluate(src.block(data.block_index_pressure),
410 integrator_flags_p.cell_evaluate);
412 velocity.reinit(cell);
413 velocity.gather_evaluate(src.block(data.block_index_velocity),
414 integrator_flags_u.cell_evaluate);
416 do_cell_integral(pressure, velocity);
418 pressure.integrate_scatter(integrator_flags_p.cell_integrate,
419 dst.block(data.block_index_pressure));
420 velocity.integrate_scatter(integrator_flags_u.cell_integrate,
421 dst.block(data.block_index_velocity));
426 face_loop(dealii::MatrixFree<dim, Number>
const & matrix_free_in,
427 BlockVectorType & dst,
428 BlockVectorType
const & src,
429 Range
const & face_range)
const
431 FaceIntegratorP pressure_m(matrix_free_in,
true, data.dof_index_pressure, data.quad_index);
432 FaceIntegratorP pressure_p(matrix_free_in,
false, data.dof_index_pressure, data.quad_index);
434 FaceIntegratorU velocity_m(matrix_free_in,
true, data.dof_index_velocity, data.quad_index);
435 FaceIntegratorU velocity_p(matrix_free_in,
false, data.dof_index_velocity, data.quad_index);
437 for(
unsigned int face = face_range.first; face < face_range.second; face++)
439 pressure_m.reinit(face);
440 pressure_m.gather_evaluate(src.block(data.block_index_pressure),
441 integrator_flags_p.face_evaluate);
442 pressure_p.reinit(face);
443 pressure_p.gather_evaluate(src.block(data.block_index_pressure),
444 integrator_flags_p.face_evaluate);
446 velocity_m.reinit(face);
447 velocity_m.gather_evaluate(src.block(data.block_index_velocity),
448 integrator_flags_u.face_evaluate);
449 velocity_p.reinit(face);
450 velocity_p.gather_evaluate(src.block(data.block_index_velocity),
451 integrator_flags_u.face_evaluate);
453 do_face_integral<true>(pressure_m, pressure_p, velocity_m, velocity_p);
455 pressure_m.integrate_scatter(integrator_flags_p.face_integrate,
456 dst.block(data.block_index_pressure));
457 pressure_p.integrate_scatter(integrator_flags_p.face_integrate,
458 dst.block(data.block_index_pressure));
460 velocity_m.integrate_scatter(integrator_flags_u.face_integrate,
461 dst.block(data.block_index_velocity));
462 velocity_p.integrate_scatter(integrator_flags_u.face_integrate,
463 dst.block(data.block_index_velocity));
468 boundary_face_loop(dealii::MatrixFree<dim, Number>
const & matrix_free_in,
469 BlockVectorType & dst,
470 BlockVectorType
const & src,
471 Range
const & face_range)
const
473 FaceIntegratorP pressure_m(matrix_free_in,
true, data.dof_index_pressure, data.quad_index);
476 FaceIntegratorU velocity_m(matrix_free_in,
true, data.dof_index_velocity, data.quad_index);
482 for(
unsigned int face = face_range.first; face < face_range.second; face++)
484 pressure_m.reinit(face);
485 pressure_m.gather_evaluate(src.block(data.block_index_pressure),
486 integrator_flags_p.face_evaluate);
487 pressure_p.reinit(face, evaluation_time);
489 velocity_m.reinit(face);
490 velocity_m.gather_evaluate(src.block(data.block_index_velocity),
491 integrator_flags_u.face_evaluate);
492 velocity_p.reinit(face, evaluation_time);
494 do_face_integral<false>(pressure_m, pressure_p, velocity_m, velocity_p);
496 pressure_m.integrate_scatter(integrator_flags_p.face_integrate,
497 dst.block(data.block_index_pressure));
498 velocity_m.integrate_scatter(integrator_flags_u.face_integrate,
499 dst.block(data.block_index_velocity));
505 mutable Number evaluation_time;
507 dealii::SmartPointer<dealii::MatrixFree<dim, Number>
const> matrix_free;