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