From: Daniel Arndt Date: Thu, 27 Oct 2016 16:06:29 +0000 (+0200) Subject: Add test for the inverse diagonal X-Git-Tag: v8.5.0-rc1~539^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F3296%2Fhead;p=dealii.git Add test for the inverse diagonal --- diff --git a/include/deal.II/matrix_free/operators.h b/include/deal.II/matrix_free/operators.h index 24e6e510b6..2762304780 100644 --- a/include/deal.II/matrix_free/operators.h +++ b/include/deal.II/matrix_free/operators.h @@ -780,27 +780,19 @@ namespace MatrixFreeOperators LinearAlgebra::distributed::Vector ones; Base::initialize_dof_vector(Base::inverse_diagonal_entries); Base::initialize_dof_vector(ones); - Base::inverse_diagonal_entries = 1.; ones = 1.; ones.update_ghost_values(); apply_add(Base::inverse_diagonal_entries, ones); - const IndexSet &locally_owned_elements - = Base::inverse_diagonal_entries.locally_owned_elements(); - locally_owned_elements.print(std::cout); - const std::vector constrained_dofs + const std::vector &constrained_dofs = Base::data->get_constrained_dofs(); - for (unsigned int i=0; i::inverse_diagonal_entries(constrained_dofs[i]) = 1.; + for (unsigned int i=0; i< constrained_dofs.size(); ++i) + Base::inverse_diagonal_entries.local_element(constrained_dofs[i]) = 1.; const unsigned int local_size = Base::inverse_diagonal_entries.local_size(); for (unsigned int i=0; i::inverse_diagonal_entries.local_element(i); - Base::inverse_diagonal_entries.local_element(i) - =1./Base::inverse_diagonal_entries.local_element(i); - } + Base::inverse_diagonal_entries.local_element(i) + =1./Base::inverse_diagonal_entries.local_element(i); Base::inverse_diagonal_entries.compress(VectorOperation::insert); Base::inverse_diagonal_entries.update_ghost_values(); diff --git a/tests/matrix_free/mass_operator_02.cc b/tests/matrix_free/mass_operator_02.cc new file mode 100644 index 0000000000..1b83bb0c1b --- /dev/null +++ b/tests/matrix_free/mass_operator_02.cc @@ -0,0 +1,194 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2016 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// this tests the correctness of matrix_diagonal_inverse of +// MatrixFreeOperators::MassOperator. + +#include "../tests.h" + +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +#include + + + + +template +void test () +{ + typedef double number; + + parallel::distributed::Triangulation tria (MPI_COMM_WORLD); + GridGenerator::hyper_cube (tria); + tria.refine_global(3); + + FE_Q fe (fe_degree); + DoFHandler dof (tria); + dof.distribute_dofs(fe); + + IndexSet owned_set = dof.locally_owned_dofs(); + IndexSet relevant_set; + DoFTools::extract_locally_relevant_dofs (dof, relevant_set); + + ConstraintMatrix constraints (relevant_set); + DoFTools::make_hanging_node_constraints(dof, constraints); + VectorTools::interpolate_boundary_values (dof, 0, ZeroFunction(), + constraints); + constraints.close(); + + deallog << "Testing " << dof.get_fe().get_name() << std::endl; + //std::cout << "Number of cells: " << tria.n_global_active_cells() << std::endl; + //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl; + //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl; + + MatrixFree mf_data; + { + const QGauss<1> quad (fe_degree+1); + typename MatrixFree::AdditionalData data; + data.mpi_communicator = MPI_COMM_WORLD; + data.tasks_parallel_scheme = + MatrixFree::AdditionalData::none; + data.tasks_block_size = 7; + mf_data.reinit (dof, constraints, quad, data); + } + + MatrixFreeOperators::MassOperator mf; + mf.initialize(mf_data); + mf.compute_diagonal(); + const LinearAlgebra::distributed::Vector diagonal + = mf.get_matrix_diagonal_inverse(); + + LinearAlgebra::distributed::Vector in, out, ref; + mf_data.initialize_dof_vector (in); + out.reinit (in); + ref.reinit (in); + + for (unsigned int i=0; i quadrature_formula(fe_degree+1); + + FEValues fe_values (dof.get_fe(), quadrature_formula, + update_values | update_gradients | + update_JxW_values); + + const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; + const unsigned int n_q_points = quadrature_formula.size(); + + FullMatrix cell_matrix (dofs_per_cell, dofs_per_cell); + std::vector local_dof_indices (dofs_per_cell); + + typename DoFHandler::active_cell_iterator + cell = dof.begin_active(), + endc = dof.end(); + for (; cell!=endc; ++cell) + if (cell->is_locally_owned()) + { + cell_matrix = 0; + fe_values.reinit (cell); + + for (unsigned int q_point=0; q_pointget_dof_indices(local_dof_indices); + constraints.distribute_local_to_global (cell_matrix, + local_dof_indices, + sparse_matrix); + } + } + sparse_matrix.compress(VectorOperation::add); + + sparse_matrix.vmult (ref, in); + + for (unsigned int i=0; i(); + test<2,2>(); + deallog.pop(); + + deallog.push("3d"); + test<3,1>(); + test<3,2>(); + deallog.pop(); +} + diff --git a/tests/matrix_free/mass_operator_02.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output b/tests/matrix_free/mass_operator_02.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..0b7375c116 --- /dev/null +++ b/tests/matrix_free/mass_operator_02.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output @@ -0,0 +1,25 @@ + +DEAL:2d::Testing FE_Q<2>(1) +DEAL:2d::Norm of difference: 4.26326e-14 +DEAL:2d::l2_norm: 504.352 +DEAL:2d::l1_norm: 3536.64 +DEAL:2d::linfty_norm: 92.1600 +DEAL:2d:: +DEAL:2d::Testing FE_Q<2>(2) +DEAL:2d::Norm of difference: 4.54747e-13 +DEAL:2d::l2_norm: 5051.63 +DEAL:2d::l1_norm: 68866.9 +DEAL:2d::linfty_norm: 576.000 +DEAL:2d:: +DEAL:3d::Testing FE_Q<3>(1) +DEAL:3d::Norm of difference: 5.68434e-13 +DEAL:3d::l2_norm: 11325.6 +DEAL:3d::l1_norm: 207861. +DEAL:3d::linfty_norm: 884.736 +DEAL:3d:: +DEAL:3d::Testing FE_Q<3>(2) +DEAL:3d::Norm of difference: 2.18279e-11 +DEAL:3d::l2_norm: 359043. +DEAL:3d::l1_norm: 1.80487e+07 +DEAL:3d::linfty_norm: 13824.0 +DEAL:3d:: diff --git a/tests/matrix_free/mass_operator_02.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=4.output b/tests/matrix_free/mass_operator_02.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=4.output new file mode 100644 index 0000000000..0b7375c116 --- /dev/null +++ b/tests/matrix_free/mass_operator_02.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=4.output @@ -0,0 +1,25 @@ + +DEAL:2d::Testing FE_Q<2>(1) +DEAL:2d::Norm of difference: 4.26326e-14 +DEAL:2d::l2_norm: 504.352 +DEAL:2d::l1_norm: 3536.64 +DEAL:2d::linfty_norm: 92.1600 +DEAL:2d:: +DEAL:2d::Testing FE_Q<2>(2) +DEAL:2d::Norm of difference: 4.54747e-13 +DEAL:2d::l2_norm: 5051.63 +DEAL:2d::l1_norm: 68866.9 +DEAL:2d::linfty_norm: 576.000 +DEAL:2d:: +DEAL:3d::Testing FE_Q<3>(1) +DEAL:3d::Norm of difference: 5.68434e-13 +DEAL:3d::l2_norm: 11325.6 +DEAL:3d::l1_norm: 207861. +DEAL:3d::linfty_norm: 884.736 +DEAL:3d:: +DEAL:3d::Testing FE_Q<3>(2) +DEAL:3d::Norm of difference: 2.18279e-11 +DEAL:3d::l2_norm: 359043. +DEAL:3d::l1_norm: 1.80487e+07 +DEAL:3d::linfty_norm: 13824.0 +DEAL:3d::