Mesh Oriented datABase  (version 5.6.0)
An array-based unstructured mesh library
CoupleMGen.cpp File Reference
+ Include dependency graph for CoupleMGen.cpp:

Go to the source code of this file.

Functions

double physField (double x, double y, double z)
 
int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
CoupleMGen.cpp.

Definition at line 78 of file CoupleMGen.cpp.

79 {
80  int proc_id = 0, size = 1;
81 
82  MPI_Init( &argc, &argv );
83  MPI_Comm_rank( MPI_COMM_WORLD, &proc_id );
84  MPI_Comm_size( MPI_COMM_WORLD, &size );
85 
86  Core mcore;
87  Interface* mb = &mcore;
88  EntityHandle fileset1, fileset2; // for 2 different meshes
90  // default options
91  opts.A = opts.B = opts.C = 1;
92  opts.M = opts.N = opts.K = 1;
93  opts.blockSize = 4;
94  opts.xsize = opts.ysize = opts.zsize = 1.;
95  opts.ui = CartVect( 1., 0, 0. );
96  opts.uj = CartVect( 0., 1., 0. );
97  opts.uk = CartVect( 0., 0., 1. );
98  opts.newMergeMethod = opts.quadratic = opts.keep_skins = opts.tetra = false;
99  opts.adjEnts = opts.parmerge = false;
100  opts.GL = 0;
101 
102  ProgOptions popts;
103 
104  popts.addOpt< int >( string( "blockSize,b" ), string( "Block size of mesh (default=4)" ), &opts.blockSize );
105  popts.addOpt< int >( string( "xproc,M" ), string( "Number of processors in x dir (default=1)" ), &opts.M );
106  popts.addOpt< int >( string( "yproc,N" ), string( "Number of processors in y dir (default=1)" ), &opts.N );
107  popts.addOpt< int >( string( "zproc,K" ), string( "Number of processors in z dir (default=1)" ), &opts.K );
108 
109  popts.addOpt< int >( string( "xblocks,A" ), string( "Number of blocks on a task in x dir (default=2)" ), &opts.A );
110  popts.addOpt< int >( string( "yblocks,B" ), string( "Number of blocks on a task in y dir (default=2)" ), &opts.B );
111  popts.addOpt< int >( string( "zblocks,C" ), string( "Number of blocks on a task in x dir (default=2)" ), &opts.C );
112 
113  popts.addOpt< double >( string( "xsize,x" ), string( "Total size in x direction (default=1.)" ), &opts.xsize );
114  popts.addOpt< double >( string( "ysize,y" ), string( "Total size in y direction (default=1.)" ), &opts.ysize );
115  popts.addOpt< double >( string( "zsize,z" ), string( "Total size in z direction (default=1.)" ), &opts.zsize );
116 
117  popts.addOpt< void >( "newMerge,w", "use new merging method", &opts.newMergeMethod );
118 
119  popts.addOpt< void >( "quadratic,q", "use hex 27 elements", &opts.quadratic );
120 
121  popts.addOpt< void >( "keep_skins,k", "keep skins with shared entities", &opts.keep_skins );
122 
123  popts.addOpt< void >( "tetrahedrons,t", "generate tetrahedrons", &opts.tetra );
124 
125  popts.addOpt< void >( "faces_edges,f", "create all faces and edges", &opts.adjEnts );
126 
127  popts.addOpt< int >( string( "ghost_layers,g" ), string( "Number of ghost layers (default=0)" ), &opts.GL );
128 
129  popts.addOpt< void >( "parallel_merge,p", "use parallel mesh merge, not vertex ID based merge", &opts.parmerge );
130 
131  Coupler::Method method = Coupler::LINEAR_FE;
132 
133  double toler = 1.e-6;
134  popts.addOpt< double >( string( "eps,e" ), string( "tolerance for coupling, used in locating points" ), &toler );
135 
136  bool writeMeshes = false;
137  popts.addOpt< void >( "print,p", "write meshes", &writeMeshes );
138 
139  popts.parseCommandLine( argc, argv );
140 
141  double start_time = MPI_Wtime();
142 
143  MB_CHK_ERR( mb->create_meshset( MESHSET_SET, fileset1 ) );
144  MB_CHK_ERR( mb->create_meshset( MESHSET_SET, fileset2 ) );
145 
146  ParallelComm* pc1 = new ParallelComm( mb, MPI_COMM_WORLD );
147  MeshGeneration* mgen1 = new MeshGeneration( mb, pc1, fileset1 );
148 
149  MB_CHK_ERR( mgen1->BrickInstance( opts ) ); // this will generate first mesh on fileset1
150 
151  double instance_time = MPI_Wtime();
152  double current = instance_time;
153  if( !proc_id ) std::cout << " instantiate first mesh " << instance_time - start_time << "\n";
154  // set an interpolation tag on source mesh, from phys field
155  std::string interpTag( "interp_tag" );
156  Tag tag;
157  MB_CHK_ERR( mb->tag_get_handle( interpTag.c_str(), 1, MB_TYPE_DOUBLE, tag, MB_TAG_CREAT | MB_TAG_DENSE ) );
158 
159  Range src_elems;
160  MB_CHK_ERR( pc1->get_part_entities( src_elems, 3 ) );
161  Range src_verts;
162  MB_CHK_ERR( mb->get_connectivity( src_elems, src_verts ) );
163  for( Range::iterator vit = src_verts.begin(); vit != src_verts.end(); ++vit )
164  {
165  EntityHandle vert = *vit; //?
166 
167  double vertPos[3];
168  mb->get_coords( &vert, 1, vertPos );
169 
170  double fieldValue = physField( vertPos[0], vertPos[1], vertPos[2] );
171 
172  MB_CHK_ERR( mb->tag_set_data( tag, &vert, 1, &fieldValue ) );
173  }
174 
175  double setTag_time = MPI_Wtime();
176  if( !proc_id ) std::cout << " set tag " << setTag_time - current;
177  current = instance_time;
178  // change some options, so it is a different mesh
179  int tmp1 = opts.K;
180  opts.K = opts.M;
181  opts.M = tmp1; // swap (opts.K, opts.M)
182  opts.tetra = !opts.tetra;
183  opts.blockSize++;
184 
185  ParallelComm* pc2 = new ParallelComm( mb, MPI_COMM_WORLD );
186  MeshGeneration* mgen2 = new MeshGeneration( mb, pc2, fileset2 );
187 
188  MB_CHK_ERR( mgen2->BrickInstance( opts ) ); // this will generate second mesh on fileset2
189 
190  double instance_second = MPI_Wtime();
191  if( !proc_id ) std::cout << " instance second mesh" << instance_second - current << "\n";
192  current = instance_second;
193 
194  // test the sets are fine
195  if( writeMeshes )
196  {
197  MB_CHK_SET_ERR( mb->write_file( "mesh1.h5m", 0, ";;PARALLEL=WRITE_PART;CPUTIME;PARALLEL_COMM=0;", &fileset1,
198  1 ),
199  "Can't write in parallel mesh 1" );
200  MB_CHK_SET_ERR( mb->write_file( "mesh2.h5m", 0, ";;PARALLEL=WRITE_PART;CPUTIME;PARALLEL_COMM=1;", &fileset2,
201  1 ),
202  "Can't write in parallel mesh 1" );
203  double write_files = MPI_Wtime();
204  if( !proc_id ) std::cout << " write files " << write_files - current << "\n";
205  current = write_files;
206  }
207 
208  // Instantiate a coupler, which also initializes the tree
209  Coupler mbc( mb, pc1, src_elems, 0 );
210 
211  double instancecoupler = MPI_Wtime();
212  if( !proc_id ) std::cout << " instance coupler " << instancecoupler - current << "\n";
213  current = instancecoupler;
214 
215  // Get points from the target mesh to interpolate
216  // We have to treat differently the case when the target is a spectral mesh
217  // In that case, the points of interest are the GL points, not the vertex nodes
218  std::vector< double > vpos; // This will have the positions we are interested in
219  int numPointsOfInterest = 0;
220 
221  Range targ_elems;
222  Range targ_verts;
223 
224  // First get all vertices adj to partition entities in target mesh
225  MB_CHK_ERR( pc2->get_part_entities( targ_elems, 3 ) );
226 
227  MB_CHK_ERR( mb->get_adjacencies( targ_elems, 0, false, targ_verts, Interface::UNION ) );
228  Range tmp_verts;
229  // Then get non-owned verts and subtract
230  MB_CHK_ERR( pc2->get_pstatus_entities( 0, PSTATUS_NOT_OWNED, tmp_verts ) );
231  targ_verts = subtract( targ_verts, tmp_verts );
232  // get position of these entities; these are the target points
233  numPointsOfInterest = (int)targ_verts.size();
234  vpos.resize( 3 * targ_verts.size() );
235  MB_CHK_ERR( mb->get_coords( targ_verts, &vpos[0] ) );
236  // Locate those points in the source mesh
237  // std::cout<<"rank "<< proc_id<< " points of interest: " << numPointsOfInterest << "\n";
238  MB_CHK_ERR( mbc.locate_points( &vpos[0], numPointsOfInterest, 0, toler ) );
239 
240  double locatetime = MPI_Wtime();
241  if( !proc_id ) std::cout << " locate points: " << locatetime - current << "\n";
242  current = locatetime;
243 
244  // Now interpolate tag onto target points
245  std::vector< double > field( numPointsOfInterest );
246 
247  MB_CHK_ERR( mbc.interpolate( method, interpTag, &field[0] ) );
248 
249  // compare with the actual phys field
250  double err_max = 0;
251  for( int i = 0; i < numPointsOfInterest; i++ )
252  {
253  double trval = physField( vpos[3 * i], vpos[3 * i + 1], vpos[3 * i + 2] );
254  double err2 = fabs( trval - field[i] );
255  if( err2 > err_max ) err_max = err2;
256  }
257 
258  double interpolateTime = MPI_Wtime();
259  if( !proc_id ) std::cout << " interpolate points: " << interpolateTime - current << "\n";
260  current = interpolateTime;
261 
262  double gerr;
263  MPI_Allreduce( &err_max, &gerr, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD );
264  if( 0 == proc_id ) std::cout << "max err " << gerr << "\n";
265 
266  delete mgen1;
267  delete mgen2;
268 
269  MPI_Finalize();
270 
271  return 0;
272 }

