diff --git a/tests/numerics/petsc_vector_test.C b/tests/numerics/petsc_vector_test.C index 72828b4476a..f9e6b369a64 100644 --- a/tests/numerics/petsc_vector_test.C +++ b/tests/numerics/petsc_vector_test.C @@ -3,11 +3,18 @@ #ifdef LIBMESH_HAVE_PETSC #include "numeric_vector_test.h" +#include +#include +#include +#include +#include +#include using namespace libMesh; class PetscVectorTest : public NumericVectorTest> { + public: PetscVectorTest() : NumericVectorTest>() { @@ -25,6 +32,8 @@ public: CPPUNIT_TEST( testPetscOperations ); + CPPUNIT_TEST( testPetscVectorGetSetThread ); + CPPUNIT_TEST_SUITE_END(); void testGetArray() @@ -129,6 +138,89 @@ public: libMesh::TOLERANCE*libMesh::TOLERANCE); } + void testPetscVectorGetSetThread() + { + LOG_UNIT_TEST; + + Mesh mesh(*TestCommWorld); + EquationSystems es(mesh); + + LinearImplicitSystem & sys = es.add_system ("test"); + + sys.add_variable("u", FEType(0, MONOMIAL)); + + MeshTools::Generation::build_line (mesh, + 50, + 0., 1., + EDGE2); + es.init(); + + *sys.solution = 1.0; + + ConstElemRange active_local_elem_range(mesh.active_local_elements_begin(), + mesh.active_local_elements_end()); + std::vector cache(mesh.parallel_n_elem()); + for (const auto & elem : active_local_elem_range) + { + std::vector indices; + sys.get_dof_map().dof_indices(elem, indices, 0); + libmesh_assert(indices.size() == 1); + cache[elem->id()] = indices[0]; + } + + ConstMonomialGetSetThread test_thread(sys, cache); + + Threads::parallel_for(active_local_elem_range, test_thread); + + sys.solution->close(); + + const libMesh::dof_id_type first = sys.solution->first_local_index(); + const libMesh::dof_id_type last = sys.solution->last_local_index(); + + for (libMesh::dof_id_type n=first; n != last; n++) + LIBMESH_ASSERT_FP_EQUAL(libMesh::libmesh_real((*sys.solution)(n)), + libMesh::Real(0.5), + libMesh::TOLERANCE*libMesh::TOLERANCE); + + } + +private: + + class ConstMonomialGetSetThread + { + public: + ConstMonomialGetSetThread (LinearImplicitSystem & system, + const std::vector & cache) + : + _system(system), + _cache(cache) + {} + + ConstMonomialGetSetThread (ConstMonomialGetSetThread & x, + Threads::split /*split*/) + : + _system(x._system), + _cache(x._cache) + {} + + void operator()(const ConstElemRange & range) const + { + PetscVector * solution = + dynamic_cast *>(_system.solution.get()); + libmesh_assert(solution); + + for (const auto & elem : range) + { + const auto modified_value = (*solution)(_cache[elem->id()]) / 2.0; + solution->set(_cache[elem->id()], modified_value); + } + } + + private: + LinearImplicitSystem & _system; + const std::vector & _cache; + }; + }; CPPUNIT_TEST_SUITE_REGISTRATION( PetscVectorTest );