Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
mbtempest.cpp File Reference
#include <iostream>
#include <iomanip>
#include <cstdlib>
#include <vector>
#include <string>
#include <sstream>
#include <cassert>
#include "moab/Core.hpp"
#include "moab/IntxMesh/IntxUtils.hpp"
#include "moab/Remapping/TempestRemapper.hpp"
#include "moab/Remapping/TempestOnlineMap.hpp"
#include "moab/ProgOptions.hpp"
#include "moab/CpuTimer.hpp"
#include "DebugOutput.hpp"
+ Include dependency graph for mbtempest.cpp:

Go to the source code of this file.

Classes

struct  ToolContext
 

Functions

static moab::ErrorCode CreateTempestMesh (ToolContext &, moab::TempestRemapper &remapper, Mesh *)
 
double sample_slow_harmonic (double dLon, double dLat)
 
double sample_fast_harmonic (double dLon, double dLat)
 
double sample_constant (double dLon, double dLat)
 
double sample_stationary_vortex (double dLon, double dLat)
 
std::string get_file_read_options (ToolContext &ctx, std::string filename)
 
int main (int argc, char *argv[])
 

Function Documentation

◆ CreateTempestMesh()

static moab::ErrorCode CreateTempestMesh ( ToolContext ctx,
moab::TempestRemapper remapper,
Mesh *  tempest_mesh 
)
static

Definition at line 1043 of file mbtempest.cpp.

1044 {
1046  int err;
1047  moab::DebugOutput& outputFormatter = ctx.outputFormatter;
1048 
1049  if( !ctx.proc_id )
1050  {
1051  outputFormatter.printf( 0, "Creating TempestRemap Mesh object ...\n" );
1052  }
1053 
1055  {
1056  ctx.timer_push( "create Tempest OverlapMesh" );
1057  // For the overlap method, choose between: "fuzzy", "exact" or "mixed"
1058  err = GenerateOverlapMesh( ctx.inFilenames[0], ctx.inFilenames[1], *tempest_mesh, ctx.outFilename, "NetCDF4",
1059  "exact", true );
1060  if( err )
1061  rval = moab::MB_FAILURE;
1062  else
1063  ctx.meshes.push_back( tempest_mesh );
1064  ctx.timer_pop();
1065  }
1067  {
1068  // Load the meshes and validate
1069  ctx.meshsets.resize( 3 );
1070  ctx.meshes.resize( 3 );
1071  ctx.meshsets[0] = remapper.GetMeshSet( moab::Remapper::SourceMesh );
1072  ctx.meshsets[1] = remapper.GetMeshSet( moab::Remapper::TargetMesh );
1073  ctx.meshsets[2] = remapper.GetMeshSet( moab::Remapper::OverlapMesh );
1074 
1075  ctx.timer_push( "load MOAB Source mesh" );
1076  // First the source
1078  ctx.meshes[0] = remapper.GetMesh( moab::Remapper::SourceMesh );
1079  ctx.timer_pop();
1080 
1081  ctx.timer_push( "load MOAB Target mesh" );
1082  // Next the target
1084  ctx.meshes[1] = remapper.GetMesh( moab::Remapper::TargetMesh );
1085  ctx.timer_pop();
1086 
1087  ctx.timer_push( "generate TempestRemap OverlapMesh" );
1088  // Now let us construct the overlap mesh, by calling TempestRemap interface directly
1089  // For the overlap method, choose between: "fuzzy", "exact" or "mixed"
1090  err = GenerateOverlapWithMeshes( *ctx.meshes[0], *ctx.meshes[1], *tempest_mesh, "" /*ctx.outFilename*/,
1091  "NetCDF4", "exact", false );
1092  ctx.timer_pop();
1093  if( err )
1094  rval = moab::MB_FAILURE;
1095  else
1096  {
1097  remapper.SetMesh( moab::Remapper::OverlapMesh, tempest_mesh );
1098  ctx.meshes[2] = remapper.GetMesh( moab::Remapper::OverlapMesh );
1099  // ctx.meshes.push_back(*tempest_mesh);
1100  }
1101  }
1103  {
1104  ctx.meshsets.resize( 3 );
1105  ctx.meshes.resize( 3 );
1106  ctx.meshsets[0] = remapper.GetMeshSet( moab::Remapper::SourceMesh );
1107  ctx.meshsets[1] = remapper.GetMeshSet( moab::Remapper::TargetMesh );
1108  ctx.meshsets[2] = remapper.GetMeshSet( moab::Remapper::OverlapMesh );
1109 
1110  const double radius_src = 1.0 /*2.0*acos(-1.0)*/;
1111  const double radius_dest = 1.0 /*2.0*acos(-1.0)*/;
1112 
1113  std::vector< int > smetadata, tmetadata;
1114  //const char* additional_read_opts = ( ctx.n_procs > 1 ? "NO_SET_CONTAINING_PARENTS;" : "" );
1115  std::string additional_read_opts_src = get_file_read_options( ctx, ctx.inFilenames[0] );
1116 
1117  ctx.timer_push( "load MOAB Source mesh" );
1118  // Load the source mesh and validate
1119  rval =
1120  remapper.LoadNativeMesh( ctx.inFilenames[0], ctx.meshsets[0], smetadata, additional_read_opts_src.c_str() );MB_CHK_ERR( rval );
1121  if( smetadata.size() ) remapper.SetMeshType( moab::Remapper::SourceMesh, smetadata );
1122  ctx.timer_pop();
1123 
1124  ctx.timer_push( "preprocess MOAB Source mesh" );
1125  // Rescale the radius of both to compute the intersection
1126  rval = moab::IntxUtils::ScaleToRadius( ctx.mbcore, ctx.meshsets[0], radius_src );MB_CHK_ERR( rval );
1127 
1128  ctx.timer_pop();
1129 
1130  // Load the target mesh and validate
1131  std::string addititional_read_opts_tgt = get_file_read_options( ctx, ctx.inFilenames[1] );
1132  // addititional_read_opts_tgt += "PARALLEL_GHOSTS=2.0.3;PARALLEL_THIN_GHOST_LAYER;";
1133 
1134  ctx.timer_push( "load MOAB Target mesh" );
1135  rval = remapper.LoadNativeMesh( ctx.inFilenames[1], ctx.meshsets[1], tmetadata,
1136  addititional_read_opts_tgt.c_str() );MB_CHK_ERR( rval );
1137  if( tmetadata.size() ) remapper.SetMeshType( moab::Remapper::TargetMesh, tmetadata );
1138  ctx.timer_pop();
1139 
1140  ctx.timer_push( "preprocess MOAB Target mesh" );
1141  rval = moab::IntxUtils::ScaleToRadius( ctx.mbcore, ctx.meshsets[1], radius_dest );MB_CHK_ERR( rval );
1142  ctx.timer_pop();
1143 
1144  if( ctx.computeWeights )
1145  {
1146  ctx.timer_push( "convert MOAB meshes to TempestRemap meshes in memory" );
1147  // convert MOAB representation to TempestRemap's Mesh for source
1148  rval = remapper.ConvertMeshToTempest( moab::Remapper::SourceMesh );MB_CHK_ERR( rval );
1149  ctx.meshes[0] = remapper.GetMesh( moab::Remapper::SourceMesh );
1150 
1151  // convert MOAB representation to TempestRemap's Mesh for target
1152  rval = remapper.ConvertMeshToTempest( moab::Remapper::TargetMesh );MB_CHK_ERR( rval );
1153  ctx.meshes[1] = remapper.GetMesh( moab::Remapper::TargetMesh );
1154  ctx.timer_pop();
1155  }
1156  }
1157  else if( ctx.meshType == moab::TempestRemapper::ICO )
1158  {
1159  ctx.timer_push( "generate ICO mesh with TempestRemap" );
1160  err = GenerateICOMesh( *tempest_mesh, ctx.blockSize, ctx.computeDual, ctx.outFilename, "NetCDF4" );
1161  ctx.timer_pop();
1162 
1163  if( err )
1164  rval = moab::MB_FAILURE;
1165  else
1166  ctx.meshes.push_back( tempest_mesh );
1167  }
1168  else if( ctx.meshType == moab::TempestRemapper::RLL )
1169  {
1170  ctx.timer_push( "generate RLL mesh with TempestRemap" );
1171  err = GenerateRLLMesh( *tempest_mesh, // Mesh& meshOut,
1172  ctx.blockSize * 2, ctx.blockSize, // int nLongitudes, int nLatitudes,
1173  0.0, 360.0, // double dLonBegin, double dLonEnd,
1174  -90.0, 90.0, // double dLatBegin, double dLatEnd,
1175  false, false, false, // bool fGlobalCap, bool fFlipLatLon, bool fForceGlobal,
1176  "" /*ctx.inFilename*/, "",
1177  "", // std::string strInputFile, std::string strInputFileLonName, std::string
1178  // strInputFileLatName,
1179  ctx.outFilename, "NetCDF4", // std::string strOutputFile, std::string strOutputFormat
1180  true // bool fVerbose
1181  );
1182  ctx.timer_pop();
1183 
1184  if( err )
1185  rval = moab::MB_FAILURE;
1186  else
1187  ctx.meshes.push_back( tempest_mesh );
1188  }
1189  else // default
1190  {
1191  ctx.timer_push( "generate CS mesh with TempestRemap" );
1192  err = GenerateCSMesh( *tempest_mesh, ctx.blockSize, ctx.outFilename, "NetCDF4" );
1193  ctx.timer_pop();
1194  if( err )
1195  rval = moab::MB_FAILURE;
1196  else
1197  ctx.meshes.push_back( tempest_mesh );
1198  }
1199 
1200  if( ctx.meshType != moab::TempestRemapper::OVERLAP_MOAB && !tempest_mesh )
1201  {
1202  std::cout << "Tempest Mesh is not a complete object; Quitting...";
1203  exit( -1 );
1204  }
1205 
1206  return rval;
1207 }

