libMesh
overlapping_coupling_test.C
Go to the documentation of this file.
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> // libmesh_make_unique
16 
17 #include "test_comm.h"
18 #include "libmesh_cppunit.h"
19 
20 
21 using namespace libMesh;
22 
23 
24 // This coupling functor is attempting to mimic in a basic way
25 // a particular use case. The idea is that the mesh has overlapping domains
26 // (we assume only two in this test) and that the elements are not
27 // topologically connected between the overlapping domains, but when the user
28 // solves their problem, some dofs on are coupled between the overlapping
29 // elements. So we want to exercise the GhostingFunctor for this scenario,
30 // in particular both algebraic ghosting (dof ghosting) and coupling
31 // ghosting (sparsity pattern). In the use case, the coupling happens at
32 // quadrature points, so we look for the overlap based on the quadrature points.
33 //
34 // Additionally, the use case also has an issue that the coupling functor
35 // cannot be added until after the System is initialized because the overlap
36 // depends on the solution. So we will exercise that case as well.
38 {
39  public:
40 
42  : GhostingFunctor(),
43  _system(system),
44  _mesh(system.get_mesh()),
45  _point_locator(system.get_mesh().sub_point_locator()),
46  _coupling_matrix(nullptr)
47  {
48  // This is a highly specalized coupling functor where we are assuming
49  // only two subdomains that are overlapping so make sure there's only two.
50  CPPUNIT_ASSERT_EQUAL(2, (int)_mesh.n_subdomains());
51 
52  // We're assuming a specific set of subdomain ids
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() );
57 
58  _point_locator->enable_out_of_mesh_mode();
59  }
60 
62 
63  void set_coupling_matrix (std::unique_ptr<CouplingMatrix> & coupling_matrix)
64  { _coupling_matrix = std::move(coupling_matrix); }
65 
66  virtual void operator()
67  (const MeshBase::const_element_iterator & range_begin,
68  const MeshBase::const_element_iterator & range_end,
70  std::unordered_map<const Elem *,const CouplingMatrix*> & coupled_elements) override
71  {
72  std::unique_ptr<PointLocatorBase> sub_point_locator = _mesh.sub_point_locator();
73 
74  for( const auto & elem : as_range(range_begin,range_end) )
75  {
76  std::set<subdomain_id_type> allowed_subdomains;
77  unsigned int var = 0;
78  if( elem->subdomain_id() == 1 )
79  {
80  var = _system.variable_number("V");
81  allowed_subdomains.insert(2);
82  }
83  else if( elem->subdomain_id() == 2 )
84  {
85  var = _system.variable_number("U");
86  allowed_subdomains.insert(1);
87  }
88  else
89  CPPUNIT_FAIL("Something bad happpend");
90 
91  unsigned int dim = 2;
92  FEType fe_type = _system.variable_type(var);
93  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
94  QGauss qrule (dim, fe_type.default_quadrature_order());
95  fe->attach_quadrature_rule (&qrule);
96 
97  const std::vector<libMesh::Point> & qpoints = fe->get_xyz();
98 
99  fe->reinit(elem);
100 
101  for ( const auto & qp : qpoints )
102  {
103  const Elem * overlapping_elem = (*sub_point_locator)( qp, &allowed_subdomains );
104 
105  if( overlapping_elem && overlapping_elem->processor_id() != p )
106  coupled_elements.insert( std::make_pair(overlapping_elem,_coupling_matrix.get()) );
107  }
108  }
109  }
110 
111 private:
112 
114 
115  const MeshBase & _mesh;
116 
117  std::unique_ptr<PointLocatorBase> _point_locator;
118 
119  std::unique_ptr<CouplingMatrix> _coupling_matrix;
120 };
121 
122 
123 // We use a custom partioner for this test to help ensure we're
124 // testing the cases we want. In particular, we're going to
125 // ensure that elements in subdomain one are on different processors
126 // from elements subdomain two so we know ghosting will be required
127 // for algebraic and coupling ghosting between the two subdomains.
129 {
130 public:
131 
132  OverlappingTestPartitioner () = default;
137  virtual ~OverlappingTestPartitioner() = default;
138 
142  virtual std::unique_ptr<Partitioner> clone () const override
143  {
144  return libmesh_make_unique<OverlappingTestPartitioner>(*this);
145  }
146 
147 protected:
148 
152  virtual void _do_partition (MeshBase & mesh,
153  const unsigned int n) override
154  {
155  // If we're on one partition, then everyone gets to be on that partition
156  if (n == 1)
157  this->single_partition_range (mesh.active_elements_begin(), mesh.active_elements_end());
158  else
159  {
160  libmesh_assert_greater (n, 0);
161 
162  // If there are more partitions than elements, then just
163  // assign one element to each partition
164  if (mesh.n_active_elem() <= n)
165  {
166  processor_id_type e = 0;
167  for (auto & elem : mesh.active_element_ptr_range() )
168  {
169  elem->processor_id() = e;
170  e++;
171  }
172  }
173  // Otherwise, we'll split up the partitions into two groups
174  // and then assign the elements of each subdomain into each of
175  // the two groups
176  else
177  {
178  unsigned int n_sub_two_elems = std::distance( mesh.active_subdomain_elements_begin(2),
180 
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;
184  else
185  n_parts_sub_two = n/2;
186 
187  const unsigned int n_parts_sub_one = n - n_parts_sub_two;
188 
189  const dof_id_type sub_two_blk_size = cast_int<dof_id_type>
191  mesh.active_subdomain_elements_end(2) )/n_parts_sub_two );
192 
193  this->assign_proc_id_subdomain( mesh, 2, sub_two_blk_size, n_parts_sub_two, 0 );
194 
195 
196  const dof_id_type sub_one_blk_size = cast_int<dof_id_type>
198  mesh.active_subdomain_elements_end(1) )/n_parts_sub_one );
199 
200  this->assign_proc_id_subdomain( mesh, 1, sub_one_blk_size, n_parts_sub_one, n_parts_sub_two );
201  }
202  }
203  }
204 
206  const subdomain_id_type sid,
207  const dof_id_type blksize,
208  const unsigned int n_parts,
209  const processor_id_type offset )
210  {
211  dof_id_type e = 0;
212  for (auto & elem : mesh.active_subdomain_elements_ptr_range(sid))
213  {
214  if ((e/blksize) < n_parts)
215  elem->processor_id() = offset + cast_int<processor_id_type>(e/blksize);
216  else
217  elem->processor_id() = offset;
218 
219  e++;
220  }
221  }
222 
223 };
224 
225 
226 // Common functionality for all the different tests we'll be running
227 // We're creating two elements that will live in subdomain 1 and one
228 // element that is lives in subdomain 2 that overlaps both elements in
229 // subdomain 1, but is not topologically connected.
231 {
232 protected:
233 
234  std::unique_ptr<MeshBase> _mesh;
235 
236  std::unique_ptr<EquationSystems> _es;
237 
238  void build_quad_mesh(unsigned int n_refinements = 0)
239  {
240  // We are making assumptions in various places about the presence
241  // of the elements on the current processor so we're restricting to
242  // ReplicatedMesh for now.
243  _mesh = libmesh_make_unique<ReplicatedMesh>(*TestCommWorld);
244 
245  _mesh->set_mesh_dimension(2);
246 
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 );
251 
252  {
253  Elem* elem = _mesh->add_elem( new Quad4 );
254  elem->set_id(0);
255  elem->subdomain_id() = 1;
256 
257  for (unsigned int n=0; n<4; n++)
258  elem->set_node(n) = _mesh->node_ptr(n);
259  }
260 
261  _mesh->add_point( Point(1.0,2.0),4 );
262  _mesh->add_point( Point(0.0,2.0),5 );
263 
264  {
265  Elem* elem = _mesh->add_elem( new Quad4 );
266  elem->set_id(1);
267  elem->subdomain_id() = 1;
268 
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);
273  }
274 
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 );
279 
280  {
281  Elem* elem = _mesh->add_elem( new Quad4 );
282  elem->set_id(2);
283  elem->subdomain_id() = 2;
284 
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);
289  }
290 
291  _mesh->partitioner() = libmesh_make_unique<OverlappingTestPartitioner>();
292 
293  _mesh->prepare_for_use();
294 
295 #ifdef LIBMESH_ENABLE_AMR
296  if (n_refinements > 0)
297  {
298  MeshRefinement refine(*_mesh);
299  refine.uniformly_refine(n_refinements);
300  }
301 #endif // LIBMESH_ENABLE_AMR
302  }
303 
305  {
306  _es = libmesh_make_unique<EquationSystems>(*_mesh);
307  LinearImplicitSystem & sys = _es->add_system<LinearImplicitSystem> ("SimpleSystem");
308 
309  std::set<subdomain_id_type> sub_one;
310  sub_one.insert(1);
311 
312  std::set<subdomain_id_type> sub_two;
313  sub_two.insert(2);
314 
315  sys.add_variable("U", FIRST, LAGRANGE, &sub_two);
316  sys.add_variable("L", FIRST, LAGRANGE, &sub_two);
317 
318  sys.add_variable("V", FIRST, LAGRANGE, &sub_one);
319  sys.add_variable("p", FIRST, LAGRANGE, &sub_one);
320 
321  _es->init();
322  }
323 
324  void clear()
325  {
326  _es.reset();
327  _mesh.reset();
328  }
329 
330  void setup_coupling_matrix( std::unique_ptr<CouplingMatrix> & coupling )
331  {
332  System & system = _es->get_system("SimpleSystem");
333 
334  coupling = libmesh_make_unique<CouplingMatrix>(system.n_vars());
335 
336  const unsigned int u_var = system.variable_number("U");
337  const unsigned int l_var = system.variable_number("L");
338  const unsigned int v_var = system.variable_number("V");
339  const unsigned int p_var = system.variable_number("p");
340 
341  // Only adding the overlapping couplings since the primary ones should
342  // be there by default.
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;
348  }
349 
350 };
351 
352 // Mainly just to sanity check the mesh construction and
353 // the assumptions made in the OverlappingCouplingFunctor
354 // as well as the custom Partitioner.
355 class OverlappingFunctorTest : public CppUnit::TestCase,
356  public OverlappingTestBase
357 {
358 public:
359 
360  CPPUNIT_TEST_SUITE( OverlappingFunctorTest );
361 
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 );
369 #endif
370 
371  CPPUNIT_TEST_SUITE_END();
372 
373 public:
374 
375  public:
376  void setUp()
377  {
378  this->build_quad_mesh();
379  this->init(*_mesh);
380  }
381 
382  void tearDown()
383  { this->clear(); }
384 
386  { this->run_coupling_functor_test(0); }
387 
389  { this->run_coupling_functor_test(1); }
390 
392  {
394  _es->reinit();
395  this->run_coupling_functor_test(0);
396  }
397 
399  {
401  _es->reinit();
402  this->run_coupling_functor_test(1);
403  }
404 
406  {
407  this->run_partitioner_test(0);
408  }
409 
411  {
412  this->run_partitioner_test(1);
413  }
414 
415 private:
416 
417  // This is basically to sanity check the coupling functor
418  // with the supplied mesh to make sure all the assumptions
419  // are kosher. This test requires AMR and some kind of a
420  // linear solver, so let's only run it if we have PETSc.
421  void run_coupling_functor_test(unsigned int n_refinements)
422  {
423 #if defined(LIBMESH_ENABLE_AMR) && defined(LIBMESH_HAVE_SOLVER)
424  if( n_refinements > 0 )
425  {
426  MeshRefinement refine(*_mesh);
427  refine.uniformly_refine(n_refinements);
428  _es->reinit();
429  }
430 #endif
431 
432  System & system = _es->get_system("SimpleSystem");
433 
434  OverlappingCouplingFunctor coupling_functor(system);
435 
436  std::unordered_map<const Elem *,const CouplingMatrix*> subdomain_one_couplings;
437  std::unordered_map<const Elem *,const CouplingMatrix*> subdomain_two_couplings;
438 
439  coupling_functor( _mesh->active_subdomain_elements_begin(1),
440  _mesh->active_subdomain_elements_end(1),
442  subdomain_one_couplings );
443 
444  coupling_functor( _mesh->active_subdomain_elements_begin(2),
445  _mesh->active_subdomain_elements_end(2),
447  subdomain_two_couplings );
448 
449  dof_id_type n_elems_subdomain_two = std::distance( _mesh->active_subdomain_elements_begin(2),
450  _mesh->active_subdomain_elements_end(2) );
451 
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() );
454  }
455 
456  void run_partitioner_test(unsigned int n_refinements)
457  {
458 #ifdef LIBMESH_ENABLE_AMR
459  if( n_refinements > 0 )
460  {
461  MeshRefinement refine(*_mesh);
462  refine.uniformly_refine(n_refinements);
463  _es->reinit();
464  }
465 #endif // LIBMESH_ENABLE_AMR
466 
467  System & system = _es->get_system("SimpleSystem");
468 
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());
472 
473  for (auto & elem : _mesh->active_subdomain_elements_ptr_range(2))
474  sub_two_proc_ids.insert(elem->processor_id());
475 
476 
477  // Everyone should be on the same processor if only 1 processor
478  if( system.n_processors() == 1 )
479  {
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 );
483  }
484  // Otherwise these sets should be disjoint
485  else
486  {
487  // Make sure no subdomain one ids in the subdomain two ids
488  for (auto & id : sub_one_proc_ids )
489  CPPUNIT_ASSERT( sub_two_proc_ids.find(id) == sub_two_proc_ids.end() );
490 
491  // Vice-versa
492  for (auto & id : sub_two_proc_ids )
493  CPPUNIT_ASSERT( sub_one_proc_ids.find(id) == sub_one_proc_ids.end() );
494  }
495  }
496 
497 };
498 
499 // In this testing, we're relying on the presence of PETSc
500 #ifdef LIBMESH_HAVE_PETSC
501 
502 
503 // This testing class the algebraic ghosting of the
504 // OverlappingCouplingFunctor.
505 class OverlappingAlgebraicGhostingTest : public CppUnit::TestCase,
506  public OverlappingTestBase
507 {
508 public:
509  CPPUNIT_TEST_SUITE( OverlappingAlgebraicGhostingTest );
510 
511  CPPUNIT_TEST( testGhostingCouplingMatrix );
512  CPPUNIT_TEST( testGhostingNullCouplingMatrix );
513  CPPUNIT_TEST( testGhostingNullCouplingMatrixUnifRef );
514 
515  CPPUNIT_TEST_SUITE_END();
516 
517 public:
518 
519  void setUp()
520  {}
521 
522  void tearDown()
523  { this->clear(); }
524 
526  {
527  this->run_ghosting_test(0, true);
528  }
529 
531  {
532  this->run_ghosting_test(0, false);
533  }
534 
536  {
537  std::unique_ptr<CouplingMatrix> coupling_matrix;
538 
539  this->run_ghosting_test(2, false);
540  }
541 
542 private:
543 
544  void run_ghosting_test(const unsigned int n_refinements, bool build_coupling_matrix)
545  {
546  this->build_quad_mesh(n_refinements);
547  this->init(*_mesh);
548 
549  std::unique_ptr<CouplingMatrix> coupling_matrix;
550  if (build_coupling_matrix)
551  this->setup_coupling_matrix(coupling_matrix);
552 
553  LinearImplicitSystem & system = _es->get_system<LinearImplicitSystem>("SimpleSystem");
554 
555  // If we don't add this coupling functor and properly recompute the
556  // sparsity pattern, then PETSc will throw a malloc error when we
557  // try to assemble into the global matrix
558  OverlappingCouplingFunctor functor(system);
559  functor.set_coupling_matrix(coupling_matrix);
560 
561  DofMap & dof_map = system.get_dof_map();
562  dof_map.add_algebraic_ghosting_functor(functor);
563  dof_map.reinit_send_list(system.get_mesh());
564 
565  // Update current local solution
567 
568  system.current_local_solution->init(system.n_dofs(), system.n_local_dofs(),
569  dof_map.get_send_list(), false,
571 
572  system.solution->localize(*(system.current_local_solution),dof_map.get_send_list());
573 
574  std::unique_ptr<PointLocatorBase> point_locator = _mesh->sub_point_locator();
575 
576  const unsigned int u_var = system.variable_number("U");
577 
579 
580  FEMContext subdomain_one_context(system);
581  FEMContext subdomain_two_context(system);
582 
583  // The use case on which this test is based, we only add terms to the residual
584  // corresponding to the dofs in the second subdomain, but that have couplings
585  // to dofs in the first subdomain.
586  for (const auto & elem : _mesh->active_local_subdomain_elements_ptr_range(2))
587  {
588  // A little extra unit testing on the range iterator
589  CPPUNIT_ASSERT_EQUAL(2, (int)elem->subdomain_id());
590 
591  const std::vector<libMesh::Point> & qpoints = subdomain_two_context.get_element_fe(u_var)->get_xyz();
592 
593  // Setup the context for the current element
594  subdomain_two_context.pre_fe_reinit(system,elem);
595  subdomain_two_context.elem_fe_reinit();
596 
597  std::set<subdomain_id_type> allowed_subdomains;
598  allowed_subdomains.insert(1);
599 
600  // Now loop over the quadrature points and find the subdomain-one element that overlaps
601  // with the current subdomain-two element and then initialize the FEMContext for the
602  // subdomain-one element. If the algebraic ghosting has not been done properly,
603  // this will error.
604  for ( const auto & qp : qpoints )
605  {
606  const Elem * overlapping_elem = (*point_locator)( qp, &allowed_subdomains );
607  CPPUNIT_ASSERT(overlapping_elem);
608 
609  // Setup the context for the overlapping element
610  subdomain_one_context.pre_fe_reinit(system,overlapping_elem);
611  subdomain_one_context.elem_fe_reinit();
612  }
613  }
614  }
615 
616 };
617 
618 
619 // This testing class now exercises testing the sparsity
620 // pattern augmented by the OverlappingCouplingFunctor
621 class OverlappingCouplingGhostingTest : public CppUnit::TestCase,
622  public OverlappingTestBase
623 {
624 public:
625  CPPUNIT_TEST_SUITE( OverlappingCouplingGhostingTest );
626 
627  CPPUNIT_TEST( testSparsityCouplingMatrix );
628  CPPUNIT_TEST( testSparsityNullCouplingMatrix );
629  CPPUNIT_TEST( testSparsityNullCouplingMatrixUnifRef );
630 
631  CPPUNIT_TEST_SUITE_END();
632 
633 public:
634 
635  void setUp()
636  {}
637 
638  void tearDown()
639  { this->clear(); }
640 
642  {
643  this->run_sparsity_pattern_test(0, true);
644  }
645 
647  {
648  this->run_sparsity_pattern_test(0, false);
649  }
650 
652  {
653  this->run_sparsity_pattern_test(1, false);
654  }
655 
656 private:
657 
658  void run_sparsity_pattern_test(const unsigned int n_refinements, bool build_coupling_matrix)
659  {
660  this->build_quad_mesh(n_refinements);
661  this->init(*_mesh);
662 
663  std::unique_ptr<CouplingMatrix> coupling_matrix;
664  if (build_coupling_matrix)
665  this->setup_coupling_matrix(coupling_matrix);
666 
667  LinearImplicitSystem & system = _es->get_system<LinearImplicitSystem>("SimpleSystem");
668 
669  // If we don't add this coupling functor and properly recompute the
670  // sparsity pattern, then PETSc will throw a malloc error when we
671  // try to assemble into the global matrix
672  OverlappingCouplingFunctor coupling_functor(system);
673  coupling_functor.set_coupling_matrix(coupling_matrix);
674 
675  DofMap & dof_map = system.get_dof_map();
676  dof_map.add_coupling_functor(coupling_functor);
677  dof_map.reinit_send_list(system.get_mesh());
678 
679  // Update current local solution
681 
682  system.current_local_solution->init(system.n_dofs(), system.n_local_dofs(),
683  dof_map.get_send_list(), false,
685 
686  system.solution->localize(*(system.current_local_solution),dof_map.get_send_list());
687 
688  // Now that we've added the coupling functor, we need
689  // to recompute the sparsity
690  dof_map.clear_sparsity();
691  dof_map.compute_sparsity(system.get_mesh());
692 
693  // Now that we've recomputed the sparsity pattern, we need
694  // to reinitialize the system matrix.
695  libMesh::SparseMatrix<libMesh::Number> & matrix = system.get_matrix("System Matrix");
696  libmesh_assert(dof_map.is_attached(matrix));
697  matrix.init();
698 
699  std::unique_ptr<PointLocatorBase> point_locator = _mesh->sub_point_locator();
700 
701  const unsigned int u_var = system.variable_number("U");
702  const unsigned int v_var = system.variable_number("V");
703 
704  DenseMatrix<Number> K12, K21;
705 
706  FEMContext subdomain_one_context(system);
707  FEMContext subdomain_two_context(system);
708 
709  // Add normally coupled parts of the matrix
710  for (const auto & elem : _mesh->active_local_subdomain_elements_ptr_range(1))
711  {
712  subdomain_one_context.pre_fe_reinit(system,elem);
713  subdomain_one_context.elem_fe_reinit();
714 
715  std::vector<dof_id_type> & rows = subdomain_one_context.get_dof_indices();
716 
717  // Fill with ones in case PETSc ignores the zeros at some point
718  std::fill( subdomain_one_context.get_elem_jacobian().get_values().begin(),
719  subdomain_one_context.get_elem_jacobian().get_values().end(),
720  1);
721 
722  // Insert the Jacobian for the dofs for this element
723  system.matrix->add_matrix( subdomain_one_context.get_elem_jacobian(), rows );
724  }
725 
726  for (const auto & elem : _mesh->active_local_subdomain_elements_ptr_range(2))
727  {
728  // A little extra unit testing on the range iterator
729  CPPUNIT_ASSERT_EQUAL(2, (int)elem->subdomain_id());
730 
731  const std::vector<libMesh::Point> & qpoints = subdomain_two_context.get_element_fe(u_var)->get_xyz();
732 
733  // Setup the context for the current element
734  subdomain_two_context.pre_fe_reinit(system,elem);
735  subdomain_two_context.elem_fe_reinit();
736 
737  // We're only assembling rows for the dofs on subdomain 2 (U,L), so
738  // the current element will have all those dof_indices.
739  std::vector<dof_id_type> & rows = subdomain_two_context.get_dof_indices();
740 
741  std::fill( subdomain_two_context.get_elem_jacobian().get_values().begin(),
742  subdomain_two_context.get_elem_jacobian().get_values().end(),
743  1);
744 
745  // Insert the Jacobian for the normally coupled dofs for this element
746  system.matrix->add_matrix( subdomain_two_context.get_elem_jacobian(), rows );
747 
748  std::set<subdomain_id_type> allowed_subdomains;
749  allowed_subdomains.insert(1);
750 
751  // Now loop over the quadrature points and find the subdomain-one element that overlaps
752  // with the current subdomain-two element and then add a local element matrix with
753  // the coupling to the global matrix to try and trip any issues with sparsity pattern
754  // construction
755  for ( const auto & qp : qpoints )
756  {
757  const Elem * overlapping_elem = (*point_locator)( qp, &allowed_subdomains );
758  CPPUNIT_ASSERT(overlapping_elem);
759 
760  // Setup the context for the overlapping element
761  subdomain_one_context.pre_fe_reinit(system,overlapping_elem);
762  subdomain_one_context.elem_fe_reinit();
763 
764  // We're only coupling to the "V" variable so only need those dof indices
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() );
768 
769  // This will also zero the matrix so we can just insert zeros for this test
770  K21.resize( rows.size(), columns.size() );
771 
772  std::fill(K21.get_values().begin(), K21.get_values().end(), 1);
773 
774  // Now adding this local matrix to the global would trip a PETSc
775  // malloc error if the sparsity pattern hasn't been correctly
776  // built to include the overlapping coupling.
777  system.matrix->add_matrix (K21, rows, columns);
778 
779  // Now add the other part of the overlapping coupling
780  K12.resize(v_indices.size(), rows.size());
781  std::fill(K12.get_values().begin(), K12.get_values().end(), 1);
782  system.matrix->add_matrix(K12,v_indices,rows);
783  }
784  } // end element loop
785 
786  // We need to make sure to close the matrix for this test. There could still
787  // be PETSc malloc errors tripped here if we didn't allocate the off-processor
788  // part of the sparsity pattern correctly.
789  system.matrix->close();
790  }
791 
792 };
793 #endif // LIBMESH_HAVE_PETSC
794 
795 
797 
798 #ifdef LIBMESH_HAVE_PETSC
801 #endif
libMesh::MeshBase::active_subdomain_elements_begin
virtual element_iterator active_subdomain_elements_begin(subdomain_id_type subdomain_id)=0
libMesh::SparseMatrix::close
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
OverlappingAlgebraicGhostingTest::testGhostingNullCouplingMatrixUnifRef
void testGhostingNullCouplingMatrixUnifRef()
Definition: overlapping_coupling_test.C:535
libMesh::DofMap::add_algebraic_ghosting_functor
void add_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor, bool to_mesh=true)
Adds a functor which can specify algebraic ghosting requirements for use with distributed vectors.
Definition: dof_map.C:1862
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
OverlappingFunctorTest
Definition: overlapping_coupling_test.C:355
OverlappingCouplingFunctor::_coupling_matrix
std::unique_ptr< CouplingMatrix > _coupling_matrix
Definition: overlapping_coupling_test.C:119
libMesh::System::n_vars
unsigned int n_vars() const
Definition: system.h:2155
libMesh::DofMap::compute_sparsity
void compute_sparsity(const MeshBase &)
Computes the sparsity pattern for the matrices corresponding to proc_id and sends that data to Linear...
Definition: dof_map.C:1768
libMesh::DofMap::is_attached
bool is_attached(SparseMatrix< Number > &matrix)
Matrices should not be attached more than once.
Definition: dof_map.C:310
OverlappingTestPartitioner::clone
virtual std::unique_ptr< Partitioner > clone() const override
Definition: overlapping_coupling_test.C:142
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::DiffContext::get_dof_indices
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:367
OverlappingCouplingGhostingTest::run_sparsity_pattern_test
void run_sparsity_pattern_test(const unsigned int n_refinements, bool build_coupling_matrix)
Definition: overlapping_coupling_test.C:658
libMesh::SparseMatrix::init
virtual void init(const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type nnz=30, const numeric_index_type noz=10, const numeric_index_type blocksize=1)=0
Initialize SparseMatrix with the specified sizes.
OverlappingFunctorTest::tearDown
void tearDown()
Definition: overlapping_coupling_test.C:382
OverlappingTestPartitioner
Definition: overlapping_coupling_test.C:128
libMesh::MeshBase::active_element_ptr_range
virtual SimpleRange< element_iterator > active_element_ptr_range()=0
libMesh::DofObject::set_id
dof_id_type & set_id()
Definition: dof_object.h:776
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
OverlappingAlgebraicGhostingTest::tearDown
void tearDown()
Definition: overlapping_coupling_test.C:522
libMesh::FEMContext::pre_fe_reinit
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1642
OverlappingTestBase
Definition: overlapping_coupling_test.C:230
OverlappingCouplingFunctor::OverlappingCouplingFunctor
OverlappingCouplingFunctor(System &system)
Definition: overlapping_coupling_test.C:41
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::DenseMatrix< Number >
libMesh::Partitioner
The Partitioner class provides a uniform interface for partitioning algorithms.
Definition: partitioner.h:50
libMesh::GHOSTED
Definition: enum_parallel_type.h:37
OverlappingCouplingGhostingTest::testSparsityNullCouplingMatrix
void testSparsityNullCouplingMatrix()
Definition: overlapping_coupling_test.C:646
OverlappingFunctorTest::checkOverlappingPartitioner
void checkOverlappingPartitioner()
Definition: overlapping_coupling_test.C:405
libMesh::DenseMatrix::get_values
std::vector< T > & get_values()
Definition: dense_matrix.h:371
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
OverlappingAlgebraicGhostingTest::run_ghosting_test
void run_ghosting_test(const unsigned int n_refinements, bool build_coupling_matrix)
Definition: overlapping_coupling_test.C:544
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
OverlappingAlgebraicGhostingTest::testGhostingNullCouplingMatrix
void testGhostingNullCouplingMatrix()
Definition: overlapping_coupling_test.C:530
libMesh::SparseMatrix
Generic sparse matrix.
Definition: vector_fe_ex5.C:45
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::NumericVector::build
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
Definition: numeric_vector.C:49
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
OverlappingCouplingGhostingTest::setUp
void setUp()
Definition: overlapping_coupling_test.C:635
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
OverlappingCouplingFunctor::set_coupling_matrix
void set_coupling_matrix(std::unique_ptr< CouplingMatrix > &coupling_matrix)
Definition: overlapping_coupling_test.C:63
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
OverlappingCouplingFunctor
Definition: overlapping_coupling_test.C:37
libMesh::libmesh_assert
libmesh_assert(ctx)
OverlappingCouplingGhostingTest
Definition: overlapping_coupling_test.C:621
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::MeshBase::active_elements_begin
virtual element_iterator active_elements_begin()=0
Active, local, and negation forms of the element iterators described above.
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::System::n_local_dofs
dof_id_type n_local_dofs() const
Definition: system.C:187
OverlappingTestBase::_mesh
std::unique_ptr< MeshBase > _mesh
Definition: overlapping_coupling_test.C:234
libMesh::GhostingFunctor
This abstract base class defines the interface by which library code and user code can report associa...
Definition: ghosting_functor.h:153
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
OverlappingCouplingGhostingTest::testSparsityNullCouplingMatrixUnifRef
void testSparsityNullCouplingMatrixUnifRef()
Definition: overlapping_coupling_test.C:651
libMesh::MeshBase::active_subdomain_elements_end
virtual element_iterator active_subdomain_elements_end(subdomain_id_type subdomain_id)=0
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::System::add_variable
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1069
libMesh::System::current_local_solution
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1551
OverlappingTestBase::_es
std::unique_ptr< EquationSystems > _es
Definition: overlapping_coupling_test.C:236
TestCommWorld
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:111
OverlappingCouplingFunctor::_system
System & _system
Definition: overlapping_coupling_test.C:113
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::MeshBase::active_subdomain_elements_ptr_range
virtual SimpleRange< element_iterator > active_subdomain_elements_ptr_range(subdomain_id_type subdomain_id)=0
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
OverlappingTestPartitioner::_do_partition
virtual void _do_partition(MeshBase &mesh, const unsigned int n) override
Partition the MeshBase into n subdomains.
Definition: overlapping_coupling_test.C:152
OverlappingCouplingFunctor::_mesh
const MeshBase & _mesh
Definition: overlapping_coupling_test.C:115
OverlappingFunctorTest::checkCouplingFunctorTriUnifRef
void checkCouplingFunctorTriUnifRef()
Definition: overlapping_coupling_test.C:398
OverlappingFunctorTest::run_coupling_functor_test
void run_coupling_functor_test(unsigned int n_refinements)
Definition: overlapping_coupling_test.C:421
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
OverlappingFunctorTest::checkOverlappingPartitionerUnifRef
void checkOverlappingPartitionerUnifRef()
Definition: overlapping_coupling_test.C:410
libMesh::Quad4
The QUAD4 is an element in 2D composed of 4 nodes.
Definition: face_quad4.h:51
OverlappingCouplingGhostingTest::tearDown
void tearDown()
Definition: overlapping_coupling_test.C:638
libMesh::DofMap::get_send_list
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:496
libMesh::MeshBase::const_element_iterator
The definition of the const_element_iterator struct.
Definition: mesh_base.h:1891
libMesh::Elem::set_node
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2059
OverlappingCouplingGhostingTest::testSparsityCouplingMatrix
void testSparsityCouplingMatrix()
Definition: overlapping_coupling_test.C:641
libMesh::DofMap::reinit_send_list
void reinit_send_list(MeshBase &mesh)
Clears the _send_list vector and then rebuilds it.
Definition: dof_map.C:1690
OverlappingTestPartitioner::assign_proc_id_subdomain
void assign_proc_id_subdomain(const MeshBase &mesh, const subdomain_id_type sid, const dof_id_type blksize, const unsigned int n_parts, const processor_id_type offset)
Definition: overlapping_coupling_test.C:205
OverlappingTestBase::setup_coupling_matrix
void setup_coupling_matrix(std::unique_ptr< CouplingMatrix > &coupling)
Definition: overlapping_coupling_test.C:330
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::FEMContext::elem_fe_reinit
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
Definition: fem_context.C:1436
libMesh::MeshTools::Modification::all_tri
void all_tri(MeshBase &mesh)
Converts the 2D quadrilateral elements of a Mesh into triangular elements.
Definition: mesh_modification.C:280
OverlappingTestBase::build_quad_mesh
void build_quad_mesh(unsigned int n_refinements=0)
Definition: overlapping_coupling_test.C:238
libMesh::FEMContext::get_element_fe
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh.
Definition: fem_context.h:275
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::SparseMatrix::add_matrix
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
OverlappingAlgebraicGhostingTest
Definition: overlapping_coupling_test.C:505
OverlappingFunctorTest::checkCouplingFunctorQuad
void checkCouplingFunctorQuad()
Definition: overlapping_coupling_test.C:385
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
CPPUNIT_TEST_SUITE_REGISTRATION
CPPUNIT_TEST_SUITE_REGISTRATION(OverlappingFunctorTest)
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libmesh_cppunit.h
OverlappingTestBase::init
void init(MeshBase &mesh)
Definition: overlapping_coupling_test.C:304
libMesh::DofMap::clear_sparsity
void clear_sparsity()
Clears the sparsity pattern.
Definition: dof_map.C:1800
libMesh::Elem::subdomain_id
subdomain_id_type subdomain_id() const
Definition: elem.h:2069
operator=
Iterator & operator=(const Iterator &rhs)
Assignment operator.
Definition: variant_filter_iterator.h:493
libMesh::DiffContext::get_elem_jacobian
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:283
libMesh::DofObject::invalid_processor_id
static const processor_id_type invalid_processor_id
An invalid processor_id to distinguish DoFs that have not been assigned to a processor.
Definition: dof_object.h:432
OverlappingCouplingFunctor::~OverlappingCouplingFunctor
virtual ~OverlappingCouplingFunctor()
Definition: overlapping_coupling_test.C:61
OverlappingCouplingFunctor::_point_locator
std::unique_ptr< PointLocatorBase > _point_locator
Definition: overlapping_coupling_test.C:117
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
OverlappingFunctorTest::setUp
void setUp()
Definition: overlapping_coupling_test.C:376
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::System::n_dofs
dof_id_type n_dofs() const
Definition: system.C:150
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
OverlappingFunctorTest::checkCouplingFunctorQuadUnifRef
void checkCouplingFunctorQuadUnifRef()
Definition: overlapping_coupling_test.C:388
libMesh::TestClass
Definition: id_types.h:33
test_comm.h
OverlappingFunctorTest::run_partitioner_test
void run_partitioner_test(unsigned int n_refinements)
Definition: overlapping_coupling_test.C:456
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::LinearImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Definition: linear_implicit_system.h:55
OverlappingFunctorTest::checkCouplingFunctorTri
void checkCouplingFunctorTri()
Definition: overlapping_coupling_test.C:391
libMesh::MeshBase::active_elements_end
virtual element_iterator active_elements_end()=0
libMesh::MeshRefinement::uniformly_refine
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
Definition: mesh_refinement.C:1678
libMesh::ImplicitSystem::get_matrix
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
Definition: implicit_system.C:262
OverlappingTestBase::clear
void clear()
Definition: overlapping_coupling_test.C:324
libMesh::FIRST
Definition: enum_order.h:42
libMesh::MeshBase::n_active_elem
virtual dof_id_type n_active_elem() const =0
OverlappingAlgebraicGhostingTest::setUp
void setUp()
Definition: overlapping_coupling_test.C:519
libMesh::DofMap::add_coupling_functor
void add_coupling_functor(GhostingFunctor &coupling_functor, bool to_mesh=true)
Adds a functor which can specify coupling requirements for creation of sparse matrices.
Definition: dof_map.C:1838
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
OverlappingAlgebraicGhostingTest::testGhostingCouplingMatrix
void testGhostingCouplingMatrix()
Definition: overlapping_coupling_test.C:525