Loading [MathJax]/extensions/tex2jax.js
Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
visuMapVtk.cpp File Reference
#include "moab/MOABConfig.h"
#include "moab/ProgOptions.hpp"
#include "moab/Core.hpp"
#include "moab/Range.hpp"
#include "netcdf.h"
#include <cmath>
#include <sstream>
#include <map>
#include <Eigen/Sparse>
+ Include dependency graph for visuMapVtk.cpp:

Go to the source code of this file.

Macros

#define ERR_NC(e)
 
#define GET_DIM1(ncdim, name, val)
 
#define GET_VAR1(name, id, dims)
 
#define GET_1D_INT_VAR1(name, id, vals)
 
#define GET_1D_DBL_VAR1(name, id, vals)
 

Functions

int main (int argc, char *argv[])
 

Variables

int ncFile1
 

Macro Definition Documentation

◆ ERR_NC

#define ERR_NC (   e)
Value:
{ \ printf( "Error: %s\n", nc_strerror( e ) ); \ exit( 2 ); \ }

Definition at line 30 of file visuMapVtk.cpp.

◆ GET_1D_DBL_VAR1

#define GET_1D_DBL_VAR1 (   name,
  id,
  vals 
)
Value:
{ \ std::vector< int > dum_dims; \ GET_VAR1( name, id, dum_dims ); \ if( -1 != ( id ) ) \ { \ size_t ntmp; \ int dvfail = nc_inq_dimlen( ncFile1, dum_dims[0], &ntmp ); \ if( NC_NOERR != dvfail ) \ { \ ERR_NC( dvfail ) \ } \ ( vals ).resize( ntmp ); \ size_t ntmp1 = 0; \ dvfail = nc_get_vara_double( ncFile1, id, &ntmp1, &ntmp, &( vals )[0] ); \ if( NC_NOERR != dvfail ) \ { \ ERR_NC( dvfail ) \ } \ } \ }

Definition at line 100 of file visuMapVtk.cpp.

◆ GET_1D_INT_VAR1

#define GET_1D_INT_VAR1 (   name,
  id,
  vals 
)
Value:
{ \ GET_VAR1( name, id, vals ); \ if( -1 != ( id ) ) \ { \ size_t ntmp; \ int ivfail = nc_inq_dimlen( ncFile1, ( vals )[0], &ntmp ); \ if( NC_NOERR != ivfail ) \ { \ ERR_NC( ivfail ) \ } \ ( vals ).resize( ntmp ); \ size_t ntmp1 = 0; \ ivfail = nc_get_vara_int( ncFile1, id, &ntmp1, &ntmp, &( vals )[0] ); \ if( NC_NOERR != ivfail ) \ { \ ERR_NC( ivfail ) \ } \ } \ }

Definition at line 79 of file visuMapVtk.cpp.

◆ GET_DIM1

#define GET_DIM1 (   ncdim,
  name,
  val 
)
Value:
{ \ int gdfail = nc_inq_dimid( ncFile1, name, &( ncdim ) ); \ if( NC_NOERR == gdfail ) \ { \ size_t tmp_val; \ gdfail = nc_inq_dimlen( ncFile1, ncdim, &tmp_val ); \ if( NC_NOERR != gdfail ) \ { \ ERR_NC( gdfail ) \ } \ else \ ( val ) = tmp_val; \ } \ else \ ( val ) = 0; \ }

Definition at line 41 of file visuMapVtk.cpp.

◆ GET_VAR1

#define GET_VAR1 (   name,
  id,
  dims 
)
Value:
{ \ ( id ) = -1; \ int gvfail = nc_inq_varid( ncFile1, name, &( id ) ); \ if( NC_NOERR == gvfail ) \ { \ int ndims; \ gvfail = nc_inq_varndims( ncFile1, id, &ndims ); \ if( NC_NOERR == gvfail ) \ { \ ( dims ).resize( ndims ); \ gvfail = nc_inq_vardimid( ncFile1, id, &( dims )[0] ); \ if( NC_NOERR != gvfail ) \ { \ ERR_NC( gvfail ) \ } \ } \ } \ }

Definition at line 59 of file visuMapVtk.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 125 of file visuMapVtk.cpp.

126 { 127  128  ProgOptions opts; 129  int dimSource = 2; // for FV meshes is 2; for SE meshes, use fine mesh, dim will be 0 130  int dimTarget = 2; // 131  int otype = 0; 132  133  std::string inputfile1, inputSource, inputTarget; 134  opts.addOpt< std::string >( "map,m", "input map ", &inputfile1 ); 135  opts.addOpt< std::string >( "source,s", "source mesh", &inputSource ); 136  opts.addOpt< std::string >( "target,t", "target mesh", &inputTarget ); 137  opts.addOpt< int >( "dimSource,d", "dimension of source ", &dimSource ); 138  opts.addOpt< int >( "dimTarget,g", "dimension of target ", &dimTarget ); 139  opts.addOpt< int >( "typeOutput,o", " output type vtk(0), h5m(1)", &otype ); 140  141  int startSourceID = -1, endSourceID = -1, startTargetID = -1, endTargetID = -1; 142  opts.addOpt< int >( "startSourceID,b", "start source id ", &startSourceID ); 143  opts.addOpt< int >( "endSourceID,e", "end source id ", &endSourceID ); 144  opts.addOpt< int >( "startTargetID,c", "start target id ", &startTargetID ); 145  opts.addOpt< int >( "endTargetID,f", "end target id ", &endTargetID ); 146  // -b startSourceID -e endSourceID -c startTargetID -f endTargetID 147  148  opts.parseCommandLine( argc, argv ); 149  150  std::string extension = ".vtk"; 151  if( 1 == otype ) extension = ".h5m"; 152  // Open netcdf/exodus file 153  int fail = nc_open( inputfile1.c_str(), 0, &ncFile1 ); 154  if( NC_NOWRITE != fail ) 155  { 156  ERR_NC( fail ) 157  } 158  159  std::cout << " opened " << inputfile1 << " for map 1 \n"; 160  161  int temp_dim; 162  int na1, nb1, ns1; 163  GET_DIM1( temp_dim, "n_a", na1 ); 164  GET_DIM1( temp_dim, "n_b", nb1 ); 165  GET_DIM1( temp_dim, "n_s", ns1 ); 166  std::cout << " n_a, n_b, n_s : " << na1 << ", " << nb1 << ", " << ns1 << " for map 1 \n"; 167  std::vector< int > col1( ns1 ), row1( ns1 ); 168  std::vector< double > val1( ns1 ); 169  int idrow1, idcol1, ids1; 170  GET_1D_INT_VAR1( "row", idrow1, row1 ); 171  GET_1D_INT_VAR1( "col", idcol1, col1 ); 172  GET_1D_DBL_VAR1( "S", ids1, val1 ); 173  // first matrix 174  typedef Eigen::Triplet< double > Triplet; 175  std::vector< Triplet > tripletList; 176  tripletList.reserve( ns1 ); 177  for( int iv = 0; iv < ns1; iv++ ) 178  { 179  // all row and column indices are 1-based in the map file. 180  // inside Eigen3, we will use 0-based; then we will have to add back 1 when 181  // we dump out the files 182  tripletList.push_back( Triplet( row1[iv] - 1, col1[iv] - 1, val1[iv] ) ); 183  } 184  Eigen::SparseMatrix< double > weight1( nb1, na1 ); 185  186  weight1.setFromTriplets( tripletList.begin(), tripletList.end() ); 187  weight1.makeCompressed(); 188  // we read the matrix; now read moab source and target 189  Core core; 190  Interface* mb = &core; 191  ErrorCode rval; 192  Tag gtag = mb->globalId_tag(); 193  EntityHandle sourceSet, targetSet; 194  // those are the maps from global ids to the moab entity handles corresponding to those global ids 195  // which are corresponding to the global DOFs 196  map< int, EntityHandle > sourceHandles; 197  map< int, EntityHandle > targetHandles; 198  rval = mb->create_meshset( MESHSET_SET, sourceSet );MB_CHK_SET_ERR( rval, "can't create source mesh set" ); 199  rval = mb->create_meshset( MESHSET_SET, targetSet );MB_CHK_SET_ERR( rval, "can't create target mesh set" ); 200  const char* readopts = ""; 201  rval = mb->load_file( inputSource.c_str(), &sourceSet, readopts );MB_CHK_SET_ERR( rval, "Failed to read" ); 202  rval = mb->load_file( inputTarget.c_str(), &targetSet, readopts );MB_CHK_SET_ERR( rval, "Failed to read" ); 203  Range sRange; 204  rval = mb->get_entities_by_dimension( sourceSet, dimSource, sRange );MB_CHK_SET_ERR( rval, "Failed to get sRange" ); 205  vector< int > sids; 206  sids.resize( sRange.size() ); 207  rval = mb->tag_get_data( gtag, sRange, &sids[0] );MB_CHK_SET_ERR( rval, "Failed to get ids for srange" ); 208  // all global ids are 1 based in the source file, and they correspond to the dofs in the map file 209  for( size_t i = 0; i < sids.size(); i++ ) 210  { 211  EntityHandle eh = sRange[i]; 212  int gid = sids[i]; 213  sourceHandles[gid] = eh; 214  } 215  Range tRange; 216  rval = mb->get_entities_by_dimension( targetSet, dimTarget, tRange );MB_CHK_SET_ERR( rval, "Failed to get tRange" ); 217  vector< int > tids; 218  tids.resize( tRange.size() ); 219  rval = mb->tag_get_data( gtag, tRange, &tids[0] );MB_CHK_SET_ERR( rval, "Failed to get ids for trange" ); 220  // all global ids are 1 based in the target file, and they correspond to the dofs in the map file 221  for( size_t i = 0; i < tids.size(); i++ ) 222  { 223  EntityHandle eh = tRange[i]; 224  int gid = tids[i]; 225  targetHandles[gid] = eh; 226  } 227  // a dense tag for weights 228  Tag wtag; 229  double defVal = 0; 230  231  std::string name_map = inputfile1; 232  // strip last 3 chars (.nc extension) 233  name_map.erase( name_map.begin() + name_map.length() - 3, name_map.end() ); 234  // if path , remove from name 235  size_t pos = name_map.rfind( '/', name_map.length() ); 236  if( pos != std::string::npos ) name_map = name_map.erase( 0, pos + 1 ); 237  238  rval = mb->tag_get_handle( "weight", 1, MB_TYPE_DOUBLE, wtag, MB_TAG_CREAT | MB_TAG_DENSE, &defVal );MB_CHK_SET_ERR( rval, "Failed to create weight" ); 239  EntityHandle partialSet; 240  rval = mb->create_meshset( MESHSET_SET, partialSet );MB_CHK_SET_ERR( rval, "can't create partial set" ); 241  // how to get a complete row in sparse matrix? Or a complete column ? 242  for( int col = startSourceID - 1; col <= endSourceID - 1; col++ ) 243  { 244  Range targetEnts; // will find entries for column col-1 in sparse matrix weight1, and its entries 245  // will assign a dense tag with values, and write out the file 246  if( col < 0 ) continue; 247  Eigen::SparseVector< double > colVect = weight1.col( col ); 248  // the row indices correspond to target cells 249  for( Eigen::SparseVector< double >::InnerIterator it( colVect ); it; ++it ) 250  { 251  double weight = it.value(); // == vec[ it.index() ] 252  // we add back the 1 that we subtract 253  int globalIdRow = it.index() + 1; 254  EntityHandle th = targetHandles[globalIdRow]; 255  targetEnts.insert( th ); 256  rval = mb->tag_set_data( wtag, &th, 1, &weight );MB_CHK_SET_ERR( rval, "Failed to set weight tag on target" ); 257  } 258  259  if( dimTarget == 0 ) 260  { 261  Range adjCells; 262  rval = mb->get_adjacencies( targetEnts, 2, false, adjCells, Interface::UNION );MB_CHK_SET_ERR( rval, " can't get adj cells " ); 263  targetEnts.merge( adjCells ); 264  } 265  266  rval = mb->add_entities( partialSet, targetEnts );MB_CHK_SET_ERR( rval, "Failed to add target entities to partial set" ); 267  // write now the set in a numbered file 268  std::stringstream fff; 269  fff << name_map << "_column" << col + 1 << extension; 270  rval = mb->write_mesh( fff.str().c_str(), &partialSet, 1 );MB_CHK_ERR( rval ); 271  // remove from partial set the entities it has 272  rval = mb->clear_meshset( &partialSet, 1 );MB_CHK_SET_ERR( rval, "Failed to empty partial set" ); 273  } 274  275  // how to get a complete row in sparse matrix? 276  for( int row = startTargetID - 1; row <= endTargetID - 1; row++ ) 277  { 278  Range sourceEnts; // will find entries for row in sparse matrix weight1, and its entries 279  // will assign a dense tag with values, and write out the file 280  if( row < 0 ) continue; 281  Eigen::SparseVector< double > rowVect = weight1.row( row ); 282  // the row indices correspond to target cells 283  for( Eigen::SparseVector< double >::InnerIterator it( rowVect ); it; ++it ) 284  { 285  double weight = it.value(); // == vec[ it.index() ] 286  int globalIdCol = it.index() + 1; 287  EntityHandle sh = sourceHandles[globalIdCol]; 288  sourceEnts.insert( sh ); 289  rval = mb->tag_set_data( wtag, &sh, 1, &weight );MB_CHK_SET_ERR( rval, "Failed to set weight tag on source" ); 290  } 291  if( dimSource == 0 ) 292  { 293  Range adjCells; 294  rval = mb->get_adjacencies( sourceEnts, 2, false, adjCells, Interface::UNION );MB_CHK_SET_ERR( rval, " can't get adj cells " ); 295  sourceEnts.merge( adjCells ); 296  } 297  rval = mb->add_entities( partialSet, sourceEnts );MB_CHK_SET_ERR( rval, "Failed to add source entities" ); 298  // write now the set in a numbered file 299  std::stringstream fff; 300  fff << name_map << "_row" << row + 1 << extension; 301  rval = mb->write_mesh( fff.str().c_str(), &partialSet, 1 );MB_CHK_ERR( rval ); 302  rval = mb->clear_meshset( &partialSet, 1 );MB_CHK_SET_ERR( rval, "Failed to empty partial set" ); 303  } 304  return 0; 305 }

References moab::Core::add_entities(), ProgOptions::addOpt(), moab::Core::clear_meshset(), moab::Core::create_meshset(), ERR_NC, ErrorCode, moab::fail(), GET_1D_DBL_VAR1, GET_1D_INT_VAR1, moab::Core::get_adjacencies(), GET_DIM1, moab::Core::get_entities_by_dimension(), moab::Core::globalId_tag(), moab::Range::index(), moab::Range::insert(), moab::Core::load_file(), mb, MB_CHK_ERR, MB_CHK_SET_ERR, MB_TAG_CREAT, MB_TAG_DENSE, MB_TYPE_DOUBLE, moab::Range::merge(), MESHSET_SET, ncFile1, ProgOptions::parseCommandLine(), moab::Range::size(), moab::Core::tag_get_data(), moab::Core::tag_get_handle(), moab::Core::tag_set_data(), moab::Interface::UNION, and moab::Core::write_mesh().

Variable Documentation

◆ ncFile1

int ncFile1

Definition at line 39 of file visuMapVtk.cpp.

Referenced by main().