References ToolContext::blockSize, ToolContext::computeDual, ToolContext::computeWeights, moab::TempestRemapper::ConvertMeshToTempest(), moab::TempestRemapper::DEFAULT, ErrorCode, get_file_read_options(), moab::TempestRemapper::GetMesh(), moab::TempestRemapper::GetMeshSet(), moab::TempestRemapper::ICO, ToolContext::inFilenames, moab::TempestRemapper::LoadMesh(), moab::Remapper::LoadNativeMesh(), MB_CHK_ERR, MB_SUCCESS, ToolContext::mbcore, ToolContext::meshes, ToolContext::meshsets, ToolContext::meshType, ToolContext::outFilename, ToolContext::outputFormatter, moab::TempestRemapper::OVERLAP_FILES, moab::TempestRemapper::OVERLAP_MEMORY, moab::TempestRemapper::OVERLAP_MOAB, moab::Remapper::OverlapMesh, moab::DebugOutput::printf(), ToolContext::proc_id, moab::TempestRemapper::RLL, moab::IntxUtils::ScaleToRadius(), moab::TempestRemapper::SetMesh(), moab::TempestRemapper::SetMeshType(), moab::Remapper::SourceMesh, moab::Remapper::TargetMesh, ToolContext::timer_pop(), and ToolContext::timer_push().

Referenced by main().

◆ get_file_read_options()

std::string get_file_read_options ( ToolContext ctx,
std::string  filename 
)

Definition at line 417 of file mbtempest.cpp.

