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;
127 double toler = 5.e-10;
129 ierr = MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
130 assert( MPI_SUCCESS == ierr );
131 ierr = MPI_Comm_rank( MPI_COMM_WORLD, &rank );
132 assert( MPI_SUCCESS == ierr );
134 result =
get_file_options( argc, argv, nprocs, rank, meshFiles, method, interpTag, gNormTag, ssNormTag, ssTagNames,
135 ssTagValues, readOpts, outFile, writeOpts, dbgFile,
help, toler );
140 ierr = MPI_Finalize();
145 if( !dbgFile.empty() )
147 std::stringstream dfname;
148 dfname << dbgFile << rank <<
".txt";
149 if( !std::freopen( dfname.str().c_str(),
"a", stdout ) )
return -1;
150 if( !std::freopen( dfname.str().c_str(),
"a", stderr ) )
return -1;
155 if( NULL == mbImpl )
return 1;
158 std::vector< ParallelComm* > pcs( meshFiles.size() );
164 for(
unsigned int i = 0; i < meshFiles.size(); i++ )
167 int index = pcs[i]->get_id();
168 std::string newReadopts;
169 std::ostringstream extraOpt;
170 extraOpt <<
";PARALLEL_COMM=" << index;
171 newReadopts = readOpts + extraOpt.str();
175 result = mbImpl->
load_file( meshFiles[i].c_str(), &roots[i], newReadopts.c_str() );
MB_CHK_ERR( result );
180 result = report_iface_ents( mbImpl, pcs,
true );
MB_CHK_ERR( result );
182 double instant_time = 0.0, pointloc_time = 0.0,
interp_time = 0.0, gnorm_time = 0.0, ssnorm_time = 0.0;
185 result = test_interpolation( mbImpl, method, interpTag, gNormTag, ssNormTag, ssTagNames, ssTagValues, roots, pcs,
188 reduceMax( instant_time );
189 reduceMax( pointloc_time );
193 printf(
"\nMax time : %g %g %g (inst loc interp -- %d procs)\n", instant_time, pointloc_time,
interp_time,
197 if( !outFile.empty() )
202 std::string newwriteOpts;
203 std::ostringstream extraOpt;
204 if( nprocs > 1 ) extraOpt <<
";PARALLEL_COMM=" << 1;
205 newwriteOpts = writeOpts + extraOpt.str();
206 result = mbImpl->
write_file( outFile.c_str(), NULL, newwriteOpts.c_str(), partSets );
MB_CHK_ERR( result );
209 std::cout <<
"Wrote " << outFile << std::endl;
210 std::cout <<
"mbcoupler_test complete." << std::endl;
216 for(
unsigned int i = 0; i < meshFiles.size(); i++ )
221 ierr = MPI_Finalize();
222 assert( MPI_SUCCESS == ierr );
226 ErrorCode report_iface_ents(
Interface* mbImpl, std::vector< ParallelComm* >& pcs,
const bool print_results )
232 for(
unsigned int p = 0; p < pcs.size(); p++ )
234 for(
int i = 0; i < 4; i++ )
236 tmp_result = pcs[p]->get_iface_entities( -1, i, iface_ents[i] );
240 std::cerr <<
"get_iface_entities returned error on proc " << pcs[p]->proc_config().proc_rank()
241 <<
"; message: " << std::endl;
242 std::string last_error;
244 if( last_error.empty() )
245 std::cerr <<
"(none)" << std::endl;
247 std::cerr << last_error << std::endl;
250 if( 0 != i ) iface_ents[4].
merge( iface_ents[i] );
258 MPI_Comm_rank( MPI_COMM_WORLD, &rank );
260 if( print_results || iface_ents[0].
size() != iface_ents[5].
size() )
262 std::cerr <<
"Proc " << rank <<
" iface entities: " << std::endl;
263 for(
int i = 0; i < 4; i++ )
264 std::cerr <<
" " << iface_ents[i].
size() <<
" " << i <<
"d iface entities." << std::endl;
265 std::cerr <<
" (" << iface_ents[5].
size() <<
" verts adj to other iface ents)" << std::endl;
286 std::vector< std::string >& meshFiles,
288 std::string& interpTag,
289 std::string& gNormTag,
290 std::string& ssNormTag,
291 std::vector< const char* >& ssTagNames,
292 std::vector< const char* >& ssTagValues,
293 std::string& readOpts,
294 std::string& outFile,
295 std::string& writeOpts,
296 std::string& dbgFile,
304 readOpts = ( nprocs > 1 ?
"PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"
305 "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=3.0.1;CPUTIME"
306 :
"PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"
307 "PARALLEL_RESOLVE_SHARED_ENTS;CPUTIME" );
309 writeOpts = ( nprocs > 1 ?
"PARALLEL=WRITE_PART;CPUTIME" :
"" );
311 std::string defaultDbgFile = argv[0];
314 bool haveMeshes =
false;
315 bool haveInterpTag =
false;
320 if( argc > 1 && argv[1] == std::string(
"-h" ) )
328 if( argv[npos] == std::string(
"-meshes" ) )
333 meshFiles.resize( numFiles );
334 for(
int i = 0; i < numFiles; i++ )
337 meshFiles[i] = argv[npos++];
340 std::cerr <<
" ERROR - missing correct number of mesh filenames" << std::endl;
347 else if( argv[npos] == std::string(
"-itag" ) )
352 interpTag = argv[npos++];
355 std::cerr <<
" ERROR - missing <interp_tag>" << std::endl;
359 haveInterpTag =
true;
361 else if( argv[npos] == std::string(
"-meth" ) )
365 if( argv[npos][0] ==
'0' )
367 else if( argv[npos][0] ==
'1' )
369 else if( argv[npos][0] ==
'2' )
371 else if( argv[npos][0] ==
'3' )
373 else if( argv[npos][0] ==
'4' )
377 std::cerr <<
" ERROR - unrecognized method number " << method << std::endl;
382 else if( argv[npos] == std::string(
"-eps" ) )
387 epsilon = atof( argv[npos++] );
390 std::cerr <<
" ERROR - missing <epsilon>" << std::endl;
394 else if( argv[npos] == std::string(
"-gnorm" ) )
399 gNormTag = argv[npos++];
402 std::cerr <<
" ERROR - missing <gnorm_tag>" << std::endl;
406 else if( argv[npos] == std::string(
"-ssnorm" ) )
411 ssNormTag = argv[npos++];
414 std::cerr <<
" ERROR - missing <ssnorm_tag>" << std::endl;
420 char* opts = argv[npos++];
421 char sep1[1] = {
';' };
422 char sep2[1] = {
'=' };
423 bool end_vals_seen =
false;
424 std::vector< char* > tmpTagOpts;
427 for(
char* i = strtok( opts, sep1 ); i; i = strtok( 0, sep1 ) )
428 tmpTagOpts.push_back( i );
431 for(
unsigned int j = 0; j < tmpTagOpts.size(); j++ )
433 char* e = strtok( tmpTagOpts[j], sep2 );
434 ssTagNames.push_back( e );
435 e = strtok( 0, sep2 );
442 std::cerr <<
" ERROR - new value seen after end of values in "
450 ssTagValues.push_back( (
const char*)valp );
455 end_vals_seen =
true;
456 ssTagValues.push_back( (
const char*)0 );
462 std::cerr <<
" ERROR - missing <ssnorm_selection>" << std::endl;
466 else if( argv[npos] == std::string(
"-ropts" ) )
471 readOpts = argv[npos++];
474 std::cerr <<
" ERROR - missing <roptions>" << std::endl;
478 else if( argv[npos] == std::string(
"-outfile" ) )
483 outFile = argv[npos++];
486 std::cerr <<
" ERROR - missing <out_file>" << std::endl;
490 else if( argv[npos] == std::string(
"-wopts" ) )
495 writeOpts = argv[npos++];
498 std::cerr <<
" ERROR - missing <woptions>" << std::endl;
502 else if( argv[npos] == std::string(
"-dbgout" ) )
508 dbgFile = argv[npos++];
510 dbgFile = defaultDbgFile;
515 std::cerr <<
" ERROR - Unrecognized parameter:" << argv[npos] << std::endl;
516 std::cerr <<
" Skipping..." << std::endl;
523 meshFiles.resize( 2 );
524 meshFiles[0] = std::string( TestDir +
"unittest/64bricks_1khex.h5m" );
525 meshFiles[1] = std::string( TestDir +
"unittest/64bricks_12ktet.h5m" );
527 std::cout <<
"Mesh files not entered; using default files " << meshFiles[0] <<
" and " << meshFiles[1]
533 interpTag =
"vertex_field";
534 std::cout <<
"Interpolation field name not given, using default of " << interpTag << std::endl;
537 #ifdef MOAB_HAVE_HDF5
540 std::stringstream dfname;
541 dfname <<
"dum" << nprocs <<
".h5m";
542 outFile = dfname.str();
543 if( 0 == rank ) std::cout <<
"No arguments given; using output file " << outFile << std::endl;
554 std::string& interpTag,
555 std::string& gNormTag,
556 std::string& ssNormTag,
557 std::vector< const char* >& ssTagNames,
558 std::vector< const char* >& ssTagValues,
560 std::vector< ParallelComm* >& pcs,
561 double& instant_time,
562 double& pointloc_time,
571 Range src_elems, targ_elems, targ_verts;
577 Coupler mbc( mbImpl, pcs[0], src_elems, 0,
false );
579 mbc.initialize_tree();
582 bool specSou =
false, specTar =
false;
583 result = mbc.initialize_spectral_elements( roots[0], roots[1], specSou, specTar );
585 instant_time = MPI_Wtime();
590 std::vector< double > vpos;
591 int numPointsOfInterest = 0;
597 result = pcs[1]->get_part_entities( targ_elems, 3 );
MB_CHK_ERR( result );
600 targ_verts = targ_elems;
606 targ_verts =
subtract( targ_verts, tmp_verts );
608 numPointsOfInterest = (int)targ_verts.
size();
609 vpos.resize( 3 * targ_verts.
size() );
612 std::cout <<
"rank " << pcs[0]->proc_config().proc_rank() <<
" points of interest: " << numPointsOfInterest
614 result = mbc.locate_points( &vpos[0], numPointsOfInterest, 0, toler );
MB_CHK_ERR( result );
621 result = pcs[1]->get_part_entities( targ_elems, 3 );
MB_CHK_ERR( result );
622 result = mbc.get_gl_points_on_elements( targ_elems, vpos, numPointsOfInterest );
MB_CHK_ERR( result );
623 std::cout <<
"rank " << pcs[0]->proc_config().proc_rank() <<
" points of interest: " << numPointsOfInterest
627 pointloc_time = MPI_Wtime();
630 std::vector< double > field( numPointsOfInterest );
632 result = mbc.interpolate( method, interpTag, &field[0] );
MB_CHK_ERR( result );
637 if( !gNormTag.empty() )
646 gnorm_time = MPI_Wtime();
650 if( !ssNormTag.empty() )
653 result = mbc.normalize_subset( roots[0], ssNormTag.c_str(), &ssTagNames[0], ssTagNames.size(), &ssTagValues[0],
656 result = mbc.normalize_subset( roots[1], ssNormTag.c_str(), &ssTagNames[0], ssTagNames.size(), &ssTagValues[0],
660 ssnorm_time = MPI_Wtime();
662 ssnorm_time -= gnorm_time;
665 pointloc_time -= instant_time;
673 std::string newtag = interpTag +
"_TAR";
694 int ntot = numPointsOfInterest / targ_elems.
size();
696 std::string newtag = interpTag +
"_TAR";