MOAB: Mesh Oriented datABase  (version 5.5.0)
gcrm_par.cpp File Reference
#include "TestUtil.hpp"
#include "moab/Core.hpp"
#include "moab/ParallelComm.hpp"
#include "moab/ProgOptions.hpp"
#include "MBParallelConventions.h"
#include "moab/ReadUtilIface.hpp"
#include "MBTagConventions.hpp"
#include <sstream>
+ Include dependency graph for gcrm_par.cpp:

Go to the source code of this file.

Functions

void test_read_onevar_trivial ()
 
void test_read_mesh_parallel_trivial ()
 
void test_gather_onevar_on_rank0 ()
 
void test_gather_onevar_on_rank1 ()
 
void test_multiple_loads_of_same_file ()
 
void read_one_cell_var (bool rcbzoltan)
 
void read_mesh_parallel (bool rcbzoltan)
 
void gather_one_cell_var (int gather_set_rank)
 
void multiple_loads_of_same_file ()
 
int main (int argc, char *argv[])
 
void test_read_onevar_rcbzoltan ()
 
void test_read_mesh_parallel_rcbzoltan ()
 

Variables

std::string example = TestDir + "unittest/io/gcrm_r3.nc"
 
std::string read_options
 
const double eps = 1e-6
 
const int layers = 3
 

Function Documentation

◆ gather_one_cell_var()

void gather_one_cell_var ( int  gather_set_rank)

Definition at line 470 of file gcrm_par.cpp.

471 {
472  Core moab;
473  Interface& mb = moab;
474 
475  EntityHandle file_set;
476  ErrorCode rval = mb.create_meshset( MESHSET_SET, file_set );CHECK_ERR( rval );
477 
478  read_options = "PARALLEL=READ_PART;PARTITION_METHOD=TRIVIAL;PARALLEL_RESOLVE_SHARED_ENTS";
479  std::ostringstream gather_set_option;
480  gather_set_option << ";GATHER_SET=" << gather_set_rank;
481  read_options += gather_set_option.str();
482 
483  rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
484 
485  ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
486  int procs = pcomm->proc_config().proc_size();
487  int rank = pcomm->proc_config().proc_rank();
488 
489  // Make sure gather_set_rank is valid
490  if( gather_set_rank < 0 || gather_set_rank >= procs ) return;
491 
492  Range cells, cells_owned;
493  rval = mb.get_entities_by_type( file_set, MBPOLYGON, cells );CHECK_ERR( rval );
494 
495  // Get local owned cells
496  rval = pcomm->filter_pstatus( cells, PSTATUS_NOT_OWNED, PSTATUS_NOT, -1, &cells_owned );CHECK_ERR( rval );
497 
498  EntityHandle gather_set = 0;
499  if( gather_set_rank == rank )
500  {
501  // Get gather set
502  ReadUtilIface* readUtilIface;
503  mb.query_interface( readUtilIface );
504  rval = readUtilIface->get_gather_set( gather_set );CHECK_ERR( rval );
505  assert( gather_set != 0 );
506  }
507 
508  Tag vorticity_tag0, gid_tag;
509  rval = mb.tag_get_handle( "vorticity0", layers, MB_TYPE_DOUBLE, vorticity_tag0, MB_TAG_DENSE );CHECK_ERR( rval );
510 
511  gid_tag = mb.globalId_tag();
512 
513  pcomm->gather_data( cells_owned, vorticity_tag0, gid_tag, gather_set, gather_set_rank );
514 
515  if( gather_set_rank == rank )
516  {
517  // Get gather set cells
518  Range gather_set_cells;
519  rval = mb.get_entities_by_type( gather_set, MBPOLYGON, gather_set_cells );CHECK_ERR( rval );
520  CHECK_EQUAL( (size_t)642, gather_set_cells.size() );
521  CHECK_EQUAL( (size_t)1, gather_set_cells.psize() );
522 
523  // Check vorticity0 tag values on 4 gather set cells: first cell, two median cells, and last
524  // cell
525  EntityHandle cell_ents[] = { gather_set_cells[0], gather_set_cells[320], gather_set_cells[321],
526  gather_set_cells[641] };
527  double vorticity0_val[4 * layers];
528  rval = mb.tag_get_data( vorticity_tag0, &cell_ents[0], 4, vorticity0_val );CHECK_ERR( rval );
529 
530  // Only check first two layers
531  // Layer 0
532  CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
533  CHECK_REAL_EQUAL( 0.131688, vorticity0_val[1 * layers], eps );
534  CHECK_REAL_EQUAL( -0.554888, vorticity0_val[2 * layers], eps );
535  CHECK_REAL_EQUAL( -0.554888, vorticity0_val[3 * layers], eps );
536  // Layer 1
537  CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
538  CHECK_REAL_EQUAL( 0.131686, vorticity0_val[1 * layers + 1], eps );
539  CHECK_REAL_EQUAL( -0.554881, vorticity0_val[2 * layers + 1], eps );
540  CHECK_REAL_EQUAL( -0.554881, vorticity0_val[3 * layers + 1], eps );
541  }
542 }

References CHECK_EQUAL, CHECK_ERR, CHECK_REAL_EQUAL, moab::Core::create_meshset(), eps, ErrorCode, example, moab::ParallelComm::filter_pstatus(), moab::ParallelComm::gather_data(), moab::Core::get_entities_by_type(), moab::ReadUtilIface::get_gather_set(), moab::ParallelComm::get_pcomm(), moab::Core::globalId_tag(), layers, moab::Core::load_file(), mb, MB_TAG_DENSE, MB_TYPE_DOUBLE, MBPOLYGON, MESHSET_SET, moab::ParallelComm::proc_config(), moab::ProcConfig::proc_rank(), moab::ProcConfig::proc_size(), moab::Range::psize(), PSTATUS_NOT, PSTATUS_NOT_OWNED, moab::Interface::query_interface(), rank, read_options, moab::Range::size(), moab::Core::tag_get_data(), and moab::Core::tag_get_handle().

Referenced by test_gather_onevar_on_rank0(), and test_gather_onevar_on_rank1().

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 40 of file gcrm_par.cpp.

41 {
42  MPI_Init( &argc, &argv );
43  int result = 0;
44 
46 #if defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_ZOLTAN )
48 #endif
49 
51 #if defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_ZOLTAN )
53 #endif
54 
57 
59 
60  MPI_Finalize();
61  return result;
62 }

References RUN_TEST, test_gather_onevar_on_rank0(), test_gather_onevar_on_rank1(), test_multiple_loads_of_same_file(), test_read_mesh_parallel_rcbzoltan(), test_read_mesh_parallel_trivial(), test_read_onevar_rcbzoltan(), and test_read_onevar_trivial().

◆ multiple_loads_of_same_file()

void multiple_loads_of_same_file ( )

Definition at line 544 of file gcrm_par.cpp.

545 {
546  Core moab;
547  Interface& mb = moab;
548 
549  // Need a file set for nomesh to work right
550  EntityHandle file_set;
551  ErrorCode rval;
552  rval = mb.create_meshset( MESHSET_SET, file_set );CHECK_ERR( rval );
553 
554  // Read first only header information, no mesh, no variable
555  read_options = "PARALLEL=READ_PART;PARTITION;NOMESH;VARIABLE=;PARTITION_METHOD=TRIVIAL";
556 
557  rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
558 
559  // Create mesh, no variable
560  read_options = "PARALLEL=READ_PART;PARTITION;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION_METHOD="
561  "TRIVIAL;VARIABLE=";
562 
563  rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
564 
565  // Read variable vorticity at timestep 0, no mesh
566  read_options = "PARALLEL=READ_PART;PARTITION;PARTITION_METHOD=TRIVIAL;NOMESH;VARIABLE="
567  "vorticity;TIMESTEP=0";
568 
569  rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
570 
571  Range local_verts;
572  rval = mb.get_entities_by_type( file_set, MBVERTEX, local_verts );CHECK_ERR( rval );
573 
574  Range local_edges;
575  rval = mb.get_entities_by_type( file_set, MBEDGE, local_edges );CHECK_ERR( rval );
576 
577  Range local_cells;
578  rval = mb.get_entities_by_type( file_set, MBPOLYGON, local_cells );CHECK_ERR( rval );
579  // No mixed elements
580  CHECK_EQUAL( (size_t)1, local_cells.psize() );
581 
582  ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
583  int procs = pcomm->proc_config().proc_size();
584  int rank = pcomm->proc_config().proc_rank();
585 
586  // Make check runs this test on two processors
587  if( 2 == procs )
588  {
589  CHECK_EQUAL( (size_t)321, local_cells.size() );
590 
591  // Check tag for cell variable vorticity at timestep 0
592  Tag vorticity_tag0;
593  rval = mb.tag_get_handle( "vorticity0", layers, MB_TYPE_DOUBLE, vorticity_tag0 );CHECK_ERR( rval );
594 
595  // Get vorticity0 tag values on 3 local cells
596  double vorticity0_val[3 * layers];
597  EntityHandle cell_ents[] = { local_cells[0], local_cells[160], local_cells[320] };
598  rval = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
599 
600  if( 0 == rank )
601  {
602  CHECK_EQUAL( (size_t)687, local_verts.size() );
603  CHECK_EQUAL( (size_t)1007, local_edges.size() );
604 
605  // Layer 0
606  CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
607  CHECK_REAL_EQUAL( -1.708188, vorticity0_val[1 * layers], eps );
608  CHECK_REAL_EQUAL( 0.131688, vorticity0_val[2 * layers], eps );
609  // Layer 1
610  CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
611  CHECK_REAL_EQUAL( -1.708164, vorticity0_val[1 * layers + 1], eps );
612  CHECK_REAL_EQUAL( 0.131686, vorticity0_val[2 * layers + 1], eps );
613  }
614  else if( 1 == rank )
615  {
616  CHECK_EQUAL( (size_t)688, local_verts.size() );
617  CHECK_EQUAL( (size_t)1008, local_edges.size() );
618 
619  // Layer 0
620  CHECK_REAL_EQUAL( -0.554888, vorticity0_val[0 * layers], eps );
621  CHECK_REAL_EQUAL( 2.434397, vorticity0_val[1 * layers], eps );
622  CHECK_REAL_EQUAL( -0.554888, vorticity0_val[2 * layers], eps );
623  // Layer 1
624  CHECK_REAL_EQUAL( -0.554881, vorticity0_val[0 * layers + 1], eps );
625  CHECK_REAL_EQUAL( 2.434363, vorticity0_val[1 * layers + 1], eps );
626  CHECK_REAL_EQUAL( -0.554881, vorticity0_val[2 * layers + 1], eps );
627  }
628  }
629 }

References CHECK_EQUAL, CHECK_ERR, CHECK_REAL_EQUAL, moab::Core::create_meshset(), eps, ErrorCode, example, moab::Core::get_entities_by_type(), moab::ParallelComm::get_pcomm(), layers, moab::Core::load_file(), mb, MB_TYPE_DOUBLE, MBEDGE, MBPOLYGON, MBVERTEX, MESHSET_SET, moab::ParallelComm::proc_config(), moab::ProcConfig::proc_rank(), moab::ProcConfig::proc_size(), moab::Range::psize(), rank, read_options, moab::Range::size(), moab::Core::tag_get_data(), and moab::Core::tag_get_handle().

Referenced by test_multiple_loads_of_same_file().

◆ read_mesh_parallel()

void read_mesh_parallel ( bool  rcbzoltan)

Definition at line 288 of file gcrm_par.cpp.