418 {
419  // if running in serial, blank option may suffice in general
420  std::string opts = "";
421  if( ctx.n_procs > 1 )
422  {
423  size_t lastindex = filename.find_last_of( "." );
424  std::string extension = filename.substr( lastindex + 1, filename.size() );
425  if( extension == "h5m" ) return "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_SHARED_ENTS;";
426  // "PARALLEL_GHOSTS=2.0.3;PARALLEL_THIN_GHOST_LAYER;";
427  else if( extension == "nc" )
428  {
429  // set default set of options
430  opts = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;";
431 
432  if( ctx.proc_id == 0 )
433  {
434  {
435  NcFile ncFile( filename.c_str(), NcFile::ReadOnly );
436  if( !ncFile.is_valid() )
437  _EXCEPTION1( "Unable to open grid file \"%s\" for reading", filename.c_str() );
438 
439  // Check for dimension names "grid_size", "grid_rank" and "grid_corners"
440  int iSCRIPFormat = 0, iMPASFormat = 0, iESMFFormat = 0;
441  for( int i = 0; i < ncFile.num_dims(); i++ )
442  {
443  std::string strDimName = ncFile.get_dim( i )->name();
444  if( strDimName == "grid_size" || strDimName == "grid_corners" || strDimName == "grid_rank" )
445  iSCRIPFormat++;
446  if( strDimName == "nodeCount" || strDimName == "elementCount" ||
447  strDimName == "maxNodePElement" )
448  iESMFFormat++;
449  if( strDimName == "nCells" || strDimName == "nEdges" || strDimName == "nVertices" ||
450  strDimName == "vertexDegree" )
451  iMPASFormat++;
452  }
453  ncFile.close();
454 
455  if( iESMFFormat == 3 ) // Input from a NetCDF ESMF file
456  {
457  opts = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;PARALLEL_RESOLVE_SHARED_ENTS;VARIABLE=;";
458  }
459  if( iSCRIPFormat == 3 ) // Input from a NetCDF SCRIP file
460  {
461  opts = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;";
462  }
463  if( iMPASFormat == 4 ) // Input from a NetCDF SCRIP file
464  {
465  opts = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;"
466  "PARALLEL_RESOLVE_SHARED_ENTS;NO_EDGES;NO_MIXED_ELEMENTS;VARIABLE=;";
467  }
468  }
469  }
470 
471 #ifdef MOAB_HAVE_MPI
472  int line_size = opts.size();
473  MPI_Bcast( &line_size, 1, MPI_INT, 0, MPI_COMM_WORLD );
474  if( ctx.proc_id != 0 ) opts.resize( line_size );
475  MPI_Bcast( const_cast< char* >( opts.data() ), line_size, MPI_CHAR, 0, MPI_COMM_WORLD );
476 #endif
477  }
478  else
479  return "PARALLEL=BCAST_DELETE;PARTITION=TRIVIAL;PARALLEL_RESOLVE_SHARED_ENTS;";
480  }
481  return opts;
482 }

References ToolContext::n_procs, ncFile, and ToolContext::proc_id.

Referenced by CreateTempestMesh().

◆ main()

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

Definition at line 484 of file mbtempest.cpp.

