Mesh Oriented datABase  (version 5.6.0)
An array-based unstructured mesh library
mbcoupler_test.cpp
Go to the documentation of this file.
1 // MOAB includes
2 #include "moab/ParallelComm.hpp"
4 #include "moab/Core.hpp"
5 #include "Coupler.hpp"
6 #include "moab_mpi.h"
7 #include "ElemUtil.hpp"
8 #include "TestUtil.hpp"
9 
10 // C++ includes
11 #include <iostream>
12 #include <iomanip>
13 #include <sstream>
14 #include <cassert>
15 
16 using namespace moab;
17 
18 /*
19  Sample usages:
20  1) P_0 interpolation: ./mbcoupler_test -meshes <src_mesh> <target_mesh> -itag <interp_tag> -meth
21  0 -outfile <output> 2) P_1 interpolation: ./mbcoupler_test -meshes <src_mesh> <target_mesh> -itag
22  <interp_tag> -meth 1 -outfile <output> 3) P_1 interpolation with epsilon control: ./mbcoupler_test
23  -meshes <src_mesh> <target_mesh> -itag <interp_tag> -meth 1 -eps <tolerance for locating points;
24  say 0.01> -outfile <output> 3) P_0 interpolation with global normalization: ./mbcoupler_test
25  -meshes <src_mesh> <target_mesh> -itag <interp_tag> -meth 0 -gnorm <gnorm_tag_name> -outfile
26  <output> 4) P_1 interpolation with subset normalization: ./mbcoupler_test -meshes <src_mesh>
27  <target_mesh> -itag <interp_tag> -meth 1 -ssnorm <snorm_tag_name> <snorm_criteria: MATERIAL_SET>
28  -outfile <output> 5) P_1 interpolation for meshes on sphere ./mbcoupler_test -meshes <src_mesh>
29  <target_mesh> -itag <interp_tag> -meth 4 -outfile <output>
30 */
31 
32 // Print usage
34 {
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>]]"
39  << std::endl;
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 "
50  "are selected"
51  << std::endl;
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;
64 }
65 
66 #ifdef MOAB_HAVE_HDF5
67 
68 ErrorCode get_file_options( int argc,
69  char** argv,
70  int nprocs,
71  int rank,
72  std::vector< std::string >& meshFiles,
73  Coupler::Method& method,
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,
80  std::string& outFile,
81  std::string& writeOpts,
82  std::string& dbgFile,
83  bool& help,
84  double& epsilon );
85 
86 ErrorCode report_iface_ents( Interface* mbImpl, std::vector< ParallelComm* >& pcs, bool print_results );
87 
88 ErrorCode test_interpolation( Interface* mbImpl,
89  Coupler::Method method,
90  std::string& interpTag,
91  std::string& gNormTag,
92  std::string& ssNormTag,
93  std::vector< const char* >& ssTagNames,
94  std::vector< const char* >& ssTagValues,
95  EntityHandle* roots,
96  std::vector< ParallelComm* >& pcs,
97  double& instant_time,
98  double& pointloc_time,
99  double& interp_time,
100  double& gnorm_time,
101  double& ssnorm_time,
102  double& toler );
103 
104 void reduceMax( double& v )
105 {
106  double buf;
107 
108  MPI_Barrier( MPI_COMM_WORLD );
109  MPI_Allreduce( &v, &buf, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD );
110 
111  v = buf;
112 }
113 
114 int main( int argc, char** argv )
115 {
116  // Need to init MPI first, to tell how many procs and rank
117  int ierr = MPI_Init( &argc, &argv );
118  assert( MPI_SUCCESS == ierr );
119 
120  std::vector< const char* > ssTagNames, ssTagValues;
121  std::vector< std::string > meshFiles;
122  std::string interpTag, gNormTag, ssNormTag, readOpts, outFile, writeOpts, dbgFile;
124 
125  bool help = false;
126  double toler = 5.e-10;
127  int nprocs, rank;
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 );
132 
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 ) );
135 
136  if( help )
137  {
138  print_usage();
139  return MPI_Finalize();
140  }
141 
142  // Redirect stdout and stderr if dbgFile is not null
143  if( !dbgFile.empty() )
144  {
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;
149  }
150 
151  // Create MOAB instance based on that
152  Interface* mbImpl = new( std::nothrow ) Core();
153  if( NULL == mbImpl ) return 1;
154 
155  // Read in mesh(es)
156  std::vector< ParallelComm* > pcs( meshFiles.size() );
157 
158  // Create root sets for each mesh using the iMesh API. Then pass these
159  // to the load_file functions to be populated.
160  EntityHandle* roots = (EntityHandle*)malloc( sizeof( EntityHandle ) * meshFiles.size() );
161 
162  for( unsigned int i = 0; i < meshFiles.size(); i++ )
163  {
164  pcs[i] = new ParallelComm( mbImpl, MPI_COMM_WORLD );
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();
170 
171  MB_CHK_ERR( mbImpl->create_meshset( MESHSET_SET, roots[i] ) );
172 
173  MB_CHK_ERR( mbImpl->load_file( meshFiles[i].c_str(), &roots[i], newReadopts.c_str() ) );
174  // MB_CHK_ERR( rps[i]->load_file(meshFiles[i].c_str(), &roots[i],
175  // FileOptions(readOpts.c_str())) );
176  }
177 
178  MB_CHK_ERR( report_iface_ents( mbImpl, pcs, true ) );
179 
180  double instant_time = 0.0, pointloc_time = 0.0, interp_time = 0.0, gnorm_time = 0.0, ssnorm_time = 0.0;
181  // Test interpolation and global normalization and subset normalization
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 ) );
184 
185  reduceMax( instant_time );
186  reduceMax( pointloc_time );
187  reduceMax( interp_time );
188 
189  if( 0 == rank )
190  printf( "\nMax time : %g %g %g (inst loc interp -- %d procs)\n", instant_time, pointloc_time, interp_time,
191  nprocs );
192 
193  // Output mesh
194  if( !outFile.empty() )
195  {
196  Range partSets;
197  // Only save the target mesh
198  partSets.insert( (EntityHandle)roots[1] );
199  std::string newwriteOpts;
200  std::ostringstream extraOpt;
201  if( nprocs > 1 ) extraOpt << ";PARALLEL_COMM=" << 1;
202  newwriteOpts = writeOpts + extraOpt.str();
203  MB_CHK_ERR( mbImpl->write_file( outFile.c_str(), NULL, newwriteOpts.c_str(), partSets ) );
204  if( 0 == rank )
205  {
206  std::cout << "Wrote " << outFile << std::endl;
207  std::cout << "mbcoupler_test complete." << std::endl;
208  }
209  }
210 
211  free( roots );
212  for( unsigned int i = 0; i < meshFiles.size(); i++ )
213  delete pcs[i];
214 
215  delete mbImpl;
216 
217  ierr = MPI_Finalize();
218  assert( MPI_SUCCESS == ierr );
219  return 0;
220 }
221 
222 ErrorCode report_iface_ents( Interface* mbImpl, std::vector< ParallelComm* >& pcs, const bool print_results )
223 {
224  Range iface_ents[6];
225 
226  // Now figure out which vertices are shared
227  for( unsigned int p = 0; p < pcs.size(); p++ )
228  {
229  for( int i = 0; i < 4; i++ )
230  {
231  if( pcs[p]->get_iface_entities( -1, i, iface_ents[i] ) != MB_SUCCESS )
232  {
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;
236  MB_CHK_ERR( mbImpl->get_last_error( last_error ) );
237  if( last_error.empty() )
238  std::cerr << "(none)" << std::endl;
239  else
240  std::cerr << last_error << std::endl;
241  }
242  if( 0 != i ) iface_ents[4].merge( iface_ents[i] );
243  }
244  }
245 
246  // Report # iface entities
247  MB_CHK_ERR( mbImpl->get_adjacencies( iface_ents[4], 0, false, iface_ents[5], Interface::UNION ) );
248 
249  int rank;
250  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
251 
252  if( print_results || iface_ents[0].size() != iface_ents[5].size() )
253  {
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;
258  }
259 
260  return moab::MB_SUCCESS;
261 }
262 
263 // Check first character for a '-'.
264 // Return true if one is found. False otherwise.
265 bool check_for_flag( const char* str )
266 {
267  if( '-' == str[0] )
268  return true;
269  else
270  return false;
271 }
272 
273 // get_file_options() function with added possibilities for mbcoupler_test.
274 ErrorCode get_file_options( int argc,
275  char** argv,
276  int nprocs,
277  int rank,
278  std::vector< std::string >& meshFiles,
279  Coupler::Method& method,
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,
289  bool& help,
290  double& epsilon )
291 {
292  // Initialize some of the outputs to null values indicating not present
293  // in the argument list.
294  gNormTag = "";
295  ssNormTag = "";
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" );
300  outFile = "";
301  writeOpts = ( nprocs > 1 ? "PARALLEL=WRITE_PART;CPUTIME" : "" );
302  dbgFile = "";
303  std::string defaultDbgFile = argv[0]; // The executable name will be the default debug output file.
304 
305  // These will indicate if we've gotten our required parameters at the end of parsing.
306  bool haveMeshes = false;
307  bool haveInterpTag = false;
308 
309  // Loop over the values in argv pulling out an parsing each one
310  int npos = 1;
311 
312  if( argc > 1 && argv[1] == std::string( "-h" ) )
313  {
314  help = true;
315  return MB_SUCCESS;
316  }
317 
318  while( npos < argc )
319  {
320  if( argv[npos] == std::string( "-meshes" ) )
321  {
322  // Parse out the mesh filenames
323  npos++;
324  int numFiles = 2;
325  meshFiles.resize( numFiles );
326  for( int i = 0; i < numFiles; i++ )
327  {
328  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
329  meshFiles[i] = argv[npos++];
330  else
331  {
332  std::cerr << " ERROR - missing correct number of mesh filenames" << std::endl;
333  return MB_FAILURE;
334  }
335  }
336 
337  haveMeshes = true;
338  }
339  else if( argv[npos] == std::string( "-itag" ) )
340  {
341  // Parse out the interpolation tag
342  npos++;
343  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
344  interpTag = argv[npos++];
345  else
346  {
347  std::cerr << " ERROR - missing <interp_tag>" << std::endl;
348  return MB_FAILURE;
349  }
350 
351  haveInterpTag = true;
352  }
353  else if( argv[npos] == std::string( "-meth" ) )
354  {
355  // Parse out the interpolation tag
356  npos++;
357  if( argv[npos][0] == '0' )
358  method = Coupler::CONSTANT;
359  else if( argv[npos][0] == '1' )
360  method = Coupler::LINEAR_FE;
361  else if( argv[npos][0] == '2' )
362  method = Coupler::QUADRATIC_FE;
363  else if( argv[npos][0] == '3' )
364  method = Coupler::SPECTRAL;
365  else if( argv[npos][0] == '4' )
366  method = Coupler::SPHERICAL;
367  else
368  {
369  std::cerr << " ERROR - unrecognized method number " << method << std::endl;
370  return MB_FAILURE;
371  }
372  npos++;
373  }
374  else if( argv[npos] == std::string( "-eps" ) )
375  {
376  // Parse out the tolerance
377  npos++;
378  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
379  epsilon = atof( argv[npos++] );
380  else
381  {
382  std::cerr << " ERROR - missing <epsilon>" << std::endl;
383  return MB_FAILURE;
384  }
385  }
386  else if( argv[npos] == std::string( "-gnorm" ) )
387  {
388  // Parse out the global normalization tag
389  npos++;
390  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
391  gNormTag = argv[npos++];
392  else
393  {
394  std::cerr << " ERROR - missing <gnorm_tag>" << std::endl;
395  return MB_FAILURE;
396  }
397  }
398  else if( argv[npos] == std::string( "-ssnorm" ) )
399  {
400  // Parse out the subset normalization tag and selection criteria
401  npos++;
402  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
403  ssNormTag = argv[npos++];
404  else
405  {
406  std::cerr << " ERROR - missing <ssnorm_tag>" << std::endl;
407  return MB_FAILURE;
408  }
409 
410  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
411  {
412  char* opts = argv[npos++];
413  char sep1[1] = { ';' };
414  char sep2[1] = { '=' };
415  bool end_vals_seen = false;
416  std::vector< char* > tmpTagOpts;
417 
418  // First get the options
419  for( char* i = strtok( opts, sep1 ); i; i = strtok( 0, sep1 ) )
420  tmpTagOpts.push_back( i );
421 
422  // Parse out the name and val or just name.
423  for( unsigned int j = 0; j < tmpTagOpts.size(); j++ )
424  {
425  char* e = strtok( tmpTagOpts[j], sep2 );
426  ssTagNames.push_back( e );
427  e = strtok( 0, sep2 );
428  if( e != NULL )
429  {
430  // We have a value
431  if( end_vals_seen )
432  {
433  // ERROR we should not have a value after none are seen
434  std::cerr << " ERROR - new value seen after end of values in "
435  "<ssnorm_selection>"
436  << std::endl;
437  return MB_FAILURE;
438  }
439  // Otherwise get the value string from e and convert it to an int
440  int* valp = new int;
441  *valp = atoi( e );
442  ssTagValues.push_back( (const char*)valp );
443  }
444  else
445  {
446  // Otherwise there is no '=' so push a null on the list
447  end_vals_seen = true;
448  ssTagValues.push_back( (const char*)0 );
449  }
450  }
451  }
452  else
453  {
454  std::cerr << " ERROR - missing <ssnorm_selection>" << std::endl;
455  return MB_FAILURE;
456  }
457  }
458  else if( argv[npos] == std::string( "-ropts" ) )
459  {
460  // Parse out the mesh file read options
461  npos++;
462  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
463  readOpts = argv[npos++];
464  else
465  {
466  std::cerr << " ERROR - missing <roptions>" << std::endl;
467  return MB_FAILURE;
468  }
469  }
470  else if( argv[npos] == std::string( "-outfile" ) )
471  {
472  // Parse out the output file name
473  npos++;
474  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
475  outFile = argv[npos++];
476  else
477  {
478  std::cerr << " ERROR - missing <out_file>" << std::endl;
479  return MB_FAILURE;
480  }
481  }
482  else if( argv[npos] == std::string( "-wopts" ) )
483  {
484  // Parse out the output file write options
485  npos++;
486  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
487  writeOpts = argv[npos++];
488  else
489  {
490  std::cerr << " ERROR - missing <woptions>" << std::endl;
491  return MB_FAILURE;
492  }
493  }
494  else if( argv[npos] == std::string( "-dbgout" ) )
495  {
496  // Parse out the debug output file name.
497  // If no name then use the default.
498  npos++;
499  if( ( npos < argc ) && ( !check_for_flag( argv[npos] ) ) )
500  dbgFile = argv[npos++];
501  else
502  dbgFile = defaultDbgFile;
503  }
504  else
505  {
506  // Unrecognized parameter. Skip it and move along.
507  std::cerr << " ERROR - Unrecognized parameter:" << argv[npos] << std::endl;
508  std::cerr << " Skipping..." << std::endl;
509  npos++;
510  }
511  }
512 
513  if( !haveMeshes )
514  {
515  meshFiles.resize( 2 );
516  meshFiles[0] = std::string( TestDir + "unittest/64bricks_1khex.h5m" );
517  meshFiles[1] = std::string( TestDir + "unittest/64bricks_12ktet.h5m" );
518  if( 0 == rank )
519  std::cout << "Mesh files not entered; using default files " << meshFiles[0] << " and " << meshFiles[1]
520  << std::endl;
521  }
522 
523  if( !haveInterpTag )
524  {
525  interpTag = "vertex_field";
526  std::cout << "Interpolation field name not given, using default of " << interpTag << std::endl;
527  }
528 
529 #ifdef MOAB_HAVE_HDF5
530  if( 1 == argc )
531  {
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;
536  }
537 #endif
538 
539  return MB_SUCCESS;
540 }
541 
542 // End get_file_options()
543 
544 ErrorCode test_interpolation( Interface* mbImpl,
545  Coupler::Method method,
546  std::string& interpTag,
547  std::string& gNormTag,
548  std::string& ssNormTag,
549  std::vector< const char* >& ssTagNames,
550  std::vector< const char* >& ssTagValues,
551  EntityHandle* roots,
552  std::vector< ParallelComm* >& pcs,
553  double& instant_time,
554  double& pointloc_time,
555  double& interp_time,
556  double& gnorm_time,
557  double& ssnorm_time,
558  double& toler )
559 {
560  assert( method >= Coupler::CONSTANT && method <= Coupler::SPHERICAL );
561 
562  // Source is 1st mesh, target is 2nd
563  Range src_elems, targ_elems, targ_verts;
564  MB_CHK_ERR( pcs[0]->get_part_entities( src_elems, 3 ) );
565 
566  double start_time = MPI_Wtime();
567 
568  // Instantiate a coupler, which does not initialize the tree yet
569  Coupler mbc( mbImpl, pcs[0], src_elems, 0, false ); // do not initialize tree yet
570  if( Coupler::SPHERICAL == method ) mbc.set_spherical();
571  mbc.initialize_tree(); // it is expensive, but do something different for spherical
572 
573  // Initialize spectral elements, if they exist
574  bool specSou = false, specTar = false;
575  MB_CHK_ERR( mbc.initialize_spectral_elements( roots[0], roots[1], specSou, specTar ) );
576 
577  instant_time = MPI_Wtime();
578 
579  // Get points from the target mesh to interpolate
580  // We have to treat differently the case when the target is a spectral mesh
581  // In that case, the points of interest are the GL points, not the vertex nodes
582  std::vector< double > vpos; // This will have the positions we are interested in
583  int numPointsOfInterest = 0;
584  if( !specTar )
585  {
586  // Usual case
587  Range tmp_verts;
588 
589  // First get all vertices adj to partition entities in target mesh
590  MB_CHK_ERR( pcs[1]->get_part_entities( targ_elems, 3 ) );
591  if( Coupler::SPHERICAL == method ) MB_CHK_ERR( pcs[1]->get_part_entities( targ_elems, 2 ) ); // get the polygons/quads on a sphere.
592  if( Coupler::CONSTANT == method )
593  targ_verts = targ_elems;
594  else
595  MB_CHK_ERR( mbImpl->get_adjacencies( targ_elems, 0, false, targ_verts, Interface::UNION ) );
596 
597  // Then get non-owned verts and subtract
598  MB_CHK_ERR( pcs[1]->get_pstatus_entities( 0, PSTATUS_NOT_OWNED, tmp_verts ) );
599  targ_verts = subtract( targ_verts, tmp_verts );
600 
601  // get position of these entities; these are the target points
602  numPointsOfInterest = (int)targ_verts.size();
603  vpos.resize( 3 * targ_verts.size() );
604  MB_CHK_ERR( mbImpl->get_coords( targ_verts, &vpos[0] ) );
605 
606  // Locate those points in the source mesh
607  std::cout << "rank " << pcs[0]->proc_config().proc_rank() << " points of interest: " << numPointsOfInterest
608  << "\n";
609  MB_CHK_ERR( mbc.locate_points( &vpos[0], numPointsOfInterest, 0, toler ) );
610  }
611  else
612  {
613  // In this case, the target mesh is spectral, we want values
614  // interpolated on the GL positions; for each element, get the GL points, and construct
615  // CartVect!!!
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
619  << "\n";
620  }
621 
622  pointloc_time = MPI_Wtime();
623 
624  // Now interpolate tag onto target points
625  std::vector< double > field( numPointsOfInterest );
626 
627  MB_CHK_ERR( mbc.interpolate( method, interpTag, &field[0] ) );
628 
629  interp_time = MPI_Wtime();
630 
631  // Do global normalization if specified
632  if( !gNormTag.empty() )
633  {
634  // Normalize the source mesh
635  MB_CHK_ERR( mbc.normalize_mesh( roots[0], gNormTag.c_str(), Coupler::VOLUME, 4 ) );
636 
637  // Normalize the target mesh
638  MB_CHK_ERR( mbc.normalize_mesh( roots[1], gNormTag.c_str(), Coupler::VOLUME, 4 ) );
639  }
640 
641  gnorm_time = MPI_Wtime();
642 
643  // Do subset normalization if specified
644 
645  if( !ssNormTag.empty() )
646  {
647  MB_CHK_ERR( mbc.normalize_subset( roots[0], ssNormTag.c_str(), &ssTagNames[0], ssTagNames.size(), &ssTagValues[0],
648  Coupler::VOLUME, 4 ) );
649 
650  MB_CHK_ERR( mbc.normalize_subset( roots[1], ssNormTag.c_str(), &ssTagNames[0], ssTagNames.size(), &ssTagValues[0],
651  Coupler::VOLUME, 4 ) );
652  }
653 
654  ssnorm_time = MPI_Wtime();
655 
656  ssnorm_time -= gnorm_time;
657  gnorm_time -= interp_time;
658  interp_time -= pointloc_time;
659  pointloc_time -= instant_time;
660  instant_time -= start_time;
661 
662  // Set field values as tag on target vertices
663  if( specSou )
664  {
665  // Create a new tag for the values on the target
666  Tag tag;
667  std::string newtag = interpTag + "_TAR";
668  MB_CHK_ERR( mbImpl->tag_get_handle( newtag.c_str(), 1, MB_TYPE_DOUBLE, tag, MB_TAG_CREAT | MB_TAG_DENSE ) );
669  MB_CHK_ERR( mbImpl->tag_set_data( tag, targ_verts, &field[0] ) );
670  }
671  else
672  {
673  if( !specTar )
674  {
675  // Use original tag
676  Tag tag;
677  MB_CHK_ERR( mbImpl->tag_get_handle( interpTag.c_str(), 1, MB_TYPE_DOUBLE, tag ) );
678  MB_CHK_ERR( mbImpl->tag_set_data( tag, targ_verts, &field[0] ) );
679  }
680  else
681  {
682  // We have the field values computed at each GL points, on each element
683  // in the target mesh
684  // We need to create a new tag, on elements, of size _ntot, to hold
685  // all those values.
686  // So it turns out we need ntot. Maybe we can compute it from the
687  // number of values computed, divided by number of elements
688  int ntot = numPointsOfInterest / targ_elems.size();
689  Tag tag;
690  std::string newtag = interpTag + "_TAR";
691  MB_CHK_ERR( mbImpl->tag_get_handle( newtag.c_str(), ntot, MB_TYPE_DOUBLE, tag, MB_TAG_CREAT | MB_TAG_DENSE ) );
692  MB_CHK_ERR( mbImpl->tag_set_data( tag, targ_elems, &field[0] ) );
693  }
694  }
695 
696  // Done
697  return MB_SUCCESS;
698 }
699 
700 #else
701 
702 int main( int /*argc*/, char** /*argv*/ )
703 {
704  print_usage();
705  return 0;
706 }
707 
708 #endif