289 {
290  Core moab;
291  Interface& mb = moab;
292 
293  read_options = "PARALLEL=READ_PART;PARTITION_METHOD=TRIVIAL;PARALLEL_RESOLVE_SHARED_ENTS;VARIABLE=";
294  if( rcbzoltan )
295  read_options = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;PARALLEL_RESOLVE_SHARED_ENTS;VARIABLE=";
296 
297  ErrorCode rval = mb.load_file( example.c_str(), NULL, read_options.c_str() );CHECK_ERR( rval );
298 
299  ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
300  int procs = pcomm->proc_config().proc_size();
301  int rank = pcomm->proc_config().proc_rank();
302 
303  rval = pcomm->check_all_shared_handles();CHECK_ERR( rval );
304 
305  // Get local vertices
306  Range local_verts;
307  rval = mb.get_entities_by_type( 0, MBVERTEX, local_verts );CHECK_ERR( rval );
308 
309  int verts_num = local_verts.size();
310  if( 2 == procs )
311  {
312  if( rcbzoltan )
313  {
314  if( 0 == rank )
315  CHECK_EQUAL( 688, verts_num );
316  else if( 1 == rank )
317  CHECK_EQUAL( 687, verts_num ); // Not owned vertices included
318  }
319  else
320  {
321  if( 0 == rank )
322  CHECK_EQUAL( 687, verts_num );
323  else if( 1 == rank )
324  CHECK_EQUAL( 688, verts_num ); // Not owned vertices included
325  }
326  }
327 
328  rval = pcomm->filter_pstatus( local_verts, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval );
329 
330  verts_num = local_verts.size();
331  if( 2 == procs )
332  {
333  if( rcbzoltan )
334  {
335  if( 0 == rank )
336  CHECK_EQUAL( 688, verts_num );
337  else if( 1 == rank )
338  CHECK_EQUAL( 592, verts_num ); // Not owned vertices excluded
339  }
340  else
341  {
342  if( 0 == rank )
343  CHECK_EQUAL( 687, verts_num );
344  else if( 1 == rank )
345  CHECK_EQUAL( 593, verts_num ); // Not owned vertices excluded
346  }
347  }
348 
349  // Get local edges
350  Range local_edges;
351  rval = mb.get_entities_by_type( 0, MBEDGE, local_edges );CHECK_ERR( rval );
352 
353  int edges_num = local_edges.size();
354  if( 2 == procs )
355  {
356  if( rcbzoltan )
357  {
358  if( 0 == rank )
359  CHECK_EQUAL( 1008, edges_num );
360  else if( 1 == rank )
361  CHECK_EQUAL( 1007, edges_num ); // Not owned edges included
362  }
363  else
364  {
365  if( 0 == rank )
366  CHECK_EQUAL( 1007, edges_num );
367  else if( 1 == rank )
368  CHECK_EQUAL( 1008, edges_num ); // Not owned edges included
369  }
370  }
371 
372  rval = pcomm->filter_pstatus( local_edges, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval );
373 
374  edges_num = local_edges.size();
375  if( 2 == procs )
376  {
377  if( rcbzoltan )
378  {
379  if( 0 == rank )
380  CHECK_EQUAL( 1008, edges_num );
381  else if( 1 == rank )
382  CHECK_EQUAL( 912, edges_num ); // Not owned edges excluded
383  }
384  else
385  {
386  if( 0 == rank )
387  CHECK_EQUAL( 1007, edges_num );
388  else if( 1 == rank )
389  CHECK_EQUAL( 913, edges_num ); // Not owned edges excluded
390  }
391  }
392 
393  // Get local cells
394  Range local_cells;
395  rval = mb.get_entities_by_type( 0, MBPOLYGON, local_cells );CHECK_ERR( rval );
396  // No mixed elements
397  CHECK_EQUAL( (size_t)1, local_cells.psize() );
398 
399  int cells_num = local_cells.size();
400  if( 2 == procs )
401  {
402  if( rcbzoltan )
403  {
404  if( 0 == rank )
405  CHECK_EQUAL( 321, cells_num );
406  else
407  CHECK_EQUAL( 321, cells_num );
408  }
409  else
410  CHECK_EQUAL( 321, cells_num );
411  }
412 
413  rval = pcomm->filter_pstatus( local_cells, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval );
414 
415  cells_num = local_cells.size();
416  if( 2 == procs )
417  {
418  if( rcbzoltan )
419  {
420  if( 0 == rank )
421  CHECK_EQUAL( 321, cells_num );
422  else
423  CHECK_EQUAL( 321, cells_num );
424  }
425  else
426  CHECK_EQUAL( 321, cells_num );
427  }
428 
429  std::cout << "proc: " << rank << " verts:" << verts_num << "\n";
430 
431  int total_verts_num;
432  MPI_Reduce( &verts_num, &total_verts_num, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() );
433  if( 0 == rank )
434  {
435  std::cout << "total vertices: " << total_verts_num << "\n";
436  CHECK_EQUAL( 1280, total_verts_num );
437  }
438 
439  std::cout << "proc: " << rank << " edges:" << edges_num << "\n";
440 
441  int total_edges_num;
442  MPI_Reduce( &edges_num, &total_edges_num, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() );
443  if( 0 == rank )
444  {
445  std::cout << "total edges: " << total_edges_num << "\n";
446  CHECK_EQUAL( 1920, total_edges_num );
447  }
448 
449  std::cout << "proc: " << rank << " cells:" << cells_num << "\n";
450 
451  int total_cells_num;
452  MPI_Reduce( &cells_num, &total_cells_num, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() );
453  if( 0 == rank )
454  {
455  std::cout << "total cells: " << total_cells_num << "\n";
456  CHECK_EQUAL( 642, total_cells_num );
457  }
458 
459 #ifdef MOAB_HAVE_HDF5_PARALLEL
460  std::string write_options( "PARALLEL=WRITE_PART;" );
461 
462  std::string output_file = "test_gcrm";
463  if( rcbzoltan ) output_file += "_rcbzoltan";
464  output_file += ".h5m";
465 
466  mb.write_file( output_file.c_str(), NULL, write_options.c_str() );
467 #endif
468 }

