MOAB: Mesh Oriented datABase  (version 5.3.1)
mbtempest.cpp
Go to the documentation of this file.
00001 /*
00002  * Usage: MOAB-Tempest tool
00003  *
00004  * Generate a Cubed-Sphere mesh: ./mbtempest -t 0 -res 25 -f cubed_sphere_mesh.exo
00005  * Generate a RLL mesh: ./mbtempest -t 1 -res 25 -f rll_mesh.exo
00006  * Generate a Icosahedral-Sphere mesh: ./mbtempest -t 2 -res 25 <-dual> -f icosa_mesh.exo
00007  *
00008  * Now you can compute the intersections between the meshes too!
00009  *
00010  * Generate the overlap mesh: ./mbtempest -t 3 -l cubed_sphere_mesh.exo -l rll_mesh.exo -f
00011  * overlap_mesh.exo
00012  *
00013  */
00014 
00015 #include <iostream>
00016 #include <iomanip>
00017 #include <cstdlib>
00018 #include <vector>
00019 #include <string>
00020 #include <sstream>
00021 #include <cassert>
00022 
00023 #include "moab/Core.hpp"
00024 #include "moab/IntxMesh/IntxUtils.hpp"
00025 #include "moab/Remapping/TempestRemapper.hpp"
00026 #include "moab/Remapping/TempestOnlineMap.hpp"
00027 #include "moab/ProgOptions.hpp"
00028 #include "moab/CpuTimer.hpp"
00029 #include "DebugOutput.hpp"
00030 
00031 //#ifndef MOAB_HAVE_MPI
00032 //    #error mbtempest tool requires MPI configuration
00033 //#endif
00034 
00035 #ifdef MOAB_HAVE_MPI
00036 // MPI includes
00037 #include "moab_mpi.h"
00038 #include "moab/ParallelComm.hpp"
00039 #include "MBParallelConventions.h"
00040 #endif
00041 
00042 struct ToolContext
00043 {
00044     moab::Interface* mbcore;
00045 #ifdef MOAB_HAVE_MPI
00046     moab::ParallelComm* pcomm;
00047 #endif
00048     const int proc_id, n_procs;
00049     moab::DebugOutput outputFormatter;
00050     int blockSize;
00051     std::vector< std::string > inFilenames;
00052     std::vector< Mesh* > meshes;
00053     std::vector< moab::EntityHandle > meshsets;
00054     std::vector< int > disc_orders;
00055     std::vector< std::string > disc_methods;
00056     std::vector< std::string > doftag_names;
00057     std::string outFilename;
00058     std::string intxFilename;
00059     std::string baselineFile;
00060     moab::TempestRemapper::TempestMeshType meshType;
00061     bool computeDual;
00062     bool computeWeights;
00063     bool verifyWeights;
00064     bool enforceConvexity;
00065     int ensureMonotonicity;
00066     bool fNoConservation;
00067     bool fVolumetric;
00068     bool rrmGrids;
00069     bool kdtreeSearch;
00070     bool fNoBubble, fInputConcave, fOutputConcave, fCheck;
00071 
00072 #ifdef MOAB_HAVE_MPI
00073     ToolContext( moab::Interface* icore, moab::ParallelComm* p_pcomm )
00074         : mbcore( icore ), pcomm( p_pcomm ), proc_id( pcomm->rank() ), n_procs( pcomm->size() ),
00075           outputFormatter( std::cout, pcomm->rank(), 0 ),
00076 #else
00077     ToolContext( moab::Interface* icore )
00078         : mbcore( icore ), proc_id( 0 ), n_procs( 1 ), outputFormatter( std::cout, 0, 0 ),
00079 #endif
00080           blockSize( 5 ), outFilename( "output.exo" ), intxFilename( "" ), baselineFile( "" ),
00081           meshType( moab::TempestRemapper::DEFAULT ), computeDual( false ), computeWeights( false ),
00082           verifyWeights( false ), enforceConvexity( false ), ensureMonotonicity( 0 ), fNoConservation( false ),
00083           fVolumetric( false ), rrmGrids( false ), kdtreeSearch( true ), fNoBubble( true ), fInputConcave( false ),
00084           fOutputConcave( false ), fCheck( n_procs > 1 ? false : true )
00085     {
00086         inFilenames.resize( 2 );
00087         doftag_names.resize( 2 );
00088         timer = new moab::CpuTimer();
00089 
00090         outputFormatter.set_prefix( "[MBTempest]: " );
00091     }
00092 
00093     ~ToolContext()
00094     {
00095         // for (unsigned i=0; i < meshes.size(); ++i) delete meshes[i];
00096         meshes.clear();
00097         inFilenames.clear();
00098         disc_orders.clear();
00099         disc_methods.clear();
00100         doftag_names.clear();
00101         outFilename.clear();
00102         intxFilename.clear();
00103         baselineFile.clear();
00104         meshsets.clear();
00105         delete timer;
00106     }
00107 
00108     void timer_push( std::string operation )
00109     {
00110         timer_ops = timer->time_since_birth();
00111         opName    = operation;
00112     }
00113 
00114     void timer_pop()
00115     {
00116         double locElapsed = timer->time_since_birth() - timer_ops, avgElapsed = 0, maxElapsed = 0;
00117 #ifdef MOAB_HAVE_MPI
00118         MPI_Reduce( &locElapsed, &maxElapsed, 1, MPI_DOUBLE, MPI_MAX, 0, pcomm->comm() );
00119         MPI_Reduce( &locElapsed, &avgElapsed, 1, MPI_DOUBLE, MPI_SUM, 0, pcomm->comm() );
00120 #else
00121         maxElapsed = locElapsed;
00122         avgElapsed = locElapsed;
00123 #endif
00124         if( !proc_id )
00125         {
00126             avgElapsed /= n_procs;
00127             std::cout << "[LOG] Time taken to " << opName.c_str() << ": max = " << maxElapsed
00128                       << ", avg = " << avgElapsed << "\n";
00129         }
00130         // std::cout << "\n[LOG" << proc_id << "] Time taken to " << opName << " = " <<
00131         // timer->time_since_birth() - timer_ops << std::endl;
00132         opName.clear();
00133     }
00134 
00135     void ParseCLOptions( int argc, char* argv[] )
00136     {
00137         ProgOptions opts;
00138         int imeshType                  = 0;
00139         std::string expectedFName      = "output.exo";
00140         std::string expectedMethod     = "fv";
00141         std::string expectedDofTagName = "GLOBAL_ID";
00142         int expectedOrder              = 1;
00143         bool fBubble                   = false;
00144 
00145         if( !proc_id )
00146         {
00147             std::cout << "Command line options provided to mbtempest:\n  ";
00148             for( int iarg = 0; iarg < argc; ++iarg )
00149                 std::cout << argv[iarg] << " ";
00150             std::cout << std::endl << std::endl;
00151         }
00152 
00153         opts.addOpt< int >( "type,t",
00154                             "Type of mesh (default=CS; Choose from [CS=0, RLL=1, ICO=2, OVERLAP_FILES=3, "
00155                             "OVERLAP_MEMORY=4, OVERLAP_MOAB=5])",
00156                             &imeshType );
00157         opts.addOpt< int >( "res,r", "Resolution of the mesh (default=5)", &blockSize );
00158         opts.addOpt< void >( "dual,d", "Output the dual of the mesh (relevant only for ICO mesh type)", &computeDual );
00159         opts.addOpt< std::string >( "file,f", "Output computed mesh or remapping weights to specified filename",
00160                                     &outFilename );
00161         opts.addOpt< std::string >(
00162             "load,l", "Input mesh filenames for source and target meshes. (relevant only when computing weights)",
00163             &expectedFName );
00164         opts.addOpt< void >( "noconserve,c",
00165                              "Do not apply conservation to the resultant weights (relevant only "
00166                              "when computing weights)",
00167                              &fNoConservation );
00168         opts.addOpt< void >( "volumetric,v",
00169                              "Apply a volumetric projection to compute the weights (relevant only "
00170                              "when computing weights)",
00171                              &fVolumetric );
00172         opts.addOpt< void >( "advfront,a",
00173                              "Use the advancing front intersection instead of the Kd-tree based algorithm "
00174                              "to compute mesh intersections. (relevant only when computing weights)" );
00175         opts.addOpt< int >( "monotonic,n", "Ensure monotonicity in the weight generation", &ensureMonotonicity );
00176         opts.addOpt< int >( "order,o", "Discretization orders for the source and target solution fields",
00177                             &expectedOrder );
00178         opts.addOpt< std::string >( "method,m", "Discretization method for the source and target solution fields",
00179                                     &expectedMethod );
00180         opts.addOpt< std::string >( "global_id,g",
00181                                     "Tag name that contains the global DoF IDs for source and target solution fields",
00182                                     &expectedDofTagName );
00183 
00184         opts.addOpt< void >( "weights,w",
00185                              "Compute and output the weights using the overlap mesh (generally "
00186                              "relevant only for OVERLAP mesh)",
00187                              &computeWeights );
00188         opts.addOpt< std::string >( "intx,i", "Output TempestRemap intersection mesh filename", &intxFilename );
00189 
00190         opts.addOpt< std::string >( "baseline", "Output baseline file", &baselineFile );
00191 
00192         opts.addOpt< void >( "rrmgrids",
00193                              "At least one of the meshes is a regionally refined grid (relevant to "
00194                              "accelerate intersection computation)",
00195                              &rrmGrids );
00196         opts.addOpt< void >( "checkmap", "Check the generated map for conservation and consistency", &fCheck );
00197         opts.addOpt< void >( "verify",
00198                              "Verify the accuracy of the maps by projecting analytical functions "
00199                              "from source to target "
00200                              "grid by applying the maps",
00201                              &verifyWeights );
00202         opts.addOpt< void >( "enforce_convexity", "check convexity of input meshes to compute mesh intersections",
00203                              &enforceConvexity );
00204         opts.addOpt< void >( "bubble", "use bubble on interior of spectral element nodes", &fBubble );
00205 
00206         opts.parseCommandLine( argc, argv );
00207 
00208         // By default - use Kd-tree based search; if user asks for advancing front, disable Kd-tree
00209         // algorithm
00210         kdtreeSearch = opts.numOptSet( "advfront,a" ) == 0;
00211 
00212         // negate the option
00213         fNoBubble = !fBubble;
00214 
00215         switch( imeshType )
00216         {
00217             case 0:
00218                 meshType = moab::TempestRemapper::CS;
00219                 break;
00220 
00221             case 1:
00222                 meshType = moab::TempestRemapper::RLL;
00223                 break;
00224 
00225             case 2:
00226                 meshType = moab::TempestRemapper::ICO;
00227                 break;
00228 
00229             case 3:
00230                 meshType = moab::TempestRemapper::OVERLAP_FILES;
00231                 break;
00232 
00233             case 4:
00234                 meshType = moab::TempestRemapper::OVERLAP_MEMORY;
00235                 break;
00236 
00237             case 5:
00238                 meshType = moab::TempestRemapper::OVERLAP_MOAB;
00239                 break;
00240 
00241             default:
00242                 meshType = moab::TempestRemapper::DEFAULT;
00243                 break;
00244         }
00245 
00246         if( meshType > moab::TempestRemapper::ICO )
00247         {
00248             opts.getOptAllArgs( "load,l", inFilenames );
00249             opts.getOptAllArgs( "order,o", disc_orders );
00250             opts.getOptAllArgs( "method,m", disc_methods );
00251             opts.getOptAllArgs( "global_id,i", doftag_names );
00252 
00253             if( disc_orders.size() == 0 )
00254             {
00255                 disc_orders.resize( 2, 1 );
00256             }
00257 
00258             if( disc_orders.size() == 1 )
00259             {
00260                 disc_orders.push_back( 1 );
00261             }
00262 
00263             if( disc_methods.size() == 0 )
00264             {
00265                 disc_methods.resize( 2, "fv" );
00266             }
00267 
00268             if( disc_methods.size() == 1 )
00269             {
00270                 disc_methods.push_back( "fv" );
00271             }
00272 
00273             if( doftag_names.size() == 0 )
00274             {
00275                 doftag_names.resize( 2, "GLOBAL_ID" );
00276             }
00277 
00278             if( doftag_names.size() == 1 )
00279             {
00280                 doftag_names.push_back( "GLOBAL_ID" );
00281             }
00282 
00283             assert( inFilenames.size() == 2 );
00284             assert( disc_orders.size() == 2 );
00285             assert( disc_methods.size() == 2 );
00286         }
00287 
00288         expectedFName.clear();
00289     }
00290 
00291   private:
00292     moab::CpuTimer* timer;
00293     double timer_ops;
00294     std::string opName;
00295 };
00296 
00297 // Forward declare some methods
00298 static moab::ErrorCode CreateTempestMesh( ToolContext&, moab::TempestRemapper& remapper, Mesh* );
00299 inline double sample_slow_harmonic( double dLon, double dLat );
00300 inline double sample_fast_harmonic( double dLon, double dLat );
00301 inline double sample_constant( double dLon, double dLat );
00302 inline double sample_stationary_vortex( double dLon, double dLat );
00303 
00304 int main( int argc, char* argv[] )
00305 {
00306     moab::ErrorCode rval;
00307     NcError error( NcError::verbose_nonfatal );
00308     std::stringstream sstr;
00309 
00310     int proc_id = 0, nprocs = 1;
00311 #ifdef MOAB_HAVE_MPI
00312     MPI_Init( &argc, &argv );
00313     MPI_Comm_rank( MPI_COMM_WORLD, &proc_id );
00314     MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
00315 #endif
00316 
00317     moab::Interface* mbCore = new( std::nothrow ) moab::Core;
00318 
00319     if( NULL == mbCore )
00320     {
00321         return 1;
00322     }
00323 
00324     ToolContext* runCtx;
00325 #ifdef MOAB_HAVE_MPI
00326     moab::ParallelComm* pcomm = new moab::ParallelComm( mbCore, MPI_COMM_WORLD, 0 );
00327 
00328     runCtx                   = new ToolContext( mbCore, pcomm );
00329     const char* writeOptions = ( nprocs > 1 ? "PARALLEL=WRITE_PART" : "" );
00330 #else
00331     runCtx                   = new ToolContext( mbCore );
00332     const char* writeOptions = "";
00333 #endif
00334     runCtx->ParseCLOptions( argc, argv );
00335 
00336     const double radius_src  = 1.0 /*2.0*acos(-1.0)*/;
00337     const double radius_dest = 1.0 /*2.0*acos(-1.0)*/;
00338 
00339     moab::DebugOutput& outputFormatter = runCtx->outputFormatter;
00340 
00341 #ifdef MOAB_HAVE_MPI
00342     moab::TempestRemapper remapper( mbCore, pcomm );
00343 #else
00344     moab::TempestRemapper remapper( mbCore );
00345 #endif
00346     remapper.meshValidate     = true;
00347     remapper.constructEdgeMap = false;
00348     remapper.initialize();
00349 
00350     // Default area_method = lHuiller; Options: Girard, GaussQuadrature (if TR is available)
00351 #ifdef MOAB_HAVE_TEMPESTREMAP
00352     moab::IntxAreaUtils areaAdaptor( moab::IntxAreaUtils::GaussQuadrature );
00353 #else
00354     moab::IntxAreaUtils areaAdaptor( moab::IntxAreaUtils::lHuiller );
00355 #endif
00356 
00357     Mesh* tempest_mesh = new Mesh();
00358     runCtx->timer_push( "create Tempest mesh" );
00359     rval = CreateTempestMesh( *runCtx, remapper, tempest_mesh );MB_CHK_ERR( rval );
00360     runCtx->timer_pop();
00361 
00362     const double epsrel = ReferenceTolerance;  // ReferenceTolerance is defined in Defines.h in tempestremap;
00363                                                // Defines.h is included in SparseMatrix.h
00364                                                // SparseMatrix.h is included in OfflineMap.h
00365                                                // OfflineMap.h is included in TempestOnlineMap.hpp
00366                                                // TempestOnlineMap.hpp is included in this file, and is part of MOAB
00367     // Some constant parameters
00368 
00369     const double boxeps = 1e-6;
00370 
00371     if( runCtx->meshType == moab::TempestRemapper::OVERLAP_MEMORY )
00372     {
00373         // Compute intersections with MOAB
00374         // For the overlap method, choose between: "fuzzy", "exact" or "mixed"
00375         assert( runCtx->meshes.size() == 3 );
00376 
00377 #ifdef MOAB_HAVE_MPI
00378         rval = pcomm->check_all_shared_handles();MB_CHK_ERR( rval );
00379 #endif
00380 
00381         // Load the meshes and validate
00382         rval = remapper.ConvertTempestMesh( moab::Remapper::SourceMesh );MB_CHK_ERR( rval );
00383         rval = remapper.ConvertTempestMesh( moab::Remapper::TargetMesh );MB_CHK_ERR( rval );
00384         remapper.SetMeshType( moab::Remapper::OverlapMesh, moab::TempestRemapper::OVERLAP_FILES );
00385         rval = remapper.ConvertTempestMesh( moab::Remapper::OverlapMesh );MB_CHK_ERR( rval );
00386         rval = mbCore->write_mesh( "tempest_intersection.h5m", &runCtx->meshsets[2], 1 );MB_CHK_ERR( rval );
00387 
00388         // print verbosely about the problem setting
00389         {
00390             moab::Range rintxverts, rintxelems;
00391             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 0, rintxverts );MB_CHK_ERR( rval );
00392             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 2, rintxelems );MB_CHK_ERR( rval );
00393             outputFormatter.printf( 0, "The source set contains %lu vertices and %lu elements \n", rintxverts.size(),
00394                                     rintxelems.size() );
00395 
00396             moab::Range bintxverts, bintxelems;
00397             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 0, bintxverts );MB_CHK_ERR( rval );
00398             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 2, bintxelems );MB_CHK_ERR( rval );
00399             outputFormatter.printf( 0, "The target set contains %lu vertices and %lu elements \n", bintxverts.size(),
00400                                     bintxelems.size() );
00401         }
00402 
00403         moab::EntityHandle intxset;  // == remapper.GetMeshSet(moab::Remapper::OverlapMesh);
00404 
00405         // Compute intersections with MOAB
00406         {
00407             // Create the intersection on the sphere object
00408             runCtx->timer_push( "setup the intersector" );
00409 
00410             moab::Intx2MeshOnSphere* mbintx = new moab::Intx2MeshOnSphere( mbCore );
00411             mbintx->set_error_tolerance( epsrel );
00412             mbintx->set_box_error( boxeps );
00413             mbintx->set_radius_source_mesh( radius_src );
00414             mbintx->set_radius_destination_mesh( radius_dest );
00415 #ifdef MOAB_HAVE_MPI
00416             mbintx->set_parallel_comm( pcomm );
00417 #endif
00418             rval = mbintx->FindMaxEdges( runCtx->meshsets[0], runCtx->meshsets[1] );MB_CHK_ERR( rval );
00419 
00420 #ifdef MOAB_HAVE_MPI
00421             moab::Range local_verts;
00422             rval = mbintx->build_processor_euler_boxes( runCtx->meshsets[1], local_verts );MB_CHK_ERR( rval );
00423 
00424             runCtx->timer_pop();
00425 
00426             moab::EntityHandle covering_set;
00427             runCtx->timer_push( "communicate the mesh" );
00428             rval = mbintx->construct_covering_set( runCtx->meshsets[0], covering_set );MB_CHK_ERR( rval );  // lots of communication if mesh is distributed very differently
00429             runCtx->timer_pop();
00430 #else
00431             moab::EntityHandle covering_set = runCtx->meshsets[0];
00432 #endif
00433             // Now let's invoke the MOAB intersection algorithm in parallel with a
00434             // source and target mesh set representing two different decompositions
00435             runCtx->timer_push( "compute intersections with MOAB" );
00436             rval = mbCore->create_meshset( moab::MESHSET_SET, intxset );MB_CHK_SET_ERR( rval, "Can't create new set" );
00437             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" );
00438             runCtx->timer_pop();
00439 
00440             // free the memory
00441             delete mbintx;
00442         }
00443 
00444         {
00445             moab::Range intxelems, intxverts;
00446             rval = mbCore->get_entities_by_dimension( intxset, 2, intxelems );MB_CHK_ERR( rval );
00447             rval = mbCore->get_entities_by_dimension( intxset, 0, intxverts, true );MB_CHK_ERR( rval );
00448             outputFormatter.printf( 0, "The intersection set contains %lu elements and %lu vertices \n",
00449                                     intxelems.size(), intxverts.size() );
00450 
00451             double initial_sarea =
00452                 areaAdaptor.area_on_sphere( mbCore, runCtx->meshsets[0],
00453                                             radius_src );  // use the target to compute the initial area
00454             double initial_tarea =
00455                 areaAdaptor.area_on_sphere( mbCore, runCtx->meshsets[1],
00456                                             radius_dest );  // use the target to compute the initial area
00457             double intx_area = areaAdaptor.area_on_sphere( mbCore, intxset, radius_src );
00458 
00459             outputFormatter.printf( 0, "mesh areas: source = %12.10f, target = %12.10f, intersection = %12.10f \n",
00460                                     initial_sarea, initial_tarea, intx_area );
00461             outputFormatter.printf( 0, "relative error w.r.t source = %12.10e, target = %12.10e \n",
00462                                     fabs( intx_area - initial_sarea ) / initial_sarea,
00463                                     fabs( intx_area - initial_tarea ) / initial_tarea );
00464         }
00465 
00466         // Write out our computed intersection file
00467         rval = mbCore->write_mesh( "moab_intersection.h5m", &intxset, 1 );MB_CHK_ERR( rval );
00468 
00469         if( runCtx->computeWeights )
00470         {
00471             runCtx->timer_push( "compute weights with the Tempest meshes" );
00472             // Call to generate an offline map with the tempest meshes
00473             OfflineMap weightMap;
00474             int err = GenerateOfflineMapWithMeshes(
00475                 weightMap, *runCtx->meshes[0], *runCtx->meshes[1], *runCtx->meshes[2], "",
00476                 "",  // std::string strInputMeta, std::string strOutputMeta,
00477                 runCtx->disc_methods[0],
00478                 runCtx->disc_methods[1],                         // std::string strInputType, std::string strOutputType,
00479                 runCtx->disc_orders[0], runCtx->disc_orders[1],  // int nPin=4, int nPout=4,
00480                 runCtx->fNoBubble, true,                         // bool fNoBubble = false, bool fCorrectAreas = false,
00481                 runCtx->ensureMonotonicity                       // int fMonotoneTypeID = 0
00482             );
00483             runCtx->timer_pop();
00484 
00485             std::map< std::string, std::string > mapAttributes;
00486             if( err )
00487             {
00488                 rval = moab::MB_FAILURE;
00489             }
00490             else
00491             {
00492                 weightMap.Write( "outWeights.nc", mapAttributes );
00493             }
00494         }
00495     }
00496     else if( runCtx->meshType == moab::TempestRemapper::OVERLAP_MOAB )
00497     {
00498         // Usage: mpiexec -n 2 tools/mbtempest -t 5 -l mycs_2.h5m -l myico_2.h5m -f myoverlap_2.h5m
00499 #ifdef MOAB_HAVE_MPI
00500         rval = pcomm->check_all_shared_handles();MB_CHK_ERR( rval );
00501 #endif
00502         // print verbosely about the problem setting
00503         {
00504             moab::Range rintxverts, rintxelems;
00505             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 0, rintxverts );MB_CHK_ERR( rval );
00506             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[0], 2, rintxelems );MB_CHK_ERR( rval );
00507             rval = moab::IntxUtils::fix_degenerate_quads( mbCore, runCtx->meshsets[0] );MB_CHK_ERR( rval );
00508             if( runCtx->enforceConvexity )
00509             {
00510                 rval = moab::IntxUtils::enforce_convexity( mbCore, runCtx->meshsets[0], proc_id );MB_CHK_ERR( rval );
00511             }
00512             rval = areaAdaptor.positive_orientation( mbCore, runCtx->meshsets[0], radius_src );MB_CHK_ERR( rval );
00513             if( !proc_id )
00514                 outputFormatter.printf( 0, "The source set contains %lu vertices and %lu elements \n",
00515                                         rintxverts.size(), rintxelems.size() );
00516 
00517             moab::Range bintxverts, bintxelems;
00518             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 0, bintxverts );MB_CHK_ERR( rval );
00519             rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 2, bintxelems );MB_CHK_ERR( rval );
00520             rval = moab::IntxUtils::fix_degenerate_quads( mbCore, runCtx->meshsets[1] );MB_CHK_ERR( rval );
00521             if( runCtx->enforceConvexity )
00522             {
00523                 rval = moab::IntxUtils::enforce_convexity( mbCore, runCtx->meshsets[1], proc_id );MB_CHK_ERR( rval );
00524             }
00525             rval = areaAdaptor.positive_orientation( mbCore, runCtx->meshsets[1], radius_dest );MB_CHK_ERR( rval );
00526             if( !proc_id )
00527                 outputFormatter.printf( 0, "The target set contains %lu vertices and %lu elements \n",
00528                                         bintxverts.size(), bintxelems.size() );
00529         }
00530 
00531         // First compute the covering set such that the target elements are fully covered by the
00532         // lcoal source grid
00533         runCtx->timer_push( "construct covering set for intersection" );
00534         rval = remapper.ConstructCoveringSet( epsrel, 1.0, 1.0, boxeps, runCtx->rrmGrids );MB_CHK_ERR( rval );
00535         runCtx->timer_pop();
00536 
00537         // Compute intersections with MOAB with either the Kd-tree or the advancing front algorithm
00538         runCtx->timer_push( "setup and compute mesh intersections" );
00539         rval = remapper.ComputeOverlapMesh( runCtx->kdtreeSearch, false );MB_CHK_ERR( rval );
00540         runCtx->timer_pop();
00541 
00542         // print some diagnostic checks to see if the overlap grid resolved the input meshes
00543         // correctly
00544         {
00545             double local_areas[3],
00546                 global_areas[3];  // Array for Initial area, and through Method 1 and Method 2
00547             // local_areas[0] = area_on_sphere_lHuiller ( mbCore, runCtx->meshsets[1], radius_src );
00548             local_areas[0] = areaAdaptor.area_on_sphere( mbCore, runCtx->meshsets[0], radius_src );
00549             local_areas[1] = areaAdaptor.area_on_sphere( mbCore, runCtx->meshsets[1], radius_dest );
00550             local_areas[2] = areaAdaptor.area_on_sphere( mbCore, runCtx->meshsets[2], radius_src );
00551 
00552 #ifdef MOAB_HAVE_MPI
00553             MPI_Allreduce( &local_areas[0], &global_areas[0], 3, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD );
00554 #else
00555             global_areas[0]                 = local_areas[0];
00556             global_areas[1]                 = local_areas[1];
00557             global_areas[2]                 = local_areas[2];
00558 #endif
00559             if( !proc_id )
00560             {
00561                 outputFormatter.printf( 0,
00562                                         "initial area: source mesh = %12.14f, target mesh = "
00563                                         "%12.14f, overlap mesh = %12.14f\n",
00564                                         global_areas[0], global_areas[1], global_areas[2] );
00565                 // outputFormatter.printf ( 0, " area with l'Huiller: %12.14f with Girard:
00566                 // %12.14f\n", global_areas[2], global_areas[3] ); outputFormatter.printf ( 0, "
00567                 // relative difference areas = %12.10e\n", fabs ( global_areas[2] - global_areas[3]
00568                 // ) / global_areas[2] );
00569                 outputFormatter.printf( 0, "relative error w.r.t source = %12.14e, and target = %12.14e\n",
00570                                         fabs( global_areas[0] - global_areas[2] ) / global_areas[0],
00571                                         fabs( global_areas[1] - global_areas[2] ) / global_areas[1] );
00572             }
00573         }
00574 
00575         if( runCtx->intxFilename.size() )
00576         {
00577             moab::EntityHandle writableOverlapSet;
00578             rval = mbCore->create_meshset( moab::MESHSET_SET, writableOverlapSet );MB_CHK_SET_ERR( rval, "Can't create new set" );
00579             moab::EntityHandle meshOverlapSet = remapper.GetMeshSet( moab::Remapper::OverlapMesh );
00580             moab::Range ovEnts;
00581             rval = mbCore->get_entities_by_dimension( meshOverlapSet, 2, ovEnts );MB_CHK_SET_ERR( rval, "Can't create new set" );
00582             rval = mbCore->get_entities_by_dimension( meshOverlapSet, 0, ovEnts );MB_CHK_SET_ERR( rval, "Can't create new set" );
00583 
00584 #ifdef MOAB_HAVE_MPI
00585             // Do not remove ghosted entities if we still haven't computed weights
00586             // Remove ghosted entities from overlap set before writing the new mesh set to file
00587             if( nprocs > 1 )
00588             {
00589                 moab::Range ghostedEnts;
00590                 rval = remapper.GetOverlapAugmentedEntities( ghostedEnts );MB_CHK_ERR( rval );
00591                 ovEnts = moab::subtract( ovEnts, ghostedEnts );
00592             }
00593 #endif
00594             rval = mbCore->add_entities( writableOverlapSet, ovEnts );MB_CHK_SET_ERR( rval, "Deleting ghosted entities failed" );
00595 
00596             size_t lastindex = runCtx->intxFilename.find_last_of( "." );
00597             sstr.str( "" );
00598             sstr << runCtx->intxFilename.substr( 0, lastindex ) << ".h5m";
00599             if( !runCtx->proc_id )
00600                 std::cout << "Writing out the MOAB intersection mesh file to " << sstr.str() << std::endl;
00601 
00602             // Write out our computed intersection file
00603             rval = mbCore->write_file( sstr.str().c_str(), NULL, writeOptions, &writableOverlapSet, 1 );MB_CHK_ERR( rval );
00604         }
00605 
00606         if( runCtx->computeWeights )
00607         {
00608             runCtx->meshes[2] = remapper.GetMesh( moab::Remapper::OverlapMesh );
00609             if( !runCtx->proc_id ) std::cout << std::endl;
00610 
00611             runCtx->timer_push( "setup computation of weights" );
00612             // Call to generate the remapping weights with the tempest meshes
00613             moab::TempestOnlineMap* weightMap = new moab::TempestOnlineMap( &remapper );
00614             runCtx->timer_pop();
00615 
00616             runCtx->timer_push( "compute weights with TempestRemap" );
00617 
00618             rval = weightMap->GenerateRemappingWeights(
00619                 runCtx->disc_methods[0],                         // std::string strInputType
00620                 runCtx->disc_methods[1],                         // std::string strOutputType,
00621                 runCtx->disc_orders[0], runCtx->disc_orders[1],  // int nPin, int nPout,
00622                 runCtx->fNoBubble,                               // bool fNoBubble
00623                 runCtx->ensureMonotonicity,                      // int fMonotoneTypeID,
00624                 runCtx->fVolumetric, runCtx->fNoConservation,    // bool fVolumetric, bool fNoConservation,
00625                 !runCtx->fCheck,                                 //  bool fNoCheck,
00626                 runCtx->doftag_names[0],                         // std::string source_tag_name
00627                 runCtx->doftag_names[1],                         // std::string target_tag_name,
00628                 runCtx->fInputConcave,                           // bool fInputConcave
00629                 runCtx->fOutputConcave                           // bool fOutputConcave
00630             );MB_CHK_ERR( rval );
00631             runCtx->timer_pop();
00632 
00633             // Invoke the CheckMap routine on the TempestRemap serial interface directly, if running
00634             // on a single process
00635             if( nprocs == 1 )
00636             {
00637                 const double dNormalTolerance = 1.0E-8;
00638                 const double dStrictTolerance = 1.0E-12;
00639                 weightMap->CheckMap( runCtx->fCheck, runCtx->fCheck, runCtx->fCheck && ( runCtx->ensureMonotonicity ),
00640                                      dNormalTolerance, dStrictTolerance );
00641             }
00642 
00643             if( runCtx->outFilename.size() )
00644             {
00645                 // Write the map file to disk in parallel using either HDF5 or SCRIP interface
00646                 rval = weightMap->WriteParallelMap( runCtx->outFilename.c_str() );MB_CHK_ERR( rval );
00647 
00648                 // Write out the metadata information for the map file
00649                 if( proc_id == 0 )
00650                 {
00651                     size_t lastindex = runCtx->outFilename.find_last_of( "." );
00652                     sstr.str( "" );
00653                     sstr << runCtx->outFilename.substr( 0, lastindex ) << ".meta";
00654 
00655                     std::ofstream metafile( sstr.str() );
00656                     metafile << "Generator = MOAB-TempestRemap (mbtempest) Offline Regridding "
00657                                 "Weight Generator"
00658                              << std::endl;
00659                     metafile << "domain_a = " << runCtx->inFilenames[0] << std::endl;
00660                     metafile << "domain_b = " << runCtx->inFilenames[1] << std::endl;
00661                     metafile << "domain_aNb = "
00662                              << ( runCtx->intxFilename.size() ? runCtx->intxFilename : "outOverlap.h5m" ) << std::endl;
00663                     metafile << "map_aPb = " << runCtx->outFilename << std::endl;
00664                     metafile << "type_src = " << runCtx->disc_methods[0] << std::endl;
00665                     metafile << "np_src = " << runCtx->disc_orders[0] << std::endl;
00666                     metafile << "concave_src = " << ( runCtx->fInputConcave ? "true" : "false" ) << std::endl;
00667                     metafile << "type_dst = " << runCtx->disc_methods[1] << std::endl;
00668                     metafile << "np_dst = " << runCtx->disc_orders[1] << std::endl;
00669                     metafile << "concave_dst = " << ( runCtx->fOutputConcave ? "true" : "false" ) << std::endl;
00670                     metafile << "mono_type = " << runCtx->ensureMonotonicity << std::endl;
00671                     metafile << "bubble = " << ( runCtx->fNoBubble ? "false" : "true" ) << std::endl;
00672                     metafile << "version = "
00673                              << "MOAB v5.1.0+" << std::endl;
00674                     metafile.close();
00675                 }
00676             }
00677 
00678             if( runCtx->verifyWeights )
00679             {
00680                 // Let us pick a sampling test function for solution evaluation
00681                 moab::TempestOnlineMap::sample_function testFunction =
00682                     &sample_stationary_vortex;  // &sample_slow_harmonic;
00683 
00684                 runCtx->timer_push( "describe a solution on source grid" );
00685                 moab::Tag srcAnalyticalFunction;
00686                 rval = weightMap->DefineAnalyticalSolution( srcAnalyticalFunction, "AnalyticalSolnSrcExact",
00687                                                             moab::Remapper::SourceMesh, testFunction );MB_CHK_ERR( rval );
00688                 runCtx->timer_pop();
00689                 // rval = mbCore->write_file ( "srcWithSolnTag.h5m", NULL, writeOptions,
00690                 // &runCtx->meshsets[0], 1 ); MB_CHK_ERR ( rval );
00691 
00692                 runCtx->timer_push( "describe a solution on target grid" );
00693                 moab::Tag tgtAnalyticalFunction;
00694                 moab::Tag tgtProjectedFunction;
00695                 rval = weightMap->DefineAnalyticalSolution( tgtAnalyticalFunction, "AnalyticalSolnTgtExact",
00696                                                             moab::Remapper::TargetMesh, testFunction,
00697                                                             &tgtProjectedFunction, "ProjectedSolnTgt" );MB_CHK_ERR( rval );
00698                 // rval = mbCore->write_file ( "tgtWithSolnTag.h5m", NULL, writeOptions,
00699                 // &runCtx->meshsets[1], 1 ); MB_CHK_ERR ( rval );
00700                 runCtx->timer_pop();
00701 
00702                 runCtx->timer_push( "compute solution projection on target grid" );
00703                 rval = weightMap->ApplyWeights( srcAnalyticalFunction, tgtProjectedFunction );MB_CHK_ERR( rval );
00704                 runCtx->timer_pop();
00705                 rval = mbCore->write_file( "tgtWithSolnTag2.h5m", NULL, writeOptions, &runCtx->meshsets[1], 1 );MB_CHK_ERR( rval );
00706 
00707                 if( nprocs == 1 && runCtx->baselineFile.size() )
00708                 {
00709                     // save the field from tgtWithSolnTag2 in a text file, and global ids for cells
00710                     moab::Range tgtCells;
00711                     rval = mbCore->get_entities_by_dimension( runCtx->meshsets[1], 2, tgtCells );MB_CHK_ERR( rval );
00712                     std::vector< int > globIds;
00713                     globIds.resize( tgtCells.size() );
00714                     std::vector< double > vals;
00715                     vals.resize( tgtCells.size() );
00716                     moab::Tag projTag;
00717                     rval = mbCore->tag_get_handle( "ProjectedSolnTgt", projTag );MB_CHK_ERR( rval );
00718                     moab::Tag gid = mbCore->globalId_tag();
00719                     rval          = mbCore->tag_get_data( gid, tgtCells, &globIds[0] );MB_CHK_ERR( rval );
00720                     rval = mbCore->tag_get_data( projTag, tgtCells, &vals[0] );MB_CHK_ERR( rval );
00721                     std::fstream fs;
00722                     fs.open( runCtx->baselineFile.c_str(), std::fstream::out );
00723                     fs << std::setprecision( 15 );  // maximum precision for doubles
00724                     for( size_t i = 0; i < tgtCells.size(); i++ )
00725                         fs << globIds[i] << " " << vals[i] << "\n";
00726                     fs.close();
00727                     // for good measure, save the source file too, with the tag AnalyticalSolnSrcExact
00728                     // it will be used later to test, along with a target file
00729                     rval = mbCore->write_file( "srcWithSolnTag.h5m", NULL, writeOptions, &runCtx->meshsets[0], 1 );MB_CHK_ERR( rval );
00730                 }
00731                 runCtx->timer_push( "compute error metrics against analytical solution on target grid" );
00732                 std::map< std::string, double > errMetrics;
00733                 rval = weightMap->ComputeMetrics( moab::Remapper::TargetMesh, tgtAnalyticalFunction,
00734                                                   tgtProjectedFunction, errMetrics, true );MB_CHK_ERR( rval );
00735                 runCtx->timer_pop();
00736             }
00737 
00738             delete weightMap;
00739         }
00740     }
00741 
00742     // Clean up
00743     remapper.clear();
00744     delete runCtx;
00745     delete mbCore;
00746 
00747 #ifdef MOAB_HAVE_MPI
00748     MPI_Finalize();
00749 #endif
00750     exit( 0 );
00751 }
00752 
00753 static moab::ErrorCode CreateTempestMesh( ToolContext& ctx, moab::TempestRemapper& remapper, Mesh* tempest_mesh )
00754 {
00755     moab::ErrorCode rval = moab::MB_SUCCESS;
00756     int err;
00757     moab::DebugOutput& outputFormatter = ctx.outputFormatter;
00758 
00759     if( !ctx.proc_id )
00760     {
00761         outputFormatter.printf( 0, "Creating TempestRemap Mesh object ...\n" );
00762     }
00763 
00764     if( ctx.meshType == moab::TempestRemapper::OVERLAP_FILES )
00765     {
00766         // For the overlap method, choose between: "fuzzy", "exact" or "mixed"
00767         err = GenerateOverlapMesh( ctx.inFilenames[0], ctx.inFilenames[1], *tempest_mesh, ctx.outFilename, "NetCDF4",
00768                                    "exact", true );
00769 
00770         if( err )
00771         {
00772             rval = moab::MB_FAILURE;
00773         }
00774         else
00775         {
00776             ctx.meshes.push_back( tempest_mesh );
00777         }
00778     }
00779     else if( ctx.meshType == moab::TempestRemapper::OVERLAP_MEMORY )
00780     {
00781         // Load the meshes and validate
00782         ctx.meshsets.resize( 3 );
00783         ctx.meshes.resize( 3 );
00784         ctx.meshsets[0] = remapper.GetMeshSet( moab::Remapper::SourceMesh );
00785         ctx.meshsets[1] = remapper.GetMeshSet( moab::Remapper::TargetMesh );
00786         ctx.meshsets[2] = remapper.GetMeshSet( moab::Remapper::OverlapMesh );
00787 
00788         // First the source
00789         rval = remapper.LoadMesh( moab::Remapper::SourceMesh, ctx.inFilenames[0], moab::TempestRemapper::DEFAULT );MB_CHK_ERR( rval );
00790         ctx.meshes[0] = remapper.GetMesh( moab::Remapper::SourceMesh );
00791 
00792         // Next the target
00793         rval = remapper.LoadMesh( moab::Remapper::TargetMesh, ctx.inFilenames[1], moab::TempestRemapper::DEFAULT );MB_CHK_ERR( rval );
00794         ctx.meshes[1] = remapper.GetMesh( moab::Remapper::TargetMesh );
00795 
00796         // Now let us construct the overlap mesh, by calling TempestRemap interface directly
00797         // For the overlap method, choose between: "fuzzy", "exact" or "mixed"
00798         err = GenerateOverlapWithMeshes( *ctx.meshes[0], *ctx.meshes[1], *tempest_mesh, "" /*ctx.outFilename*/,
00799                                          "NetCDF4", "exact", false );
00800 
00801         if( err )
00802         {
00803             rval = moab::MB_FAILURE;
00804         }
00805         else
00806         {
00807             remapper.SetMesh( moab::Remapper::OverlapMesh, tempest_mesh );
00808             ctx.meshes[2] = remapper.GetMesh( moab::Remapper::OverlapMesh );
00809             // ctx.meshes.push_back(*tempest_mesh);
00810         }
00811     }
00812     else if( ctx.meshType == moab::TempestRemapper::OVERLAP_MOAB )
00813     {
00814         ctx.meshsets.resize( 3 );
00815         ctx.meshes.resize( 3 );
00816         ctx.meshsets[0] = remapper.GetMeshSet( moab::Remapper::SourceMesh );
00817         ctx.meshsets[1] = remapper.GetMeshSet( moab::Remapper::TargetMesh );
00818         ctx.meshsets[2] = remapper.GetMeshSet( moab::Remapper::OverlapMesh );
00819 
00820         const double radius_src  = 1.0 /*2.0*acos(-1.0)*/;
00821         const double radius_dest = 1.0 /*2.0*acos(-1.0)*/;
00822 
00823         const char* additional_read_opts = ( ctx.n_procs > 1 ? "NO_SET_CONTAINING_PARENTS;" : "" );
00824         std::vector< int > smetadata, tmetadata;
00825         // Load the source mesh and validate
00826         rval = remapper.LoadNativeMesh( ctx.inFilenames[0], ctx.meshsets[0], smetadata, additional_read_opts );MB_CHK_ERR( rval );
00827         if( smetadata.size() )
00828         {
00829             // remapper.SetMeshType( moab::Remapper::SourceMesh,
00830             //                       static_cast< moab::TempestRemapper::TempestMeshType >( smetadata[0] ), &smetadata
00831             //                       );
00832             remapper.SetMeshType( moab::Remapper::SourceMesh,
00833                                   static_cast< moab::TempestRemapper::TempestMeshType >( smetadata[0] ) );
00834         }
00835         // Rescale the radius of both to compute the intersection
00836         rval = moab::IntxUtils::ScaleToRadius( ctx.mbcore, ctx.meshsets[0], radius_src );MB_CHK_ERR( rval );
00837         rval = remapper.ConvertMeshToTempest( moab::Remapper::SourceMesh );MB_CHK_ERR( rval );
00838         ctx.meshes[0] = remapper.GetMesh( moab::Remapper::SourceMesh );
00839 
00840         // Load the target mesh and validate
00841         rval = remapper.LoadNativeMesh( ctx.inFilenames[1], ctx.meshsets[1], tmetadata, additional_read_opts );MB_CHK_ERR( rval );
00842         if( tmetadata.size() )
00843         {
00844             // remapper.SetMeshType( moab::Remapper::TargetMesh,
00845             //                       static_cast< moab::TempestRemapper::TempestMeshType >( tmetadata[0] ), &tmetadata
00846             //                       );
00847             remapper.SetMeshType( moab::Remapper::TargetMesh,
00848                                   static_cast< moab::TempestRemapper::TempestMeshType >( tmetadata[0] ) );
00849         }
00850         rval = moab::IntxUtils::ScaleToRadius( ctx.mbcore, ctx.meshsets[1], radius_dest );MB_CHK_ERR( rval );
00851         rval = remapper.ConvertMeshToTempest( moab::Remapper::TargetMesh );MB_CHK_ERR( rval );
00852         ctx.meshes[1] = remapper.GetMesh( moab::Remapper::TargetMesh );
00853     }
00854     else if( ctx.meshType == moab::TempestRemapper::ICO )
00855     {
00856         err = GenerateICOMesh( *tempest_mesh, ctx.blockSize, ctx.computeDual, ctx.outFilename, "NetCDF4" );
00857 
00858         if( err )
00859         {
00860             rval = moab::MB_FAILURE;
00861         }
00862         else
00863         {
00864             ctx.meshes.push_back( tempest_mesh );
00865         }
00866     }
00867     else if( ctx.meshType == moab::TempestRemapper::RLL )
00868     {
00869         err = GenerateRLLMesh( *tempest_mesh,                     // Mesh& meshOut,
00870                                ctx.blockSize * 2, ctx.blockSize,  // int nLongitudes, int nLatitudes,
00871                                0.0, 360.0,                        // double dLonBegin, double dLonEnd,
00872                                -90.0, 90.0,                       // double dLatBegin, double dLatEnd,
00873                                false, false, false,  // bool fGlobalCap, bool fFlipLatLon, bool fForceGlobal,
00874                                "" /*ctx.inFilename*/, "",
00875                                "",  // std::string strInputFile, std::string strInputFileLonName, std::string
00876                                     // strInputFileLatName,
00877                                ctx.outFilename, "NetCDF4",  // std::string strOutputFile, std::string strOutputFormat
00878                                true                         // bool fVerbose
00879         );
00880 
00881         if( err )
00882         {
00883             rval = moab::MB_FAILURE;
00884         }
00885         else
00886         {
00887             ctx.meshes.push_back( tempest_mesh );
00888         }
00889     }
00890     else  // default
00891     {
00892         err = GenerateCSMesh( *tempest_mesh, ctx.blockSize, ctx.outFilename, "NetCDF4" );
00893 
00894         if( err )
00895         {
00896             rval = moab::MB_FAILURE;
00897         }
00898         else
00899         {
00900             ctx.meshes.push_back( tempest_mesh );
00901         }
00902     }
00903 
00904     if( ctx.meshType != moab::TempestRemapper::OVERLAP_MOAB && !tempest_mesh )
00905     {
00906         std::cout << "Tempest Mesh is not a complete object; Quitting...";
00907         exit( -1 );
00908     }
00909 
00910     return rval;
00911 }
00912 
00913 ///////////////////////////////////////////////
00914 //         Test functions
00915 
00916 double sample_slow_harmonic( double dLon, double dLat )
00917 {
00918     return ( 2.0 + cos( dLat ) * cos( dLat ) * cos( 2.0 * dLon ) );
00919 }
00920 
00921 double sample_fast_harmonic( double dLon, double dLat )
00922 {
00923     return ( 2.0 + pow( sin( 2.0 * dLat ), 16.0 ) * cos( 16.0 * dLon ) );
00924     // return (2.0 + pow(cos(2.0 * dLat), 16.0) * cos(16.0 * dLon));
00925 }
00926 
00927 double sample_constant( double /*dLon*/, double /*dLat*/ )
00928 {
00929     return 1.0;
00930 }
00931 
00932 double sample_stationary_vortex( double dLon, double dLat )
00933 {
00934     const double dLon0 = 0.0;
00935     const double dLat0 = 0.6;
00936     const double dR0   = 3.0;
00937     const double dD    = 5.0;
00938     const double dT    = 6.0;
00939 
00940     ///     Find the rotated longitude and latitude of a point on a sphere
00941     ///     with pole at (dLonC, dLatC).
00942     {
00943         double dSinC = sin( dLat0 );
00944         double dCosC = cos( dLat0 );
00945         double dCosT = cos( dLat );
00946         double dSinT = sin( dLat );
00947 
00948         double dTrm = dCosT * cos( dLon - dLon0 );
00949         double dX   = dSinC * dTrm - dCosC * dSinT;
00950         double dY   = dCosT * sin( dLon - dLon0 );
00951         double dZ   = dSinC * dSinT + dCosC * dTrm;
00952 
00953         dLon = atan2( dY, dX );
00954         if( dLon < 0.0 )
00955         {
00956             dLon += 2.0 * M_PI;
00957         }
00958         dLat = asin( dZ );
00959     }
00960 
00961     double dRho = dR0 * cos( dLat );
00962     double dVt  = 3.0 * sqrt( 3.0 ) / 2.0 / cosh( dRho ) / cosh( dRho ) * tanh( dRho );
00963 
00964     double dOmega;
00965     if( dRho == 0.0 )
00966     {
00967         dOmega = 0.0;
00968     }
00969     else
00970     {
00971         dOmega = dVt / dRho;
00972     }
00973 
00974     return ( 1.0 - tanh( dRho / dD * sin( dLon - dOmega * dT ) ) );
00975 }
00976 
00977 ///////////////////////////////////////////////
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines