MOAB: Mesh Oriented datABase  (version 5.4.1)
par_hex_smooth_laplace.cpp
Go to the documentation of this file.
00001 /* *****************************************************************
00002     MESQUITE -- The Mesh Quality Improvement Toolkit
00003 
00004     Copyright 2004 Sandia Corporation and Argonne National
00005     Laboratory.  Under the terms of Contract DE-AC04-94AL85000
00006     with Sandia Corporation, the U.S. Government retains certain
00007     rights in this software.
00008 
00009     This library is free software; you can redistribute it and/or
00010     modify it under the terms of the GNU Lesser General Public
00011     License as published by the Free Software Foundation; either
00012     version 2.1 of the License, or (at your option) any later version.
00013 
00014     This library is distributed in the hope that it will be useful,
00015     but WITHOUT ANY WARRANTY; without even the implied warranty of
00016     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00017     Lesser General Public License for more details.
00018 
00019     You should have received a copy of the GNU Lesser General Public License
00020     (lgpl.txt) along with this library; if not, write to the Free Software
00021     Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00022 
00023     [email protected], [email protected], [email protected],
00024     [email protected], [email protected], [email protected]
00025 
00026   ***************************************************************** */
00027 // -*- Mode : c++; tab-width: 2; c-tab-always-indent: t; indent-tabs-mode: nil; c-basic-offset: 2
00028 // -*-
00029 //
00030 //   SUMMARY:
00031 //     USAGE:
00032 //
00033 // ORIG-DATE: 24-Jan-12
00034 //  LAST-MOD: 26-Jan-1 by Stephen Kennon
00035 //
00036 //
00037 // DESCRIPTION:
00038 // ============
00039 /*! \file par_hex.cpp
00040 
00041 A test of MBMesquite's parallel capabilities.  Reads a split vtk file, smooths in parallel using
00042 Laplace smoothing, writes out the result (which can be compared with the "gold" copy of the same
00043 name in the meshFiles VTK directory).
00044 
00045 See the MBMesquite User's Guide, section "Using MBMesquite in Parallel" - this code is very similar
00046 to the example code shown therein.
00047 
00048  */
00049 // DESCRIP-END.
00050 //
00051 
00052 #include "TestUtil.hpp"
00053 #include "MeshImpl.hpp"
00054 #include "MeshUtil.hpp"
00055 #include "MsqTimer.hpp"
00056 #include "MsqDebug.hpp"
00057 #include "Mesquite.hpp"
00058 #include "MsqError.hpp"
00059 #include "Vector3D.hpp"
00060 #include "InstructionQueue.hpp"
00061 #include "LaplaceWrapper.hpp"
00062 #include "PatchData.hpp"
00063 #include "TerminationCriterion.hpp"
00064 #include "QualityAssessor.hpp"
00065 
00066 /* MBMesquite includes */
00067 #include "ParallelMeshImpl.hpp"
00068 #include "ParallelHelper.hpp"
00069 
00070 // algorithms
00071 #include "Randomize.hpp"
00072 #include "ConditionNumberQualityMetric.hpp"
00073 #include "UntangleBetaQualityMetric.hpp"
00074 #include "LPtoPTemplate.hpp"
00075 #include "LInfTemplate.hpp"
00076 #include "SteepestDescent.hpp"
00077 #include "ConjugateGradient.hpp"
00078 #include "PlanarDomain.hpp"
00079 
00080 #include <iostream>
00081 using std::cout;
00082 using std::endl;
00083 #include <cstdlib>
00084 
00085 #include <mpi.h>
00086 #include <sstream>
00087 
00088 using namespace MBMesquite;
00089 
00090 std::string VTK_3D_DIR = TestDir + "unittest/mesquite/3D/vtk/hexes/tangled/";
00091 
00092 using namespace std;
00093 
00094 int main( int argc, char* argv[] )
00095 {
00096     /* init MPI */
00097     int rank, nprocs;
00098     if( MPI_SUCCESS != MPI_Init( &argc, &argv ) )
00099     {
00100         cerr << "MPI_Init failed." << endl;
00101         return 2;
00102     }
00103     MPI_Comm_rank( MPI_COMM_WORLD, &rank );
00104     MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
00105 
00106     if( nprocs > 2 )
00107     {
00108         cerr << "parallel_laplace_smooth test can only be run with 1 or 2 processors" << std::endl;
00109         return 0;
00110     }
00111 
00112     const int debug = 0;  // 1,2,3 for more debug info
00113     if( debug )
00114     {
00115         MsqDebug::enable( 1 );
00116         if( debug > 1 ) MsqDebug::enable( 2 );
00117         if( debug > 2 ) MsqDebug::enable( 3 );
00118     }
00119 
00120     /* create processor-specific file names */
00121     ostringstream in_name, out_name, gold_name;
00122     in_name << VTK_3D_DIR << "par_original_hex_mesh." << nprocs << "." << rank << ".vtk";
00123     gold_name << VTK_3D_DIR << "par_smoothed_hex_mesh." << nprocs << "." << rank << ".vtk";
00124     out_name << "par_smoothed_hex_mesh." << nprocs << "." << rank << ".vtk";
00125 
00126     /* load different mesh files on each processor */
00127     MBMesquite::MsqError err;
00128     MBMesquite::MeshImpl mesh;
00129     mesh.read_vtk( in_name.str().c_str(), err );
00130     if( err )
00131     {
00132         cerr << err << endl;
00133         return 1;
00134     }
00135 
00136     /* create parallel mesh instance, specifying tags
00137      * containing parallel data */
00138     MBMesquite::ParallelMeshImpl parallel_mesh( &mesh, "GLOBAL_ID", "PROCESSOR_ID" );
00139     MBMesquite::ParallelHelperImpl helper;
00140     helper.set_communicator( MPI_COMM_WORLD );
00141     helper.set_parallel_mesh( &parallel_mesh );
00142     parallel_mesh.set_parallel_helper( &helper );
00143 
00144     /* do Laplacian smooth */
00145     LaplaceWrapper optimizer;
00146     optimizer.set_vertex_movement_limit_factor( 1.e-10 );
00147     optimizer.set_iteration_limit( 2000 );
00148     optimizer.enable_culling( false );
00149     optimizer.run_instructions( &parallel_mesh, err );
00150     if( err )
00151     {
00152         cerr << err << endl;
00153         return 1;
00154     }
00155 
00156     /* write mesh */
00157     mesh.write_vtk( out_name.str().c_str(), err );
00158     if( err )
00159     {
00160         cerr << err << endl;
00161         return 1;
00162     }
00163 
00164     /* compare mesh with gold copy */
00165     MeshImpl gold;
00166     gold.read_vtk( gold_name.str().c_str(), err );
00167     if( err )
00168     {
00169         cerr << err << endl;
00170         return 1;
00171     }
00172 
00173     bool do_print = true;
00174     double tol    = 1.e-4;
00175     bool diff     = MeshUtil::meshes_are_different( mesh, gold, err, tol, do_print );
00176     if( err )
00177     {
00178         cerr << err << endl;
00179         return 1;
00180     }
00181 
00182     if( diff )
00183     {
00184         cerr << "Error, computed mesh different from gold copy" << std::endl;
00185         return 1;
00186     }
00187 
00188     print_timing_diagnostics( cout );
00189 
00190     MPI_Finalize();
00191     return 0;
00192 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines