Mesh Oriented datABase  (version 5.6.0)
An array-based unstructured mesh library
visuMap.cpp File Reference
#include "moab/MOABConfig.h"
#include "moab/ProgOptions.hpp"
#include "moab/Core.hpp"
#include "moab/Range.hpp"
#include "moab/IntxMesh/IntxUtils.hpp"
#include "moab/ReadUtilIface.hpp"
#include "netcdf.h"
#include <cmath>
#include <sstream>
#include <map>
#include <Eigen/Sparse>
+ Include dependency graph for visuMap.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 51 of file visuMap.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 121 of file visuMap.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 100 of file visuMap.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 62 of file visuMap.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 80 of file visuMap.cpp.

Function Documentation

◆ main()

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

Definition at line 146 of file visuMap.cpp.

147 {
148 
149  ProgOptions opts;
150  int dimSource = 2; // for FV meshes is 2; for SE meshes, use fine mesh, dim will be 0
151  int dimTarget = 2; //
152  int otype = 2;
153  double fraction = 0.05;
154  std::string inputfile1, inputSource, inputTarget;
155  opts.addOpt< std::string >( "map,m", "input map ", &inputfile1 );
156  opts.addOpt< std::string >( "source,s", "source mesh", &inputSource );
157  opts.addOpt< std::string >( "target,t", "target mesh", &inputTarget );
158  opts.addOpt< int >( "dimSource,d", "dimension of source ", &dimSource );
159  opts.addOpt< int >( "dimTarget,g", "dimension of target ", &dimTarget );
160  opts.addOpt< int >( "typeOutput,o", " output type vtk(0), h5m(1), view(default = 2) ", &otype );
161 
162  int startSourceID = -1, endSourceID = -1, startTargetID = -1, endTargetID = -1;
163  opts.addOpt< int >( "startSourceID,b", "start source id ", &startSourceID );
164  opts.addOpt< int >( "endSourceID,e", "end source id ", &endSourceID );
165  opts.addOpt< int >( "startTargetID,c", "start target id ", &startTargetID );
166  opts.addOpt< int >( "endTargetID,f", "end target id ", &endTargetID );
167  opts.addOpt< double >( "raiseFraction,r", "fraction for raising source points height (default 0.05)", &fraction );
168  // -b startSourceID -e endSourceID -c startTargetID -f endTargetID
169 
170  opts.parseCommandLine( argc, argv );
171 
172  std::string extension = ".vtk";
173  if( 1 <= otype ) extension = ".h5m";
174 
175  // Open netcdf map file
176  int fail = nc_open( inputfile1.c_str(), 0, &ncFile1 );
177  if( NC_NOWRITE != fail )
178  {
179  ERR_NC( fail )
180  }
181 
182  std::cout << " opened " << inputfile1 << " for map 1 \n";
183 
184  int temp_dim;
185  int na1, nb1, ns1;
186  GET_DIM1( temp_dim, "n_a", na1 );
187  GET_DIM1( temp_dim, "n_b", nb1 );
188  GET_DIM1( temp_dim, "n_s", ns1 );
189  std::cout << " n_a, n_b, n_s : " << na1 << ", " << nb1 << ", " << ns1 << " for map 1 \n";
190  std::vector< int > col1( ns1 ), row1( ns1 );
191  std::vector< double > val1( ns1 );
192  int idrow1, idcol1, ids1;
193  GET_1D_INT_VAR1( "row", idrow1, row1 );
194  GET_1D_INT_VAR1( "col", idcol1, col1 );
195  GET_1D_DBL_VAR1( "S", ids1, val1 );
196 
197  // we read the matrix; now read moab source and target
198  Core core;
199  Interface* mb = &core;
200  Tag gtag = mb->globalId_tag();
201 
202  // a dense tag for weights
203  Tag wtag;
204  double defVal = 0;
205 
206  std::string name_map = inputfile1;
207  // strip last 3 chars (.nc extension)
208  name_map.erase( name_map.begin() + name_map.length() - 3, name_map.end() );
209  // if path , remove from name
210  size_t pos = name_map.rfind( '/', name_map.length() );
211  if( pos != std::string::npos ) name_map = name_map.erase( 0, pos + 1 );
212 
213  MB_CHK_SET_ERR( mb->tag_get_handle( "weight", 1, MB_TYPE_DOUBLE, wtag, MB_TAG_CREAT | MB_TAG_DENSE, &defVal ),
214  "Failed to create weight" );
215 
216  if( 2 == otype )
217  {
218  // create a view of the full map, in which each weight is shown on an edge that starts at the source cell center
219  // and ends at the target cell center
220  // source cell centers are raised a little, let's say a fraction 0.05 * radius, which is 1
221  // each edge gets the associated weight as a tag
222  // first read the cell centers for source and target meshes, directly from the map file
223 
224  std::vector< double > xc_a( na1 ), yc_a( na1 );
225  std::vector< double > xc_b( nb1 ), yc_b( nb1 );
226  int idxc_a, idxc_b, idyc_a, idyc_b;
227  GET_1D_DBL_VAR1( "xc_a", idxc_a, xc_a );
228  GET_1D_DBL_VAR1( "xc_b", idxc_b, xc_b );
229  GET_1D_DBL_VAR1( "yc_a", idyc_a, yc_a );
230  GET_1D_DBL_VAR1( "yc_b", idyc_b, yc_b );
231  // create source vertices, and target vertices
232  std::vector< double > vertex_coords_src( 3 * na1 );
233  // xc_a and yc_a are in degrees, usually
234  for( int i = 0; i < na1; i++ )
235  {
237  sph.R = 1 + fraction; // slightly higher
238  sph.lon = xc_a[i] * M_PI / 180;
239  sph.lat = yc_a[i] * M_PI / 180;
240  CartVect pos = IntxUtils::spherical_to_cart( sph );
241  vertex_coords_src[3 * i] = pos[0];
242  vertex_coords_src[3 * i + 1] = pos[1];
243  vertex_coords_src[3 * i + 2] = pos[2];
244  }
245  Range source_verts;
246  MB_CHK_SET_ERR( mb->create_vertices( &vertex_coords_src[0], na1, source_verts ),
247  "can't create source vertices" );
248  // create a set with source vertices
249  EntityHandle srcSet;
250  MB_CHK_SET_ERR( mb->create_meshset( MESHSET_SET, srcSet ), "can't create source set for vertices" );
251  MB_CHK_SET_ERR( mb->add_entities( srcSet, source_verts ), "can't add vertices" );
252  std::vector< int > vgid( na1 );
253  for( int i = 0; i < na1; i++ )
254  vgid[i] = i + 1;
255  MB_CHK_SET_ERR( mb->tag_set_data( gtag, source_verts, &vgid[0] ), "can't set global id on source verts" );
256 
257  std::vector< double > vertex_coords_tgt( 3 * nb1 );
258  // xc_a and yc_a are in degrees, usually
259  for( int i = 0; i < nb1; i++ )
260  {
262  sph.R = 1; //
263  sph.lon = xc_b[i] * M_PI / 180;
264  sph.lat = yc_b[i] * M_PI / 180;
265  CartVect pos = IntxUtils::spherical_to_cart( sph );
266  vertex_coords_tgt[3 * i] = pos[0];
267  vertex_coords_tgt[3 * i + 1] = pos[1];
268  vertex_coords_tgt[3 * i + 2] = pos[2];
269  }
270  Range target_verts;
271  MB_CHK_SET_ERR( mb->create_vertices( &vertex_coords_tgt[0], nb1, target_verts ),
272  "can't create target vertices" );
273  EntityHandle tgtSet;
274  MB_CHK_SET_ERR( mb->create_meshset( MESHSET_SET, tgtSet ), "can't create target set for vertices" );
275  MB_CHK_SET_ERR( mb->add_entities( tgtSet, target_verts ), "can't add vertices" );
276  vgid.resize( nb1 );
277  for( int i = 0; i < nb1; i++ )
278  vgid[i] = i + 1;
279  MB_CHK_SET_ERR( mb->tag_set_data( gtag, target_verts, &vgid[0] ), "can't set global id on target verts" );
280  // create ns1 edges
281 
282  ReadUtilIface* read_iface;
283  MB_CHK_ERR( mb->query_interface( read_iface ) );
284 
285  EntityHandle actual_start_handle;
286  EntityHandle* array = nullptr;
287  MB_CHK_ERR( read_iface->get_element_connect( ns1, 2, MBEDGE, 1, actual_start_handle, array ) );
288 
289  for( int i = 0; i < ns1; i++ )
290  {
291  array[2 * i] = source_verts[col1[i] - 1]; // 1 based to 0 based index
292  array[2 * i + 1] = target_verts[row1[i] - 1]; // 1 based to 0 based index
293  }
294  Range edges( actual_start_handle, actual_start_handle + ns1 - 1 );
295 
296  MB_CHK_SET_ERR( mb->tag_set_data( wtag, edges, &val1[0] ), "can't set tag on edges" );
297 
298  vgid.resize( ns1 );
299  for( int i = 0; i < ns1; i++ )
300  vgid[i] = i + 1;
301  MB_CHK_SET_ERR( mb->tag_set_data( gtag, edges, &vgid[0] ), "can't set global id on edges" );
302 
303  std::string name_file = name_map + extension;
304  MB_CHK_ERR( mb->write_mesh( name_file.c_str() ) );
305  std::cout << " wrote view map file " << name_file << " with source fraction height: " << fraction << "\n";
306 
307  return 0; // do not bother with other files created, just one file with weights on edges
308  }
309  // first matrix
310  typedef Eigen::Triplet< double > Triplet;
311  std::vector< Triplet > tripletList;
312  tripletList.reserve( ns1 );
313  for( int iv = 0; iv < ns1; iv++ )
314  {
315  // all row and column indices are 1-based in the map file.
316  // inside Eigen3, we will use 0-based; then we will have to add back 1 when
317  // we dump out the files
318  tripletList.push_back( Triplet( row1[iv] - 1, col1[iv] - 1, val1[iv] ) );
319  }
320  Eigen::SparseMatrix< double > weight1( nb1, na1 );
321 
322  weight1.setFromTriplets( tripletList.begin(), tripletList.end() );
323  weight1.makeCompressed();
324  EntityHandle sourceSet, targetSet;
325  // those are the maps from global ids to the moab entity handles corresponding to those global ids
326  // which are corresponding to the global DOFs
327  map< int, EntityHandle > sourceHandles;
328  map< int, EntityHandle > targetHandles;
329  MB_CHK_SET_ERR( mb->create_meshset( MESHSET_SET, sourceSet ), "can't create source mesh set" );
330  MB_CHK_SET_ERR( mb->create_meshset( MESHSET_SET, targetSet ), "can't create target mesh set" );
331  const char* readopts = "";
332  MB_CHK_SET_ERR( mb->load_file( inputSource.c_str(), &sourceSet, readopts ), "Failed to read" );
333  MB_CHK_SET_ERR( mb->load_file( inputTarget.c_str(), &targetSet, readopts ), "Failed to read" );
334  Range sRange;
335  MB_CHK_SET_ERR( mb->get_entities_by_dimension( sourceSet, dimSource, sRange ), "Failed to get sRange" );
336  vector< int > sids;
337  sids.resize( sRange.size() );
338  MB_CHK_SET_ERR( mb->tag_get_data( gtag, sRange, &sids[0] ), "Failed to get ids for srange" );
339  // all global ids are 1 based in the source file, and they correspond to the dofs in the map file
340  for( size_t i = 0; i < sids.size(); i++ )
341  {
342  EntityHandle eh = sRange[i];
343  int gid = sids[i];
344  sourceHandles[gid] = eh;
345  }
346  Range tRange;
347  MB_CHK_SET_ERR( mb->get_entities_by_dimension( targetSet, dimTarget, tRange ), "Failed to get tRange" );
348  vector< int > tids;
349  tids.resize( tRange.size() );
350  MB_CHK_SET_ERR( mb->tag_get_data( gtag, tRange, &tids[0] ), "Failed to get ids for trange" );
351  // all global ids are 1 based in the target file, and they correspond to the dofs in the map file
352  for( size_t i = 0; i < tids.size(); i++ )
353  {
354  EntityHandle eh = tRange[i];
355  int gid = tids[i];
356  targetHandles[gid] = eh;
357  }
358  EntityHandle partialSet;
359  MB_CHK_SET_ERR( mb->create_meshset( MESHSET_SET, partialSet ), "can't create partial set" );
360  // how to get a complete row in sparse matrix? Or a complete column ?
361  for( int col = startSourceID - 1; col <= endSourceID - 1; col++ )
362  {
363  Range targetEnts; // will find entries for column col-1 in sparse matrix weight1, and its entries
364  // will assign a dense tag with values, and write out the file
365  if( col < 0 ) continue;
366  Eigen::SparseVector< double > colVect = weight1.col( col );
367  // the row indices correspond to target cells
368  for( Eigen::SparseVector< double >::InnerIterator it( colVect ); it; ++it )
369  {
370  double weight = it.value(); // == vec[ it.index() ]
371  // we add back the 1 that we subtract
372  int globalIdRow = it.index() + 1;
373  EntityHandle th = targetHandles[globalIdRow];
374  targetEnts.insert( th );
375  MB_CHK_SET_ERR( mb->tag_set_data( wtag, &th, 1, &weight ), "Failed to set weight tag on target" );
376  }
377 
378  if( dimTarget == 0 )
379  {
380  Range adjCells;
381  MB_CHK_SET_ERR( mb->get_adjacencies( targetEnts, 2, false, adjCells, Interface::UNION ),
382  " can't get adj cells " );
383  targetEnts.merge( adjCells );
384  }
385 
386  MB_CHK_SET_ERR( mb->add_entities( partialSet, targetEnts ), "Failed to add target entities to partial set" );
387  // write now the set in a numbered file
388  std::stringstream fff;
389  fff << name_map << "_column" << col + 1 << extension;
390  MB_CHK_ERR( mb->write_mesh( fff.str().c_str(), &partialSet, 1 ) );
391  // remove from partial set the entities it has
392  MB_CHK_SET_ERR( mb->clear_meshset( &partialSet, 1 ), "Failed to empty partial set" );
393  }
394 
395  // how to get a complete row in sparse matrix?
396  for( int row = startTargetID - 1; row <= endTargetID - 1; row++ )
397  {
398  Range sourceEnts; // will find entries for row in sparse matrix weight1, and its entries
399  // will assign a dense tag with values, and write out the file
400  if( row < 0 ) continue;
401  Eigen::SparseVector< double > rowVect = weight1.row( row );
402  // the row indices correspond to target cells
403  for( Eigen::SparseVector< double >::InnerIterator it( rowVect ); it; ++it )
404  {
405  double weight = it.value(); // == vec[ it.index() ]
406  int globalIdCol = it.index() + 1;
407  EntityHandle sh = sourceHandles[globalIdCol];
408  sourceEnts.insert( sh );
409  MB_CHK_SET_ERR( mb->tag_set_data( wtag, &sh, 1, &weight ), "Failed to set weight tag on source" );
410  }
411  if( dimSource == 0 )
412  {
413  Range adjCells;
414  MB_CHK_SET_ERR( mb->get_adjacencies( sourceEnts, 2, false, adjCells, Interface::UNION ),
415  " can't get adj cells " );
416  sourceEnts.merge( adjCells );
417  }
418  MB_CHK_SET_ERR( mb->add_entities( partialSet, sourceEnts ), "Failed to add source entities" );
419  // write now the set in a numbered file
420  std::stringstream fff;
421  fff << name_map << "_row" << row + 1 << extension;
422  MB_CHK_ERR( mb->write_mesh( fff.str().c_str(), &partialSet, 1 ) );
423  MB_CHK_SET_ERR( mb->clear_meshset( &partialSet, 1 ), "Failed to empty partial set" );
424  }
425  return 0;
426 }

References moab::Core::add_entities(), ProgOptions::addOpt(), moab::Core::clear_meshset(), moab::Core::create_meshset(), moab::Core::create_vertices(), ERR_NC, moab::fail(), GET_1D_DBL_VAR1, GET_1D_INT_VAR1, moab::Core::get_adjacencies(), GET_DIM1, moab::ReadUtilIface::get_element_connect(), moab::Core::get_entities_by_dimension(), moab::Core::globalId_tag(), moab::Range::index(), moab::Range::insert(), moab::IntxUtils::SphereCoords::lat, moab::Core::load_file(), moab::IntxUtils::SphereCoords::lon, mb, MB_CHK_ERR, MB_CHK_SET_ERR, MB_TAG_CREAT, MB_TAG_DENSE, MB_TYPE_DOUBLE, MBEDGE, moab::Range::merge(), MESHSET_SET, ncFile1, ProgOptions::parseCommandLine(), moab::Interface::query_interface(), moab::IntxUtils::SphereCoords::R, moab::Range::size(), moab::IntxUtils::spherical_to_cart(), 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 60 of file visuMap.cpp.

Referenced by main().