1 #include <libmesh/equation_systems.h>
2 #include <libmesh/replicated_mesh.h>
3 #include <libmesh/numeric_vector.h>
4 #include <libmesh/dof_map.h>
5 #include <libmesh/elem.h>
6 #include <libmesh/face_quad4.h>
7 #include <libmesh/ghosting_functor.h>
8 #include <libmesh/coupling_matrix.h>
9 #include <libmesh/quadrature_gauss.h>
10 #include <libmesh/linear_implicit_system.h>
11 #include <libmesh/mesh_refinement.h>
12 #include <libmesh/mesh_modification.h>
13 #include <libmesh/partitioner.h>
14 #include <libmesh/sparse_matrix.h>
15 #include <libmesh/auto_ptr.h>
44 _mesh(system.get_mesh()),
45 _point_locator(system.get_mesh().sub_point_locator()),
46 _coupling_matrix(nullptr)
50 CPPUNIT_ASSERT_EQUAL(2, (
int)_mesh.n_subdomains());
53 std::set<subdomain_id_type> ids;
54 _mesh.subdomain_ids(ids);
55 CPPUNIT_ASSERT( ids.find(1) != ids.end() );
56 CPPUNIT_ASSERT( ids.find(2) != ids.end() );
58 _point_locator->enable_out_of_mesh_mode();
64 { _coupling_matrix = std::move(coupling_matrix); }
66 virtual void operator()
70 std::unordered_map<const Elem *,const CouplingMatrix*> & coupled_elements)
override
72 std::unique_ptr<PointLocatorBase> sub_point_locator = _mesh.sub_point_locator();
74 for(
const auto & elem :
as_range(range_begin,range_end) )
76 std::set<subdomain_id_type> allowed_subdomains;
78 if( elem->subdomain_id() == 1 )
80 var = _system.variable_number(
"V");
81 allowed_subdomains.insert(2);
83 else if( elem->subdomain_id() == 2 )
85 var = _system.variable_number(
"U");
86 allowed_subdomains.insert(1);
89 CPPUNIT_FAIL(
"Something bad happpend");
92 FEType fe_type = _system.variable_type(var);
95 fe->attach_quadrature_rule (&qrule);
97 const std::vector<libMesh::Point> & qpoints = fe->get_xyz();
101 for (
const auto & qp : qpoints )
103 const Elem * overlapping_elem = (*sub_point_locator)( qp, &allowed_subdomains );
105 if( overlapping_elem && overlapping_elem->
processor_id() != p )
106 coupled_elements.insert( std::make_pair(overlapping_elem,_coupling_matrix.get()) );
142 virtual std::unique_ptr<Partitioner>
clone ()
const override
144 return libmesh_make_unique<OverlappingTestPartitioner>(*
this);
153 const unsigned int n)
override
160 libmesh_assert_greater (n, 0);
169 elem->processor_id() = e;
181 unsigned int n_parts_sub_two = 0;
182 if( n_sub_two_elems < n/2 )
183 n_parts_sub_two = n_sub_two_elems;
185 n_parts_sub_two = n/2;
187 const unsigned int n_parts_sub_one = n - n_parts_sub_two;
189 const dof_id_type sub_two_blk_size = cast_int<dof_id_type>
193 this->assign_proc_id_subdomain(
mesh, 2, sub_two_blk_size, n_parts_sub_two, 0 );
196 const dof_id_type sub_one_blk_size = cast_int<dof_id_type>
200 this->assign_proc_id_subdomain(
mesh, 1, sub_one_blk_size, n_parts_sub_one, n_parts_sub_two );
208 const unsigned int n_parts,
214 if ((e/blksize) < n_parts)
215 elem->processor_id() = offset + cast_int<processor_id_type>(e/blksize);
217 elem->processor_id() = offset;
236 std::unique_ptr<EquationSystems>
_es;
245 _mesh->set_mesh_dimension(2);
247 _mesh->add_point(
Point(0.0,0.0),0 );
248 _mesh->add_point(
Point(1.0,0.0),1 );
249 _mesh->add_point(
Point(1.0,1.0),2 );
250 _mesh->add_point(
Point(0.0,1.0),3 );
253 Elem* elem = _mesh->add_elem(
new Quad4 );
257 for (
unsigned int n=0; n<4; n++)
258 elem->
set_node(n) = _mesh->node_ptr(n);
261 _mesh->add_point(
Point(1.0,2.0),4 );
262 _mesh->add_point(
Point(0.0,2.0),5 );
265 Elem* elem = _mesh->add_elem(
new Quad4 );
269 elem->
set_node(0) = _mesh->node_ptr(3);
270 elem->
set_node(1) = _mesh->node_ptr(2);
271 elem->
set_node(2) = _mesh->node_ptr(4);
272 elem->
set_node(3) = _mesh->node_ptr(5);
275 _mesh->add_point(
Point(0.0,0.0),6 );
276 _mesh->add_point(
Point(1.0,0.0),7 );
277 _mesh->add_point(
Point(1.0,2.0),8 );
278 _mesh->add_point(
Point(0.0,2.0),9 );
281 Elem* elem = _mesh->add_elem(
new Quad4 );
285 elem->
set_node(0) = _mesh->node_ptr(6);
286 elem->
set_node(1) = _mesh->node_ptr(7);
287 elem->
set_node(2) = _mesh->node_ptr(8);
288 elem->
set_node(3) = _mesh->node_ptr(9);
291 _mesh->partitioner() = libmesh_make_unique<OverlappingTestPartitioner>();
293 _mesh->prepare_for_use();
295 #ifdef LIBMESH_ENABLE_AMR
296 if (n_refinements > 0)
301 #endif // LIBMESH_ENABLE_AMR
306 _es = libmesh_make_unique<EquationSystems>(*_mesh);
309 std::set<subdomain_id_type> sub_one;
312 std::set<subdomain_id_type> sub_two;
332 System & system = _es->get_system(
"SimpleSystem");
334 coupling = libmesh_make_unique<CouplingMatrix>(system.
n_vars());
343 (*coupling)(u_var,v_var) =
true;
344 (*coupling)(l_var,v_var) =
true;
345 (*coupling)(l_var,p_var) =
true;
346 (*coupling)(v_var,u_var) =
true;
347 (*coupling)(v_var,l_var) =
true;
362 #ifdef LIBMESH_HAVE_SOLVER
363 CPPUNIT_TEST( checkCouplingFunctorQuad );
364 CPPUNIT_TEST( checkCouplingFunctorQuadUnifRef );
365 CPPUNIT_TEST( checkCouplingFunctorTri );
366 CPPUNIT_TEST( checkCouplingFunctorTriUnifRef );
367 CPPUNIT_TEST( checkOverlappingPartitioner );
368 CPPUNIT_TEST( checkOverlappingPartitionerUnifRef );
371 CPPUNIT_TEST_SUITE_END();
378 this->build_quad_mesh();
386 { this->run_coupling_functor_test(0); }
389 { this->run_coupling_functor_test(1); }
395 this->run_coupling_functor_test(0);
402 this->run_coupling_functor_test(1);
407 this->run_partitioner_test(0);
412 this->run_partitioner_test(1);
423 #if defined(LIBMESH_ENABLE_AMR) && defined(LIBMESH_HAVE_SOLVER)
424 if( n_refinements > 0 )
432 System & system = _es->get_system(
"SimpleSystem");
436 std::unordered_map<const Elem *,const CouplingMatrix*> subdomain_one_couplings;
437 std::unordered_map<const Elem *,const CouplingMatrix*> subdomain_two_couplings;
439 coupling_functor( _mesh->active_subdomain_elements_begin(1),
440 _mesh->active_subdomain_elements_end(1),
442 subdomain_one_couplings );
444 coupling_functor( _mesh->active_subdomain_elements_begin(2),
445 _mesh->active_subdomain_elements_end(2),
447 subdomain_two_couplings );
450 _mesh->active_subdomain_elements_end(2) );
452 CPPUNIT_ASSERT_EQUAL( n_elems_subdomain_two, (
dof_id_type) subdomain_one_couplings.size() );
453 CPPUNIT_ASSERT_EQUAL(
dof_id_type(2*n_elems_subdomain_two), (
dof_id_type) subdomain_two_couplings.size() );
458 #ifdef LIBMESH_ENABLE_AMR
459 if( n_refinements > 0 )
465 #endif // LIBMESH_ENABLE_AMR
467 System & system = _es->get_system(
"SimpleSystem");
469 std::set<processor_id_type> sub_one_proc_ids, sub_two_proc_ids;
470 for (
auto & elem : _mesh->active_subdomain_elements_ptr_range(1))
471 sub_one_proc_ids.insert(elem->processor_id());
473 for (
auto & elem : _mesh->active_subdomain_elements_ptr_range(2))
474 sub_two_proc_ids.insert(elem->processor_id());
480 CPPUNIT_ASSERT_EQUAL( 1, (
int)sub_one_proc_ids.size() );
481 CPPUNIT_ASSERT_EQUAL( 1, (
int)sub_two_proc_ids.size() );
482 CPPUNIT_ASSERT( sub_one_proc_ids == sub_two_proc_ids );
488 for (
auto &
id : sub_one_proc_ids )
489 CPPUNIT_ASSERT( sub_two_proc_ids.find(
id) == sub_two_proc_ids.end() );
492 for (
auto &
id : sub_two_proc_ids )
493 CPPUNIT_ASSERT( sub_one_proc_ids.find(
id) == sub_one_proc_ids.end() );
500 #ifdef LIBMESH_HAVE_PETSC
511 CPPUNIT_TEST( testGhostingCouplingMatrix );
512 CPPUNIT_TEST( testGhostingNullCouplingMatrix );
513 CPPUNIT_TEST( testGhostingNullCouplingMatrixUnifRef );
515 CPPUNIT_TEST_SUITE_END();
527 this->run_ghosting_test(0,
true);
532 this->run_ghosting_test(0,
false);
537 std::unique_ptr<CouplingMatrix> coupling_matrix;
539 this->run_ghosting_test(2,
false);
546 this->build_quad_mesh(n_refinements);
549 std::unique_ptr<CouplingMatrix> coupling_matrix;
550 if (build_coupling_matrix)
551 this->setup_coupling_matrix(coupling_matrix);
559 functor.set_coupling_matrix(coupling_matrix);
574 std::unique_ptr<PointLocatorBase> point_locator = _mesh->sub_point_locator();
586 for (
const auto & elem : _mesh->active_local_subdomain_elements_ptr_range(2))
589 CPPUNIT_ASSERT_EQUAL(2, (
int)elem->subdomain_id());
591 const std::vector<libMesh::Point> & qpoints = subdomain_two_context.
get_element_fe(u_var)->get_xyz();
597 std::set<subdomain_id_type> allowed_subdomains;
598 allowed_subdomains.insert(1);
604 for (
const auto & qp : qpoints )
606 const Elem * overlapping_elem = (*point_locator)( qp, &allowed_subdomains );
607 CPPUNIT_ASSERT(overlapping_elem);
610 subdomain_one_context.
pre_fe_reinit(system,overlapping_elem);
627 CPPUNIT_TEST( testSparsityCouplingMatrix );
628 CPPUNIT_TEST( testSparsityNullCouplingMatrix );
629 CPPUNIT_TEST( testSparsityNullCouplingMatrixUnifRef );
631 CPPUNIT_TEST_SUITE_END();
643 this->run_sparsity_pattern_test(0,
true);
648 this->run_sparsity_pattern_test(0,
false);
653 this->run_sparsity_pattern_test(1,
false);
660 this->build_quad_mesh(n_refinements);
663 std::unique_ptr<CouplingMatrix> coupling_matrix;
664 if (build_coupling_matrix)
665 this->setup_coupling_matrix(coupling_matrix);
673 coupling_functor.set_coupling_matrix(coupling_matrix);
699 std::unique_ptr<PointLocatorBase> point_locator = _mesh->sub_point_locator();
710 for (
const auto & elem : _mesh->active_local_subdomain_elements_ptr_range(1))
715 std::vector<dof_id_type> & rows = subdomain_one_context.
get_dof_indices();
726 for (
const auto & elem : _mesh->active_local_subdomain_elements_ptr_range(2))
729 CPPUNIT_ASSERT_EQUAL(2, (
int)elem->subdomain_id());
731 const std::vector<libMesh::Point> & qpoints = subdomain_two_context.
get_element_fe(u_var)->get_xyz();
739 std::vector<dof_id_type> & rows = subdomain_two_context.
get_dof_indices();
748 std::set<subdomain_id_type> allowed_subdomains;
749 allowed_subdomains.insert(1);
755 for (
const auto & qp : qpoints )
757 const Elem * overlapping_elem = (*point_locator)( qp, &allowed_subdomains );
758 CPPUNIT_ASSERT(overlapping_elem);
761 subdomain_one_context.
pre_fe_reinit(system,overlapping_elem);
765 std::vector<dof_id_type> & v_indices = subdomain_one_context.
get_dof_indices(v_var);
766 std::vector<dof_id_type> columns(rows);
767 columns.insert( columns.end(), v_indices.begin(), v_indices.end() );
770 K21.
resize( rows.size(), columns.size() );
780 K12.
resize(v_indices.size(), rows.size());
793 #endif // LIBMESH_HAVE_PETSC
798 #ifdef LIBMESH_HAVE_PETSC