MOAB: Mesh Oriented datABase
(version 5.4.1)
|
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: 19-Feb-02 at 10:57:52 00034 // LAST-MOD: 23-Jul-03 at 18:04:37 by Thomas Leurent 00035 // 00036 // 00037 // DESCRIPTION: 00038 // ============ 00039 /*! \file main.cpp 00040 00041 describe main.cpp here 00042 00043 */ 00044 // DESCRIP-END. 00045 // 00046 00047 #include <iostream> 00048 using std::cout; 00049 using std::endl; 00050 #include <cstdlib> 00051 00052 #include "Mesquite.hpp" 00053 #include "MsqError.hpp" 00054 #include "MeshImpl.hpp" 00055 #include "Vector3D.hpp" 00056 #include "InstructionQueue.hpp" 00057 #include "PatchData.hpp" 00058 #include "TerminationCriterion.hpp" 00059 #include "QualityAssessor.hpp" 00060 #include "PlanarDomain.hpp" 00061 #include "MsqTimer.hpp" 00062 #include "TestUtil.hpp" 00063 00064 // algorythms 00065 #include "ConditionNumberQualityMetric.hpp" 00066 #include "LInfTemplate.hpp" 00067 #include "SteepestDescent.hpp" 00068 #include "LaplacianSmoother.hpp" 00069 #include "EdgeLengthQualityMetric.hpp" 00070 using namespace MBMesquite; 00071 00072 std::string DEFAULT_INPUT = TestDir + "unittest/mesquite/2D/vtk/quads/untangled/square_quad_2.vtk"; 00073 00074 void help( const char* argv0 ) 00075 { 00076 std::cerr << "Usage: " << argv0 << " [<input_file>] [<output_file>]" << std::endl 00077 << " default input file is: " << DEFAULT_INPUT << std::endl 00078 << " defualt is no output file" << std::endl 00079 << " Warning: input mesh is assumed to lie in Z=5 plane" << std::endl; 00080 exit( 1 ); 00081 } 00082 00083 int main( int argc, char* argv[] ) 00084 { 00085 const char* input_file = DEFAULT_INPUT.c_str(); 00086 const char* output_file = NULL; 00087 switch( argc ) 00088 { 00089 default: 00090 help( argv[0] ); 00091 case 3: 00092 if( !strcmp( argv[2], "-h" ) ) help( argv[0] ); 00093 output_file = argv[2]; 00094 case 2: 00095 if( !strcmp( argv[1], "-h" ) ) help( argv[0] ); 00096 input_file = argv[1]; 00097 case 1:; 00098 } 00099 00100 /* Read a VTK Mesh file */ 00101 MsqPrintError err( cout ); 00102 MBMesquite::MeshImpl mesh; 00103 mesh.read_vtk( input_file, err ); 00104 if( err ) return 1; 00105 00106 // creates an intruction queue 00107 InstructionQueue queue1; 00108 00109 // creates a mean ratio quality metric ... 00110 ConditionNumberQualityMetric shape_metric; 00111 EdgeLengthQualityMetric lapl_met; 00112 lapl_met.set_averaging_method( QualityMetric::RMS ); 00113 00114 // creates the laplacian smoother procedures 00115 LaplacianSmoother lapl1; 00116 QualityAssessor stop_qa = QualityAssessor( &shape_metric ); 00117 stop_qa.add_quality_assessment( &lapl_met ); 00118 00119 //**************Set stopping criterion**************** 00120 TerminationCriterion sc2; 00121 sc2.add_iteration_limit( 10 ); 00122 if( err ) return 1; 00123 lapl1.set_outer_termination_criterion( &sc2 ); 00124 00125 // adds 1 pass of pass1 to mesh_set1 00126 queue1.add_quality_assessor( &stop_qa, err ); 00127 if( err ) return 1; 00128 queue1.set_master_quality_improver( &lapl1, err ); 00129 if( err ) return 1; 00130 queue1.add_quality_assessor( &stop_qa, err ); 00131 if( err ) return 1; 00132 // adds 1 passes of pass2 to mesh_set1 00133 // mesh_set1.add_quality_pass(pass2); 00134 00135 // writeVtkMesh("original_mesh", mesh, err); MSQ_CHKERR(err); 00136 00137 PlanarDomain plane( Vector3D( 0, 0, 1 ), Vector3D( 0, 0, 5 ) ); 00138 00139 // launches optimization on mesh_set1 00140 MeshDomainAssoc mesh_and_domain = MeshDomainAssoc( &mesh, &plane ); 00141 Timer t; 00142 queue1.run_instructions( &mesh_and_domain, err ); 00143 if( err ) return 1; 00144 double secs = t.since_birth(); 00145 std::cout << "Optimization completed in " << secs << " seconds" << std::endl; 00146 00147 if( output_file ) 00148 { 00149 mesh.write_vtk( output_file, err ); 00150 if( err ) return 1; 00151 std::cout << "Wrote file: " << output_file << std::endl; 00152 } 00153 00154 // check that smoother is working: 00155 // the one free vertex must be at the origin 00156 if( !DEFAULT_INPUT.compare( input_file ) ) 00157 { 00158 std::vector< Mesh::VertexHandle > vertices; 00159 mesh.get_all_vertices( vertices, err ); 00160 if( err ) return 1; 00161 00162 std::vector< bool > fixed_flags; 00163 mesh.vertices_get_fixed_flag( arrptr( vertices ), fixed_flags, vertices.size(), err ); 00164 if( err ) return 1; 00165 00166 // find one free vertex 00167 int idx = -1; 00168 for( unsigned i = 0; i < vertices.size(); ++i ) 00169 { 00170 if( fixed_flags[i] == true ) continue; 00171 if( idx != -1 ) 00172 { 00173 std::cerr << "Multiple free vertices in mesh." << std::endl; 00174 return 1; 00175 } 00176 idx = i; 00177 } 00178 00179 if( idx == -1 ) 00180 { 00181 std::cerr << "No free vertex in mesh!!!!!" << std::endl; 00182 return 1; 00183 } 00184 00185 Mesh::VertexHandle vertex = vertices[idx]; 00186 MsqVertex coords; 00187 mesh.vertices_get_coordinates( &vertex, &coords, 1, err ); 00188 if( err ) return 1; 00189 00190 // calculate distance from origin 00191 double dist = sqrt( coords[0] * coords[0] + coords[1] * coords[1] ); 00192 if( dist > 1e-8 ) 00193 { 00194 std::cerr << "Free vertex not at origin after Laplace smooth." << std::endl 00195 << "Expected location: (0,0)" << std::endl 00196 << "Actual location: (" << coords[0] << "," << coords[1] << ")" << std::endl; 00197 return 2; 00198 } 00199 } 00200 00201 return 0; 00202 }