References moab::ParallelComm::check_all_shared_handles(), CHECK_EQUAL, CHECK_ERR, ErrorCode, example, moab::ParallelComm::filter_pstatus(), moab::Core::get_entities_by_type(), moab::ParallelComm::get_pcomm(), moab::Core::load_file(), mb, MBEDGE, MBPOLYGON, MBVERTEX, moab::ProcConfig::proc_comm(), moab::ParallelComm::proc_config(), moab::ProcConfig::proc_rank(), moab::ProcConfig::proc_size(), moab::Range::psize(), PSTATUS_NOT, PSTATUS_NOT_OWNED, rank, read_options, moab::Range::size(), and moab::Core::write_file().

Referenced by test_read_mesh_parallel_rcbzoltan(), and test_read_mesh_parallel_trivial().

◆ read_one_cell_var()

void read_one_cell_var ( bool  rcbzoltan)

Definition at line 100 of file gcrm_par.cpp.

101 {
102  Core moab;
103  Interface& mb = moab;
104 
105  read_options = "PARALLEL=READ_PART;PARTITION_METHOD=TRIVIAL;NO_EDGES;VARIABLE=vorticity";
106  if( rcbzoltan )
107  read_options = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;NO_EDGES;VARIABLE=vorticity;DEBUG_IO=1";
108 
109  ErrorCode rval = mb.load_file( example.c_str(), NULL, read_options.c_str() );CHECK_ERR( rval );
110 
111  // Get local edges
112  Range local_edges;
113  rval = mb.get_entities_by_type( 0, MBEDGE, local_edges );CHECK_ERR( rval );
114  CHECK_EQUAL( (size_t)0, local_edges.size() );
115 
116  // Get local cells
117  Range local_cells;
118  rval = mb.get_entities_by_type( 0, MBPOLYGON, local_cells );CHECK_ERR( rval );
119  // No mixed elements
120  CHECK_EQUAL( (size_t)1, local_cells.psize() );
121 
122  Tag gid_tag = mb.globalId_tag();
123 
124  std::vector< int > gids( local_cells.size() );
125  rval = mb.tag_get_data( gid_tag, local_cells, &gids[0] );
126  Range local_cell_gids;
127  std::copy( gids.rbegin(), gids.rend(), range_inserter( local_cell_gids ) );
128 
129  ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
130  int procs = pcomm->proc_config().proc_size();
131  int rank = pcomm->proc_config().proc_rank();
132 
133  // Make check runs this test on two processors
134  if( 2 == procs )
135  {
136  // Check tag for cell variable vorticity at timestep 0
137  Tag vorticity_tag0;
138  rval = mb.tag_get_handle( "vorticity0", layers, MB_TYPE_DOUBLE, vorticity_tag0 );CHECK_ERR( rval );
139 
140  // Check tag for cell variable vorticity at timestep 1
141  Tag vorticity_tag1;
142  rval = mb.tag_get_handle( "vorticity1", layers, MB_TYPE_DOUBLE, vorticity_tag1 );CHECK_ERR( rval );
143 
144  // Get vorticity0 and vorticity1 tag values on 3 local cells
145  double vorticity0_val[3 * layers];
146  double vorticity1_val[3 * layers];
147 
148  if( rcbzoltan )
149  {
150  CHECK_EQUAL( (size_t)14, local_cell_gids.psize() );
151 
152  if( 0 == rank )
153  {
154  CHECK_EQUAL( (size_t)321, local_cells.size() );
155  CHECK_EQUAL( (size_t)321, local_cell_gids.size() );
156 
157  CHECK_EQUAL( 3, (int)local_cell_gids[0] );
158  CHECK_EQUAL( 162, (int)local_cell_gids[159] );
159  CHECK_EQUAL( 640, (int)local_cell_gids[318] );
160 
161  EntityHandle cell_ents[] = { local_cells[0], local_cells[159], local_cells[318] };
162  rval = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
163 
164  // Timestep 0
165  // Layer 0
166  CHECK_REAL_EQUAL( -0.725999, vorticity0_val[0 * layers], eps );
167  CHECK_REAL_EQUAL( -1.814997, vorticity0_val[1 * layers], eps );
168  CHECK_REAL_EQUAL( 0.131688, vorticity0_val[2 * layers], eps );
169  // Layer 1
170  CHECK_REAL_EQUAL( -0.725989, vorticity0_val[0 * layers + 1], eps );
171  CHECK_REAL_EQUAL( -1.814972, vorticity0_val[1 * layers + 1], eps );
172  CHECK_REAL_EQUAL( 0.131686, vorticity0_val[2 * layers + 1], eps );
173 
174  rval = mb.tag_get_data( vorticity_tag1, cell_ents, 3, vorticity1_val );CHECK_ERR( rval );
175 
176  // Timestep 1
177  // Layer 0
178  CHECK_REAL_EQUAL( -0.706871, vorticity1_val[0 * layers], eps );
179  CHECK_REAL_EQUAL( -1.767178, vorticity1_val[1 * layers], eps );
180  CHECK_REAL_EQUAL( 0.128218, vorticity1_val[2 * layers], eps );
181  // Layer 1
182  CHECK_REAL_EQUAL( -0.706861, vorticity1_val[0 * layers + 1], eps );
183  CHECK_REAL_EQUAL( -1.767153, vorticity1_val[1 * layers + 1], eps );
184  CHECK_REAL_EQUAL( 0.128216, vorticity1_val[2 * layers + 1], eps );
185  }
186  else if( 1 == rank )
187  {
188  CHECK_EQUAL( (size_t)321, local_cells.size() );
189  CHECK_EQUAL( (size_t)321, local_cell_gids.size() );
190 
191  CHECK_EQUAL( 1, (int)local_cell_gids[0] );
192  CHECK_EQUAL( 366, (int)local_cell_gids[161] );
193  CHECK_EQUAL( 1, (int)local_cell_gids[322] );
194 
195  EntityHandle cell_ents[] = { local_cells[0], local_cells[161], local_cells[322] };
196  rval = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
197 
198  // Timestep 0
199  // Layer 0
200  CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
201  CHECK_REAL_EQUAL( -1.512985, vorticity0_val[1 * layers], eps );
202  CHECK_REAL_EQUAL( 3.629994, vorticity0_val[2 * layers], eps );
203  // Layer 1
204  CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
205  CHECK_REAL_EQUAL( -1.512964, vorticity0_val[1 * layers + 1], eps );
206  CHECK_REAL_EQUAL( 3.629944, vorticity0_val[2 * layers + 1], eps );
207 
208  rval = mb.tag_get_data( vorticity_tag1, cell_ents, 3, vorticity1_val );CHECK_ERR( rval );
209 
210  // Timestep 1
211  // Layer 0
212  CHECK_REAL_EQUAL( 3.534355, vorticity1_val[0 * layers], eps );
213  CHECK_REAL_EQUAL( -1.473122, vorticity1_val[1 * layers], eps );
214  CHECK_REAL_EQUAL( 3.534355, vorticity1_val[2 * layers], eps );
215  // Layer 1
216  CHECK_REAL_EQUAL( 3.534306, vorticity1_val[0 * layers + 1], eps );
217  CHECK_REAL_EQUAL( -1.473102, vorticity1_val[1 * layers + 1], eps );
218  CHECK_REAL_EQUAL( 3.534306, vorticity1_val[2 * layers + 1], eps );
219  }
220  }
221  else
222  {
223  CHECK_EQUAL( (size_t)321, local_cells.size() );
224  CHECK_EQUAL( (size_t)321, local_cell_gids.size() );
225  CHECK_EQUAL( (size_t)1, local_cell_gids.psize() );
226 
227  EntityHandle cell_ents[] = { local_cells[0], local_cells[160], local_cells[320] };
228  rval = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
229 
230  rval = mb.tag_get_data( vorticity_tag1, cell_ents, 3, vorticity1_val );CHECK_ERR( rval );
231 
232  if( 0 == rank )
233  {
234  CHECK_EQUAL( 1, (int)local_cell_gids[0] );
235  CHECK_EQUAL( 161, (int)local_cell_gids[160] );
236  CHECK_EQUAL( 321, (int)local_cell_gids[320] );
237 
238  // Timestep 0
239  // Layer 0
240  CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
241  CHECK_REAL_EQUAL( -1.708188, vorticity0_val[1 * layers], eps );
242  CHECK_REAL_EQUAL( 0.131688, vorticity0_val[2 * layers], eps );
243  // Layer 1
244  CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
245  CHECK_REAL_EQUAL( -1.708164, vorticity0_val[1 * layers + 1], eps );
246  CHECK_REAL_EQUAL( 0.131686, vorticity0_val[2 * layers + 1], eps );
247 
248  // Timestep 1
249  // Layer 0
250  CHECK_REAL_EQUAL( 3.534355, vorticity1_val[0 * layers], eps );
251  CHECK_REAL_EQUAL( -1.663182, vorticity1_val[1 * layers], eps );
252  CHECK_REAL_EQUAL( 0.128218, vorticity1_val[2 * layers], eps );
253  // Layer 1
254  CHECK_REAL_EQUAL( 3.534306, vorticity1_val[0 * layers + 1], eps );
255  CHECK_REAL_EQUAL( -1.663160, vorticity1_val[1 * layers + 1], eps );
256  CHECK_REAL_EQUAL( 0.128216, vorticity1_val[2 * layers + 1], eps );
257  }
258  else if( 1 == rank )
259  {
260  CHECK_EQUAL( 322, (int)local_cell_gids[0] );
261  CHECK_EQUAL( 482, (int)local_cell_gids[160] );
262  CHECK_EQUAL( 642, (int)local_cell_gids[320] );
263 
264  // Timestep 0
265  // Layer 0
266  CHECK_REAL_EQUAL( -0.554888, vorticity0_val[0 * layers], eps );
267  CHECK_REAL_EQUAL( 2.434397, vorticity0_val[1 * layers], eps );
268  CHECK_REAL_EQUAL( -0.554888, vorticity0_val[2 * layers], eps );
269  // Layer 1
270  CHECK_REAL_EQUAL( -0.554881, vorticity0_val[0 * layers + 1], eps );
271  CHECK_REAL_EQUAL( 2.434363, vorticity0_val[1 * layers + 1], eps );
272  CHECK_REAL_EQUAL( -0.554881, vorticity0_val[2 * layers + 1], eps );
273 
274  // Timestep 1
275  // Layer 0
276  CHECK_REAL_EQUAL( -0.540269, vorticity1_val[0 * layers], eps );
277  CHECK_REAL_EQUAL( 2.370258, vorticity1_val[1 * layers], eps );
278  CHECK_REAL_EQUAL( -0.540269, vorticity1_val[2 * layers], eps );
279  // Layer 1
280  CHECK_REAL_EQUAL( -0.540262, vorticity1_val[0 * layers + 1], eps );
281  CHECK_REAL_EQUAL( 2.370226, vorticity1_val[1 * layers + 1], eps );
282  CHECK_REAL_EQUAL( -0.540262, vorticity1_val[2 * layers + 1], eps );
283  }
284  }
285  }
286 }

