MOAB: Mesh Oriented datABase  (version 5.5.0)
augment_with_ghosts.cpp
Go to the documentation of this file.
1 #include "moab/ParallelComm.hpp"
2 #include "moab/Core.hpp"
3 #include "moab_mpi.h"
4 #include "TestUtil.hpp"
5 #include "MBTagConventions.hpp"
6 #include <iostream>
7 #include <sstream>
8 
9 std::string filename = TestDir + "unittest/io/p8ex1.h5m";
10 
11 using namespace moab;
12 void report_sets( moab::Core* mb, int rank, int nproc )
13 {
14  // check neumann and material sets, and see if their number of quads / hexes in them
15  const char* const shared_set_tag_names[] = { MATERIAL_SET_TAG_NAME, DIRICHLET_SET_TAG_NAME, NEUMANN_SET_TAG_NAME,
17 
18  int num_tags = sizeof( shared_set_tag_names ) / sizeof( shared_set_tag_names[0] );
19 
20  for( int p = 0; p < nproc; p++ )
21  {
22  if( rank == p )
23  {
24  std::cout << " Task no:" << rank << "\n";
25  for( int i = 0; i < num_tags; i++ )
26  {
27  Tag tag;
28  ErrorCode rval = mb->tag_get_handle( shared_set_tag_names[i], 1, MB_TYPE_INTEGER, tag, MB_TAG_ANY );CHECK_ERR( rval );
29  Range sets;
30  rval = mb->get_entities_by_type_and_tag( 0, MBENTITYSET, &tag, 0, 1, sets, Interface::UNION );CHECK_ERR( rval );
31 
32  std::vector< int > vals( sets.size() );
33  rval = mb->tag_get_data( tag, sets, &vals[0] );CHECK_ERR( rval );
34  std::cout << " sets: " << shared_set_tag_names[i] << "\n";
35  int j = 0;
36  for( Range::iterator it = sets.begin(); it != sets.end(); it++, j++ )
37  {
38  Range ents;
39  rval = mb->get_entities_by_handle( *it, ents );CHECK_ERR( rval );
40  std::cout << " set " << mb->id_from_handle( *it ) << " with tagval=" << vals[j] << " has "
41  << ents.size() << " entities\n";
42  }
43  }
44  }
45  MPI_Barrier( MPI_COMM_WORLD );
46  }
47 }
49 {
50  int nproc, rank;
51  MPI_Comm_size( MPI_COMM_WORLD, &nproc );
52  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
53  moab::Core* mb = new moab::Core();
54 
55  ErrorCode rval = MB_SUCCESS;
56 
57  char read_opts[] = "PARALLEL=READ_PART;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=3.0.1.3;"
58  "PARTITION=PARALLEL_PARTITION";
59  rval = mb->load_file( filename.c_str(), 0, read_opts );CHECK_ERR( rval );
60 
61  report_sets( mb, rank, nproc );
62  // write in serial the database , on each rank
63  std::ostringstream outfile;
64  outfile << "testReadGhost_n" << nproc << "." << rank << ".h5m";
65  // the mesh contains ghosts too, but they were not part of mat/neumann set
66  // write in serial the file, to see what tags are missing / or not
67  rval = mb->write_file( outfile.str().c_str() ); // everything on root
68  CHECK_ERR( rval );
69  delete mb;
70 }
71 
73 {
74  int nproc, rank;
75  MPI_Comm_size( MPI_COMM_WORLD, &nproc );
76  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
77  moab::Core* mb = new moab::Core();
79  ErrorCode rval = MB_SUCCESS;
80 
81  // first read in parallel, then ghost, then augment
82  char read_opts[] = "PARALLEL=READ_PART;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION=PARALLEL_PARTITION";
83  rval = mb->load_file( filename.c_str(), 0, read_opts );CHECK_ERR( rval );
84 
85  int ghost_dim = 3, bridge = 0, layers = 1, addl_ents = 3;
86  rval = pc->exchange_ghost_cells( ghost_dim, bridge, layers, addl_ents, true, true );CHECK_ERR( rval );
87 
88  //
89  pc->set_debug_verbosity( 1 );
90  rval = pc->augment_default_sets_with_ghosts( 0 );CHECK_ERR( rval );
91 
92  report_sets( mb, rank, nproc );
93 
94  // write in serial the database , on each rank
95  std::ostringstream outfile;
96  outfile << "TaskMesh_n" << nproc << "." << rank << ".h5m";
97  // the mesh contains ghosts too, but they were not part of mat/neumann set
98  // write in serial the file, to see what tags are missing / or not
99  rval = mb->write_file( outfile.str().c_str() ); // everything on root
100  CHECK_ERR( rval );
101 
102  delete pc;
103  delete mb;
104 }
105 
107 {
108  int nproc, rank;
109  MPI_Comm_size( MPI_COMM_WORLD, &nproc );
110  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
111  moab::Core* mb = new moab::Core();
112 
113  ErrorCode rval = MB_SUCCESS;
114 
115  char read_opts[] = "PARALLEL=READ_PART;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=3.0.1.3;"
116  "PARTITION=PARALLEL_PARTITION;SKIP_AUGMENT_WITH_GHOSTS";
117  rval = mb->load_file( filename.c_str(), 0, read_opts );CHECK_ERR( rval );
118 
119  report_sets( mb, rank, nproc );
120  // write in serial the database , on each rank
121  std::ostringstream outfile;
122  outfile << "testReadGhost_n" << nproc << "." << rank << ".h5m";
123  // the mesh contains ghosts too, but they were not part of mat/neumann set
124  // write in serial the file, to see what tags are missing / or not
125  rval = mb->write_file( outfile.str().c_str() ); // everything on root
126  CHECK_ERR( rval );
127  delete mb;
128 }
129 int main( int argc, char* argv[] )
130 {
131  MPI_Init( &argc, &argv );
132 
133  int result = 0;
134  if( argc >= 2 ) filename = argv[1];
135 
136  result += RUN_TEST( test_read_with_ghost );
137  result += RUN_TEST( test_read_and_ghost_after );
139 
140  MPI_Finalize();
141  return 0;
142 }