From e8fdad65fb1aac6522cc7b97fd5b796a35a9ff39 Mon Sep 17 00:00:00 2001 From: bangerth Date: Mon, 8 Jul 2013 02:23:31 +0000 Subject: [PATCH] A much reduced testcase. Should produce the same number every time, but doesn't. git-svn-id: https://svn.dealii.org/trunk@29945 0785d39b-7218-0410-832d-ea1e28bc413d --- tests/base/work_stream_03.cc | 271 +++++++++++++++++++++++++++++------ 1 file changed, 231 insertions(+), 40 deletions(-) diff --git a/tests/base/work_stream_03.cc b/tests/base/work_stream_03.cc index 954c4b4ef3..123f21f9f5 100644 --- a/tests/base/work_stream_03.cc +++ b/tests/base/work_stream_03.cc @@ -19,7 +19,8 @@ // the failing tests (deal.II/project_q_01) multiple times and // verifying that it indeed works -//#include "../tests.h" +#include "../tests.h" +#include #include #include #include @@ -43,9 +44,23 @@ #include #include #include -#include -#include -#include + +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include #include #include @@ -83,36 +98,218 @@ class F : public Function }; +Triangulation<2> triangulation; +FE_Q<2> fe(3); +DoFHandler<2> dof_handler(triangulation); +ConstraintMatrix constraints; + +namespace +{ + template + struct Scratch + { + Scratch (const ::dealii::hp::FECollection &fe, + const UpdateFlags update_flags, + const Function *coefficient, + const Function *rhs_function, + const ::dealii::hp::QCollection &quadrature, + const ::dealii::hp::MappingCollection &mapping) + : + fe_collection (fe), + quadrature_collection (quadrature), + mapping_collection (mapping), + x_fe_values (mapping_collection, + fe_collection, + quadrature_collection, + update_flags), + coefficient_values(quadrature_collection.max_n_quadrature_points()), + coefficient_vector_values (quadrature_collection.max_n_quadrature_points(), + dealii::Vector (fe_collection.n_components())), + rhs_values(quadrature_collection.max_n_quadrature_points()), + rhs_vector_values(quadrature_collection.max_n_quadrature_points(), + dealii::Vector (fe_collection.n_components())), + coefficient (coefficient), + rhs_function (rhs_function), + update_flags (update_flags) + {} + + Scratch (const Scratch &data) + : + fe_collection (data.fe_collection), + quadrature_collection (data.quadrature_collection), + mapping_collection (data.mapping_collection), + x_fe_values (mapping_collection, + fe_collection, + quadrature_collection, + data.update_flags), + coefficient_values (data.coefficient_values), + coefficient_vector_values (data.coefficient_vector_values), + rhs_values (data.rhs_values), + rhs_vector_values (data.rhs_vector_values), + coefficient (data.coefficient), + rhs_function (data.rhs_function), + update_flags (data.update_flags) + {} + + const ::dealii::hp::FECollection &fe_collection; + const ::dealii::hp::QCollection &quadrature_collection; + const ::dealii::hp::MappingCollection &mapping_collection; + + ::dealii::hp::FEValues x_fe_values; + + std::vector coefficient_values; + std::vector > coefficient_vector_values; + std::vector rhs_values; + std::vector > rhs_vector_values; + + std::vector old_JxW; + + const Function *coefficient; + const Function *rhs_function; + + const UpdateFlags update_flags; + }; + + + struct CopyData + { + std::vector dof_indices; + FullMatrix cell_matrix; + dealii::Vector cell_rhs; + const ConstraintMatrix *constraints; + }; +} + +template + void + mass_assembler(const CellIterator &cell, Scratch &data, + CopyData ©_data) + { + data.x_fe_values.reinit(cell); + const FEValues &fe_values = + data.x_fe_values.get_present_fe_values(); + + const unsigned int dofs_per_cell = fe_values.dofs_per_cell, n_q_points = + fe_values.n_quadrature_points; + const FiniteElement &fe = fe_values.get_fe(); + const unsigned int n_components = fe.n_components(); + + { + copy_data.cell_matrix.reinit(dofs_per_cell, dofs_per_cell); + copy_data.cell_rhs.reinit(dofs_per_cell); + + copy_data.dof_indices.resize(dofs_per_cell); + } + + Assert (copy_data.cell_matrix.frobenius_norm() == 0, ExcInternalError()); + cell->get_dof_indices(copy_data.dof_indices); + + data.rhs_values.resize(n_q_points); + data.rhs_function->value_list(fe_values.get_quadrature_points(), + data.rhs_values); + + for (unsigned int i = 0; i < dofs_per_cell; ++i) + { + const unsigned int component_i = fe.system_to_component_index(i).first; + const double *phi_i = &fe_values.shape_value(i, 0); + + double add_data = 0; + for (unsigned int point = 0; point < n_q_points; ++point) + add_data += phi_i[point] * fe_values.JxW(point) + * data.rhs_values[point]; + copy_data.cell_rhs(i) = add_data; + } + } + + +template +void copy_local_to_global (const CopyData &data, + MatrixType *matrix, + VectorType *right_hand_side) +{ + const unsigned int dofs_per_cell = data.dof_indices.size(); + Assert (data.cell_matrix.frobenius_norm() == 0, ExcInternalError()); + + Assert (matrix->frobenius_norm() == 0, ExcInternalError()); + data.constraints->distribute_local_to_global(data.cell_matrix, + data.cell_rhs, + data.dof_indices, + *matrix, *right_hand_side); + Assert (matrix->frobenius_norm() == 0, ExcInternalError()); +//Q: why does this write anything into the matrix??? +} + + + +template +void create_mass_matrix (const Mapping &mapping, + const DoFHandler &dof, + const Quadrature &q, + SparseMatrix &matrix, + const Function &rhs, + Vector &rhs_vector, + const Function *const coefficient, + const ConstraintMatrix &constraints) +{ + Assert (matrix.m() == dof.n_dofs(), + ExcDimensionMismatch (matrix.m(), dof.n_dofs())); + Assert (matrix.n() == dof.n_dofs(), + ExcDimensionMismatch (matrix.n(), dof.n_dofs())); + + hp::FECollection fe_collection (dof.get_fe()); + hp::QCollection q_collection (q); + hp::MappingCollection mapping_collection (mapping); + Scratch + assembler_data (fe_collection, + update_values | + update_JxW_values | update_quadrature_points, + coefficient, &rhs, + q_collection, mapping_collection); + CopyData copy_data; + copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(), + assembler_data.fe_collection.max_dofs_per_cell()); + copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell()); + copy_data.dof_indices.resize (assembler_data.fe_collection.max_dofs_per_cell()); + copy_data.constraints = &constraints; + + dealii::WorkStream::run (dof.begin_active(), + static_cast::active_cell_iterator>(dof.end()), + &mass_assembler::active_cell_iterator>, + std_cxx1x::bind(& + copy_local_to_global, Vector >, + std_cxx1x::_1, &matrix, &rhs_vector), + assembler_data, + copy_data, + 2*multithread_info.n_default_threads, + 1); +} + + template void do_project (const unsigned int p) { - Triangulation triangulation; - GridGenerator::hyper_cube (triangulation); - triangulation.refine_global (3); + SparsityPattern sparsity; + { + CompressedSimpleSparsityPattern csp (dof_handler.n_dofs(), dof_handler.n_dofs()); + DoFTools::make_sparsity_pattern (dof_handler, csp, constraints, + false); - std::cout << "Start: " << __PRETTY_FUNCTION__ << ' ' << p << std::endl; - FE_Q fe(p); - DoFHandler dof_handler(triangulation); - dof_handler.distribute_dofs (fe); + sparsity.copy_from (csp); + } - ConstraintMatrix constraints; - constraints.close (); + SparseMatrix mass_matrix (sparsity); + Vector tmp (mass_matrix.n()); - Vector projection (dof_handler.n_dofs()); - Vector error (triangulation.n_active_cells()); - for (unsigned int q=0; q<4; ++q) - { - // project the function - VectorTools::project (dof_handler, - constraints, - QGauss(p+2), - F (), - projection); - Assert (std::fabs(projection.l1_norm() - 3750.000000000079) < 1e-10, - ExcInternalError()); - } - std::cout << "Done: " << __PRETTY_FUNCTION__ << ' ' << p - << std::endl; + const Function<2>* dummy = 0; + create_mass_matrix(StaticMappingQ1 < dim > ::mapping, + dof_handler, QGauss < dim > (5), mass_matrix, F(), tmp, dummy, + constraints); + std::ostringstream x; + x.precision(18); + x << "Check1: " << mass_matrix.frobenius_norm() << " " << tmp.l1_norm() << std::endl; + std::cout << x.str(); } @@ -121,10 +318,8 @@ void do_project (const unsigned int p) template void test () { - std::cout.precision(16); - Threads::TaskGroup<> g; - for (unsigned int p=1; p<12; ++p) + for (unsigned int p=1; p<48*3; ++p) g += Threads::new_task (&do_project, 3); g.join_all (); } @@ -139,16 +334,12 @@ int main () deallog.depth_console(0); deallog.threshold_double(1.e-10); - test<2>(); + GridGenerator::hyper_cube (triangulation, 0, 1); + triangulation.refine_global (1); + dof_handler.distribute_dofs (fe); + constraints.close (); - Threads::TaskGroup<> g; - for (unsigned int i=0; i<2; ++i) - { - g += Threads::new_task (&test<1>); - g += Threads::new_task (&test<2>); - g += Threads::new_task (&test<3>); - } - g.join_all (); + test<2>(); deallog << "OK" << std::endl; } -- 2.39.5