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 diachin2@llnl.gov, djmelan@sandia.gov, mbrewer@sandia.gov, 00024 pknupp@sandia.gov, tleurent@mcs.anl.gov, tmunson@mcs.anl.gov 00025 00026 ***************************************************************** */ 00027 /*! 00028 \file ParallelHelper.hpp 00029 \brief 00030 00031 Implements ParallelHelper Class 00032 00033 \author Martin Isenburg 00034 \date 2008-03-04 00035 */ 00036 00037 #ifndef Mesquite_ParallelHelper_hpp 00038 #define Mesquite_ParallelHelper_hpp 00039 00040 #include "ParallelHelperInterface.hpp" 00041 00042 #include <vector> 00043 #include <map> 00044 00045 namespace MBMesquite 00046 { 00047 typedef struct VertexIdMapKey 00048 { 00049 size_t glob_id; 00050 int proc_id; 00051 } VertexIdMapKey; 00052 00053 struct VertexIdLessFunc 00054 { 00055 bool operator()( const VertexIdMapKey& that1, const VertexIdMapKey& that2 ) const 00056 { 00057 return ( ( that1.proc_id < that2.proc_id ) || 00058 ( ( that1.proc_id == that2.proc_id ) && ( that1.glob_id < that2.glob_id ) ) ); 00059 } 00060 }; 00061 00062 typedef std::map< VertexIdMapKey, int, VertexIdLessFunc > VertexIdMap; 00063 00064 int get_parallel_rank(); 00065 int get_parallel_size(); 00066 double reduce_parallel_max( double value ); 00067 void parallel_barrier(); 00068 00069 class ParallelHelperImpl : public ParallelHelper 00070 { 00071 public: 00072 enum CommunicationModel 00073 { 00074 TrulyNonBlocking = 0, 00075 TrulyNonBlockingAvoidAllReduce, 00076 NonBlocking, 00077 NonBlockingAvoidAllReduce, 00078 Blocking, 00079 BlockingAvoidAllReduce 00080 }; 00081 00082 ParallelHelperImpl(); 00083 ~ParallelHelperImpl(); 00084 00085 // function called by application during set-up 00086 void set_parallel_mesh( ParallelMesh* mesh ); 00087 void set_communicator( size_t comm ); 00088 void set_communicator( const void* comm ) 00089 { 00090 set_communicator( reinterpret_cast< size_t >( comm ) ); 00091 } 00092 void set_communication_model( int model, MsqError& err ); 00093 void set_generate_random_numbers( int grn, MsqError& err ); 00094 00095 protected: 00096 friend class VertexMover; 00097 // functions called by VertexMover::loop_over_mesh() 00098 void smoothing_init( MsqError& err ); 00099 void compute_first_independent_set( std::vector< Mesh::VertexHandle >& fixed_vertices ); 00100 void communicate_first_independent_set( MsqError& err ); 00101 bool compute_next_independent_set(); 00102 bool get_next_partition_boundary_vertex( MBMesquite::Mesh::VertexHandle& vertex_handle ); 00103 void communicate_next_independent_set( MsqError& err ); 00104 void smoothing_close( MsqError& err ); 00105 00106 protected: 00107 friend class QualityAssessor; 00108 // functions called by QualityAssessor::loop_over_mesh() 00109 int get_rank() const; 00110 int get_nprocs() const; 00111 bool is_our_element( Mesh::ElementHandle element_handle, MsqError& err ) const; 00112 bool is_our_vertex( Mesh::VertexHandle vertex_handle, MsqError& err ) const; 00113 void communicate_min_max_to_all( double* minimum, double* maximum, MsqError& ) const; 00114 void communicate_min_max_to_zero( double* minimum, double* maximum, MsqError& ) const; 00115 void communicate_sums_to_zero( size_t* freeElementCount, 00116 int* invertedElementCount, 00117 size_t* elementCount, 00118 int* invertedSampleCount, 00119 size_t* sampleCount, 00120 long unsigned int* count, 00121 long unsigned int* invalid, 00122 double* sum, 00123 double* sqrSum, 00124 MsqError& ) const; 00125 void communicate_power_sum_to_zero( double* pMean, MsqError& ) const; 00126 void communicate_histogram_to_zero( std::vector< int >& histogram, MsqError& ) const; 00127 void communicate_any_true( bool& value, MsqError& err ) const; 00128 void communicate_all_true( bool& value, MsqError& err ) const; 00129 00130 private: 00131 ParallelMesh* mesh; 00132 00133 size_t communicator; 00134 int communication_model; 00135 00136 int rank; 00137 int nprocs; 00138 00139 // variables for VertexMover::loop_over_mesh() 00140 int generate_random_numbers; 00141 std::vector< MBMesquite::Mesh::VertexHandle > vertices; 00142 int num_vertex; 00143 std::vector< char > vtx_in_partition_boundary; 00144 int num_vtx_partition_boundary; 00145 int num_vtx_partition_boundary_local; 00146 int num_vtx_partition_boundary_remote; 00147 std::vector< MBMesquite::Mesh::VertexHandle > part_vertices; 00148 std::vector< int > part_proc_owner; 00149 std::vector< size_t > part_gid; 00150 std::vector< int > part_smoothed_flag; 00151 std::vector< double > part_rand_number; 00152 int num_exportVtx; 00153 std::vector< size_t > exportVtxGIDs; 00154 std::vector< int > exportVtxLIDs; 00155 std::vector< int > exportProc; 00156 std::vector< bool > in_independent_set; 00157 VertexIdMap vid_map; 00158 int total_num_vertices_to_smooth; 00159 int total_num_vertices_to_recv; 00160 std::vector< int > neighbourProcSend; 00161 std::vector< int > neighbourProcRecv; 00162 std::vector< int > neighbourProcSendRemain; 00163 std::vector< int > neighbourProcRecvRemain; 00164 int num_already_smoothed_vertices; 00165 int num_already_recv_vertices; 00166 std::vector< std::vector< int > > vtx_off_proc_list; 00167 std::vector< int > neighbourProc; 00168 int iteration; 00169 int global_work_remains; 00170 int next_vtx_partition_boundary; 00171 /* for exchanging unused ghost node information */ 00172 int unghost_num_vtx; 00173 std::vector< MBMesquite::Mesh::VertexHandle > unghost_vertices; 00174 int unghost_num_procs; 00175 std::vector< int > unghost_procs; 00176 std::vector< int > unghost_procs_num_vtx; 00177 std::vector< int > unghost_procs_offset; 00178 int update_num_vtx; 00179 std::vector< size_t > update_gid; 00180 int update_num_procs; 00181 std::vector< int > update_procs; 00182 std::vector< int > update_procs_num_vtx; 00183 std::vector< int > update_procs_offset; 00184 00185 // functions for VertexMover::loop_over_mesh() 00186 void compute_independent_set(); 00187 int comm_smoothed_vtx_b( MsqError& err ); 00188 int comm_smoothed_vtx_b_no_all( MsqError& err ); 00189 int comm_smoothed_vtx_nb( MsqError& err ); 00190 int comm_smoothed_vtx_nb_no_all( MsqError& err ); 00191 int comm_smoothed_vtx_tnb( MsqError& err ); 00192 int comm_smoothed_vtx_tnb_no_all( MsqError& err ); 00193 }; 00194 00195 } // namespace MBMesquite 00196 #endif // Mesquite_ParallelHelper_hpp