Mesh Oriented datABase  (version 5.5.0)
An array-based unstructured mesh library
TempestOnlineMap.hpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 ///
3 /// \file TempestOnlineMap.h
4 /// \author Vijay Mahadevan
5 /// \version November 20, 2017
6 ///
7 
8 #ifndef _TEMPESTONLINEMAP_H_
9 #define _TEMPESTONLINEMAP_H_
10 
11 #include "moab/MOABConfig.h"
12 
13 // Tempest includes
14 #ifdef MOAB_HAVE_TEMPESTREMAP
16 #include "OfflineMap.h"
17 #else
18 #error Re-configure with TempestRemap
19 #endif
20 
21 #include <string>
22 #include <vector>
23 
24 #ifdef MOAB_HAVE_EIGEN3
25 #include <Eigen/Sparse>
26 #endif
27 
28 ///////////////////////////////////////////////////////////////////////////////
29 
30 #if !defined( RECTANGULAR_TRUNCATION ) && !defined( TRIANGULAR_TRUNCATION )
31 #define RECTANGULAR_TRUNCATION
32 // #define TRIANGULAR_TRUNCATION
33 #endif
34 
35 ///////////////////////////////////////////////////////////////////////////////
36 
37 // Forward declarations
38 class Mesh;
39 
40 ///////////////////////////////////////////////////////////////////////////////
41 
42 namespace moab
43 {
44 
45 /// <summary>
46 /// An offline map between two Meshes.
47 /// </summary>
48 class TempestOnlineMap : public OfflineMap
49 {
50 
51  public:
52  /// <summary>
53  /// Generate the metadata associated with the offline map.
54  /// </summary>
56 
57  /// <summary>
58  /// Define a virtual destructor.
59  /// </summary>
60  virtual ~TempestOnlineMap();
61 
62  public:
63  // Input / Output types
65  {
70  };
71 
72  /// <summary>
73  /// Generate the offline map, given the source and target mesh and discretization details.
74  /// This method generates the mapping between the two meshes based on the overlap and stores
75  /// the result in the SparseMatrix.
76  /// </summary>
77  moab::ErrorCode GenerateRemappingWeights( std::string strInputType,
78  std::string strOutputType,
79  const GenerateOfflineMapAlgorithmOptions& mapOptions,
80  const std::string& srcDofTagName = "GLOBAL_ID",
81  const std::string& tgtDofTagName = "GLOBAL_ID" );
82 
83  /// <summary>
84  /// Generate the metadata associated with the offline map.
85  /// </summary>
86  // moab::ErrorCode GenerateMetaData();
87 
88  /// <summary>
89  /// Read the OfflineMap from a NetCDF file.
90  /// </summary>
91  moab::ErrorCode ReadParallelMap( const char* strSource,
92  const std::vector< int >& owned_dof_ids,
93  bool row_major_ownership = true );
94 
95  /// <summary>
96  /// Write the TempestOnlineMap to a parallel NetCDF file.
97  /// </summary>
98  moab::ErrorCode WriteParallelMap( const std::string& strTarget,
99  const std::map< std::string, std::string >& attrMap );
100 
101  /// <summary>
102  /// Determine if the map is first-order accurate.
103  /// </summary>
104  virtual int IsConsistent( double dTolerance );
105 
106  /// <summary>
107  /// Determine if the map is conservative.
108  /// </summary>
109  virtual int IsConservative( double dTolerance );
110 
111  /// <summary>
112  /// Determine if the map is monotone.
113  /// </summary>
114  virtual int IsMonotone( double dTolerance );
115 
116  /// <summary>
117  /// If we computed the reduction, get the vector representing the source areas for all entities
118  /// in the mesh
119  /// </summary>
120  const DataArray1D< double >& GetGlobalSourceAreas() const;
121 
122  /// <summary>
123  /// If we computed the reduction, get the vector representing the target areas for all entities
124  /// in the mesh
125  /// </summary>
126  const DataArray1D< double >& GetGlobalTargetAreas() const;
127 
128  private:
129  /// <summary>
130  /// Compute the remapping weights as a permutation matrix that relates DoFs on the source
131  /// mesh
132  /// to DoFs on the target mesh.
133  /// </summary>
134  moab::ErrorCode LinearRemapNN_MOAB( bool use_GID_matching = false, bool strict_check = false );
135 
136  /// <summary>
137  /// Compute the remapping weights for a FV field defined on the source to a
138  /// FV field defined on the target mesh.
139  /// </summary>
140  void LinearRemapFVtoFV_Tempest_MOAB( int nOrder );
141 
142  /// <summary>
143  /// Generate the OfflineMap for linear conserative element-average
144  /// spectral element to element average remapping.
145  /// </summary>
146  void LinearRemapSE0_Tempest_MOAB( const DataArray3D< int >& dataGLLNodes,
147  const DataArray3D< double >& dataGLLJacobian );
148 
149  /// <summary>
150  /// Generate the OfflineMap for cubic conserative element-average
151  /// spectral element to element average remapping.
152  /// </summary>
153  void LinearRemapSE4_Tempest_MOAB( const DataArray3D< int >& dataGLLNodes,
154  const DataArray3D< double >& dataGLLJacobian,
155  int nMonotoneType,
156  bool fContinuousIn,
157  bool fNoConservation );
158 
159  /// <summary>
160  /// Generate the OfflineMap for remapping from finite volumes to finite
161  /// elements.
162  /// </summary>
163  void LinearRemapFVtoGLL_MOAB( const DataArray3D< int >& dataGLLNodes,
164  const DataArray3D< double >& dataGLLJacobian,
165  const DataArray1D< double >& dataGLLNodalArea,
166  int nOrder,
167  int nMonotoneType,
168  bool fContinuous,
169  bool fNoConservation );
170 
171  /// <summary>
172  /// Generate the OfflineMap for remapping from finite elements to finite
173  /// elements.
174  /// </summary>
175  void LinearRemapGLLtoGLL2_MOAB( const DataArray3D< int >& dataGLLNodesIn,
176  const DataArray3D< double >& dataGLLJacobianIn,
177  const DataArray3D< int >& dataGLLNodesOut,
178  const DataArray3D< double >& dataGLLJacobianOut,
179  const DataArray1D< double >& dataNodalAreaOut,
180  int nPin,
181  int nPout,
182  int nMonotoneType,
183  bool fContinuousIn,
184  bool fContinuousOut,
185  bool fNoConservation );
186 
187  /// <summary>
188  /// Generate the OfflineMap for remapping from finite elements to finite
189  /// elements (pointwise interpolation).
190  /// </summary>
191  void LinearRemapGLLtoGLL2_Pointwise_MOAB( const DataArray3D< int >& dataGLLNodesIn,
192  const DataArray3D< double >& dataGLLJacobianIn,
193  const DataArray3D< int >& dataGLLNodesOut,
194  const DataArray3D< double >& dataGLLJacobianOut,
195  const DataArray1D< double >& dataNodalAreaOut,
196  int nPin,
197  int nPout,
198  int nMonotoneType,
199  bool fContinuousIn,
200  bool fContinuousOut );
201 
202  /// <summary>
203  /// Copy the local matrix from Tempest SparseMatrix representation (ELL)
204  /// to the parallel CSR Eigen Matrix for scalable application of matvec
205  /// needed for projections.
206  /// </summary>
207 #ifdef MOAB_HAVE_EIGEN3
208  void copy_tempest_sparsemat_to_eigen3();
209 #endif
210 
211  /// <summary>
212  /// Parallel I/O with HDF5 to write out the remapping weights from multiple processors.
213  /// </summary>
214  moab::ErrorCode WriteSCRIPMapFile( const std::string& strOutputFile,
215  const std::map< std::string, std::string >& attrMap );
216 
217  /// <summary>
218  /// Parallel I/O with NetCDF to write out the SCRIP file from multiple processors.
219  /// </summary>
220  moab::ErrorCode WriteHDF5MapFile( const std::string& filename );
221 
222  public:
223  /// <summary>
224  /// Store the tag names associated with global DoF ids for source and target meshes
225  /// </summary>
226  moab::ErrorCode SetDOFmapTags( const std::string srcDofTagName, const std::string tgtDofTagName );
227 
228  /// <summary>
229  /// Compute the association between the solution tag global DoF numbering and
230  /// the local matrix numbering so that matvec operations can be performed
231  /// consistently.
232  /// </summary>
234  bool isSrcContinuous,
235  DataArray3D< int >* srcdataGLLNodes,
236  DataArray3D< int >* srcdataGLLNodesSrc,
237  DiscretizationType destType,
238  bool isDestContinuous,
239  DataArray3D< int >* tgtdataGLLNodes );
240 
241 #ifdef MOAB_HAVE_EIGEN3
242 
243  typedef Eigen::Matrix< double, 1, Eigen::Dynamic > WeightDRowVector;
244  typedef Eigen::Matrix< double, Eigen::Dynamic, 1 > WeightDColVector;
245  typedef Eigen::SparseVector< double > WeightSVector;
246  typedef Eigen::SparseMatrix< double, Eigen::RowMajor > WeightRMatrix;
247  typedef Eigen::SparseMatrix< double, Eigen::ColMajor > WeightCMatrix;
248 
249  typedef WeightDRowVector WeightRowVector;
250  typedef WeightDColVector WeightColVector;
251  typedef WeightRMatrix WeightMatrix;
252 
253  /// <summary>
254  /// Get the raw reference to the Eigen weight matrix representing the projection from source to
255  /// destination mesh.
256  /// </summary>
257  WeightMatrix& GetWeightMatrix();
258 
259  /// <summary>
260  /// Get the row vector that is amenable for application of A*x operation.
261  /// </summary>
262  WeightRowVector& GetRowVector();
263 
264  /// <summary>
265  /// Get the column vector that is amenable for application of A^T*x operation.
266  /// </summary>
267  WeightColVector& GetColVector();
268 
269 #endif
270 
271  /// <summary>
272  /// Get the number of total Degrees-Of-Freedom defined on the source mesh.
273  /// </summary>
275 
276  /// <summary>
277  /// Get the number of total Degrees-Of-Freedom defined on the destination mesh.
278  /// </summary>
280 
281  /// <summary>
282  /// Get the number of local Degrees-Of-Freedom defined on the source mesh.
283  /// </summary>
285 
286  /// <summary>
287  /// Get the number of local Degrees-Of-Freedom defined on the destination mesh.
288  /// </summary>
290 
291  /// <summary>
292  /// Get the number of Degrees-Of-Freedom per element on the source mesh.
293  /// </summary>
295 
296  /// <summary>
297  /// Get the number of Degrees-Of-Freedom per element on the destination mesh.
298  /// </summary>
300 
301  /// <summary>
302  /// Set the number of Degrees-Of-Freedom per element on the source mesh.
303  /// </summary>
304  void SetSourceNDofsPerElement( int ns );
305 
306  /// <summary>
307  /// Get the number of Degrees-Of-Freedom per element on the destination mesh.
308  /// </summary>
309  void SetDestinationNDofsPerElement( int nt );
310 
311  /// <summary>
312  /// Get the global Degrees-Of-Freedom ID on the destination mesh.
313  /// </summary>
314  int GetRowGlobalDoF( int localID ) const;
315 
316  /// <summary>
317  /// Get the index of globaRowDoF.
318  /// </summary>
319  inline int GetIndexOfRowGlobalDoF( int globalRowDoF ) const;
320 
321  /// <summary>
322  /// Get the global Degrees-Of-Freedom ID on the source mesh.
323  /// </summary>
324  int GetColGlobalDoF( int localID ) const;
325 
326  /// <summary>
327  /// Get the index of globaColDoF.
328  /// </summary>
329  inline int GetIndexOfColGlobalDoF( int globalColDoF ) const;
330 
331  /// <summary>
332  /// Apply the weight matrix onto the source vector provided as input, and return the column
333  /// vector (solution projection) after the map application
334  /// Compute: \p tgtVals = A(S->T) * \srcVals, or
335  /// if (transpose) \p tgtVals = [A(T->S)]^T * \srcVals
336  /// </summary>
337  moab::ErrorCode ApplyWeights( std::vector< double >& srcVals,
338  std::vector< double >& tgtVals,
339  bool transpose = false );
340 
341  /// <summary>
342  /// Apply the weight matrix onto the source vector (tag) provided as input, and return the
343  /// column vector (solution projection) in a tag, after the map application
344  /// Compute: \p tgtVals = A(S->T) * \srcVals, or
345  /// if (transpose) \p tgtVals = [A(T->S)]^T * \srcVals
346  /// </summary>
347  moab::ErrorCode ApplyWeights( moab::Tag srcSolutionTag, moab::Tag tgtSolutionTag, bool transpose = false );
348 
349  typedef double ( *sample_function )( double, double );
350 
351  /// <summary>
352  /// Define an analytical solution over the given (source or target) mesh, as specificed in
353  /// the context. This routine will define a tag that is compatible with the specified
354  /// discretization method type and order and sample the solution exactly using the
355  /// analytical function provided by the user.
356  /// </summary>
358  const std::string& solnName,
360  sample_function testFunction,
361  moab::Tag* clonedSolnTag = NULL,
362  std::string cloneSolnName = "" );
363 
364  /// <summary>
365  /// Compute the error between a sampled (exact) solution and a projected solution in various
366  /// error norms.
367  /// </summary>
369  moab::Tag& exactTag,
370  moab::Tag& approxTag,
371  std::map< std::string, double >& metrics,
372  bool verbose = true );
373 
374  moab::ErrorCode fill_row_ids( std::vector< int >& ids_of_interest )
375  {
376  ids_of_interest.reserve( row_gdofmap.size() );
377  // need to add 1
378  for( auto it = row_gdofmap.begin(); it != row_gdofmap.end(); it++ )
379  ids_of_interest.push_back( *it + 1 );
380 
381  return moab::MB_SUCCESS;
382  }
383 
384  moab::ErrorCode fill_col_ids( std::vector< int >& ids_of_interest )
385  {
386  ids_of_interest.reserve( col_gdofmap.size() );
387  // need to add 1
388  for( auto it = col_gdofmap.begin(); it != col_gdofmap.end(); it++ )
389  ids_of_interest.push_back( *it + 1 );
390  return moab::MB_SUCCESS;
391  }
392 
393  moab::ErrorCode set_col_dc_dofs( std::vector< int >& values_entities );
394 
395  moab::ErrorCode set_row_dc_dofs( std::vector< int >& values_entities );
396 
397  private:
398  void setup_sizes_dimensions();
399 
400 #ifdef MOAB_HAVE_MPI
401  int rearrange_arrays_by_dofs( const std::vector< unsigned int >& gdofmap,
402  DataArray1D< double >& vecFaceArea,
403  DataArray1D< double >& dCenterLon,
404  DataArray1D< double >& dCenterLat,
405  DataArray2D< double >& dVertexLat,
406  DataArray2D< double >& dVertexLon,
407  std::vector< int >& masks,
408  unsigned& N, // this will be output too now
409  int nv,
410  int& maxdof );
411 #endif
412 
413  /// <summary>
414  /// The fundamental remapping operator object.
415  /// </summary>
417 
418 #ifdef MOAB_HAVE_EIGEN3
419 
420  WeightMatrix m_weightMatrix;
421  WeightRowVector m_rowVector;
422  WeightColVector m_colVector;
423 
424 #endif
425 
426  /// <summary>
427  /// The reference to the moab::Core object that contains source/target and overlap sets.
428  /// </summary>
430 
431 #ifdef MOAB_HAVE_MPI
432  /// <summary>
433  /// The reference to the parallel communicator object used by the Core object.
434  /// </summary>
435  moab::ParallelComm* m_pcomm;
436 #endif
437 
438  /// <summary>
439  /// The original tag data and local to global DoF mapping to associate matrix values to
440  /// solution <summary>
442  std::vector< unsigned > row_gdofmap, col_gdofmap, srccol_gdofmap;
443 
444  // make it int, because it can be -1 in new logic
446 
447  std::map< int, int > rowMap, colMap;
448 
452 
453  // Key details about the current map
458 
459  Mesh* m_meshInput;
463 
465  int rank, size;
466 };
467 
468 ///////////////////////////////////////////////////////////////////////////////
469 
470 inline int moab::TempestOnlineMap::GetRowGlobalDoF( int localRowID ) const
471 {
472  return row_gdofmap[localRowID];
473 }
474 
475 inline int moab::TempestOnlineMap::GetIndexOfRowGlobalDoF( int globalRowDoF ) const /* 0 based */
476 {
477  return globalRowDoF + 1;
478 }
479 ///////////////////////////////////////////////////////////////////////////////
480 
481 inline int moab::TempestOnlineMap::GetColGlobalDoF( int localColID ) const
482 {
483  return col_gdofmap[localColID];
484 }
485 
486 inline int moab::TempestOnlineMap::GetIndexOfColGlobalDoF( int globalColDoF ) const /* 0 based */
487 {
488  return globalColDoF + 1; // temporary
489 }
490 ///////////////////////////////////////////////////////////////////////////////
491 
493 {
494  return m_nDofsPEl_Src;
495 }
496 
497 ///////////////////////////////////////////////////////////////////////////////
498 
500 {
501  return m_nDofsPEl_Dest;
502 }
503 
504 // setters for m_nDofsPEl_Src, m_nDofsPEl_Dest
506 {
507  m_nDofsPEl_Src = ns;
508 }
509 
511 {
512  m_nDofsPEl_Dest = nt;
513 }
514 
515 ///////////////////////////////////////////////////////////////////////////////
516 #ifdef MOAB_HAVE_EIGEN3
517 
519 {
520  return m_weightMatrix.cols(); // return the global number of rows from the weight matrix
521 }
522 
523 // ///////////////////////////////////////////////////////////////////////////////
524 
526 {
527  return m_weightMatrix.rows(); // return the global number of columns from the weight matrix
528 }
529 
530 ///////////////////////////////////////////////////////////////////////////////
531 
533 {
534  return m_weightMatrix.cols(); // return the local number of rows from the weight matrix
535 }
536 
537 ///////////////////////////////////////////////////////////////////////////////
538 
540 {
541  return m_weightMatrix.rows(); // return the local number of columns from the weight matrix
542 }
543 
544 ///////////////////////////////////////////////////////////////////////////////
545 
546 inline moab::TempestOnlineMap::WeightMatrix& moab::TempestOnlineMap::GetWeightMatrix()
547 {
548  assert( m_weightMatrix.rows() != 0 && m_weightMatrix.cols() != 0 );
549  return m_weightMatrix;
550 }
551 
552 ///////////////////////////////////////////////////////////////////////////////
553 
554 inline moab::TempestOnlineMap::WeightRowVector& moab::TempestOnlineMap::GetRowVector()
555 {
556  assert( m_rowVector.size() != 0 );
557  return m_rowVector;
558 }
559 
560 ///////////////////////////////////////////////////////////////////////////////
561 
562 inline moab::TempestOnlineMap::WeightColVector& moab::TempestOnlineMap::GetColVector()
563 {
564  assert( m_colVector.size() != 0 );
565  return m_colVector;
566 }
567 
568 ///////////////////////////////////////////////////////////////////////////////
569 
570 #endif
571 
572 } // namespace moab
573 
574 #endif