485 {
486  moab::ErrorCode rval;
487  NcError error( NcError::verbose_nonfatal );
488  std::stringstream sstr;
489  std::string historyStr;
490 
491  int proc_id = 0, nprocs = 1;
492 #ifdef MOAB_HAVE_MPI
493  MPI_Init( &argc, &argv );
494  MPI_Comm_rank( MPI_COMM_WORLD, &proc_id );
495  MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
496 #endif
497 
498  moab::Core* mbCore = new( std::nothrow ) moab::Core;
499 
500  if( nullptr == mbCore )
501  {
502  return 1;
503  }
504 
505  // Build the history string
506  for( int ia = 0; ia < argc; ++ia )
507  historyStr += std::string( argv[ia] ) + " ";
508 
509  ToolContext* runCtx;
510 #ifdef MOAB_HAVE_MPI
511  moab::ParallelComm* pcomm = new moab::ParallelComm( mbCore, MPI_COMM_WORLD, 0 );
512 
513  runCtx = new ToolContext( mbCore, pcomm );
514  const char* writeOptions = ( nprocs > 1 ? "PARALLEL=WRITE_PART" : "" );
515 #else
516  runCtx = new ToolContext( mbCore );
517  const char* writeOptions = "";
518 #endif
519  runCtx->ParseCLOptions( argc, argv );
520 
521  const double radius_src = 1.0 /*2.0*acos(-1.0)*/;
522  const double radius_dest = 1.0 /*2.0*acos(-1.0)*/;
523 
524  moab::DebugOutput& outputFormatter = runCtx->outputFormatter;
525 
526 #ifdef MOAB_HAVE_MPI
527  moab::TempestRemapper remapper( mbCore, pcomm );
528 #else
529  moab::TempestRemapper remapper( mbCore );
530 #endif
531  remapper.meshValidate = true;
532  remapper.constructEdgeMap = true;
533  remapper.initialize();
534 
535  // Default area_method = lHuiller; Options: Girard, lHuiller, GaussQuadrature (if TR is available)
537 
538  Mesh* tempest_mesh = new Mesh();
539  rval = CreateTempestMesh( *runCtx, remapper, tempest_mesh );MB_CHK_ERR( rval );
540 
542  {
543  // Compute intersections with MOAB
544  // For the overlap method, choose between: "fuzzy", "exact" or "mixed"
545  assert( runCtx->meshes.size() == 3 );
546 
547 #ifdef MOAB_HAVE_MPI
548  rval = pcomm->check_all_shared_handles();MB_CHK_ERR( rval );
549 #endif
550 
551  // Load the meshes and validate
552  rval = remapper.ConvertTempestMesh( moab::Remapper::SourceMesh );MB_CHK_ERR( rval );
553  rval = remapper.ConvertTempestMesh( moab::Remapper::TargetMesh );MB_CHK_ERR( rval );
554  rval = remapper.ConvertTempestMesh( moab::Remapper::OverlapMesh );MB_CHK_ERR( rval );
555  if( !runCtx->skip_io )
556  {
557  rval = mbCore->write_mesh( "tempest_intersection.h5m", &runCtx->meshsets[2], 1 );MB_CHK_ERR( rval );
558  }
559 
560  // print verbosely about the problem setting
561  size_t velist[6], gvelist[6];
562  {
563  moab::Range rintxverts, rintxelems;
564  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 0, rintxverts );MB_CHK_ERR( rval );
565  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 2, rintxelems );MB_CHK_ERR( rval );
566  velist[0] = rintxverts.size();
567  velist[1] = rintxelems.size();
568 
569  moab::Range bintxverts, bintxelems;
570  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 0, bintxverts );MB_CHK_ERR( rval );
571  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 2, bintxelems );MB_CHK_ERR( rval );
572  velist[2] = bintxverts.size();
573  velist[3] = bintxelems.size();
574  }
575 
576  moab::EntityHandle intxset; // == remapper.GetMeshSet(moab::Remapper::OverlapMesh);
577 
578  // Compute intersections with MOAB
579  {
580  // Create the intersection on the sphere object
581  runCtx->timer_push( "setup the intersector" );
582 
583  moab::Intx2MeshOnSphere* mbintx = new moab::Intx2MeshOnSphere( mbCore );
584  mbintx->set_error_tolerance( runCtx->epsrel );
585  mbintx->set_box_error( runCtx->boxeps );
586  mbintx->set_radius_source_mesh( radius_src );
587  mbintx->set_radius_destination_mesh( radius_dest );
588 #ifdef MOAB_HAVE_MPI
589  mbintx->set_parallel_comm( pcomm );
590 #endif
591  rval = mbintx->FindMaxEdges( runCtx->meshsets[0], runCtx->meshsets[1] );MB_CHK_ERR( rval );
592 
593 #ifdef MOAB_HAVE_MPI
594  moab::Range local_verts;
595  rval = mbintx->build_processor_euler_boxes( runCtx->meshsets[1], local_verts );MB_CHK_ERR( rval );
596 
597  runCtx->timer_pop();
598 
599  moab::EntityHandle covering_set;
600  runCtx->timer_push( "communicate the mesh" );
601  // TODO:: needs clarification
602  // we compute just intersection here, no need for extra ghost layers anyway
603  // ghost layers are needed in coverage for bilinear map, which does not actually need intersection
604  // this will be fixed in the future, bilinear map needs just coverage, not intersection
605  // so I am not passing the ghost layer here, even though there is an option in runCtx for a ghost layer
606  rval = mbintx->construct_covering_set( runCtx->meshsets[0], covering_set );MB_CHK_ERR( rval ); // lots of communication if mesh is distributed very differently
607  runCtx->timer_pop();
608 
609  // print verbosely about the problem setting
610  {
611  moab::Range cintxverts, cintxelems;
612  rval = mbCore->get_entities_by_dimension( covering_set, 0, cintxverts );MB_CHK_ERR( rval );
613  rval = mbCore->get_entities_by_dimension( covering_set, 2, cintxelems );MB_CHK_ERR( rval );
614  velist[4] = cintxverts.size();
615  velist[5] = cintxelems.size();
616  }
617 
618  MPI_Reduce( velist, gvelist, 6, MPI_UINT64_T, MPI_SUM, 0, MPI_COMM_WORLD );
619 
620 #else
621  moab::EntityHandle covering_set = runCtx->meshsets[0];
622  for( int i = 0; i < 6; i++ )
623  gvelist[i] = velist[i];
624 #endif
625 
626  if( !proc_id )
627  {
628  outputFormatter.printf( 0, "The source set contains %lu vertices and %lu elements \n", gvelist[0],
629  gvelist[0] );
630  outputFormatter.printf( 0, "The covering set contains %lu vertices and %lu elements \n", gvelist[2],
631  gvelist[2] );
632  outputFormatter.printf( 0, "The target set contains %lu vertices and %lu elements \n", gvelist[1],
633  gvelist[1] );
634  }
635 
636  // Now let's invoke the MOAB intersection algorithm in parallel with a
637  // source and target mesh set representing two different decompositions
638  runCtx->timer_push( "compute intersections with MOAB" );
639  rval = mbCore->create_meshset( moab::MESHSET_SET, intxset );MB_CHK_SET_ERR( rval, "Can't create new set" );
640  rval = mbintx->intersect_meshes( covering_set, runCtx->meshsets[1], intxset );MB_CHK_SET_ERR( rval, "Can't compute the intersection of meshes on the sphere" );
641  runCtx->timer_pop();
642 
643  // free the memory
644  delete mbintx;
645  }
646 
647  {
648  moab::Range intxelems, intxverts;
649  rval = mbCore->get_entities_by_dimension( intxset, 2, intxelems );MB_CHK_ERR( rval );
650  rval = mbCore->get_entities_by_dimension( intxset, 0, intxverts, true );MB_CHK_ERR( rval );
651  outputFormatter.printf( 0, "The intersection set contains %lu elements and %lu vertices \n",
652  intxelems.size(), intxverts.size() );
653 
654  moab::IntxAreaUtils areaAdaptorHuiller( moab::IntxAreaUtils::lHuiller ); // lHuiller, GaussQuadrature
655  double initial_sarea =
656  areaAdaptorHuiller.area_on_sphere( mbCore, runCtx->meshsets[0],
657  radius_src ); // use the target to compute the initial area
658  double initial_tarea =
659  areaAdaptorHuiller.area_on_sphere( mbCore, runCtx->meshsets[1],
660  radius_dest ); // use the target to compute the initial area
661  double intx_area = areaAdaptorHuiller.area_on_sphere( mbCore, intxset, radius_src );
662 
663  outputFormatter.printf( 0, "mesh areas: source = %12.10f, target = %12.10f, intersection = %12.10f \n",
664  initial_sarea, initial_tarea, intx_area );
665  outputFormatter.printf( 0, "relative error w.r.t source = %12.10e, target = %12.10e \n",
666  fabs( intx_area - initial_sarea ) / initial_sarea,
667  fabs( intx_area - initial_tarea ) / initial_tarea );
668  }
669 
670  // Write out our computed intersection file
671  if( !runCtx->skip_io )
672  {
673  rval = mbCore->write_mesh( "moab_intersection.h5m", &intxset, 1 );MB_CHK_ERR( rval );
674  }
675 
676  if( runCtx->computeWeights )
677  {
678  runCtx->timer_push( "compute weights with the Tempest meshes" );
679  // Call to generate an offline map with the tempest meshes
680  OfflineMap weightMap;
681  int err = GenerateOfflineMapWithMeshes( *runCtx->meshes[0], *runCtx->meshes[1], *runCtx->meshes[2],
682  runCtx->disc_methods[0], // std::string strInputType
683  runCtx->disc_methods[1], // std::string strOutputType,
684  runCtx->mapOptions, weightMap );
685  runCtx->timer_pop();
686 
687  std::map< std::string, std::string > mapAttributes;
688  if( err )
689  {
690  rval = moab::MB_FAILURE;
691  }
692  else
693  {
694  if( !runCtx->skip_io ) weightMap.Write( "outWeights.nc", mapAttributes );
695  }
696  }
697  }
698  else if( runCtx->meshType == moab::TempestRemapper::OVERLAP_MOAB )
699  {
700  // Usage: mpiexec -n 2 tools/mbtempest -t 5 -l mycs_2.h5m -l myico_2.h5m -f myoverlap_2.h5m
701 #ifdef MOAB_HAVE_MPI
702  rval = pcomm->check_all_shared_handles();MB_CHK_ERR( rval );
703 #endif
704 
705  // print verbosely about the problem setting
706  size_t velist[4] = {}, gvelist[4] = {};
707  {
708  moab::Range srcverts, srcelems;
709  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 0, srcverts );MB_CHK_ERR( rval );
710  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 2, srcelems );MB_CHK_ERR( rval );
711  rval = moab::IntxUtils::fix_degenerate_quads( mbCore, runCtx->meshsets[0] );MB_CHK_ERR( rval );
712  if( runCtx->enforceConvexity )
713  {
714  rval = moab::IntxUtils::enforce_convexity( mbCore, runCtx->meshsets[0], proc_id );MB_CHK_ERR( rval );
715  }
716  rval = areaAdaptor.positive_orientation( mbCore, runCtx->meshsets[0], radius_src );MB_CHK_ERR( rval );
717  // if( !proc_id )
718  // outputFormatter.printf( 0, "The source set contains %lu vertices and %lu elements \n", srcverts.size(),
719  // srcelems.size() );
720  velist[0] = srcverts.size();
721  velist[1] = srcelems.size();
722 
723  moab::Range tgtverts, tgtelems;
724  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 0, tgtverts );MB_CHK_ERR( rval );
725  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 2, tgtelems );MB_CHK_ERR( rval );
726  rval = moab::IntxUtils::fix_degenerate_quads( mbCore, runCtx->meshsets[1] );MB_CHK_ERR( rval );
727  if( runCtx->enforceConvexity )
728  {
729  rval = moab::IntxUtils::enforce_convexity( mbCore, runCtx->meshsets[1], proc_id );MB_CHK_ERR( rval );
730  }
731  rval = areaAdaptor.positive_orientation( mbCore, runCtx->meshsets[1], radius_dest );MB_CHK_ERR( rval );
732  // if( !proc_id )
733  // outputFormatter.printf( 0, "The target set contains %lu vertices and %lu elements \n", tgtverts.size(),
734  // tgtelems.size() );
735  velist[2] = tgtverts.size();
736  velist[3] = tgtelems.size();
737  }
738  //rval = mbCore->write_file( "source_mesh.h5m", nullptr, writeOptions, &runCtx->meshsets[0], 1 );MB_CHK_ERR( rval );
739  //rval = mbCore->write_file( "target_mesh.h5m", nullptr, writeOptions, &runCtx->meshsets[1], 1 );MB_CHK_ERR( rval );
740 
741  // if( runCtx->nlayers && nprocs > 1 )
742  // {
743  // remapper.ResetMeshSet( moab::Remapper::SourceMesh, runCtx->meshsets[3] );
744  // runCtx->meshes[0] = remapper.GetMesh( moab::Remapper::SourceMesh ); // ?
745  // }
746 
747  // First compute the covering set such that the target elements are fully covered by the
748  // local source grid
749  runCtx->timer_push( "construct covering set for intersection" );
750  // if ghosting, no gnomonic
751  if (runCtx->nlayers >=1)
752  runCtx->useGnomonicProjection = false;
753  rval = remapper.ConstructCoveringSet( runCtx->epsrel, 1.0, 1.0, runCtx->boxeps, runCtx->rrmGrids,
754  runCtx->useGnomonicProjection, runCtx->nlayers );MB_CHK_ERR( rval );
755  runCtx->timer_pop();
756 
757 #ifdef MOAB_HAVE_MPI
758  MPI_Reduce( velist, gvelist, 4, MPI_UINT64_T, MPI_SUM, 0, MPI_COMM_WORLD );
759 #else
760  for( int i = 0; i < 4; i++ )
761  gvelist[i] = velist[i];
762 #endif
763 
764  if( !proc_id && runCtx->print_diagnostics )
765  {
766  outputFormatter.printf( 0, "The source set contains %lu vertices and %lu elements \n", gvelist[0],
767  gvelist[1] );
768  outputFormatter.printf( 0, "The target set contains %lu vertices and %lu elements \n", gvelist[2],
769  gvelist[3] );
770  }
771 
772  // Compute intersections with MOAB with either the Kd-tree or the advancing front algorithm
773  runCtx->timer_push( "setup and compute mesh intersections" );
774  rval = remapper.ComputeOverlapMesh( runCtx->kdtreeSearch, false );MB_CHK_ERR( rval );
775  runCtx->timer_pop();
776 
777  // print some diagnostic checks to see if the overlap grid resolved the input meshes
778  // correctly
779  double dTotalOverlapArea = 0.0;
780  if( runCtx->print_diagnostics )
781  {
782  moab::IntxAreaUtils areaAdaptorHuiller(
783  moab::IntxAreaUtils::GaussQuadrature ); // lHuiller, GaussQuadrature
784  double local_areas[3],
785  global_areas[3]; // Array for Initial area, and through Method 1 and Method 2
786  // local_areas[0] = area_on_sphere_lHuiller ( mbCore, runCtx->meshsets[1], radius_src );
787  local_areas[0] = areaAdaptorHuiller.area_on_sphere( mbCore, runCtx->meshsets[0], radius_src );
788  local_areas[1] = areaAdaptorHuiller.area_on_sphere( mbCore, runCtx->meshsets[1], radius_dest );
789  local_areas[2] = areaAdaptorHuiller.area_on_sphere( mbCore, runCtx->meshsets[2], radius_src );
790 
791 #ifdef MOAB_HAVE_MPI
792  MPI_Allreduce( &local_areas[0], &global_areas[0], 3, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD );
793 #else
794  global_areas[0] = local_areas[0];
795  global_areas[1] = local_areas[1];
796  global_areas[2] = local_areas[2];
797 #endif
798  if( !proc_id )
799  {
800  outputFormatter.printf( 0,
801  "initial area: source mesh = %12.14f, target mesh = "
802  "%12.14f, overlap mesh = %12.14f\n",
803  global_areas[0], global_areas[1], global_areas[2] );
804  // outputFormatter.printf ( 0, " area with l'Huiller: %12.14f with Girard:
805  // %12.14f\n", global_areas[2], global_areas[3] ); outputFormatter.printf ( 0, "
806  // relative difference areas = %12.10e\n", fabs ( global_areas[2] - global_areas[3]
807  // ) / global_areas[2] );
808  outputFormatter.printf( 0, "relative error w.r.t source = %12.14e, and target = %12.14e\n",
809  fabs( global_areas[0] - global_areas[2] ) / global_areas[0],
810  fabs( global_areas[1] - global_areas[2] ) / global_areas[1] );
811  }
812  dTotalOverlapArea = global_areas[2];
813  }
814 
815  if( runCtx->intxFilename.size() )
816  {
817  moab::EntityHandle writableOverlapSet;
818  rval = mbCore->create_meshset( moab::MESHSET_SET, writableOverlapSet );MB_CHK_SET_ERR( rval, "Can't create new set" );
819  moab::EntityHandle meshOverlapSet = remapper.GetMeshSet( moab::Remapper::OverlapMesh );
820  moab::Range ovEnts;
821  rval = mbCore->get_entities_by_dimension( meshOverlapSet, 2, ovEnts );MB_CHK_SET_ERR( rval, "Can't create new set" );
822  rval = mbCore->get_entities_by_dimension( meshOverlapSet, 0, ovEnts );MB_CHK_SET_ERR( rval, "Can't create new set" );
823 
824 #ifdef MOAB_HAVE_MPI
825  // Do not remove ghosted entities if we still haven't computed weights
826  // Remove ghosted entities from overlap set before writing the new mesh set to file
827  if( nprocs > 1 )
828  {
829  moab::Range ghostedEnts;
830  rval = remapper.GetOverlapAugmentedEntities( ghostedEnts );MB_CHK_ERR( rval );
831  ovEnts = moab::subtract( ovEnts, ghostedEnts );
832 #ifdef MOAB_DBG
833  if( !runCtx->skip_io )
834  {
835  std::stringstream filename;
836  filename << "aug_overlap" << runCtx->pcomm->rank() << ".h5m";
837  rval = runCtx->mbcore->write_file( filename.str().c_str(), 0, 0, &meshOverlapSet, 1 );MB_CHK_ERR( rval );
838  }
839 #endif
840  }
841 #endif
842  rval = mbCore->add_entities( writableOverlapSet, ovEnts );MB_CHK_SET_ERR( rval, "adding local intx cells failed" );
843 
844 #ifdef MOAB_HAVE_MPI
845 #ifdef MOAB_DBG
846  if( nprocs > 1 && !runCtx->skip_io )
847  {
848  std::stringstream filename;
849  filename << "writable_intx_" << runCtx->pcomm->rank() << ".h5m";
850  rval = runCtx->mbcore->write_file( filename.str().c_str(), 0, 0, &writableOverlapSet, 1 );MB_CHK_ERR( rval );
851  }
852 #endif
853 #endif
854 
855  size_t lastindex = runCtx->intxFilename.find_last_of( "." );
856  sstr.str( "" );
857  sstr << runCtx->intxFilename.substr( 0, lastindex ) << ".h5m";
858  if( !runCtx->proc_id )
859  std::cout << "Writing out the MOAB intersection mesh file to " << sstr.str() << std::endl;
860 
861  // Write out our computed intersection file
862  if( !runCtx->skip_io )
863  {
864  rval = mbCore->write_file( sstr.str().c_str(), nullptr, writeOptions, &writableOverlapSet, 1 );MB_CHK_ERR( rval );
865  }
866  }
867 
868  if( runCtx->computeWeights )
869  {
870  runCtx->meshes[2] = remapper.GetMesh( moab::Remapper::OverlapMesh );
871  if( !runCtx->proc_id ) std::cout << std::endl;
872 
873  runCtx->timer_push( "setup computation of weights" );
874  // Call to generate the remapping weights with the tempest meshes
875  moab::TempestOnlineMap* weightMap = new moab::TempestOnlineMap( &remapper );
876  runCtx->timer_pop();
877 
878  runCtx->timer_push( "compute weights with TempestRemap" );
879  rval = weightMap->GenerateRemappingWeights(
880  runCtx->disc_methods[0], // std::string strInputType
881  runCtx->disc_methods[1], // std::string strOutputType,
882  runCtx->mapOptions, // const GenerateOfflineMapAlgorithmOptions& options
883  runCtx->doftag_names[0], // const std::string& source_tag_name
884  runCtx->doftag_names[1] // const std::string& target_tag_name
885  );MB_CHK_ERR( rval );
886  runCtx->timer_pop();
887 
888  weightMap->PrintMapStatistics();
889 
890  // Invoke the CheckMap routine on the TempestRemap serial interface directly, if running
891  // on a single process
892  if( runCtx->fCheck )
893  {
894  const double dNormalTolerance = 1.0E-8;
895  const double dStrictTolerance = 1.0E-12;
896  weightMap->CheckMap( runCtx->fCheck, runCtx->fCheck, runCtx->fCheck && ( runCtx->ensureMonotonicity ),
897  dNormalTolerance, dStrictTolerance, dTotalOverlapArea );
898  }
899 
900  if( runCtx->outFilename.size() )
901  {
902  std::map< std::string, std::string > attrMap;
903  attrMap["MOABversion"] = std::string( MOAB_VERSION );
904  attrMap["Title"] = "MOAB-TempestRemap (mbtempest) Offline Regridding Weight Generator";
905  attrMap["normalization"] = "ovarea";
906  attrMap["remap_options"] = runCtx->mapOptions.strMethod;
907  attrMap["domain_a"] = runCtx->inFilenames[0];
908  attrMap["domain_b"] = runCtx->inFilenames[1];
909  if( runCtx->intxFilename.size() ) attrMap["domain_aUb"] = runCtx->intxFilename;
910  attrMap["map_aPb"] = runCtx->outFilename;
911  attrMap["methodorder_a"] = runCtx->disc_methods[0] + ":" + std::to_string( runCtx->disc_orders[0] ) +
912  ":" + std::string( runCtx->doftag_names[0] );
913  attrMap["concave_a"] = runCtx->mapOptions.fSourceConcave ? "true" : "false";
914  attrMap["methodorder_b"] = runCtx->disc_methods[1] + ":" + std::to_string( runCtx->disc_orders[1] ) +
915  ":" + std::string( runCtx->doftag_names[1] );
916  attrMap["concave_b"] = runCtx->mapOptions.fTargetConcave ? "true" : "false";
917  attrMap["bubble"] = runCtx->mapOptions.fNoBubble ? "false" : "true";
918  attrMap["history"] = historyStr;
919 
920  // Write the map file to disk in parallel using either HDF5 or SCRIP interface
921  // in extra case; maybe need a better solution, just create it with the right meshset
922  // from the beginning;
923  if( !runCtx->skip_io )
924  {
925  rval = weightMap->WriteParallelMap( runCtx->outFilename.c_str(), attrMap );MB_CHK_ERR( rval );
926  }
927  }
928 
929  if( runCtx->verifyWeights )
930  {
931  // Let us pick a sampling test function for solution evaluation
932  // SH, SV, FH, USERVAR
933  bool userVariable = false;
935  if( !runCtx->variableToVerify.compare( "SH" ) )
936  testFunction = &sample_slow_harmonic;
937  else if( !runCtx->variableToVerify.compare( "FH" ) )
938  testFunction = &sample_fast_harmonic;
939  else if( !runCtx->variableToVerify.compare( "SV" ) )
940  testFunction = &sample_stationary_vortex;
941  else
942  {
943  userVariable = runCtx->variableToVerify.size() ? true : false;
944  testFunction = runCtx->variableToVerify.size() ? nullptr : sample_stationary_vortex;
945  }
946 
947  moab::Tag srcAnalyticalFunction;
948  moab::Tag tgtAnalyticalFunction;
949  moab::Tag tgtProjectedFunction;
950  if( testFunction )
951  {
952  runCtx->timer_push( "describe a solution on source grid" );
953  rval = weightMap->DefineAnalyticalSolution( srcAnalyticalFunction, "AnalyticalSolnSrcExact",
954  moab::Remapper::SourceMesh, testFunction );MB_CHK_ERR( rval );
955  runCtx->timer_pop();
956 
957  runCtx->timer_push( "describe a solution on target grid" );
958 
959  rval = weightMap->DefineAnalyticalSolution( tgtAnalyticalFunction, "AnalyticalSolnTgtExact",
960  moab::Remapper::TargetMesh, testFunction,
961  &tgtProjectedFunction, "ProjectedSolnTgt" );MB_CHK_ERR( rval );
962  // rval = mbCore->write_file ( "tgtWithSolnTag.h5m", nullptr, writeOptions,
963  // &runCtx->meshsets[1], 1 ); MB_CHK_ERR ( rval );
964  runCtx->timer_pop();
965  }
966  else
967  {
968  rval = mbCore->tag_get_handle( runCtx->variableToVerify.c_str(), srcAnalyticalFunction );MB_CHK_ERR( rval );
969 
970  rval = mbCore->tag_get_handle( "ProjectedSolnTgt", 1, moab::MB_TYPE_DOUBLE, tgtProjectedFunction,
972  }
973 
974  if( !runCtx->skip_io )
975  {
976  rval = mbCore->write_file( "srcWithSolnTag.h5m", nullptr, writeOptions, &runCtx->meshsets[0], 1 );MB_CHK_ERR( rval );
977  }
978 
979  runCtx->timer_push( "compute solution projection on target grid" );
980  rval = weightMap->ApplyWeights( srcAnalyticalFunction, tgtProjectedFunction, false, runCtx->cassType );MB_CHK_ERR( rval );
981  runCtx->timer_pop();
982 
983  if( !runCtx->skip_io )
984  {
985  rval = mbCore->write_file( "tgtWithSolnTag2.h5m", nullptr, writeOptions, &runCtx->meshsets[1], 1 );MB_CHK_ERR( rval );
986  }
987 
988  if( nprocs == 1 && runCtx->baselineFile.size() )
989  {
990  // save the field from tgtWithSolnTag2 in a text file, and global ids for cells
991  moab::Range tgtCells;
992  rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 2, tgtCells );MB_CHK_ERR( rval );
993  std::vector< int > globIds;
994  globIds.resize( tgtCells.size() );
995  std::vector< double > vals;
996  vals.resize( tgtCells.size() );
997  moab::Tag projTag;
998  rval = mbCore->tag_get_handle( "ProjectedSolnTgt", projTag );MB_CHK_ERR( rval );
999  moab::Tag gid = mbCore->globalId_tag();
1000  rval = mbCore->tag_get_data( gid, tgtCells, &globIds[0] );MB_CHK_ERR( rval );
1001  rval = mbCore->tag_get_data( projTag, tgtCells, &vals[0] );MB_CHK_ERR( rval );
1002  std::fstream fs;
1003  fs.open( runCtx->baselineFile.c_str(), std::fstream::out );
1004  fs << std::setprecision( 15 ); // maximum precision for doubles
1005  for( size_t i = 0; i < tgtCells.size(); i++ )
1006  fs << globIds[i] << " " << vals[i] << "\n";
1007  fs.close();
1008  // for good measure, save the source file too, with the tag AnalyticalSolnSrcExact
1009  // it will be used later to test, along with a target file
1010  if( !runCtx->skip_io )
1011  {
1012  rval =
1013  mbCore->write_file( "srcWithSolnTag.h5m", nullptr, writeOptions, &runCtx->meshsets[0], 1 );MB_CHK_ERR( rval );
1014  }
1015  }
1016 
1017  // compute error metrics if it is a known analytical functional
1018  if( !userVariable )
1019  {
1020  runCtx->timer_push( "compute error metrics against analytical solution on target grid" );
1021  std::map< std::string, double > errMetrics;
1022  rval = weightMap->ComputeMetrics( moab::Remapper::TargetMesh, tgtAnalyticalFunction,
1023  tgtProjectedFunction, errMetrics, true );MB_CHK_ERR( rval );
1024  runCtx->timer_pop();
1025  }
1026  }
1027 
1028  delete weightMap;
1029  }
1030  }
1031 
1032  // Clean up
1033  remapper.clear();
1034  delete runCtx;
1035  delete mbCore;
1036 
1037 #ifdef MOAB_HAVE_MPI
1038  MPI_Finalize();
1039 #endif
1040  exit( 0 );
1041 }

References moab::Core::add_entities(), moab::TempestOnlineMap::ApplyWeights(), moab::IntxAreaUtils::area_on_sphere(), ToolContext::baselineFile, ToolContext::boxeps, ToolContext::cassType, moab::ParallelComm::check_all_shared_handles(), moab::TempestRemapper::clear(), moab::TempestOnlineMap::ComputeMetrics(), moab::TempestRemapper::ComputeOverlapMesh(), ToolContext::computeWeights, moab::TempestRemapper::ConstructCoveringSet(), moab::TempestRemapper::constructEdgeMap, moab::TempestRemapper::ConvertTempestMesh(), moab::Core::create_meshset(), CreateTempestMesh(), moab::TempestOnlineMap::DefineAnalyticalSolution(), ToolContext::disc_methods, ToolContext::disc_orders, ToolContext::doftag_names, moab::IntxUtils::enforce_convexity(), ToolContext::enforceConvexity, ToolContext::ensureMonotonicity, ToolContext::epsrel, error(), ErrorCode, ToolContext::fCheck, moab::Intx2Mesh::FindMaxEdges(), moab::IntxUtils::fix_degenerate_quads(), moab::IntxAreaUtils::GaussQuadrature, moab::TempestOnlineMap::GenerateRemappingWeights(), moab::Core::get_entities_by_dimension(), moab::TempestRemapper::GetMesh(), moab::TempestRemapper::GetMeshSet(), moab::TempestRemapper::GetOverlapAugmentedEntities(), moab::Core::globalId_tag(), ToolContext::inFilenames, moab::TempestRemapper::initialize(), moab::Intx2Mesh::intersect_meshes(), ToolContext::intxFilename, ToolContext::kdtreeSearch, moab::IntxAreaUtils::lHuiller, ToolContext::mapOptions, MB_CHK_ERR, MB_CHK_SET_ERR, MB_TAG_CREAT, MB_TAG_DENSE, MB_TYPE_DOUBLE, ToolContext::mbcore, ToolContext::meshes, MESHSET_SET, ToolContext::meshsets, ToolContext::meshType, moab::TempestRemapper::meshValidate, MOAB_VERSION, ToolContext::nlayers, ToolContext::outFilename, ToolContext::outputFormatter, moab::TempestRemapper::OVERLAP_MEMORY, moab::TempestRemapper::OVERLAP_MOAB, moab::Remapper::OverlapMesh, ToolContext::ParseCLOptions(), moab::IntxAreaUtils::positive_orientation(), ToolContext::print_diagnostics, moab::DebugOutput::printf(), moab::TempestOnlineMap::PrintMapStatistics(), ToolContext::proc_id, ToolContext::rrmGrids, sample_fast_harmonic(), sample_slow_harmonic(), sample_stationary_vortex(), moab::Intx2Mesh::set_box_error(), moab::Intx2Mesh::set_error_tolerance(), moab::Intx2MeshOnSphere::set_radius_destination_mesh(), moab::Intx2MeshOnSphere::set_radius_source_mesh(), moab::Range::size(), ToolContext::skip_io, moab::Remapper::SourceMesh, moab::subtract(), moab::Core::tag_get_data(), moab::Core::tag_get_handle(), moab::Remapper::TargetMesh, ToolContext::timer_pop(), ToolContext::timer_push(), ToolContext::useGnomonicProjection, ToolContext::variableToVerify, ToolContext::verifyWeights, moab::Core::write_file(), moab::Core::write_mesh(), and moab::TempestOnlineMap::WriteParallelMap().

