Mesh Oriented datABase  (version 5.5.1)
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 #ifndef MOAB_HAVE_TEMPESTREMAP
13 #error Re-configure with TempestRemap
14 #endif
15 
17 
18 // Tempest includes
19 #pragma GCC diagnostic push
20 #pragma GCC diagnostic ignored "-Wunused-but-set-variable"
21 #pragma GCC diagnostic ignored "-Wunused-variable"
22 #include "OfflineMap.h"
23 
24 #ifdef MOAB_HAVE_EIGEN3
25 // #ifdef MOAB_HAVE_BLAS
26 // #define EIGEN_USE_BLAS
27 // #endif
28 // #define EIGEN_DONT_PARALLELIZE
29 // #define EIGEN_DONT_VECTORIZE
30 // #define EIGEN_STRONG_INLINE
31 #define EIGEN_DEFAULT_TO_ROW_MAJOR
32 #define EIGEN_RUNTIME_NO_MALLOC
33 // #define EIGEN_UNROLLING_LIMIT 0
34 #include <Eigen/Sparse>
35 #endif
36 
37 #include <unordered_set>
38 
39 #pragma GCC diagnostic pop
40 
41 ///////////////////////////////////////////////////////////////////////////////
42 
43 #if !defined( RECTANGULAR_TRUNCATION ) && !defined( TRIANGULAR_TRUNCATION )
44 #define RECTANGULAR_TRUNCATION
45 // #define TRIANGULAR_TRUNCATION
46 #endif
47 
48 ///////////////////////////////////////////////////////////////////////////////
49 
50 // Forward declarations
51 class Mesh;
52 
53 ///////////////////////////////////////////////////////////////////////////////
54 
55 namespace moab
56 {
57 
58 /// <summary>
59 /// An offline map between two Meshes.
60 /// </summary>
61 class TempestOnlineMap : public OfflineMap
62 {
63 
64  public:
65  /// <summary>
66  /// Generate the metadata associated with the offline map.
67  /// </summary>
69 
70  /// <summary>
71  /// Define a virtual destructor.
72  /// </summary>
73  virtual ~TempestOnlineMap();
74 
75  public:
76  // Input / Output types
78  {
83  };
84 
85  // Type of limiter
86  enum CAASType
87  {
88  CAAS_NONE = 0,
92  CAAS_QLT = 4
93  };
94 
95  /// <summary>
96  /// Generate the offline map, given the source and target mesh and discretization details.
97  /// This method generates the mapping between the two meshes based on the overlap and stores
98  /// the result in the SparseMatrix.
99  /// </summary>
100  moab::ErrorCode GenerateRemappingWeights( std::string strInputType,
101  std::string strOutputType,
102  const GenerateOfflineMapAlgorithmOptions& mapOptions,
103  const std::string& srcDofTagName = "GLOBAL_ID",
104  const std::string& tgtDofTagName = "GLOBAL_ID" );
105 
106  /// <summary>
107  /// Generate the metadata associated with the offline map.
108  /// </summary>
109  // moab::ErrorCode GenerateMetaData();
110 
111  /// <summary>
112  /// Read the OfflineMap from a NetCDF file.
113  /// </summary>
114  moab::ErrorCode ReadParallelMap( const char* strSource,
115  const std::vector< int >& tgt_dof_ids,
116  int arearead,
117  std::vector< double >& areaA,
118  int& nA,
119  std::vector< double >& areaB,
120  int& nB );
121 
122  /// <summary>
123  /// Write the TempestOnlineMap to a parallel NetCDF file.
124  /// </summary>
125  moab::ErrorCode WriteParallelMap( const std::string& strTarget,
126  const std::map< std::string, std::string >& attrMap );
127 
128  /// <summary>
129  /// Determine if the map is first-order accurate.
130  /// </summary>
131  virtual int IsConsistent( double dTolerance );
132 
133  /// <summary>
134  /// Determine if the map is conservative.
135  /// </summary>
136  virtual int IsConservative( double dTolerance );
137 
138  /// <summary>
139  /// Determine if the map is monotone.
140  /// </summary>
141  virtual int IsMonotone( double dTolerance );
142 
143  /// <summary>
144  /// If we computed the reduction, get the vector representing the source areas for all entities
145  /// in the mesh
146  /// </summary>
147  const DataArray1D< double >& GetGlobalSourceAreas() const;
148 
149  /// <summary>
150  /// If we computed the reduction, get the vector representing the target areas for all entities
151  /// in the mesh
152  /// </summary>
153  const DataArray1D< double >& GetGlobalTargetAreas() const;
154 
155  /// <summary>
156  /// Print information and metadata about the remapping weights.
157  /// </summary>
158  void PrintMapStatistics();
159 
160  private:
161  /// <summary>
162  /// Compute the remapping weights as a permutation matrix that relates DoFs on the source
163  /// mesh
164  /// to DoFs on the target mesh.
165  /// </summary>
166  moab::ErrorCode LinearRemapNN_MOAB( bool use_GID_matching = false, bool strict_check = false );
167 
168  /// <summary>
169  /// Compute the remapping weights for a FV field defined on the source to a
170  /// FV field defined on the target mesh.
171  /// </summary>
172  void LinearRemapFVtoFV_Tempest_MOAB( int nOrder );
173 
174  /// <summary>
175  /// Generate the OfflineMap for linear conserative element-average
176  /// spectral element to element average remapping.
177  /// </summary>
178  void LinearRemapSE0_Tempest_MOAB( const DataArray3D< int >& dataGLLNodes,
179  const DataArray3D< double >& dataGLLJacobian );
180 
181  /// <summary>
182  /// Generate the OfflineMap for cubic conserative element-average
183  /// spectral element to element average remapping.
184  /// </summary>
185  void LinearRemapSE4_Tempest_MOAB( const DataArray3D< int >& dataGLLNodes,
186  const DataArray3D< double >& dataGLLJacobian,
187  int nMonotoneType,
188  bool fContinuousIn,
189  bool fNoConservation );
190 
191  /// <summary>
192  /// Generate the OfflineMap for remapping from finite volumes to finite
193  /// elements.
194  /// </summary>
195  void LinearRemapFVtoGLL_MOAB( const DataArray3D< int >& dataGLLNodes,
196  const DataArray3D< double >& dataGLLJacobian,
197  const DataArray1D< double >& dataGLLNodalArea,
198  int nOrder,
199  int nMonotoneType,
200  bool fContinuous,
201  bool fNoConservation );
202 
203  /// <summary>
204  /// Generate the OfflineMap for remapping from finite elements to finite
205  /// elements.
206  /// </summary>
207  void LinearRemapGLLtoGLL2_MOAB( const DataArray3D< int >& dataGLLNodesIn,
208  const DataArray3D< double >& dataGLLJacobianIn,
209  const DataArray3D< int >& dataGLLNodesOut,
210  const DataArray3D< double >& dataGLLJacobianOut,
211  const DataArray1D< double >& dataNodalAreaOut,
212  int nPin,
213  int nPout,
214  int nMonotoneType,
215  bool fContinuousIn,
216  bool fContinuousOut,
217  bool fNoConservation );
218 
219  /// <summary>
220  /// Generate the OfflineMap for remapping from finite elements to finite
221  /// elements (pointwise interpolation).
222  /// </summary>
223  void LinearRemapGLLtoGLL2_Pointwise_MOAB( const DataArray3D< int >& dataGLLNodesIn,
224  const DataArray3D< double >& dataGLLJacobianIn,
225  const DataArray3D< int >& dataGLLNodesOut,
226  const DataArray3D< double >& dataGLLJacobianOut,
227  const DataArray1D< double >& dataNodalAreaOut,
228  int nPin,
229  int nPout,
230  int nMonotoneType,
231  bool fContinuousIn,
232  bool fContinuousOut );
233 
234  /// <summary>
235  /// Copy the local matrix from Tempest SparseMatrix representation (ELL)
236  /// to the parallel CSR Eigen Matrix for scalable application of matvec
237  /// needed for projections.
238  /// </summary>
239 #ifdef MOAB_HAVE_EIGEN3
240  void copy_tempest_sparsemat_to_eigen3();
241 #endif
242 
243  /// <summary>
244  /// Parallel I/O with HDF5 to write out the remapping weights from multiple processors.
245  /// </summary>
246  moab::ErrorCode WriteSCRIPMapFile( const std::string& strOutputFile,
247  const std::map< std::string, std::string >& attrMap );
248 
249  /// <summary>
250  /// Parallel I/O with NetCDF to write out the SCRIP file from multiple processors.
251  /// </summary>
252  moab::ErrorCode WriteHDF5MapFile( const std::string& filename );
253 
254  public:
255  /// <summary>
256  /// Store the tag names associated with global DoF ids for source and target meshes to be used for mapping.
257  /// <param name="srcDofTagName">The tag name associated with global DoF ids for the source mesh</param>
258  /// <param name="tgtDofTagName">The tag name associated with global DoF ids for the target mesh</param>
259  /// </summary>
260  moab::ErrorCode SetDOFmapTags( const std::string srcDofTagName, const std::string tgtDofTagName );
261 
262  /// <summary>
263  /// @brief Compute the association between the solution tag global DoF numbering and
264  /// the local matrix numbering so that matvec operations can be performed
265  /// consistently.
266  /// <param name="srcType">The discretization type of the source mesh</param>
267  /// <param name="srcOrder">The order of the discretization on the source mesh</param>
268  /// <param name="isSrcContinuous">The continuity of the discretization on the source mesh</param>
269  /// <param name="srcdataGLLNodes">The GLL nodes on the source mesh</param>
270  /// <param name="srcdataGLLNodesSrc">The GLL nodes on the source mesh</param>
271  /// <param name="destType">The discretization type of the destination mesh</param>
272  /// <param name="destOrder">The order of the discretization on the destination mesh</param>
273  /// <param name="isTgtContinuous">The continuity of the discretization on the destination mesh</param>
274  /// <param name="tgtdataGLLNodes">The GLL nodes on the destination mesh</param>
275  /// </summary>
277  int srcOrder,
278  bool isSrcContinuous,
279  DataArray3D< int >* srcdataGLLNodes,
280  DataArray3D< int >* srcdataGLLNodesSrc,
281  DiscretizationType destType,
282  int destOrder,
283  bool isTgtContinuous,
284  DataArray3D< int >* tgtdataGLLNodes );
285 
286  /// <summary>
287  /// @brief ApplyBoundsLimiting - Apply bounds limiting to the data field
288  /// @param dataInDouble - input data field
289  /// @param dataOutDouble - output data field
290  /// @param caasType - type of limiter
291  /// @param caasIteration - iteration number of limiter
292  /// @return - pair of mass defect pre and post limiter application
293  /// </summary>
294  std::pair< double, double > ApplyBoundsLimiting( std::vector< double >& dataInDouble,
295  std::vector< double >& dataOutDouble,
296  CAASType caasType = CAAS_GLOBAL,
297  int caasIteration = 0,
298  double mismatch = 0.0 );
299 
300  /// @brief
301  /// @param vecAdjFaces
302  /// @param nrings
303  /// @param entities
304  /// @param useMOABAdjacencies
305  /// @param trMesh
306  /// @return
307  void ComputeAdjacencyRelations( std::vector< std::unordered_set< int > >& vecAdjFaces,
308  int nrings,
309  const Range& entities,
310  bool useMOABAdjacencies = true,
311  Mesh* trMesh = nullptr );
312 
313 #ifdef MOAB_HAVE_EIGEN3
314 
315  typedef Eigen::Matrix< double, 1, Eigen::Dynamic > WeightDRowVector;
316  typedef Eigen::Matrix< double, Eigen::Dynamic, 1 > WeightDColVector;
317  typedef Eigen::SparseVector< double > WeightSVector;
318  typedef Eigen::SparseMatrix< double, Eigen::RowMajor > WeightRMatrix;
319  typedef Eigen::SparseMatrix< double, Eigen::ColMajor > WeightCMatrix;
320 
321  typedef WeightDRowVector WeightRowVector;
322  typedef WeightDColVector WeightColVector;
323  typedef WeightRMatrix WeightMatrix;
324 
325  /// <summary>
326  /// Get the raw reference to the Eigen weight matrix representing the projection from source to
327  /// destination mesh.
328  /// </summary>
329  WeightMatrix& GetWeightMatrix();
330 
331  /// <summary>
332  /// Get the row vector that is amenable for application of A*x operation.
333  /// </summary>
334  WeightRowVector& GetRowVector();
335 
336  /// <summary>
337  /// Get the column vector that is amenable for application of A^T*x operation.
338  /// </summary>
339  WeightColVector& GetColVector();
340 
341 #endif
342 
343  /// <summary>
344  /// Get the number of total Degrees-Of-Freedom defined on the source mesh.
345  /// </summary>
347 
348  /// <summary>
349  /// Get the number of total Degrees-Of-Freedom defined on the destination mesh.
350  /// </summary>
352 
353  /// <summary>
354  /// Get the number of local Degrees-Of-Freedom defined on the source mesh.
355  /// </summary>
357 
358  /// <summary>
359  /// Get the number of local Degrees-Of-Freedom defined on the destination mesh.
360  /// </summary>
362 
363  /// <summary>
364  /// Get the number of Degrees-Of-Freedom per element on the source mesh.
365  /// </summary>
367 
368  /// <summary>
369  /// Get the number of Degrees-Of-Freedom per element on the destination mesh.
370  /// </summary>
372 
373  /// <summary>
374  /// Set the number of Degrees-Of-Freedom per element on the source mesh.
375  /// </summary>
376  void SetSourceNDofsPerElement( int ns );
377 
378  /// <summary>
379  /// Get the number of Degrees-Of-Freedom per element on the destination mesh.
380  /// </summary>
381  void SetDestinationNDofsPerElement( int nt );
382 
383  /// <summary>
384  /// Get the global Degrees-Of-Freedom ID on the destination mesh.
385  /// </summary>
386  int GetRowGlobalDoF( int localID ) const;
387 
388  /// <summary>
389  /// Get the index of globaRowDoF.
390  /// </summary>
391  inline int GetIndexOfRowGlobalDoF( int globalRowDoF ) const;
392 
393  /// <summary>
394  /// Get the global Degrees-Of-Freedom ID on the source mesh.
395  /// </summary>
396  int GetColGlobalDoF( int localID ) const;
397 
398  /// <summary>
399  /// Get the index of globaColDoF.
400  /// </summary>
401  inline int GetIndexOfColGlobalDoF( int globalColDoF ) const;
402 
403  /// <summary>
404  /// Apply the weight matrix onto the source vector (tag) provided as input, and return the
405  /// column vector (solution projection) in a tag, after the map application
406  /// Compute: \p tgtVals = A(S->T) * \srcVals, or
407  /// if (transpose) \p tgtVals = [A(T->S)]^T * \srcVals
408  /// </summary>
409  moab::ErrorCode ApplyWeights( moab::Tag srcSolutionTag,
410  moab::Tag tgtSolutionTag,
411  bool transpose = false,
412  CAASType caasType = CAAS_NONE,
413  double default_projection = 0.0 );
414 
415  typedef double ( *sample_function )( double, double );
416 
417  /// <summary>
418  /// Define an analytical solution over the given (source or target) mesh, as specificed in
419  /// the context. This routine will define a tag that is compatible with the specified
420  /// discretization method type and order and sample the solution exactly using the
421  /// analytical function provided by the user.
422  /// </summary>
424  const std::string& solnName,
426  sample_function testFunction,
427  moab::Tag* clonedSolnTag = NULL,
428  std::string cloneSolnName = "" );
429 
430  /// <summary>
431  /// Compute the error between a sampled (exact) solution and a projected solution in various
432  /// error norms.
433  /// </summary>
435  moab::Tag& exactTag,
436  moab::Tag& approxTag,
437  std::map< std::string, double >& metrics,
438  bool verbose = true );
439 
440  moab::ErrorCode fill_col_ids( std::vector< int >& ids_of_interest )
441  {
442  ids_of_interest.reserve( col_gdofmap.size() );
443  // need to add 1
444  for( auto it = col_gdofmap.begin(); it != col_gdofmap.end(); it++ )
445  ids_of_interest.push_back( *it + 1 );
446  return moab::MB_SUCCESS;
447  }
448 
449  moab::ErrorCode set_col_dc_dofs( std::vector< int >& values_entities );
450 
451  moab::ErrorCode set_row_dc_dofs( std::vector< int >& values_entities );
452 
453  // hack
454  void SetMeshInput( Mesh* imesh )
455  {
456  m_meshInput = imesh;
457  };
458 
459  private:
460  template < typename SparseMatrixType >
461  void serializeSparseMatrix( const SparseMatrixType& mat, const std::string& filename );
462 
463  void setup_sizes_dimensions();
464 
465  void CAASLimiter( std::vector< double >& dataCorrectedField,
466  std::vector< double >& dataLowerBound,
467  std::vector< double >& dataUpperBound,
468  double& dMass );
469  double QLTLimiter( int caasIteration,
470  std::vector< double >& dataCorrectedField,
471  std::vector< double >& dataLowerBound,
472  std::vector< double >& dataUpperBound,
473  std::vector< double >& dMassDefect );
474 
475  /// <summary>
476  /// Apply the weight matrix onto the source vector provided as input, and return the column
477  /// vector (solution projection) after the map application
478  /// Compute: \p tgtVals = A(S->T) * \srcVals, or
479  /// if (transpose) \p tgtVals = [A(T->S)]^T * \srcVals
480  /// </summary>
481  moab::ErrorCode ApplyWeights( std::vector< double >& srcVals,
482  std::vector< double >& tgtVals,
483  bool transpose = false );
484 
485 #ifdef MOAB_HAVE_MPI
486  int rearrange_arrays_by_dofs( const std::vector< unsigned int >& gdofmap,
487  DataArray1D< double >& vecFaceArea,
488  DataArray1D< double >& dCenterLon,
489  DataArray1D< double >& dCenterLat,
490  DataArray2D< double >& dVertexLat,
491  DataArray2D< double >& dVertexLon,
492  std::vector< int >& masks,
493  unsigned& N, // this will be output too now
494  int nv,
495  int& maxdof );
496 #endif
497 
498  /// <summary>
499  /// The fundamental remapping operator object.
500  /// </summary>
502 
503 #ifdef MOAB_HAVE_EIGEN3
504 
505  int num_rows, num_cols;
506  WeightMatrix m_weightMatrix;
507  WeightRowVector m_rowVector;
508  WeightColVector m_colVector;
509 
510 #endif
511 
512  /// <summary>
513  /// The reference to the moab::Core object that contains source/target and overlap sets.
514  /// </summary>
516 
517 #ifdef MOAB_HAVE_MPI
518  /// <summary>
519  /// The reference to the parallel communicator object used by the Core object.
520  /// </summary>
521  moab::ParallelComm* m_pcomm;
522 #endif
523 
524  /// <summary>
525  /// The original tag data and local to global DoF mapping to associate matrix values to
526  /// solution <summary>
528  std::vector< unsigned > row_gdofmap, col_gdofmap, srccol_gdofmap;
529 
530  // make it int, because it can be -1 in new logic
532 
533  std::map< int, int > rowMap, colMap;
535 
539 
540  // Key details about the current map
545 
546  Mesh* m_meshInput;
550 
552  int rank, size;
553 };
554 
555 ///////////////////////////////////////////////////////////////////////////////
556 
557 inline int moab::TempestOnlineMap::GetRowGlobalDoF( int localRowID ) const
558 {
559  return row_gdofmap[localRowID];
560 }
561 
562 inline int moab::TempestOnlineMap::GetIndexOfRowGlobalDoF( int globalRowDoF ) const /* 0 based */
563 {
564  return globalRowDoF + 1;
565 }
566 ///////////////////////////////////////////////////////////////////////////////
567 
568 inline int moab::TempestOnlineMap::GetColGlobalDoF( int localColID ) const
569 {
570  return col_gdofmap[localColID];
571 }
572 
573 inline int moab::TempestOnlineMap::GetIndexOfColGlobalDoF( int globalColDoF ) const /* 0 based */
574 {
575  return globalColDoF + 1; // temporary
576 }
577 ///////////////////////////////////////////////////////////////////////////////
578 
580 {
581  return m_nDofsPEl_Src;
582 }
583 
584 ///////////////////////////////////////////////////////////////////////////////
585 
587 {
588  return m_nDofsPEl_Dest;
589 }
590 
591 // setters for m_nDofsPEl_Src, m_nDofsPEl_Dest
593 {
594  m_nDofsPEl_Src = ns;
595 }
596 
598 {
599  m_nDofsPEl_Dest = nt;
600 }
601 
602 ///////////////////////////////////////////////////////////////////////////////
603 #ifdef MOAB_HAVE_EIGEN3
604 
606 {
607  return m_weightMatrix.cols(); // return the global number of rows from the weight matrix
608 }
609 
610 // ///////////////////////////////////////////////////////////////////////////////
611 
613 {
614  return m_weightMatrix.rows(); // return the global number of columns from the weight matrix
615 }
616 
617 ///////////////////////////////////////////////////////////////////////////////
618 
620 {
621  return m_weightMatrix.cols(); // return the local number of rows from the weight matrix
622 }
623 
624 ///////////////////////////////////////////////////////////////////////////////
625 
627 {
628  return m_weightMatrix.rows(); // return the local number of columns from the weight matrix
629 }
630 
631 ///////////////////////////////////////////////////////////////////////////////
632 
633 inline moab::TempestOnlineMap::WeightMatrix& moab::TempestOnlineMap::GetWeightMatrix()
634 {
635  assert( m_weightMatrix.rows() != 0 && m_weightMatrix.cols() != 0 );
636  return m_weightMatrix;
637 }
638 
639 ///////////////////////////////////////////////////////////////////////////////
640 
641 inline moab::TempestOnlineMap::WeightRowVector& moab::TempestOnlineMap::GetRowVector()
642 {
643  assert( m_rowVector.size() != 0 );
644  return m_rowVector;
645 }
646 
647 ///////////////////////////////////////////////////////////////////////////////
648 
649 inline moab::TempestOnlineMap::WeightColVector& moab::TempestOnlineMap::GetColVector()
650 {
651  assert( m_colVector.size() != 0 );
652  return m_colVector;
653 }
654 
655 ///////////////////////////////////////////////////////////////////////////////
656 
657 #endif
658 
659 } // namespace moab
660 
661 #endif