References CHECK_EQUAL, CHECK_ERR, CHECK_REAL_EQUAL, eps, ErrorCode, example, moab::Core::get_entities_by_type(), moab::ParallelComm::get_pcomm(), moab::Core::globalId_tag(), layers, moab::Core::load_file(), mb, MB_TYPE_DOUBLE, MBEDGE, MBPOLYGON, moab::ParallelComm::proc_config(), moab::ProcConfig::proc_rank(), moab::ProcConfig::proc_size(), moab::Range::psize(), rank, read_options, moab::Range::size(), moab::Core::tag_get_data(), and moab::Core::tag_get_handle().

Referenced by test_read_onevar_rcbzoltan(), and test_read_onevar_trivial().

◆ test_gather_onevar_on_rank0()

void test_gather_onevar_on_rank0 ( )

Definition at line 84 of file gcrm_par.cpp.

85 {
87 }

References gather_one_cell_var().

Referenced by main().

◆ test_gather_onevar_on_rank1()

void test_gather_onevar_on_rank1 ( )

Definition at line 89 of file gcrm_par.cpp.

90 {
92 }

References gather_one_cell_var().

Referenced by main().

◆ test_multiple_loads_of_same_file()

void test_multiple_loads_of_same_file ( )

Definition at line 94 of file gcrm_par.cpp.