◆ sample_constant()

double sample_constant ( double  dLon,
double  dLat 
)
inline

Definition at line 1224 of file mbtempest.cpp.

1225 {
1226  return 1.0;
1227 }

◆ sample_fast_harmonic()

double sample_fast_harmonic ( double  dLon,
double  dLat 
)
inline

Definition at line 1218 of file mbtempest.cpp.

1219 {
1220  return ( 2.0 + pow( sin( 2.0 * dLat ), 16.0 ) * cos( 16.0 * dLon ) );
1221  // return (2.0 + pow(cos(2.0 * dLat), 16.0) * cos(16.0 * dLon));
1222 }

Referenced by main().

◆ sample_slow_harmonic()

double sample_slow_harmonic ( double  dLon,
double  dLat 
)
inline

Definition at line 1213 of file mbtempest.cpp.

1214 {
1215  return ( 2.0 + cos( dLat ) * cos( dLat ) * cos( 2.0 * dLon ) );
1216 }

Referenced by main().

◆ sample_stationary_vortex()

double sample_stationary_vortex ( double  dLon,
double  dLat 
)
inline

Find the rotated longitude and latitude of a point on a sphere with pole at (dLonC, dLatC).

Definition at line 1229 of file mbtempest.cpp.

1230 {
1231  const double dLon0 = 0.0;
1232  const double dLat0 = 0.6;
1233  const double dR0 = 3.0;
1234  const double dD = 5.0;
1235  const double dT = 6.0;
1236 
1237  /// Find the rotated longitude and latitude of a point on a sphere
1238  /// with pole at (dLonC, dLatC).
1239  {
1240  double dSinC = sin( dLat0 );
1241  double dCosC = cos( dLat0 );
1242  double dCosT = cos( dLat );
1243  double dSinT = sin( dLat );
1244 
1245  double dTrm = dCosT * cos( dLon - dLon0 );
1246  double dX = dSinC * dTrm - dCosC * dSinT;
1247  double dY = dCosT * sin( dLon - dLon0 );
1248  double dZ = dSinC * dSinT + dCosC * dTrm;
1249 
1250  dLon = atan2( dY, dX );
1251  if( dLon < 0.0 )
1252  {
1253  dLon += 2.0 * M_PI;
1254  }
1255  dLat = asin( dZ );
1256  }
1257 
1258  double dRho = dR0 * cos( dLat );
1259  double dVt = 3.0 * sqrt( 3.0 ) / 2.0 / cosh( dRho ) / cosh( dRho ) * tanh( dRho );
1260 
1261  double dOmega;
1262  if( dRho == 0.0 )
1263  {
1264  dOmega = 0.0;
1265  }
1266  else
1267  {
1268  dOmega = dVt / dRho;
1269  }
1270 
1271  return ( 1.0 - tanh( dRho / dD * sin( dLon - dOmega * dT ) ) );
1272 }

Referenced by main().