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