MOAB: Mesh Oriented datABase  (version 5.4.1)
uref_hirec_test_parallel.cpp
Go to the documentation of this file.
00001 /*This unit test is for the convergence study of high order reconstruction under uniform
00002  * refinement*/
00003 #include <iostream>
00004 #include <string>
00005 #include <sstream>
00006 #if defined( __MINGW32__ )
00007 #include <sys/time.h>
00008 #else
00009 #include <ctime>
00010 #endif
00011 
00012 #include <vector>
00013 #include <algorithm>
00014 #include "moab/Core.hpp"
00015 #include "moab/Range.hpp"
00016 #include "moab/MeshTopoUtil.hpp"
00017 #include "moab/NestedRefine.hpp"
00018 #include "moab/DiscreteGeometry/HiReconstruction.hpp"
00019 #include "TestUtil.hpp"
00020 #include "geomObject.cpp"
00021 #include <cmath>
00022 #include <cstdlib>
00023 #include <cassert>
00024 
00025 #ifdef MOAB_HAVE_MPI
00026 #include "moab/ParallelComm.hpp"
00027 #include "MBParallelConventions.h"
00028 #include "ReadParallel.hpp"
00029 #include "moab/FileOptions.hpp"
00030 #include "MBTagConventions.hpp"
00031 #include "moab_mpi.h"
00032 #endif
00033 
00034 using namespace moab;
00035 
00036 #define nsamples 10
00037 
00038 #ifdef MOAB_HAVE_MPI
00039 std::string read_options;
00040 #endif
00041 
00042 ErrorCode load_meshset_hirec( const char* infile,
00043                               Interface* mbimpl,
00044                               EntityHandle& meshset,
00045                               ParallelComm*& pc,
00046                               const int degree,
00047                               const int dim )
00048 {
00049     ErrorCode error;
00050     error = mbimpl->create_meshset( moab::MESHSET_SET, meshset );MB_CHK_ERR( error );
00051 #ifdef MOAB_HAVE_MPI
00052     int nprocs, rank;
00053     MPI_Comm comm = MPI_COMM_WORLD;
00054     MPI_Comm_size( comm, &nprocs );
00055     MPI_Comm_rank( comm, &rank );
00056     EntityHandle partnset;
00057     error = mbimpl->create_meshset( moab::MESHSET_SET, partnset );MB_CHK_ERR( error );
00058 
00059     if( nprocs > 1 )
00060     {
00061         pc = moab::ParallelComm::get_pcomm( mbimpl, partnset, &comm );
00062     }
00063 
00064     if( nprocs > 1 )
00065     {
00066         int nghlayers = degree > 0 ? HiReconstruction::estimate_num_ghost_layers( degree, true ) : 0;
00067         assert( nghlayers < 10 );
00068 
00069         if( nghlayers )
00070         {
00071             // get ghost layers
00072             if( dim == 2 )
00073             {
00074                 read_options = "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_"
00075                                "SHARED_ENTS;PARALLEL_GHOSTS=2.0.";
00076             }
00077             else if( dim == 1 )
00078             {
00079                 read_options = "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_"
00080                                "SHARED_ENTS;PARALLEL_GHOSTS=1.0.";
00081             }
00082             else
00083             {
00084                 MB_SET_ERR( MB_FAILURE, "unsupported dimension" );
00085             }
00086 
00087             read_options += (char)( '0' + nghlayers );
00088             // std::cout << "On processor " << rank << " Degree=" << degree << " "<< read_options <<
00089             // std::endl;
00090         }
00091         else
00092         {
00093             read_options = "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_SHARED_ENTS;";
00094         }
00095 
00096         /*for debug*/
00097         /*std::string outfile = std::string(infile); std::size_t dotpos = outfile.find_last_of(".");
00098         std::ostringstream convert; convert << rank;
00099         std::string localfile = outfile.substr(0,dotpos)+convert.str()+".vtk";
00100         //write local mesh
00101         EntityHandle local;
00102         error = mbimpl->create_meshset(moab::MESHSET_SET,local);CHECK_ERR(error);
00103         std::string local_options =
00104         "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_SHARED_ENTS;"; error =
00105         mbimpl->load_file(infile,&local,local_options.c_str()); MB_CHK_ERR(error); error =
00106         mbimpl->write_file(localfile.c_str(),0,0,&local,1);*/
00107 
00108         error = mbimpl->load_file( infile, &meshset, read_options.c_str() );MB_CHK_ERR( error );
00109         /*for debug
00110         //write local mesh with ghost layers
00111         localfile = outfile.substr(0,dotpos)+convert.str()+"_ghost.vtk";
00112         error = mbimpl->write_file(localfile.c_str(),0,0,&meshset,1);*/
00113     }
00114     else
00115     {
00116         error = mbimpl->load_file( infile, &meshset );MB_CHK_ERR( error );
00117     }
00118 
00119 #else
00120     assert( !pc && degree && dim );
00121     error = mbimpl->load_file( infile, &meshset );MB_CHK_ERR( error );
00122 #endif
00123     return error;
00124 }
00125 
00126 ErrorCode closedsurface_uref_hirec_convergence_study( const char* infile,
00127                                                       std::vector< int >& degs2fit,
00128                                                       bool interp,
00129                                                       int dim,
00130                                                       geomObject* obj,
00131                                                       int& ntestverts,
00132                                                       std::vector< double >& geoml1errs,
00133                                                       std::vector< double >& geoml2errs,
00134                                                       std::vector< double >& geomlinferrs )
00135 {
00136     Core moab;
00137     Interface* mbImpl = &moab;
00138     ParallelComm* pc  = NULL;
00139     EntityHandle meshset;
00140 
00141 #ifdef MOAB_HAVE_MPI
00142     int nprocs, rank;
00143     MPI_Comm comm = MPI_COMM_WORLD;
00144     MPI_Comm_size( comm, &nprocs );
00145     MPI_Comm_rank( comm, &rank );
00146 #endif
00147 
00148     ErrorCode error;
00149     // mesh will be loaded and communicator pc will be updated
00150     int mxdeg = 1;
00151     for( size_t i = 0; i < degs2fit.size(); ++i )
00152     {
00153         mxdeg = std::max( degs2fit[i], mxdeg );
00154     }
00155     error = load_meshset_hirec( infile, mbImpl, meshset, pc, mxdeg, dim );MB_CHK_ERR( error );
00156 
00157     Range elems, elems_owned;
00158     error = mbImpl->get_entities_by_dimension( meshset, dim, elems );MB_CHK_ERR( error );
00159 
00160 #ifdef MOAB_HAVE_MPI
00161 
00162     if( pc )
00163     {
00164         error = pc->filter_pstatus( elems, PSTATUS_GHOST, PSTATUS_NOT, -1, &elems_owned );MB_CHK_ERR( error );
00165     }
00166     else
00167     {
00168         elems_owned = elems;
00169     }
00170 
00171 #endif
00172 
00173 #ifdef MOAB_HAVE_MPI
00174     // std::cout << "Mesh has " << nelems << " elements on Processor " << rank << " in total;";
00175     // std::cout << elems_owned.size() << " of which are locally owned elements" << std::endl;
00176 #else
00177     // std::cout << "Mesh has " << nelems << " elements" << std::endl;
00178 #endif
00179 
00180     /************************
00181      *   convergence study   *
00182      *************************/
00183     // project onto exact geometry since each level with uref has only linear coordinates
00184     Range verts;
00185     error = mbImpl->get_entities_by_dimension( meshset, 0, verts );MB_CHK_ERR( error );
00186 
00187     for( Range::iterator ivert = verts.begin(); ivert != verts.end(); ++ivert )
00188     {
00189         EntityHandle currvert = *ivert;
00190         double currcoords[3], exactcoords[3];
00191         error = mbImpl->get_coords( &currvert, 1, currcoords );MB_CHK_ERR( error );
00192         obj->project_points2geom( 3, currcoords, exactcoords, NULL );
00193 
00194         error = mbImpl->set_coords( &currvert, 1, exactcoords );MB_CHK_ERR( error );
00195         // for debug
00196         /*error = mbImpl->get_coords(&currvert,1,currcoords); MB_CHK_ERR(error);
00197         assert(currcoords[0]==exactcoords[0]&&currcoords[1]==exactcoords[1]&&currcoords[2]==exactcoords[2]);*/
00198     }
00199 
00200     // test ahf on mesh with ghost layers
00201     /*Range verts_owned;
00202     #ifdef MOAB_HAVE_MPI
00203         if(pc){
00204             error =
00205     pc->filter_pstatus(verts,PSTATUS_GHOST,PSTATUS_NOT,-1,&verts_owned);MB_CHK_ERR(error); }else{
00206     verts_owned = verts;
00207         }
00208 
00209         HalfFacetRep *ahf = new HalfFacetRep(&moab,pc,meshset);
00210         error = ahf->initialize(); MB_CHK_ERR(error);
00211 
00212         for(int i=0;i<nprocs;++i){
00213             if(rank==i){
00214                 std::cout << "Processor " << rank << " Local elements: \n";
00215                 for(Range::iterator ielem=elems.begin();ielem!=elems.end();++ielem){
00216                     EntityHandle currelem = *ielem;
00217                     std::vector<EntityHandle> conn;
00218                     error = mbImpl->get_connectivity(&currelem,1,conn); MB_CHK_ERR(error);
00219                     std::cout << *ielem << ": ";
00220                     for(size_t k=0;k<conn.size();++k) std::cout << conn[k] << " ";
00221                     std::cout << std::endl;
00222                 }
00223 
00224                 std::cout << "Processor " << rank << " Local owned elements: \n";
00225                 for(Range::iterator ielem=elems_owned.begin();ielem!=elems_owned.end();++ielem){
00226                     EntityHandle currelem = *ielem;
00227                     std::vector<EntityHandle> conn;
00228                     error = mbImpl->get_connectivity(&currelem,1,conn); MB_CHK_ERR(error);
00229                     std::cout << *ielem << ": ";
00230                     for(size_t k=0;k<conn.size();++k) std::cout << conn[k] << " ";
00231                     std::cout << std::endl;
00232                 }
00233 
00234 
00235                 std::cout << "Processor " << rank << " Local vertices: ";
00236                 for(Range::iterator ivert=verts.begin();ivert!=verts.end();++ivert){
00237                     std::cout << *ivert << " ";
00238                 }
00239                 std::cout << std::endl;
00240                 std::cout << "Processor " << rank << " Local owned vertices: ";
00241                 for(Range::iterator ivert=verts_owned.begin();ivert!=verts_owned.end();++ivert){
00242                     std::cout << *ivert << " ";
00243                 }
00244                 std::cout << std::endl;
00245 
00246                 for(Range::iterator ivert=verts_owned.begin();ivert!=verts_owned.end();++ivert){
00247                     EntityHandle currvid = *ivert;
00248                     std::cout << "Processor " << rank << " local verts: " << *ivert << " has
00249     adjfaces: "; std::vector<EntityHandle> adjfaces;
00250                     //error = ahf->get_up_adjacencies(currvid,2,adjfaces); MB_CHK_ERR(error);
00251                     error = mbImpl->get_adjacencies(&currvid,1,2,false,adjfaces); MB_CHK_ERR(error);
00252                     for(size_t k=0;k<adjfaces.size();++k){
00253                         std::cout << adjfaces[k] << " ";
00254                     }
00255                     std::cout << std::endl;
00256                 }
00257 
00258                 for(Range::iterator ivert=verts.begin();ivert!=verts.end();++ivert){
00259                     EntityHandle currvid = *ivert;
00260                     std::cout << "Processor " << rank << " all verts: " << *ivert << " has adjfaces:
00261     "; std::vector<EntityHandle> adjfaces;
00262                     //error = ahf->get_up_adjacencies(*ivert,2,adjfaces); MB_CHK_ERR(error);
00263                     error = mbImpl->get_adjacencies(&currvid,1,2,false,adjfaces); MB_CHK_ERR(error);
00264                     for(size_t k=0;k<adjfaces.size();++k){
00265                         std::cout << adjfaces[k] << " ";
00266                     }
00267                     std::cout << std::endl;
00268                 }
00269 
00270             }else{
00271                 MPI_Barrier(MPI_COMM_WORLD);
00272             }
00273         }
00274         exit(0);
00275     #endif*/
00276 
00277     // generate random points on each elements, assument 3D coordinates
00278     int nvpe = TYPE_FROM_HANDLE( *elems.begin() ) == MBTRI ? 3 : 4;
00279     std::vector< double > testpnts, testnaturalcoords;
00280     ntestverts = elems_owned.size() * nsamples;
00281     testpnts.reserve( 3 * elems_owned.size() * nsamples );
00282     testnaturalcoords.reserve( nvpe * elems_owned.size() * nsamples );
00283 
00284     for( Range::iterator ielem = elems_owned.begin(); ielem != elems_owned.end(); ++ielem )
00285     {
00286         EntityHandle currelem = *ielem;
00287         std::vector< EntityHandle > conn;
00288         error = mbImpl->get_connectivity( &currelem, 1, conn );MB_CHK_ERR( error );
00289         std::vector< double > elemcoords( 3 * conn.size() );
00290         error = mbImpl->get_coords( &( conn[0] ), conn.size(), &( elemcoords[0] ) );MB_CHK_ERR( error );
00291         EntityType type = TYPE_FROM_HANDLE( currelem );
00292 
00293         for( int s = 0; s < nsamples; ++s )
00294         {
00295             if( type == MBTRI )
00296             {
00297                 double a = (double)rand() / RAND_MAX, b = (double)rand() / RAND_MAX, c = (double)rand() / RAND_MAX, sum;
00298                 sum = a + b + c;
00299 
00300                 if( sum < 1e-12 )
00301                 {
00302                     --s;
00303                     continue;
00304                 }
00305                 else
00306                 {
00307                     a /= sum, b /= sum, c /= sum;
00308                 }
00309 
00310                 assert( a > 0 && b > 0 && c > 0 && fabs( a + b + c - 1 ) < 1e-12 );
00311                 testpnts.push_back( a * elemcoords[0] + b * elemcoords[3] + c * elemcoords[6] );
00312                 testpnts.push_back( a * elemcoords[1] + b * elemcoords[4] + c * elemcoords[7] );
00313                 testpnts.push_back( a * elemcoords[2] + b * elemcoords[5] + c * elemcoords[8] );
00314                 testnaturalcoords.push_back( a ), testnaturalcoords.push_back( b ), testnaturalcoords.push_back( c );
00315             }
00316             else if( type == MBQUAD )
00317             {
00318                 double xi = (double)rand() / RAND_MAX, eta = (double)rand() / RAND_MAX, N[4];
00319                 xi   = 2 * xi - 1;
00320                 eta  = 2 * eta - 1;
00321                 N[0] = ( 1 - xi ) * ( 1 - eta ) / 4, N[1] = ( 1 + xi ) * ( 1 - eta ) / 4,
00322                 N[2] = ( 1 + xi ) * ( 1 + eta ) / 4, N[3] = ( 1 - xi ) * ( 1 + eta ) / 4;
00323                 testpnts.push_back( N[0] * elemcoords[0] + N[1] * elemcoords[3] + N[2] * elemcoords[6] +
00324                                     N[3] * elemcoords[9] );
00325                 testpnts.push_back( N[0] * elemcoords[1] + N[1] * elemcoords[4] + N[2] * elemcoords[7] +
00326                                     N[3] * elemcoords[10] );
00327                 testpnts.push_back( N[0] * elemcoords[2] + N[1] * elemcoords[5] + N[2] * elemcoords[8] +
00328                                     N[3] * elemcoords[11] );
00329                 testnaturalcoords.push_back( N[0] ), testnaturalcoords.push_back( N[1] ),
00330                     testnaturalcoords.push_back( N[2] ), testnaturalcoords.push_back( N[3] );
00331             }
00332             else
00333             {
00334                 throw std::invalid_argument( "NOT SUPPORTED ELEMENT TYPE" );
00335             }
00336         }
00337     }
00338 
00339     // Compute error of linear interpolation in PARALLEL
00340     double l1err, l2err, linferr;
00341     geoml1errs.clear();
00342     geoml2errs.clear();
00343     geomlinferrs.clear();
00344     obj->compute_projecterror( 3, elems_owned.size() * nsamples, &( testpnts[0] ), l1err, l2err, linferr );
00345     geoml1errs.push_back( l1err );
00346     geoml2errs.push_back( l2err );
00347     geomlinferrs.push_back( linferr );
00348     /*Perform high order projection and compute error*/
00349     // initialize
00350     HiReconstruction hirec( &moab, pc, meshset );
00351 
00352     for( size_t ideg = 0; ideg < degs2fit.size(); ++ideg )
00353     {
00354         // High order reconstruction
00355         error = hirec.reconstruct3D_surf_geom( degs2fit[ideg], interp, true, true );MB_CHK_ERR( error );
00356 
00357         int index = 0;
00358         // for debug
00359         /*double maxlinferr=0,elemlinferr=0,eleml1err=0,eleml2err;
00360         EntityHandle elem_maxerr;*/
00361 
00362         for( Range::iterator ielem = elems_owned.begin(); ielem != elems_owned.end(); ++ielem, ++index )
00363         {
00364             // Projection
00365             error =
00366                 hirec.hiproj_walf_in_element( *ielem, nvpe, nsamples, &( testnaturalcoords[nvpe * nsamples * index] ),
00367                                               &( testpnts[3 * nsamples * index] ) );MB_CHK_ERR( error );
00368             // for debug
00369             /*obj->compute_projecterror(3,nsamples,&(testpnts[3*nsamples*index]),eleml1err,eleml2err,elemlinferr);
00370             if(elemlinferr>maxlinferr){
00371                 maxlinferr = elemlinferr; elem_maxerr = *ielem;
00372             }*/
00373         }
00374 
00375         assert( (size_t)index == elems_owned.size() );
00376         // Estimate error
00377         obj->compute_projecterror( 3, elems_owned.size() * nsamples, &( testpnts[0] ), l1err, l2err, linferr );
00378         geoml1errs.push_back( l1err );
00379         geoml2errs.push_back( l2err );
00380         geomlinferrs.push_back( linferr );
00381 
00382         // for debug
00383         /*std::cout << "Degree " << degs2fit[ideg] << " has L-inf error " << maxlinferr << " at
00384         element " << elems_owned.index(elem_maxerr) << ":"; std::vector<EntityHandle> conn; error =
00385         mbImpl->get_connectivity(&elem_maxerr,1,conn); MB_CHK_ERR(error); for(size_t
00386         ii=0;ii<conn.size();++ii) std::cout << verts.index(conn[ii]) << " "; std::cout << std::endl;
00387 
00388         for(size_t ii=0;ii<conn.size();++ii){
00389             //assume triangle
00390             std::vector<double> vertexcoords(3,0);
00391             error = mbImpl->get_coords(&(conn[ii]),1,&(vertexcoords[0])); MB_CHK_ERR(error);
00392             std::cout << verts.index(conn[ii]) << ":" << vertexcoords[0] << " "<< vertexcoords[1] <<
00393         " "<< vertexcoords[2] << "\n"; GEOMTYPE geomtype; std::vector<double> local_coords_system,
00394         local_coeffs; int deg_out; bool local_interp; bool hasfit =
00395         hirec.get_fittings_data(conn[ii],geomtype,local_coords_system,deg_out,local_coeffs,local_interp);
00396         assert(hasfit); std::cout << "Local fitting result: " << std::endl; std::cout << "Geom Type:
00397         " << geomtype << std::endl; std::cout << "Local coordinates system: " << std::endl;
00398             std::cout << local_coords_system[0] << "\t" << local_coords_system[1] << "\t" <<
00399         local_coords_system[2] << std::endl; std::cout << local_coords_system[3] << "\t" <<
00400         local_coords_system[4] << "\t" << local_coords_system[5] << std::endl; std::cout <<
00401         local_coords_system[6] << "\t" << local_coords_system[7] << "\t" << local_coords_system[8]
00402         << std::endl; std::cout << "Fitting degree is " << deg_out << " interp is " << local_interp
00403         << std::endl; std::cout << "Fitting coefficients are :" << std::endl; for(size_t
00404         kk=0;kk<local_coeffs.size();++kk) std::cout << local_coeffs[kk] << "\t"; std::cout <<
00405         std::endl;
00406         }*/
00407     }
00408 
00409     return error;
00410 }
00411 
00412 void usage()
00413 {
00414     std::cout << "usage: mpirun -np <number of processors> ./uref_hirec_test_parallel <prefix of "
00415                  "files> -str <first number> -end <lastnumber> -suffix <suffix> -interp <0=least "
00416                  "square, 1=interpolation> -dim <mesh dimension> -geom <s=sphere,t=torus>"
00417               << std::endl;
00418 }
00419 
00420 int main( int argc, char* argv[] )
00421 {
00422 #ifdef MOAB_HAVE_MPI
00423     MPI_Init( &argc, &argv );
00424     int nprocs, rank;
00425     MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
00426     MPI_Comm_rank( MPI_COMM_WORLD, &rank );
00427 #endif
00428     std::string prefix, suffix, istr, iend;
00429     int dim = 2, geom = 0;
00430     bool interp = false;
00431     ErrorCode error;
00432 
00433     if( argc == 1 )
00434     {
00435         usage();
00436         prefix = TestDir + "unittest/sphere_quads_5_ML_";
00437         suffix = ".h5m";
00438         istr   = "0";
00439         iend   = "1";
00440         std::cout << "Using defaults: MeshFiles/unittest/sphere_quads_5_ML_*.h5m -str 0 -end 1 "
00441                      "-suffix \".h5m\" -interp 0 -dim 2 -geom s"
00442                   << std::endl;
00443     }
00444     else
00445     {
00446         prefix      = std::string( argv[1] );
00447         bool hasdim = false;
00448 
00449         for( int i = 2; i < argc; ++i )
00450         {
00451             if( i + 1 != argc )
00452             {
00453                 if( std::string( argv[i] ) == "-suffix" )
00454                 {
00455                     suffix = std::string( argv[++i] );
00456                 }
00457                 else if( std::string( argv[i] ) == "-str" )
00458                 {
00459                     istr = std::string( argv[++i] );
00460                 }
00461                 else if( std::string( argv[i] ) == "-end" )
00462                 {
00463                     iend = std::string( argv[++i] );
00464                 }
00465                 else if( std::string( argv[i] ) == "-interp" )
00466                 {
00467                     interp = atoi( argv[++i] );
00468                 }
00469                 else if( std::string( argv[i] ) == "-dim" )
00470                 {
00471                     dim    = atoi( argv[++i] );
00472                     hasdim = true;
00473                 }
00474                 else if( std::string( argv[i] ) == "-geom" )
00475                 {
00476                     geom = std::string( argv[++i] ) == "s" ? 0 : 1;
00477                 }
00478                 else
00479                 {
00480 #ifdef MOAB_HAVE_MPI
00481 
00482                     if( 0 == rank )
00483                     {
00484                         usage();
00485                     }
00486 
00487                     MPI_Finalize();
00488 #else
00489                     usage();
00490 #endif
00491                     return 0;
00492                 }
00493             }
00494         }
00495 
00496         if( !hasdim )
00497         {
00498 #ifdef MOAB_HAVE_MPI
00499 
00500             if( 0 == rank )
00501             {
00502                 std::cout << "Dimension of input mesh should be provided, positive and less than 3" << std::endl;
00503             }
00504 
00505 #else
00506             std::cout << "Dimension of input mesh should be provided, positive and less than 3" << std::endl;
00507 #endif
00508             return 0;
00509         }
00510 
00511         if( dim > 2 || dim <= 0 )
00512         {
00513 #ifdef MOAB_HAVE_MPI
00514 
00515             if( 0 == rank )
00516             {
00517                 std::cout << "Input dimesion should be positive and less than 3;" << std::endl;
00518             }
00519 
00520 #else
00521             std::cout << "Input dimesion should be positive and less than 3;" << std::endl;
00522 #endif
00523             return 0;
00524         }
00525 
00526 #ifdef MOAB_HAVE_MPI
00527 
00528         if( 0 == rank )
00529         {
00530             std::cout << "Testing on " << prefix + istr + "-" + iend + suffix << " with dimension " << dim << "\n";
00531             std::string opts = interp ? "interpolation" : "least square fitting";
00532             std::cout << "High order reconstruction in " << opts << std::endl;
00533         }
00534 
00535 #else
00536         std::cout << "Testing on " << prefix + istr + "-" + iend + suffix << " with dimension " << dim << "\n";
00537         std::string opts = interp ? "interpolation" : "least square fitting";
00538         std::cout << "High order reconstruction in " << opts << std::endl;
00539 #endif
00540     }
00541 
00542     geomObject* obj;
00543 
00544     if( geom == 0 )
00545     {
00546         obj = new sphere();
00547     }
00548     else
00549     {
00550         obj = new torus();
00551     }
00552 
00553     std::vector< int > degs2fit;
00554     for( int d = 1; d <= 6; ++d )
00555     {
00556         degs2fit.push_back( d );
00557     }
00558 
00559     int begin = atoi( istr.c_str() ), end = atoi( iend.c_str() ) + 1;
00560 #ifdef MOAB_HAVE_MPI
00561     std::vector< std::vector< double > > geoml1errs_global, geoml2errs_global, geomlinferrs_global;
00562 
00563     if( rank == 0 )
00564     {
00565         geoml1errs_global =
00566             std::vector< std::vector< double > >( 1 + degs2fit.size(), std::vector< double >( end - begin, 0 ) );
00567         geoml2errs_global =
00568             std::vector< std::vector< double > >( 1 + degs2fit.size(), std::vector< double >( end - begin, 0 ) );
00569         geomlinferrs_global =
00570             std::vector< std::vector< double > >( 1 + degs2fit.size(), std::vector< double >( end - begin, 0 ) );
00571     }
00572 
00573 #else
00574     std::vector< std::vector< double > > geoml1errs_global( 1 + degs2fit.size(),
00575                                                             std::vector< double >( end - begin, 0 ) );
00576     std::vector< std::vector< double > > geoml2errs_global( 1 + degs2fit.size(),
00577                                                             std::vector< double >( end - begin, 0 ) );
00578     std::vector< std::vector< double > > geomlinferrs_global( 1 + degs2fit.size(),
00579                                                               std::vector< double >( end - begin, 0 ) );
00580 #endif
00581 
00582     for( int i = begin; i < end; ++i )
00583     {
00584         std::ostringstream convert;
00585         convert << i;
00586         std::string infile = prefix + convert.str() + suffix;
00587         int ntestverts;
00588         std::vector< double > geoml1errs, geoml2errs, geomlinferrs;
00589 
00590 #ifdef MOAB_HAVE_MPI
00591         std::cout << "Processor " << rank << " is working on file " << infile << std::endl;
00592 #endif
00593         error = closedsurface_uref_hirec_convergence_study( infile.c_str(), degs2fit, interp, dim, obj, ntestverts,
00594                                                             geoml1errs, geoml2errs, geomlinferrs );MB_CHK_ERR( error );
00595         assert( geoml1errs.size() == 1 + degs2fit.size() && geoml2errs.size() == 1 + degs2fit.size() &&
00596                 geomlinferrs.size() == 1 + degs2fit.size() );
00597 #ifdef MOAB_HAVE_MPI
00598 
00599         if( nprocs > 1 )
00600         {
00601             int ntestverts_global = 0;
00602             MPI_Reduce( &ntestverts, &ntestverts_global, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );
00603 
00604             for( size_t d = 0; d < degs2fit.size() + 1; ++d )
00605             {
00606                 double local_l1err = ntestverts * geoml1errs[d],
00607                        local_l2err = geoml2errs[d] * ( geoml2errs[d] * ntestverts ), local_linferr = geomlinferrs[d];
00608                 std::cout << "On Processor " << rank << " with mesh " << i
00609                           << " Degree = " << ( d == 0 ? 0 : degs2fit[d - 1] ) << " L1:" << geoml1errs[d]
00610                           << " L2:" << geoml2errs[d] << " Li:" << geomlinferrs[d] << std::endl;
00611                 double global_l1err = 0, global_l2err = 0, global_linferr = 0;
00612                 MPI_Reduce( &local_l1err, &global_l1err, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD );
00613                 MPI_Reduce( &local_l2err, &global_l2err, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD );
00614                 MPI_Reduce( &local_linferr, &global_linferr, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD );
00615 
00616                 if( rank == 0 )
00617                 {
00618                     geoml1errs_global[d][i - begin]   = global_l1err / ntestverts_global;
00619                     geoml2errs_global[d][i - begin]   = sqrt( global_l2err / ntestverts_global );
00620                     geomlinferrs_global[d][i - begin] = global_linferr;
00621                 }
00622             }
00623         }
00624         else
00625         {
00626             for( size_t d = 0; d < degs2fit.size() + 1; ++d )
00627             {
00628                 geoml1errs_global[d][i - begin]   = geoml1errs[d];
00629                 geoml2errs_global[d][i - begin]   = geoml2errs[d];
00630                 geomlinferrs_global[d][i - begin] = geomlinferrs[d];
00631             }
00632         }
00633 
00634 #else
00635 
00636         for( size_t d = 0; d < degs2fit.size() + 1; ++d )
00637         {
00638             geoml1errs_global[d][i - begin]   = geoml1errs[d];
00639             geoml2errs_global[d][i - begin]   = geoml2errs[d];
00640             geomlinferrs_global[d][i - begin] = geomlinferrs[d];
00641         }
00642 
00643 #endif
00644     }
00645 
00646 #ifdef MOAB_HAVE_MPI
00647 
00648     if( rank == 0 )
00649     {
00650         std::cout << "Degrees: 0 ";
00651 
00652         for( size_t ideg = 0; ideg < degs2fit.size(); ++ideg )
00653         {
00654             std::cout << degs2fit[ideg] << " ";
00655         }
00656 
00657         std::cout << std::endl;
00658         std::cout << "L1-norm error: \n";
00659 
00660         for( size_t i = 0; i < geoml1errs_global.size(); ++i )
00661         {
00662             for( size_t j = 0; j < geoml1errs_global[i].size(); ++j )
00663             {
00664                 std::cout << geoml1errs_global[i][j] << " ";
00665             }
00666 
00667             std::cout << std::endl;
00668         }
00669 
00670         std::cout << "L2-norm error: \n";
00671 
00672         for( size_t i = 0; i < geoml2errs_global.size(); ++i )
00673         {
00674             for( size_t j = 0; j < geoml2errs_global[i].size(); ++j )
00675             {
00676                 std::cout << geoml2errs_global[i][j] << " ";
00677             }
00678 
00679             std::cout << std::endl;
00680         }
00681 
00682         std::cout << "Linf-norm error: \n";
00683 
00684         for( size_t i = 0; i < geomlinferrs_global.size(); ++i )
00685         {
00686             for( size_t j = 0; j < geomlinferrs_global[i].size(); ++j )
00687             {
00688                 std::cout << geomlinferrs_global[i][j] << " ";
00689             }
00690 
00691             std::cout << std::endl;
00692         }
00693     }
00694 
00695 #else
00696     std::cout << "Degrees: 0 ";
00697 
00698     for( size_t ideg = 0; ideg < degs2fit.size(); ++ideg )
00699     {
00700         std::cout << degs2fit[ideg] << " ";
00701     }
00702 
00703     std::cout << std::endl;
00704     std::cout << "L1-norm error: \n";
00705 
00706     for( size_t i = 0; i < geoml1errs_global.size(); ++i )
00707     {
00708         for( size_t j = 0; j < geoml1errs_global[i].size(); ++j )
00709         {
00710             std::cout << geoml1errs_global[i][j] << " ";
00711         }
00712 
00713         std::cout << std::endl;
00714     }
00715 
00716     std::cout << "L2-norm error: \n";
00717 
00718     for( size_t i = 0; i < geoml2errs_global.size(); ++i )
00719     {
00720         for( size_t j = 0; j < geoml2errs_global[i].size(); ++j )
00721         {
00722             std::cout << geoml2errs_global[i][j] << " ";
00723         }
00724 
00725         std::cout << std::endl;
00726     }
00727 
00728     std::cout << "Linf-norm error: \n";
00729 
00730     for( size_t i = 0; i < geomlinferrs_global.size(); ++i )
00731     {
00732         for( size_t j = 0; j < geomlinferrs_global[i].size(); ++j )
00733         {
00734             std::cout << geomlinferrs_global[i][j] << " ";
00735         }
00736 
00737         std::cout << std::endl;
00738     }
00739 
00740 #endif
00741 #ifdef MOAB_HAVE_MPI
00742     MPI_Finalize();
00743 #endif
00744 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines