8 #include "TestUtil.hpp"
35 std::cerr <<
"Usage: ";
36 std::cerr <<
"mbcoupler_test -meshes <source_mesh> <target_mesh> -itag <interp_tag> [-gnorm "
37 "<gnorm_tag>] [-ssnorm <ssnorm_tag> <ssnorm_selection>] [-ropts <roptions>] "
38 "[-outfile <out_file> [-wopts <woptions>]] [-dbgout [<dbg_file>]]"
40 std::cerr <<
" -meshes" << std::endl;
41 std::cerr <<
" Read in mesh files <source_mesh> and <target_mesh>." << std::endl;
42 std::cerr <<
" -itag" << std::endl;
43 std::cerr <<
" Interpolate tag <interp_tag> from source mesh to target mesh." << std::endl;
44 std::cerr <<
" -gnorm" << std::endl;
45 std::cerr <<
" Normalize the value of tag <gnorm_tag> over then entire mesh and save to" << std::endl;
46 std::cerr <<
" tag \"<gnorm_tag>_normf\" on the mesh set. Do this for all meshes." << std::endl;
47 std::cerr <<
" -ssnorm" << std::endl;
48 std::cerr <<
" Normalize the value of tag <ssnorm_tag> over subsets of a mesh and save to" << std::endl;
49 std::cerr <<
" tag \"<ssnorm_tag>_normf\" on the Entity Set for each subset. Subsets "
52 std::cerr <<
" using criteria in <ssnorm_selection>. Do this for all meshes." << std::endl;
53 std::cerr <<
" -ropts" << std::endl;
54 std::cerr <<
" Read in the mesh files using options in <roptions>." << std::endl;
55 std::cerr <<
" -outfile" << std::endl;
56 std::cerr <<
" Write out target mesh to <out_file>." << std::endl;
57 std::cerr <<
" -wopts" << std::endl;
58 std::cerr <<
" Write out mesh files using options in <woptions>." << std::endl;
59 std::cerr <<
" -dbgout" << std::endl;
60 std::cerr <<
" Write stdout and stderr streams to the file \'<dbg_file>.txt\'." << std::endl;
61 std::cerr <<
" -eps" << std::endl;
62 std::cerr <<
" epsilon" << std::endl;
63 std::cerr <<
" -meth <method> (0=CONSTANT, 1=LINEAR_FE, 2=QUADRATIC_FE, 3=SPECTRAL, 4=SPHERICAL)" << std::endl;
72 std::vector< std::string >& meshFiles,
74 std::string& interpTag,
75 std::string& gNormTag,
76 std::string& ssNormTag,
77 std::vector< const char* >& ssTagNames,
78 std::vector< const char* >& ssTagValues,
79 std::string& readOpts,
81 std::string& writeOpts,
86 ErrorCode report_iface_ents(
Interface* mbImpl, std::vector< ParallelComm* >& pcs,
bool print_results );
90 std::string& interpTag,
91 std::string& gNormTag,
92 std::string& ssNormTag,
93 std::vector< const char* >& ssTagNames,
94 std::vector< const char* >& ssTagValues,
96 std::vector< ParallelComm* >& pcs,
98 double& pointloc_time,
104 void reduceMax(
double& v )
108 MPI_Barrier( MPI_COMM_WORLD );
109 MPI_Allreduce( &v, &buf, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD );
114 int main(
int argc,
char** argv )
117 int ierr = MPI_Init( &argc, &argv );
118 assert( MPI_SUCCESS == ierr );
120 std::vector< const char* > ssTagNames, ssTagValues;
121 std::vector< std::string > meshFiles;
122 std::string interpTag, gNormTag, ssNormTag, readOpts, outFile, writeOpts, dbgFile;
126 double toler = 5.e-10;
128 ierr = MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
129 assert( MPI_SUCCESS == ierr );
130 ierr = MPI_Comm_rank( MPI_COMM_WORLD, &rank );
131 assert( MPI_SUCCESS == ierr );
133 MB_CHK_ERR(
get_file_options( argc, argv, nprocs, rank, meshFiles, method, interpTag, gNormTag, ssNormTag, ssTagNames,
134 ssTagValues, readOpts, outFile, writeOpts, dbgFile,
help, toler ) );
139 return MPI_Finalize();
143 if( !dbgFile.empty() )
145 std::stringstream dfname;
146 dfname << dbgFile << rank <<
".txt";
147 if( !std::freopen( dfname.str().c_str(),
"a", stdout ) )
return -1;
148 if( !std::freopen( dfname.str().c_str(),
"a", stderr ) )
return -1;
153 if( NULL == mbImpl )
return 1;
156 std::vector< ParallelComm* > pcs( meshFiles.size() );
162 for(
unsigned int i = 0; i < meshFiles.size(); i++ )
165 int index = pcs[i]->get_id();
166 std::string newReadopts;
167 std::ostringstream extraOpt;
168 extraOpt <<
";PARALLEL_COMM=" <<
index;
169 newReadopts = readOpts + extraOpt.str();
178 MB_CHK_ERR( report_iface_ents( mbImpl, pcs,
true ) );
180 double instant_time = 0.0, pointloc_time = 0.0,
interp_time = 0.0, gnorm_time = 0.0, ssnorm_time = 0.0;
182 MB_CHK_ERR( test_interpolation( mbImpl, method, interpTag, gNormTag, ssNormTag, ssTagNames, ssTagValues, roots, pcs,
183 instant_time, pointloc_time,
interp_time, gnorm_time, ssnorm_time, toler ) );
185 reduceMax( instant_time );
186 reduceMax( pointloc_time );
190 printf(
"\nMax time : %g %g %g (inst loc interp -- %d procs)\n", instant_time, pointloc_time,
interp_time,
194 if( !outFile.empty() )
199 std::string newwriteOpts;
200 std::ostringstream extraOpt;
201 if( nprocs > 1 ) extraOpt <<
";PARALLEL_COMM=" << 1;
202 newwriteOpts = writeOpts + extraOpt.str();
206 std::cout <<
"Wrote " << outFile << std::endl;
207 std::cout <<
"mbcoupler_test complete." << std::endl;
212 for(
unsigned int i = 0; i < meshFiles.size(); i++ )
217 ierr = MPI_Finalize();
218 assert( MPI_SUCCESS == ierr );
222 ErrorCode report_iface_ents(
Interface* mbImpl, std::vector< ParallelComm* >& pcs,
const bool print_results )
227 for(
unsigned int p = 0; p < pcs.size(); p++ )
229 for(
int i = 0; i < 4; i++ )
231 if( pcs[p]->get_iface_entities( -1, i, iface_ents[i] ) !=
MB_SUCCESS )
233 std::cerr <<
"get_iface_entities returned error on proc " << pcs[p]->proc_config().proc_rank()
234 <<
"; message: " << std::endl;
235 std::string last_error;
237 if( last_error.empty() )
238 std::cerr <<
"(none)" << std::endl;
240 std::cerr << last_error << std::endl;
242 if( 0 != i ) iface_ents[4].
merge( iface_ents[i] );
250 MPI_Comm_rank( MPI_COMM_WORLD, &rank );
252 if( print_results || iface_ents[0].size() != iface_ents[5].size() )
254 std::cerr <<
"Proc " << rank <<
" iface entities: " << std::endl;
255 for(
int i = 0; i < 4; i++ )
256 std::cerr <<
" " << iface_ents[i].size() <<
" " << i <<
"d iface entities." << std::endl;
257 std::cerr <<
" (" << iface_ents[5].
size() <<
" verts adj to other iface ents)" << std::endl;
278 std::vector< std::string >& meshFiles,
280 std::string& interpTag,
281 std::string& gNormTag,
282 std::string& ssNormTag,
283 std::vector< const char* >& ssTagNames,
284 std::vector< const char* >& ssTagValues,
285 std::string& readOpts,
286 std::string& outFile,
287 std::string& writeOpts,
288 std::string& dbgFile,
296 readOpts = ( nprocs > 1 ?
"PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"
297 "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=3.0.1;CPUTIME"
298 :
"PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"
299 "PARALLEL_RESOLVE_SHARED_ENTS;CPUTIME" );
301 writeOpts = ( nprocs > 1 ?
"PARALLEL=WRITE_PART;CPUTIME" :
"" );
303 std::string defaultDbgFile = argv[0];
306 bool haveMeshes =
false;
307 bool haveInterpTag =
false;
312 if( argc > 1 && argv[1] == std::string(
"-h" ) )
320 if( argv[npos] == std::string(
"-meshes" ) )
325 meshFiles.resize( numFiles );
326 for(
int i = 0; i < numFiles; i++ )
329 meshFiles[i] = argv[npos++];
332 std::cerr <<
" ERROR - missing correct number of mesh filenames" << std::endl;
339 else if( argv[npos] == std::string(
"-itag" ) )
344 interpTag = argv[npos++];
347 std::cerr <<
" ERROR - missing <interp_tag>" << std::endl;
351 haveInterpTag =
true;
353 else if( argv[npos] == std::string(
"-meth" ) )
357 if( argv[npos][0] ==
'0' )
359 else if( argv[npos][0] ==
'1' )
361 else if( argv[npos][0] ==
'2' )
363 else if( argv[npos][0] ==
'3' )
365 else if( argv[npos][0] ==
'4' )
369 std::cerr <<
" ERROR - unrecognized method number " << method << std::endl;
374 else if( argv[npos] == std::string(
"-eps" ) )
379 epsilon = atof( argv[npos++] );
382 std::cerr <<
" ERROR - missing <epsilon>" << std::endl;
386 else if( argv[npos] == std::string(
"-gnorm" ) )
391 gNormTag = argv[npos++];
394 std::cerr <<
" ERROR - missing <gnorm_tag>" << std::endl;
398 else if( argv[npos] == std::string(
"-ssnorm" ) )
403 ssNormTag = argv[npos++];
406 std::cerr <<
" ERROR - missing <ssnorm_tag>" << std::endl;
412 char* opts = argv[npos++];
413 char sep1[1] = {
';' };
414 char sep2[1] = {
'=' };
415 bool end_vals_seen =
false;
416 std::vector< char* > tmpTagOpts;
419 for(
char* i = strtok( opts, sep1 ); i; i = strtok( 0, sep1 ) )
420 tmpTagOpts.push_back( i );
423 for(
unsigned int j = 0; j < tmpTagOpts.size(); j++ )
425 char* e = strtok( tmpTagOpts[j], sep2 );
426 ssTagNames.push_back( e );
427 e = strtok( 0, sep2 );
434 std::cerr <<
" ERROR - new value seen after end of values in "
442 ssTagValues.push_back( (
const char*)valp );
447 end_vals_seen =
true;
448 ssTagValues.push_back( (
const char*)0 );
454 std::cerr <<
" ERROR - missing <ssnorm_selection>" << std::endl;
458 else if( argv[npos] == std::string(
"-ropts" ) )
463 readOpts = argv[npos++];
466 std::cerr <<
" ERROR - missing <roptions>" << std::endl;
470 else if( argv[npos] == std::string(
"-outfile" ) )
475 outFile = argv[npos++];
478 std::cerr <<
" ERROR - missing <out_file>" << std::endl;
482 else if( argv[npos] == std::string(
"-wopts" ) )
487 writeOpts = argv[npos++];
490 std::cerr <<
" ERROR - missing <woptions>" << std::endl;
494 else if( argv[npos] == std::string(
"-dbgout" ) )
500 dbgFile = argv[npos++];
502 dbgFile = defaultDbgFile;
507 std::cerr <<
" ERROR - Unrecognized parameter:" << argv[npos] << std::endl;
508 std::cerr <<
" Skipping..." << std::endl;
515 meshFiles.resize( 2 );
516 meshFiles[0] = std::string( TestDir +
"unittest/64bricks_1khex.h5m" );
517 meshFiles[1] = std::string( TestDir +
"unittest/64bricks_12ktet.h5m" );
519 std::cout <<
"Mesh files not entered; using default files " << meshFiles[0] <<
" and " << meshFiles[1]
525 interpTag =
"vertex_field";
526 std::cout <<
"Interpolation field name not given, using default of " << interpTag << std::endl;
529 #ifdef MOAB_HAVE_HDF5
532 std::stringstream dfname;
533 dfname <<
"dum" << nprocs <<
".h5m";
534 outFile = dfname.str();
535 if( 0 == rank ) std::cout <<
"No arguments given; using output file " << outFile << std::endl;
546 std::string& interpTag,
547 std::string& gNormTag,
548 std::string& ssNormTag,
549 std::vector< const char* >& ssTagNames,
550 std::vector< const char* >& ssTagValues,
552 std::vector< ParallelComm* >& pcs,
553 double& instant_time,
554 double& pointloc_time,
563 Range src_elems, targ_elems, targ_verts;
564 MB_CHK_ERR( pcs[0]->get_part_entities( src_elems, 3 ) );
569 Coupler mbc( mbImpl, pcs[0], src_elems, 0,
false );
571 mbc.initialize_tree();
574 bool specSou =
false, specTar =
false;
575 MB_CHK_ERR( mbc.initialize_spectral_elements( roots[0], roots[1], specSou, specTar ) );
577 instant_time = MPI_Wtime();
582 std::vector< double > vpos;
583 int numPointsOfInterest = 0;
590 MB_CHK_ERR( pcs[1]->get_part_entities( targ_elems, 3 ) );
593 targ_verts = targ_elems;
599 targ_verts =
subtract( targ_verts, tmp_verts );
602 numPointsOfInterest = (int)targ_verts.
size();
603 vpos.resize( 3 * targ_verts.
size() );
607 std::cout <<
"rank " << pcs[0]->proc_config().proc_rank() <<
" points of interest: " << numPointsOfInterest
609 MB_CHK_ERR( mbc.locate_points( &vpos[0], numPointsOfInterest, 0, toler ) );
616 MB_CHK_ERR( pcs[1]->get_part_entities( targ_elems, 3 ) );
617 MB_CHK_ERR( mbc.get_gl_points_on_elements( targ_elems, vpos, numPointsOfInterest ) );
618 std::cout <<
"rank " << pcs[0]->proc_config().proc_rank() <<
" points of interest: " << numPointsOfInterest
622 pointloc_time = MPI_Wtime();
625 std::vector< double > field( numPointsOfInterest );
627 MB_CHK_ERR( mbc.interpolate( method, interpTag, &field[0] ) );
632 if( !gNormTag.empty() )
641 gnorm_time = MPI_Wtime();
645 if( !ssNormTag.empty() )
647 MB_CHK_ERR( mbc.normalize_subset( roots[0], ssNormTag.c_str(), &ssTagNames[0], ssTagNames.size(), &ssTagValues[0],
650 MB_CHK_ERR( mbc.normalize_subset( roots[1], ssNormTag.c_str(), &ssTagNames[0], ssTagNames.size(), &ssTagValues[0],
654 ssnorm_time = MPI_Wtime();
656 ssnorm_time -= gnorm_time;
659 pointloc_time -= instant_time;
667 std::string newtag = interpTag +
"_TAR";
688 int ntot = numPointsOfInterest / targ_elems.
size();
690 std::string newtag = interpTag +
"_TAR";