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