From 51e2b3b60ea91ed916dd41efade72ec8e00f7442 Mon Sep 17 00:00:00 2001 From: Peter Munch Date: Tue, 30 Jul 2024 15:33:42 +0200 Subject: [PATCH] Add coupling terms in prototype (#1212) Former-commit-id: a091baf2a75c0128e09edbdc07a174278b85f9bc --- ...e_mortar_poisson_consistent_integration.cc | 181 +++++++++++++++++- 1 file changed, 177 insertions(+), 4 deletions(-) diff --git a/prototypes/matrix_free_mortar_poisson_consistent_integration/matrix_free_mortar_poisson_consistent_integration.cc b/prototypes/matrix_free_mortar_poisson_consistent_integration/matrix_free_mortar_poisson_consistent_integration.cc index ef357e5fec..72419f461c 100644 --- a/prototypes/matrix_free_mortar_poisson_consistent_integration/matrix_free_mortar_poisson_consistent_integration.cc +++ b/prototypes/matrix_free_mortar_poisson_consistent_integration/matrix_free_mortar_poisson_consistent_integration.cc @@ -12,6 +12,7 @@ #include #include #include +#include #include #include @@ -328,9 +329,6 @@ class PoissonOperator : public Subscriptor FEPointEvaluation<1, dim, dim, double> phi_p( mapping, dof_handler.get_fe(), update_values | update_gradients); - Vector buffer_0; - Vector buffer_1; - for (const auto &[JxWs, cell_0, points_0, @@ -346,8 +344,11 @@ class PoissonOperator : public Subscriptor phi_p.reinit(cell_1, points_1); // gather + Vector buffer_0; buffer_0.reinit(cell_dof_0->get_fe().n_dofs_per_cell()); cell_dof_0->get_dof_values(src, buffer_0); + + Vector buffer_1; buffer_1.reinit(cell_dof_1->get_fe().n_dofs_per_cell()); cell_dof_1->get_dof_values(src, buffer_1); @@ -425,6 +426,36 @@ class PoissonOperator : public Subscriptor else DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints); + // for coupling term + for (const auto &[JxWs, + cell_0, + points_0, + cell_1, + points_1, + normals, + penalty_parameter] : all_intersections) + { + const auto cell_dof_0 = + cell_0->as_dof_handler_iterator(dof_handler); + const auto cell_dof_1 = + cell_1->as_dof_handler_iterator(dof_handler); + + std::vector dof_indices_0( + cell_dof_0->get_fe().n_dofs_per_cell()); + std::vector dof_indices_1( + cell_dof_1->get_fe().n_dofs_per_cell()); + + cell_dof_0->get_dof_indices(dof_indices_0); + cell_dof_1->get_dof_indices(dof_indices_1); + + constraints.add_entries_local_to_global(dof_indices_0, + dof_indices_1, + dsp); + constraints.add_entries_local_to_global(dof_indices_1, + dof_indices_0, + dsp); + } + dsp.compress(); system_matrix.reinit(dsp); @@ -442,6 +473,147 @@ class PoissonOperator : public Subscriptor do_vmult_cell_single, this); + // for coupling term + FEPointEvaluation<1, dim, dim, double> phi_m( + mapping, dof_handler.get_fe(), update_values | update_gradients); + FEPointEvaluation<1, dim, dim, double> phi_p( + mapping, dof_handler.get_fe(), update_values | update_gradients); + + for (const auto &[JxWs, + cell_0, + points_0, + cell_1, + points_1, + normals, + penalty_parameter] : all_intersections) + { + const auto cell_dof_0 = + cell_0->as_dof_handler_iterator(dof_handler); + const auto cell_dof_1 = + cell_1->as_dof_handler_iterator(dof_handler); + + phi_m.reinit(cell_0, points_0); + phi_p.reinit(cell_1, points_1); + + const unsigned int n_dofs_per_cell_0 = + cell_dof_0->get_fe().n_dofs_per_cell(); + const unsigned int n_dofs_per_cell_1 = + cell_dof_1->get_fe().n_dofs_per_cell(); + + Vector buffer_0; + buffer_0.reinit(n_dofs_per_cell_0); + + Vector buffer_1; + buffer_1.reinit(n_dofs_per_cell_1); + + std::vector dof_indices_0( + cell_dof_0->get_fe().n_dofs_per_cell()); + std::vector dof_indices_1( + cell_dof_0->get_fe().n_dofs_per_cell()); + + cell_dof_0->get_dof_indices(dof_indices_0); + cell_dof_1->get_dof_indices(dof_indices_1); + + // evaluate + phi_m.evaluate(buffer_0, EvaluationFlags::values); + phi_p.evaluate(buffer_1, EvaluationFlags::values); + + for (unsigned int b = 0; b < 2; ++b) + { + FullMatrix matrix_0; + FullMatrix matrix_1; + + if (b == 0) + { + matrix_0.reinit(n_dofs_per_cell_0, n_dofs_per_cell_0); + matrix_1.reinit(n_dofs_per_cell_1, n_dofs_per_cell_0); + } + else + { + matrix_0.reinit(n_dofs_per_cell_0, n_dofs_per_cell_1); + matrix_1.reinit(n_dofs_per_cell_1, n_dofs_per_cell_1); + } + + for (unsigned int j = 0; + j < (b == 0 ? n_dofs_per_cell_0 : n_dofs_per_cell_1); + ++j) + { + buffer_0 = 0.0; + buffer_1 = 0.0; + + if (b == 0) + buffer_0[j] = 1.0; + else + buffer_1[j] = 1.0; + + phi_m.evaluate(buffer_0, EvaluationFlags::values); + phi_p.evaluate(buffer_1, EvaluationFlags::values); + + // quadrature loop + for (const auto q : phi_m.quadrature_point_indices()) + { + const auto normal = normals[q]; + const auto JxW = JxWs[q]; + + const auto value_m = phi_m.get_value(q); + const auto value_p = phi_p.get_value(q); + + const auto gradient_m = phi_m.get_gradient(q); + const auto gradient_p = phi_p.get_gradient(q); + + const auto jump_value = (value_m - value_p) * 0.5 * JxW; + const auto avg_gradient = + normal * (gradient_m + gradient_p) * 0.5 * JxW; + + const double sigma = penalty_parameter * panalty_factor; + + phi_m.submit_gradient(-jump_value * normal, q); + phi_m.submit_value(+jump_value * sigma * 2.0 - + avg_gradient, + q); + + phi_p.submit_gradient(-jump_value * normal, q); + phi_p.submit_value(-jump_value * sigma * 2.0 + + avg_gradient, + q); + } + + // integrate + phi_m.test_and_sum(buffer_0, EvaluationFlags::values); + phi_p.test_and_sum(buffer_1, EvaluationFlags::values); + + for (unsigned int i = 0; i < n_dofs_per_cell_0; ++i) + matrix_0[i][j] = buffer_0[i]; + for (unsigned int i = 0; i < n_dofs_per_cell_1; ++i) + matrix_1[i][j] = buffer_1[i]; + } + + + if (b == 0) + { + constraints.distribute_local_to_global(matrix_0, + dof_indices_0, + system_matrix); + constraints.distribute_local_to_global(matrix_1, + dof_indices_1, + dof_indices_0, + system_matrix); + } + else + { + constraints.distribute_local_to_global(matrix_0, + dof_indices_0, + dof_indices_1, + system_matrix); + constraints.distribute_local_to_global(matrix_1, + dof_indices_1, + system_matrix); + } + } + } + + system_matrix.compress(VectorOperation::add); + this->valid_system = true; } } @@ -585,7 +757,8 @@ main(int argc, char **argv) ReductionControl reduction_control(10000, 1e-20, 1e-4); - TrilinosWrappers::PreconditionAMG preconditioner; + TrilinosWrappers::PreconditionAMG + preconditioner; // TrilinosWrappers::SolverDirect preconditioner.initialize(op.get_system_matrix()); SolverGMRES solver(reduction_control);