95 {
97 }

References multiple_loads_of_same_file().

Referenced by main().

◆ test_read_mesh_parallel_rcbzoltan()

void test_read_mesh_parallel_rcbzoltan ( )

Definition at line 79 of file gcrm_par.cpp.

80 {
81  read_mesh_parallel( true );
82 }

References read_mesh_parallel().

Referenced by main().

◆ test_read_mesh_parallel_trivial()

void test_read_mesh_parallel_trivial ( )

Definition at line 74 of file gcrm_par.cpp.

75 {
76  read_mesh_parallel( false );
77 }

References read_mesh_parallel().

Referenced by main().

◆ test_read_onevar_rcbzoltan()

void test_read_onevar_rcbzoltan ( )

Definition at line 69 of file gcrm_par.cpp.

70 {
71  read_one_cell_var( true );
72 }

References read_one_cell_var().

Referenced by main().

◆ test_read_onevar_trivial()

void test_read_onevar_trivial ( )

Definition at line 64 of file gcrm_par.cpp.

65 {
66  read_one_cell_var( false );
67 }

References read_one_cell_var().

Referenced by main().

Variable Documentation

◆ eps

const double eps = 1e-6

◆ example

std::string example = TestDir + "unittest/io/gcrm_r3.nc"

◆ layers

const int layers = 3

◆ read_options