25 #include "grins_config.h"
27 #ifdef GRINS_HAVE_CPPUNIT
29 #include <libmesh/ignore_warnings.h>
30 #include <cppunit/extensions/HelperMacros.h>
31 #include <cppunit/TestCase.h>
32 #include <libmesh/restore_warnings.h>
35 #include "grins_test_paths.h"
44 #include "libmesh/elem.h"
45 #include "libmesh/getpot.h"
46 #include "libmesh/face_quad4.h"
47 #include "libmesh/face_quad9.h"
48 #include "libmesh/fe_interface.h"
49 #include "libmesh/mesh_refinement.h"
50 #include "libmesh/serial_mesh.h"
53 #include <libmesh/ignore_warnings.h>
114 std::string filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mesh_quad4_100elem_2D.in";
115 GetPot input_quad4(filename);
116 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = this->
build_mesh(input_quad4);
119 filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mesh_quad9_100elem_2D.in";
120 GetPot input_quad9(filename);
128 std::string filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mesh_quad4_100elem_2D.in";
129 GetPot input(filename);
130 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = this->
build_mesh(input);
132 libMesh::Point origin(0.0,6.5);
133 libMesh::Point end_point(10.0,0.25);
134 libMesh::Real theta =
calc_theta(origin,end_point);
137 rayfire->init(*mesh);
140 mesh->elem_ptr(0)->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
142 libMesh::MeshRefinement mr(*mesh);
143 mr.refine_elements();
145 rayfire->reinit(*mesh);
148 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->id()) ) );
151 for (
unsigned int i=0; i<4; i++)
152 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem( mesh->elem_ptr(0)->child_ptr(i)->id() ) ) );
168 std::string filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mesh_quad4_100elem_2D.in";
169 GetPot input_quad4(filename);
170 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = this->
build_mesh(input_quad4);
173 filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mesh_quad9_100elem_2D.in";
174 GetPot input_quad9(filename);
184 libMesh::Point origin(0.0,0.1);
185 libMesh::Point end_point(1.0,0.1);
186 libMesh::Real theta =
calc_theta(origin,end_point);
189 rayfire->init(*mesh);
191 libMesh::MeshRefinement mr(*mesh);
192 mr.uniformly_refine();
193 rayfire->reinit(*mesh);
194 mr.uniformly_refine();
195 rayfire->reinit(*mesh);
198 mesh->elem_ptr(0)->child_ptr(1)->child_ptr(1)->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
200 for (
unsigned int c=0; c<mesh->elem_ptr(0)->child_ptr(0)->n_children(); c++)
201 mesh->elem_ptr(0)->child_ptr(0)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
203 mr.refine_and_coarsen_elements();
204 rayfire->reinit(*mesh);
206 CPPUNIT_ASSERT( mesh->elem_ptr(0)->child_ptr(0)->active() );
208 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->id()) );
209 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(0)->id()) );
210 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(1)->child_ptr(0)->id()) );
211 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(1)->child_ptr(1)->id()) );
219 std::string filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mixed_quad_tri.in";
220 GetPot input(filename);
221 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = this->
build_mesh(input);
223 libMesh::Point origin(0.0,0.25);
224 libMesh::Point end_point(1.0,0.25);
225 libMesh::Real theta =
calc_theta(origin,end_point);
228 rayfire->init(*mesh);
231 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(0) );
232 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(6) );
233 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(9) );
234 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(8) );
237 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(2)) );
238 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(3)) );
239 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(4)) );
240 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(5)) );
241 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(1)) );
242 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(7)) );
245 libMesh::MeshRefinement mr(*mesh);
246 mr.uniformly_refine();
248 rayfire->reinit(*mesh);
251 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->id()) );
252 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->id()) );
253 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(6)->child_ptr(0)->id()) );
254 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(9)->child_ptr(1)->id()) );
255 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(8)->child_ptr(2)->id()) );
258 for (
unsigned int c=0; c<mesh->elem_ptr(6)->n_children(); c++)
260 mesh->elem_ptr(6)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
261 mesh->elem_ptr(9)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
264 mr.coarsen_elements();
265 rayfire->reinit(*mesh);
267 CPPUNIT_ASSERT( mesh->elem_ptr(6)->active() );
268 CPPUNIT_ASSERT( mesh->elem_ptr(9)->active() );
271 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->id()) );
272 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->id()) );
273 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(8)->child_ptr(2)->id()) );
274 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(6) );
275 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(9) );
277 for (
unsigned int c=0; c<mesh->elem_ptr(6)->n_children(); c++)
279 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(6)->child_ptr(c)->id())) );
280 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(9)->child_ptr(c)->id())) );
290 libMesh::MeshRefinement mr(*mesh);
291 mr.uniformly_refine();
294 libMesh::Point origin(0.0,0.1);
295 libMesh::Point end_point(1.0,0.1);
296 libMesh::Real theta =
calc_theta(origin,end_point);
299 rayfire->init(*mesh);
301 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->id()) );
302 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->id()) );
304 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(2)->id())) );
305 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(3)->id())) );
308 mr.uniformly_refine();
309 rayfire->reinit(*mesh);
311 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->child_ptr(0)->id()) );
312 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->child_ptr(1)->id()) );
313 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(0)->id()) );
314 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(1)->id()) );
317 mr.uniformly_coarsen();
318 rayfire->reinit(*mesh);
319 mr.uniformly_coarsen();
320 rayfire->reinit(*mesh);
323 CPPUNIT_ASSERT( mesh->elem_ptr(0)->active() );
326 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(0) );
329 for (
unsigned int e=1; e<mesh->n_elem(); e++)
330 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(e)) );
335 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh =
new libMesh::SerialMesh(*
TestCommWorld);
337 mesh->set_mesh_dimension(2);
339 mesh->add_point( libMesh::Point(0.0,0.0),0 );
340 mesh->add_point( libMesh::Point(1.0,0.0),1 );
341 mesh->add_point( libMesh::Point(1.0,1.0),2 );
342 mesh->add_point( libMesh::Point(0.0,1.0),3 );
343 mesh->add_point( libMesh::Point(1.0,2.0),4 );
344 mesh->add_point( libMesh::Point(0.0,2.0),5 );
346 libMesh::Elem* elem0 = mesh->add_elem(
new libMesh::Quad4 );
347 elem0->set_node(0) = mesh->node_ptr(0);
348 elem0->set_node(1) = mesh->node_ptr(1);
349 elem0->set_node(2) = mesh->node_ptr(2);
350 elem0->set_node(3) = mesh->node_ptr(3);
352 libMesh::Elem* elem1 = mesh->add_elem(
new libMesh::Quad4 );
353 elem1->set_node(0) = mesh->node_ptr(3);
354 elem1->set_node(1) = mesh->node_ptr(2);
355 elem1->set_node(2) = mesh->node_ptr(4);
356 elem1->set_node(3) = mesh->node_ptr(5);
358 mesh->prepare_for_use();
360 libMesh::Point origin(0.0,1.0);
361 libMesh::Real theta = 0.0;
364 rayfire->init(*mesh);
366 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(elem0->id()) );
367 CPPUNIT_ASSERT( !rayfire->map_to_rayfire_elem(elem1->id()) );
369 elem0->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
371 libMesh::MeshRefinement mr(*mesh);
372 mr.refine_elements();
374 rayfire->reinit(*mesh);
376 CPPUNIT_ASSERT( !rayfire->map_to_rayfire_elem(elem1->id()) );
378 CPPUNIT_ASSERT( !rayfire->map_to_rayfire_elem(elem0->child_ptr(0)->id()) );
379 CPPUNIT_ASSERT( !rayfire->map_to_rayfire_elem(elem0->child_ptr(1)->id()) );
380 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(elem0->child_ptr(2)->id()) );
381 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(elem0->child_ptr(3)->id()) );
387 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh =
new libMesh::SerialMesh(*
TestCommWorld);
389 mesh->set_mesh_dimension(2);
391 mesh->add_point( libMesh::Point(0.0,0.0),0 );
392 mesh->add_point( libMesh::Point(1.0,0.0),1 );
393 mesh->add_point( libMesh::Point(1.1,1.1),2 );
394 mesh->add_point( libMesh::Point(0.2,1.0),3 );
396 libMesh::Elem* e = mesh->add_elem(
new libMesh::Quad4 );
397 for (
unsigned int n=0; n<4; n++)
398 e->set_node(n) = mesh->node_ptr(n);
400 mesh->prepare_for_use();
402 CPPUNIT_ASSERT_EQUAL( mesh->n_elem(), (libMesh::dof_id_type)1 );
404 libMesh::Real y = 0.5250005;
405 libMesh::Real x = y/5.0;
407 libMesh::Point origin(x,y);
408 libMesh::Real theta = 0.0;
410 std::vector<unsigned int> children_in_rayfire;
411 children_in_rayfire.push_back(1);
412 children_in_rayfire.push_back(2);
413 children_in_rayfire.push_back(3);
415 std::vector<unsigned int> children_not_in_rayfire;
416 children_not_in_rayfire.push_back(0);
418 this->
test_deformed_elem(mesh,origin,theta,children_in_rayfire,children_not_in_rayfire);
426 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh =
new libMesh::SerialMesh(*
TestCommWorld);
428 mesh->set_mesh_dimension(2);
430 mesh->add_point( libMesh::Point(0.26591876082146976,0.016285303166482301),2 );
431 mesh->add_point( libMesh::Point(0.26539426622418877,0.016285726631637167),3 );
432 mesh->add_point( libMesh::Point(0.26538091506882988,0.015714299294800886),0 );
433 mesh->add_point( libMesh::Point(0.26590538328042834,0.015713833483130532),1 );
435 libMesh::Elem* e = mesh->add_elem(
new libMesh::Quad4 );
436 for (
unsigned int n=0; n<4; n++)
437 e->set_node(n) = mesh->node_ptr(n);
439 mesh->prepare_for_use();
441 CPPUNIT_ASSERT_EQUAL( mesh->n_elem(), (libMesh::dof_id_type)1 );
443 libMesh::Real x = 0.26538759034362924;
444 libMesh::Real y = 0.01600000000000000;
446 libMesh::Point origin(x,y);
447 libMesh::Real theta = 0.0;
449 std::vector<unsigned int> children_in_rayfire;
450 children_in_rayfire.push_back(0);
451 children_in_rayfire.push_back(2);
452 children_in_rayfire.push_back(3);
454 std::vector<unsigned int> children_not_in_rayfire;
455 children_not_in_rayfire.push_back(1);
457 this->
test_deformed_elem(mesh,origin,theta,children_in_rayfire,children_not_in_rayfire);
465 std::string filename = std::string(GRINS_TEST_UNIT_INPUT_SRCDIR)+
"/mesh_quad4_9elem_2D.in";
466 GetPot input_quad4(filename);
467 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = this->
build_mesh(input_quad4);
470 libMesh::MeshRefinement mr(*mesh);
471 mesh->elem_ptr(1)->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
472 mr.refine_elements();
474 CPPUNIT_ASSERT( !(mesh->elem_ptr(1)->active()) );
475 CPPUNIT_ASSERT( mesh->elem_ptr(1)->has_children() );
477 libMesh::Point origin(0.0,0.1);
478 libMesh::Real theta = 0.0;
482 rayfire->init(*mesh);
485 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(0) );
486 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(2) );
489 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(1)) );
492 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ref(1).child_ptr(0)->id()) );
493 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ref(1).child_ptr(1)->id()) );
498 libMesh::Real
calc_theta(libMesh::Point& start, libMesh::Point end)
500 return std::atan2( (end(1)-start(1)), (end(0)-start(0)) );
503 GRINS::SharedPtr<libMesh::UnstructuredMesh>
build_mesh(
const GetPot& input )
512 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh =
new libMesh::SerialMesh(*
TestCommWorld);
514 mesh->set_mesh_dimension(2);
516 mesh->add_point( libMesh::Point(0.0,0.0),0 );
517 mesh->add_point( libMesh::Point(1.0,0.0),1 );
518 mesh->add_point( libMesh::Point(1.0,1.0),2 );
519 mesh->add_point( libMesh::Point(0.0,1.0),3 );
521 libMesh::Elem* elem = mesh->add_elem(
new libMesh::Quad4 );
522 for (
unsigned int n=0; n<4; n++)
523 elem->set_node(n) = mesh->node_ptr(n);
525 mesh->prepare_for_use();
533 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh =
new libMesh::SerialMesh(*
TestCommWorld);
535 mesh->set_mesh_dimension(2);
537 mesh->add_point( libMesh::Point(0.0,0.0),0 );
538 mesh->add_point( libMesh::Point(1.0,0.0),1 );
539 mesh->add_point( libMesh::Point(1.0,1.0),2 );
540 mesh->add_point( libMesh::Point(0.0,1.0),3 );
541 mesh->add_point( libMesh::Point(0.5,0.0),4 );
542 mesh->add_point( libMesh::Point(1.0,0.5),5 );
543 mesh->add_point( libMesh::Point(0.5,1.0),6 );
544 mesh->add_point( libMesh::Point(0.0,0.5),7 );
545 mesh->add_point( libMesh::Point(0.5,0.5),8 );
547 libMesh::Elem* elem = mesh->add_elem(
new libMesh::Quad9 );
548 for (
unsigned int n=0; n<9; n++)
549 elem->set_node(n) = mesh->node_ptr(n);
551 mesh->prepare_for_use();
559 CPPUNIT_ASSERT_EQUAL( mesh->n_elem(), (libMesh::dof_id_type)1 );
561 libMesh::Point origin(0.0,0.1);
562 libMesh::Point end_point(1.0,0.2);
563 libMesh::Real theta =
calc_theta(origin,end_point);
566 rayfire->init(*mesh);
568 libMesh::Elem* elem = mesh->elem_ptr(0);
569 CPPUNIT_ASSERT(elem);
571 elem->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
573 libMesh::MeshRefinement mr(*mesh);
574 mr.refine_elements();
576 rayfire->reinit(*mesh);
578 const libMesh::Elem* rayfire_elem0 = rayfire->map_to_rayfire_elem( elem->child_ptr(0)->id() );
579 CPPUNIT_ASSERT(rayfire_elem0);
581 const libMesh::Elem* rayfire_elem1 = rayfire->map_to_rayfire_elem( elem->child_ptr(1)->id() );
582 CPPUNIT_ASSERT(rayfire_elem1);
584 CPPUNIT_ASSERT( (rayfire_elem0->node_ptr(0))->absolute_fuzzy_equals(origin) );
585 CPPUNIT_ASSERT( (rayfire_elem0->node_ptr(1))->absolute_fuzzy_equals(libMesh::Point(0.5,0.15)) );
586 CPPUNIT_ASSERT( (rayfire_elem1->node_ptr(0))->absolute_fuzzy_equals(libMesh::Point(0.5,0.15)) );
587 CPPUNIT_ASSERT( (rayfire_elem1->node_ptr(1))->absolute_fuzzy_equals(end_point) );
592 libMesh::Point origin(0.0,0.0);
593 libMesh::Point end_point(1.0,1.0);
594 libMesh::Real theta =
calc_theta(origin,end_point);
597 rayfire->init(*mesh);
599 libMesh::Elem* elem = mesh->elem_ptr(0);
600 CPPUNIT_ASSERT(elem);
602 elem->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
604 libMesh::MeshRefinement mr(*mesh);
605 mr.refine_elements();
607 rayfire->reinit(*mesh);
609 const libMesh::Elem* rayfire_elem0 = rayfire->map_to_rayfire_elem( elem->child_ptr(0)->id() );
610 CPPUNIT_ASSERT(rayfire_elem0);
612 const libMesh::Elem* rayfire_elem1 = rayfire->map_to_rayfire_elem( elem->child_ptr(3)->id() );
613 CPPUNIT_ASSERT(rayfire_elem1);
615 CPPUNIT_ASSERT( (rayfire_elem0->node_ptr(0))->absolute_fuzzy_equals(origin) );
616 CPPUNIT_ASSERT( (rayfire_elem0->node_ptr(1))->absolute_fuzzy_equals(libMesh::Point(0.5,0.5)) );
617 CPPUNIT_ASSERT( (rayfire_elem1->node_ptr(0))->absolute_fuzzy_equals(libMesh::Point(0.5,0.5)) );
618 CPPUNIT_ASSERT( (rayfire_elem1->node_ptr(1))->absolute_fuzzy_equals(end_point) );
623 libMesh::Point origin(0.4999,0.0);
624 libMesh::Point end_point(1.0,1.0);
625 libMesh::Real theta =
calc_theta(origin,end_point);
628 rayfire->init(*mesh);
630 libMesh::Elem* elem = mesh->elem_ptr(0);
631 CPPUNIT_ASSERT(elem);
633 elem->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
635 libMesh::MeshRefinement mr(*mesh);
636 mr.refine_elements();
638 rayfire->reinit(*mesh);
640 const libMesh::Elem* rayfire_elem0 = rayfire->map_to_rayfire_elem( elem->child_ptr(0)->id() );
641 CPPUNIT_ASSERT(rayfire_elem0);
643 const libMesh::Elem* rayfire_elem1 = rayfire->map_to_rayfire_elem( elem->child_ptr(1)->id() );
644 CPPUNIT_ASSERT(rayfire_elem1);
646 const libMesh::Elem* rayfire_elem2 = rayfire->map_to_rayfire_elem( elem->child_ptr(3)->id() );
647 CPPUNIT_ASSERT(rayfire_elem2);
650 const libMesh::Elem* non_rayfire_elem = rayfire->map_to_rayfire_elem( elem->child_ptr(2)->id() );
651 CPPUNIT_ASSERT(!non_rayfire_elem);
657 libMesh::Point origin(0.0,6.5);
658 libMesh::Point end_point(10.0,0.25);
659 libMesh::Real theta =
calc_theta(origin,end_point);
662 rayfire->init(*mesh);
667 std::map<unsigned int,std::vector<unsigned int> > refine_elems;
668 refine_elems[60] = std::vector<unsigned int>();
669 refine_elems[60].push_back(0);
670 refine_elems[60].push_back(1);
672 refine_elems[52] = std::vector<unsigned int>();
673 refine_elems[52].push_back(0);
675 refine_elems[34] = std::vector<unsigned int>();
676 refine_elems[34].push_back(1);
677 refine_elems[34].push_back(2);
678 refine_elems[34].push_back(3);
680 refine_elems[25] = std::vector<unsigned int>();
681 refine_elems[25].push_back(3);
683 refine_elems[26] = std::vector<unsigned int>();
684 refine_elems[26].push_back(0);
685 refine_elems[26].push_back(1);
686 refine_elems[26].push_back(2);
688 refine_elems[17] = std::vector<unsigned int>();
689 refine_elems[17].push_back(2);
690 refine_elems[17].push_back(3);
692 refine_elems[9] = std::vector<unsigned int>();
693 refine_elems[9].push_back(1);
694 refine_elems[9].push_back(2);
695 refine_elems[9].push_back(3);
699 std::map<unsigned int,std::vector<unsigned int> >::iterator it = refine_elems.begin();
700 for(; it != refine_elems.end(); it++)
702 libMesh::Elem* elem = mesh->elem_ptr( it->first );
703 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(elem->id()) );
704 elem->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
708 libMesh::MeshRefinement mr(*mesh);
709 mr.refine_elements();
712 rayfire->reinit(*mesh);
715 it = refine_elems.begin();
716 for(; it != refine_elems.end(); it++)
718 libMesh::Elem* elem = mesh->elem_ptr( it->first );
721 CPPUNIT_ASSERT_EQUAL( elem->refinement_flag(), libMesh::Elem::RefinementState::INACTIVE );
724 CPPUNIT_ASSERT( !( rayfire->map_to_rayfire_elem(elem->id()) ) );
728 std::vector<unsigned int> children = it->second;
729 unsigned int index = 0;
730 for(
unsigned int c=0; c<elem->n_children(); c++)
732 if (index<children.size())
734 if (c==children[index])
737 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem( elem->child_ptr(c)->id() ) );
741 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem( elem->child_ptr(c)->id() )) );
749 libMesh::Point origin(0.0,0.0);
750 libMesh::Point end_point(1.0,0.25);
751 libMesh::Real theta =
calc_theta(origin,end_point);
754 rayfire->init(*mesh);
756 libMesh::MeshRefinement mr(*mesh);
759 for (
unsigned int i=0; i<3; i++)
761 mr.uniformly_refine();
762 rayfire->reinit(*mesh);
765 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->child_ptr(0)->child_ptr(0)->id()) );
766 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->child_ptr(0)->child_ptr(1)->id()) );
767 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->child_ptr(1)->child_ptr(0)->id()) );
768 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(0)->child_ptr(1)->child_ptr(1)->id()) );
770 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(0)->child_ptr(2)->id()) );
771 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(0)->child_ptr(3)->id()) );
772 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(1)->child_ptr(2)->id()) );
773 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(1)->child_ptr(1)->child_ptr(3)->id()) );
780 libMesh::Point origin(0.0,6.5);
781 libMesh::Point end_point(10.0,0.25);
782 libMesh::Real theta =
calc_theta(origin,end_point);
785 rayfire->init(*mesh);
788 libMesh::MeshRefinement mr(*mesh);
789 mr.uniformly_refine();
790 rayfire->reinit(*mesh);
793 for(
unsigned int c=0; c<4; c++)
795 mesh->elem_ptr(60)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
796 mesh->elem_ptr(25)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
797 mesh->elem_ptr(26)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
798 mesh->elem_ptr(9)->child_ptr(c)->set_refinement_flag(libMesh::Elem::RefinementState::COARSEN);
801 mr.coarsen_elements();
802 rayfire->reinit(*mesh);
805 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(60) );
806 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(25) );
807 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(26) );
808 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(9) );
811 for(
unsigned int c=0; c<4; c++)
813 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(60)->child_ptr(c)->id())) );
814 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(25)->child_ptr(c)->id())) );
815 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(26)->child_ptr(c)->id())) );
816 CPPUNIT_ASSERT( !(rayfire->map_to_rayfire_elem(mesh->elem_ptr(9)->child_ptr(c)->id())) );
821 libMesh::Point & origin, libMesh::Real theta,
822 std::vector<unsigned int> & children_in_rayfire,
823 std::vector<unsigned int> & children_not_in_rayfire )
826 rayfire->init(*mesh);
828 libMesh::Elem* elem = mesh->elem_ptr(0);
829 CPPUNIT_ASSERT(elem);
831 elem->set_refinement_flag(libMesh::Elem::RefinementState::REFINE);
833 libMesh::MeshRefinement mr(*mesh);
834 mr.refine_elements();
836 rayfire->reinit(*mesh);
838 for (
unsigned int c=0; c<children_in_rayfire.size(); c++)
839 CPPUNIT_ASSERT( rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(children_in_rayfire[c])->id()) );
841 for (
unsigned int c=0; c<children_not_in_rayfire.size(); c++)
842 CPPUNIT_ASSERT( !rayfire->map_to_rayfire_elem(mesh->elem_ptr(0)->child_ptr(children_not_in_rayfire[c])->id()) );
852 #endif // GRINS_HAVE_CPPUNIT
void test_coarsen(GRINS::SharedPtr< libMesh::UnstructuredMesh > mesh)
Test coarsening specific elements on a large mesh.
void single_elems()
Refine a single element.
void test_large_mesh(GRINS::SharedPtr< libMesh::UnstructuredMesh > mesh)
Refine specific elements on a large mesh.
CPPUNIT_TEST_SUITE_REGISTRATION(AntiochAirNASA9ThermoTest)
void refine_and_coarsen()
Refine and coarsen before calling reinit()
void through_vertex_postrefinment()
After refinement, the rayfire will travel through a vertex.
GRINS::SharedPtr< libMesh::UnstructuredMesh > build_quad9_elem()
Build a single square QUAD9.
void mixed_type_mesh()
Mesh contains 2 QUAD4 and 2 TRI3.
libMesh::Parallel::Communicator * TestCommWorld
void refine_elem_not_on_rayfire()
Make sure refined elements not on the rayfire do not get added.
GRINS::SharedPtr< libMesh::UnstructuredMesh > build_quad4_elem()
Build a single square QUAD4.
CPPUNIT_TEST_SUITE(RayfireTestAMR)
void large_2D_mesh()
A 10x10 mesh with selectively refined elements.
void refined_above_unrefined()
2 QUAD4 elems, one is refined, one if not, rayfire travels boundary between them
void init_on_refined_mesh()
void test_multiple_refinements(GRINS::SharedPtr< libMesh::UnstructuredMesh > mesh)
Refine the mesh 3 times.
void test_deformed_elem(GRINS::SharedPtr< libMesh::UnstructuredMesh > &mesh, libMesh::Point &origin, libMesh::Real theta, std::vector< unsigned int > &children_in_rayfire, std::vector< unsigned int > &children_not_in_rayfire)
SharedPtr< libMesh::UnstructuredMesh > build(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
Builds the libMesh::Mesh according to input options.
void multiple_refinements()
Test for multiple successive refinements.
CPPUNIT_TEST(single_elems)
GRINS::SharedPtr< libMesh::UnstructuredMesh > build_mesh(const GetPot &input)
libMesh::Real calc_theta(libMesh::Point &start, libMesh::Point end)
void near_vertex_postrefinement()
After refinement, the rayfire will travel very near (but not through) a vertex.
void test_near_vertex(GRINS::SharedPtr< libMesh::UnstructuredMesh > mesh)
void refine_deformed_elem_near_tolerance()
void test_single_elem(GRINS::SharedPtr< libMesh::UnstructuredMesh > mesh)
Runs the test on single square QUADs of unit length and width.
void coarsen_elements()
Test of the coarsening functionality.
void test_through_vertex(GRINS::SharedPtr< libMesh::UnstructuredMesh > mesh)
void start_with_refined_mesh()
Mesh given to init() is already refined.
void refine_deformed_elem()
QUAD4 elem is deformed, rayfire goes within libMesh::TOLERANCE central node.