References moab::MeshGeneration::BrickOpts::A, ProgOptions::addOpt(), moab::MeshGeneration::BrickOpts::adjEnts, moab::MeshGeneration::BrickOpts::B, moab::Range::begin(), moab::MeshGeneration::BrickOpts::blockSize, moab::MeshGeneration::BrickInstance(), moab::MeshGeneration::BrickOpts::C, moab::Core::create_meshset(), moab::Range::end(), moab::Core::get_adjacencies(), moab::Core::get_connectivity(), moab::Core::get_coords(), moab::ParallelComm::get_part_entities(), moab::ParallelComm::get_pstatus_entities(), moab::MeshGeneration::BrickOpts::GL, moab::Coupler::interpolate(), moab::MeshGeneration::BrickOpts::K, moab::MeshGeneration::BrickOpts::keep_skins, moab::Coupler::LINEAR_FE, moab::Coupler::locate_points(), moab::MeshGeneration::BrickOpts::M, mb, MB_CHK_ERR, MB_CHK_SET_ERR, MB_TAG_CREAT, MB_TAG_DENSE, MB_TYPE_DOUBLE, MESHSET_SET, moab::MeshGeneration::BrickOpts::N, moab::MeshGeneration::BrickOpts::newMergeMethod, moab::MeshGeneration::BrickOpts::parmerge, ProgOptions::parseCommandLine(), moab::physField(), PSTATUS_NOT_OWNED, moab::MeshGeneration::BrickOpts::quadratic, moab::Range::size(), start_time, moab::subtract(), moab::Core::tag_get_handle(), moab::Core::tag_set_data(), moab::MeshGeneration::BrickOpts::tetra, moab::MeshGeneration::BrickOpts::ui, moab::MeshGeneration::BrickOpts::uj, moab::MeshGeneration::BrickOpts::uk, moab::Interface::UNION, moab::Core::write_file(), moab::MeshGeneration::BrickOpts::xsize, moab::MeshGeneration::BrickOpts::ysize, and moab::MeshGeneration::BrickOpts::zsize.

◆ physField()

double physField ( double  x,
double  y,
double  z 
)

Definition at line 70 of file CoupleMGen.cpp.

71 {
72 
73  double out = sin( M_PI * x ) * cos( M_PI * y ) * sin( M_PI * z );
74 
75  return out;
76 }