Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
TempestOnlineMap.cpp
Go to the documentation of this file.
1 /*
2  * =====================================================================================
3  *
4  * Filename: TempestOnlineMap.hpp
5  *
6  * Description: Interface to the TempestRemap library to compute the consistent,
7  * and accurate high-order conservative remapping weights for overlap
8  * grids on the sphere in climate simulations.
9  *
10  * Author: Vijay S. Mahadevan (vijaysm), [email protected]
11  *
12  * =====================================================================================
13  */
14 
15 #include "Announce.h"
16 #include "DataArray3D.h"
17 #include "FiniteVolumeTools.h"
18 #include "FiniteElementTools.h"
19 #include "TriangularQuadrature.h"
20 #include "GaussQuadrature.h"
21 #include "GaussLobattoQuadrature.h"
22 #include "SparseMatrix.h"
23 #include "STLStringHelper.h"
24 #include "LinearRemapFV.h"
25 
26 #include "LinearRemapSE0.h"
27 #include "LinearRemapFV.h"
28 
30 #include "DebugOutput.hpp"
31 #include "moab/TupleList.hpp"
32 #include "moab/MeshTopoUtil.hpp"
33 
34 #include <fstream>
35 #include <cmath>
36 #include <cstdlib>
37 #include <numeric>
38 #include <algorithm>
39 
40 #ifdef MOAB_HAVE_NETCDFPAR
41 #include "netcdfcpp_par.hpp"
42 #else
43 #include "netcdfcpp.h"
44 #endif
45 
46 // #define USE_NATIVE_TEMPESTREMAP_ROUTINES
47 
48 ///////////////////////////////////////////////////////////////////////////////
49 
50 // #define VERBOSE
51 // #define VVERBOSE
52 // #define CHECK_INCREASING_DOF
53 
54 ///////////////////////////////////////////////////////////////////////////////
55 
56 #define MPI_CHK_ERR( err ) \
57  if( err ) \
58  { \
59  std::cout << "MPI Failure. ErrorCode (" << ( err ) << ") "; \
60  std::cout << "\nMPI Aborting... \n"; \
61  return moab::MB_FAILURE; \
62  }
63 
64 moab::TempestOnlineMap::TempestOnlineMap( moab::TempestRemapper* remapper ) : OfflineMap(), m_remapper( remapper )
65 {
66  // Get the references for the MOAB core objects
68 #ifdef MOAB_HAVE_MPI
69  m_pcomm = m_remapper->get_parallel_communicator();
70 #endif
71 
72  // Update the references to the meshes
74  m_meshInputCov = remapper->GetCoveringMesh();
77 
78  is_parallel = remapper->is_parallel;
79  is_root = remapper->is_root;
80  rank = remapper->rank;
81  size = remapper->size;
82 
83  // set default order
85 
86  // Initialize dimension information from file
87  this->setup_sizes_dimensions();
88 }
89 
91 {
92  if( m_meshInputCov )
93  {
94  std::vector< std::string > dimNames;
95  std::vector< int > dimSizes;
96  dimNames.push_back( "num_elem" );
97  dimSizes.push_back( m_meshInputCov->faces.size() );
98 
99  this->InitializeSourceDimensions( dimNames, dimSizes );
100  }
101 
102  if( m_meshOutput )
103  {
104  std::vector< std::string > dimNames;
105  std::vector< int > dimSizes;
106  dimNames.push_back( "num_elem" );
107  dimSizes.push_back( m_meshOutput->faces.size() );
108 
109  this->InitializeTargetDimensions( dimNames, dimSizes );
110  }
111 }
112 
113 ///////////////////////////////////////////////////////////////////////////////
114 
116 {
117  m_interface = nullptr;
118 #ifdef MOAB_HAVE_MPI
119  m_pcomm = nullptr;
120 #endif
121  m_meshInput = nullptr;
122  m_meshOutput = nullptr;
123  m_meshOverlap = nullptr;
124 }
125 
126 ///////////////////////////////////////////////////////////////////////////////
127 
129  const std::string tgtDofTagName )
130 {
131  moab::ErrorCode rval;
132 
133  int tagSize = 0;
134  tagSize = ( m_eInputType == DiscretizationType_FV ? 1 : m_nDofsPEl_Src * m_nDofsPEl_Src );
135  rval =
136  m_interface->tag_get_handle( srcDofTagName.c_str(), tagSize, MB_TYPE_INTEGER, this->m_dofTagSrc, MB_TAG_ANY );
137 
138  if( rval == moab::MB_TAG_NOT_FOUND && m_eInputType != DiscretizationType_FV )
139  {
140  MB_CHK_SET_ERR( MB_FAILURE, "DoF tag is not set correctly for source mesh." );
141  }
142  else
143  MB_CHK_ERR( rval );
144 
145  tagSize = ( m_eOutputType == DiscretizationType_FV ? 1 : m_nDofsPEl_Dest * m_nDofsPEl_Dest );
146  rval =
147  m_interface->tag_get_handle( tgtDofTagName.c_str(), tagSize, MB_TYPE_INTEGER, this->m_dofTagDest, MB_TAG_ANY );
148  if( rval == moab::MB_TAG_NOT_FOUND && m_eOutputType != DiscretizationType_FV )
149  {
150  MB_CHK_SET_ERR( MB_FAILURE, "DoF tag is not set correctly for target mesh." );
151  }
152  else
153  MB_CHK_ERR( rval );
154 
155  return moab::MB_SUCCESS;
156 }
157 
158 ///////////////////////////////////////////////////////////////////////////////
159 
161  int srcOrder,
162  bool isSrcContinuous,
163  DataArray3D< int >* srcdataGLLNodes,
164  DataArray3D< int >* srcdataGLLNodesSrc,
165  DiscretizationType destType,
166  int destOrder,
167  bool isTgtContinuous,
168  DataArray3D< int >* tgtdataGLLNodes )
169 {
170  moab::ErrorCode rval;
171  std::vector< bool > dgll_cgll_row_ldofmap, dgll_cgll_col_ldofmap, dgll_cgll_covcol_ldofmap;
172  std::vector< int > src_soln_gdofs, locsrc_soln_gdofs, tgt_soln_gdofs;
173 
174  // We are assuming that these are element based tags that are sized: np * np
175  m_srcDiscType = srcType;
176  m_destDiscType = destType;
177  m_input_order = srcOrder;
178  m_output_order = destOrder;
179 
180  bool vprint = is_root && false;
181 
182  // Compute and store the total number of source and target DoFs corresponding
183  // to number of rows and columns in the mapping.
184 #ifdef VVERBOSE
185  {
186  src_soln_gdofs.resize( m_remapper->m_covering_source_entities.size() * m_nDofsPEl_Src * m_nDofsPEl_Src, -1 );
187  rval = m_interface->tag_get_data( m_dofTagSrc, m_remapper->m_covering_source_entities, &src_soln_gdofs[0] );MB_CHK_ERR( rval );
188  locsrc_soln_gdofs.resize( m_remapper->m_source_entities.size() * m_nDofsPEl_Src * m_nDofsPEl_Src );
189  rval = m_interface->tag_get_data( m_dofTagSrc, m_remapper->m_source_entities, &locsrc_soln_gdofs[0] );MB_CHK_ERR( rval );
190  tgt_soln_gdofs.resize( m_remapper->m_target_entities.size() * m_nDofsPEl_Dest * m_nDofsPEl_Dest );
191  rval = m_interface->tag_get_data( m_dofTagDest, m_remapper->m_target_entities, &tgt_soln_gdofs[0] );MB_CHK_ERR( rval );
192 
193  if( is_root )
194  {
195  {
196  std::ofstream output_file( "sourcecov-gids-0.txt" );
197  output_file << "I, GDOF\n";
198  for( unsigned i = 0; i < src_soln_gdofs.size(); ++i )
199  output_file << i << ", " << src_soln_gdofs[i] << "\n";
200 
201  output_file << "ELEMID, IDOF, LDOF, GDOF, NDOF\n";
202  m_nTotDofs_SrcCov = 0;
203  if( isSrcContinuous )
204  dgll_cgll_covcol_ldofmap.resize(
205  m_remapper->m_covering_source_entities.size() * m_nDofsPEl_Src * m_nDofsPEl_Src, false );
206  for( unsigned j = 0; j < m_remapper->m_covering_source_entities.size(); j++ )
207  {
208  for( int p = 0; p < m_nDofsPEl_Src; p++ )
209  {
210  for( int q = 0; q < m_nDofsPEl_Src; q++ )
211  {
212  const int localDOF = ( *srcdataGLLNodes )[p][q][j] - 1;
213  const int offsetDOF = j * m_nDofsPEl_Src * m_nDofsPEl_Src + p * m_nDofsPEl_Src + q;
214  if( isSrcContinuous && !dgll_cgll_covcol_ldofmap[localDOF] )
215  {
216  m_nTotDofs_SrcCov++;
217  dgll_cgll_covcol_ldofmap[localDOF] = true;
218  }
219  output_file << m_remapper->lid_to_gid_covsrc[j] << ", " << offsetDOF << ", " << localDOF
220  << ", " << src_soln_gdofs[offsetDOF] << ", " << m_nTotDofs_SrcCov << "\n";
221  }
222  }
223  }
224  output_file.flush(); // required here
225  output_file.close();
226  dgll_cgll_covcol_ldofmap.clear();
227  }
228 
229  {
230  std::ofstream output_file( "source-gids-0.txt" );
231  output_file << "I, GDOF\n";
232  for( unsigned i = 0; i < locsrc_soln_gdofs.size(); ++i )
233  output_file << i << ", " << locsrc_soln_gdofs[i] << "\n";
234 
235  output_file << "ELEMID, IDOF, LDOF, GDOF, NDOF\n";
236  m_nTotDofs_Src = 0;
237  if( isSrcContinuous )
238  dgll_cgll_col_ldofmap.resize(
239  m_remapper->m_source_entities.size() * m_nDofsPEl_Src * m_nDofsPEl_Src, false );
240  for( unsigned j = 0; j < m_remapper->m_source_entities.size(); j++ )
241  {
242  for( int p = 0; p < m_nDofsPEl_Src; p++ )
243  {
244  for( int q = 0; q < m_nDofsPEl_Src; q++ )
245  {
246  const int localDOF = ( *srcdataGLLNodesSrc )[p][q][j] - 1;
247  const int offsetDOF = j * m_nDofsPEl_Src * m_nDofsPEl_Src + p * m_nDofsPEl_Src + q;
248  if( isSrcContinuous && !dgll_cgll_col_ldofmap[localDOF] )
249  {
250  m_nTotDofs_Src++;
251  dgll_cgll_col_ldofmap[localDOF] = true;
252  }
253  output_file << m_remapper->lid_to_gid_src[j] << ", " << offsetDOF << ", " << localDOF
254  << ", " << locsrc_soln_gdofs[offsetDOF] << ", " << m_nTotDofs_Src << "\n";
255  }
256  }
257  }
258  output_file.flush(); // required here
259  output_file.close();
260  dgll_cgll_col_ldofmap.clear();
261  }
262 
263  {
264  std::ofstream output_file( "target-gids-0.txt" );
265  output_file << "I, GDOF\n";
266  for( unsigned i = 0; i < tgt_soln_gdofs.size(); ++i )
267  output_file << i << ", " << tgt_soln_gdofs[i] << "\n";
268 
269  output_file << "ELEMID, IDOF, GDOF, NDOF\n";
270  m_nTotDofs_Dest = 0;
271 
272  for( unsigned i = 0; i < tgt_soln_gdofs.size(); ++i )
273  {
274  output_file << m_remapper->lid_to_gid_tgt[i] << ", " << i << ", " << tgt_soln_gdofs[i] << ", "
275  << m_nTotDofs_Dest << "\n";
276  m_nTotDofs_Dest++;
277  }
278 
279  output_file.flush(); // required here
280  output_file.close();
281  }
282  }
283  else
284  {
285  {
286  std::ofstream output_file( "sourcecov-gids-1.txt" );
287  output_file << "I, GDOF\n";
288  for( unsigned i = 0; i < src_soln_gdofs.size(); ++i )
289  output_file << i << ", " << src_soln_gdofs[i] << "\n";
290 
291  output_file << "ELEMID, IDOF, LDOF, GDOF, NDOF\n";
292  m_nTotDofs_SrcCov = 0;
293  if( isSrcContinuous )
294  dgll_cgll_covcol_ldofmap.resize(
295  m_remapper->m_covering_source_entities.size() * m_nDofsPEl_Src * m_nDofsPEl_Src, false );
296  for( unsigned j = 0; j < m_remapper->m_covering_source_entities.size(); j++ )
297  {
298  for( int p = 0; p < m_nDofsPEl_Src; p++ )
299  {
300  for( int q = 0; q < m_nDofsPEl_Src; q++ )
301  {
302  const int localDOF = ( *srcdataGLLNodes )[p][q][j] - 1;
303  const int offsetDOF = j * m_nDofsPEl_Src * m_nDofsPEl_Src + p * m_nDofsPEl_Src + q;
304  if( isSrcContinuous && !dgll_cgll_covcol_ldofmap[localDOF] )
305  {
306  m_nTotDofs_SrcCov++;
307  dgll_cgll_covcol_ldofmap[localDOF] = true;
308  }
309  output_file << m_remapper->lid_to_gid_covsrc[j] << ", " << offsetDOF << ", " << localDOF
310  << ", " << src_soln_gdofs[offsetDOF] << ", " << m_nTotDofs_SrcCov << "\n";
311  }
312  }
313  }
314  output_file.flush(); // required here
315  output_file.close();
316  dgll_cgll_covcol_ldofmap.clear();
317  }
318 
319  {
320  std::ofstream output_file( "source-gids-1.txt" );
321  output_file << "I, GDOF\n";
322  for( unsigned i = 0; i < locsrc_soln_gdofs.size(); ++i )
323  output_file << i << ", " << locsrc_soln_gdofs[i] << "\n";
324 
325  output_file << "ELEMID, IDOF, LDOF, GDOF, NDOF\n";
326  m_nTotDofs_Src = 0;
327  if( isSrcContinuous )
328  dgll_cgll_col_ldofmap.resize(
329  m_remapper->m_source_entities.size() * m_nDofsPEl_Src * m_nDofsPEl_Src, false );
330  for( unsigned j = 0; j < m_remapper->m_source_entities.size(); j++ )
331  {
332  for( int p = 0; p < m_nDofsPEl_Src; p++ )
333  {
334  for( int q = 0; q < m_nDofsPEl_Src; q++ )
335  {
336  const int localDOF = ( *srcdataGLLNodesSrc )[p][q][j] - 1;
337  const int offsetDOF = j * m_nDofsPEl_Src * m_nDofsPEl_Src + p * m_nDofsPEl_Src + q;
338  if( isSrcContinuous && !dgll_cgll_col_ldofmap[localDOF] )
339  {
340  m_nTotDofs_Src++;
341  dgll_cgll_col_ldofmap[localDOF] = true;
342  }
343  output_file << m_remapper->lid_to_gid_src[j] << ", " << offsetDOF << ", " << localDOF
344  << ", " << locsrc_soln_gdofs[offsetDOF] << ", " << m_nTotDofs_Src << "\n";
345  }
346  }
347  }
348  output_file.flush(); // required here
349  output_file.close();
350  dgll_cgll_col_ldofmap.clear();
351  }
352 
353  {
354  std::ofstream output_file( "target-gids-1.txt" );
355  output_file << "I, GDOF\n";
356  for( unsigned i = 0; i < tgt_soln_gdofs.size(); ++i )
357  output_file << i << ", " << tgt_soln_gdofs[i] << "\n";
358 
359  output_file << "ELEMID, IDOF, GDOF, NDOF\n";
360  m_nTotDofs_Dest = 0;
361 
362  for( unsigned i = 0; i < tgt_soln_gdofs.size(); ++i )
363  {
364  output_file << m_remapper->lid_to_gid_tgt[i] << ", " << i << ", " << tgt_soln_gdofs[i] << ", "
365  << m_nTotDofs_Dest << "\n";
366  m_nTotDofs_Dest++;
367  }
368 
369  output_file.flush(); // required here
370  output_file.close();
371  }
372  }
373  }
374 #endif
375 
376  // Now compute the mapping and store it for the covering mesh
377  int srcTagSize = ( m_eInputType == DiscretizationType_FV ? 1 : m_nDofsPEl_Src * m_nDofsPEl_Src );
378  if( m_remapper->point_cloud_source )
379  {
380  assert( m_nDofsPEl_Src == 1 );
381  col_gdofmap.resize( m_remapper->m_covering_source_vertices.size(), UINT_MAX );
382  col_dtoc_dofmap.resize( m_remapper->m_covering_source_vertices.size(), UINT_MAX );
383  src_soln_gdofs.resize( m_remapper->m_covering_source_vertices.size(), UINT_MAX );
384  rval = m_interface->tag_get_data( m_dofTagSrc, m_remapper->m_covering_source_vertices, &src_soln_gdofs[0] );MB_CHK_ERR( rval );
385  srcTagSize = 1;
386  }
387  else
388  {
389  col_gdofmap.resize( m_remapper->m_covering_source_entities.size() * srcTagSize, UINT_MAX );
390  col_dtoc_dofmap.resize( m_remapper->m_covering_source_entities.size() * srcTagSize, UINT_MAX );
391  src_soln_gdofs.resize( m_remapper->m_covering_source_entities.size() * srcTagSize, UINT_MAX );
392  rval = m_interface->tag_get_data( m_dofTagSrc, m_remapper->m_covering_source_entities, &src_soln_gdofs[0] );MB_CHK_ERR( rval );
393  }
394 
395  // std::cout << "TOnlineMap: Process: " << rank << " and covering entities = [" <<
396  // col_dofmap.size() << ", " << src_soln_gdofs.size() << "]\n"; MPI_Barrier(MPI_COMM_WORLD);
397 
398 #ifdef ALTERNATE_NUMBERING_IMPLEMENTATION
399  unsigned maxSrcIndx = 0;
400 
401  // for ( unsigned j = 0; j < m_covering_source_entities.size(); j++ )
402  std::vector< int > locdofs( srcTagSize );
403  std::map< Node, moab::EntityHandle > mapLocalMBNodes;
404  double elcoords[3];
405  for( unsigned iel = 0; iel < m_remapper->m_covering_source_entities.size(); ++iel )
406  {
407  EntityHandle eh = m_remapper->m_covering_source_entities[iel];
408  rval = m_interface->get_coords( &eh, 1, elcoords );MB_CHK_ERR( rval );
409  Node elCentroid( elcoords[0], elcoords[1], elcoords[2] );
410  mapLocalMBNodes.insert( std::pair< Node, moab::EntityHandle >( elCentroid, eh ) );
411  }
412 
413  const NodeVector& nodes = m_remapper->m_covering_source->nodes;
414  for( unsigned j = 0; j < m_remapper->m_covering_source->faces.size(); j++ )
415  {
416  const Face& face = m_remapper->m_covering_source->faces[j];
417 
418  Node centroid;
419  centroid.x = centroid.y = centroid.z = 0.0;
420  for( unsigned l = 0; l < face.edges.size(); ++l )
421  {
422  centroid.x += nodes[face[l]].x;
423  centroid.y += nodes[face[l]].y;
424  centroid.z += nodes[face[l]].z;
425  }
426  const double factor = 1.0 / face.edges.size();
427  centroid.x *= factor;
428  centroid.y *= factor;
429  centroid.z *= factor;
430 
431  EntityHandle current_eh;
432  if( mapLocalMBNodes.find( centroid ) != mapLocalMBNodes.end() )
433  {
434  current_eh = mapLocalMBNodes[centroid];
435  }
436 
437  rval = m_interface->tag_get_data( m_dofTagSrc, &current_eh, 1, &locdofs[0] );MB_CHK_ERR( rval );
438  for( int p = 0; p < m_nDofsPEl_Src; p++ )
439  {
440  for( int q = 0; q < m_nDofsPEl_Src; q++ )
441  {
442  const int localDOF = ( *srcdataGLLNodes )[p][q][j] - 1;
443  const int offsetDOF = p * m_nDofsPEl_Src + q;
444  maxSrcIndx = ( localDOF > maxSrcIndx ? localDOF : maxSrcIndx );
445  std::cout << "Col: " << current_eh << ", " << m_remapper->lid_to_gid_covsrc[j] << ", " << offsetDOF
446  << ", " << localDOF << ", " << locdofs[offsetDOF] - 1 << ", " << maxSrcIndx << "\n";
447  }
448  }
449  }
450 #endif
451 
452  m_nTotDofs_SrcCov = 0;
453  if( srcdataGLLNodes == nullptr )
454  { /* we only have a mapping for elements as DoFs */
455  for( unsigned i = 0; i < col_gdofmap.size(); ++i )
456  {
457  assert( src_soln_gdofs[i] > 0 );
458  col_gdofmap[i] = src_soln_gdofs[i] - 1;
459  col_dtoc_dofmap[i] = i;
460  if( vprint ) std::cout << "Col: " << i << ", " << col_gdofmap[i] << "\n";
461  m_nTotDofs_SrcCov++;
462  }
463  }
464  else
465  {
466  if( isSrcContinuous )
467  dgll_cgll_covcol_ldofmap.resize( m_remapper->m_covering_source_entities.size() * srcTagSize, false );
468  // Put these remap coefficients into the SparseMatrix map
469  for( unsigned j = 0; j < m_remapper->m_covering_source_entities.size(); j++ )
470  {
471  for( int p = 0; p < m_nDofsPEl_Src; p++ )
472  {
473  for( int q = 0; q < m_nDofsPEl_Src; q++ )
474  {
475  const int localDOF = ( *srcdataGLLNodes )[p][q][j] - 1;
476  const int offsetDOF = j * srcTagSize + p * m_nDofsPEl_Src + q;
477  if( isSrcContinuous && !dgll_cgll_covcol_ldofmap[localDOF] )
478  {
479  m_nTotDofs_SrcCov++;
480  dgll_cgll_covcol_ldofmap[localDOF] = true;
481  }
482  if( !isSrcContinuous ) m_nTotDofs_SrcCov++;
483  assert( src_soln_gdofs[offsetDOF] > 0 );
484  col_gdofmap[localDOF] = src_soln_gdofs[offsetDOF] - 1;
485  col_dtoc_dofmap[offsetDOF] = localDOF;
486  if( vprint )
487  std::cout << "Col: " << m_remapper->lid_to_gid_covsrc[j] << ", " << offsetDOF << ", "
488  << localDOF << ", " << col_gdofmap[offsetDOF] << ", " << m_nTotDofs_SrcCov << "\n";
489  }
490  }
491  }
492  }
493 
494  if( m_remapper->point_cloud_source )
495  {
496  assert( m_nDofsPEl_Src == 1 );
497  srccol_gdofmap.resize( m_remapper->m_source_vertices.size(), UINT_MAX );
498  srccol_dtoc_dofmap.resize( m_remapper->m_covering_source_vertices.size(), UINT_MAX );
499  locsrc_soln_gdofs.resize( m_remapper->m_source_vertices.size(), UINT_MAX );
500  rval = m_interface->tag_get_data( m_dofTagSrc, m_remapper->m_source_vertices, &locsrc_soln_gdofs[0] );MB_CHK_ERR( rval );
501  }
502  else
503  {
504  srccol_gdofmap.resize( m_remapper->m_source_entities.size() * srcTagSize, UINT_MAX );
505  srccol_dtoc_dofmap.resize( m_remapper->m_source_entities.size() * srcTagSize, UINT_MAX );
506  locsrc_soln_gdofs.resize( m_remapper->m_source_entities.size() * srcTagSize, UINT_MAX );
507  rval = m_interface->tag_get_data( m_dofTagSrc, m_remapper->m_source_entities, &locsrc_soln_gdofs[0] );MB_CHK_ERR( rval );
508  }
509 
510  // Now compute the mapping and store it for the original source mesh
511  m_nTotDofs_Src = 0;
512  if( srcdataGLLNodesSrc == nullptr )
513  { /* we only have a mapping for elements as DoFs */
514  for( unsigned i = 0; i < srccol_gdofmap.size(); ++i )
515  {
516  assert( locsrc_soln_gdofs[i] > 0 );
517  srccol_gdofmap[i] = locsrc_soln_gdofs[i] - 1;
518  srccol_dtoc_dofmap[i] = i;
519  m_nTotDofs_Src++;
520  }
521  }
522  else
523  {
524  if( isSrcContinuous ) dgll_cgll_col_ldofmap.resize( m_remapper->m_source_entities.size() * srcTagSize, false );
525  // Put these remap coefficients into the SparseMatrix map
526  for( unsigned j = 0; j < m_remapper->m_source_entities.size(); j++ )
527  {
528  for( int p = 0; p < m_nDofsPEl_Src; p++ )
529  {
530  for( int q = 0; q < m_nDofsPEl_Src; q++ )
531  {
532  const int localDOF = ( *srcdataGLLNodesSrc )[p][q][j] - 1;
533  const int offsetDOF = j * srcTagSize + p * m_nDofsPEl_Src + q;
534  if( isSrcContinuous && !dgll_cgll_col_ldofmap[localDOF] )
535  {
536  m_nTotDofs_Src++;
537  dgll_cgll_col_ldofmap[localDOF] = true;
538  }
539  if( !isSrcContinuous ) m_nTotDofs_Src++;
540  assert( locsrc_soln_gdofs[offsetDOF] > 0 );
541  srccol_gdofmap[localDOF] = locsrc_soln_gdofs[offsetDOF] - 1;
542  srccol_dtoc_dofmap[offsetDOF] = localDOF;
543  }
544  }
545  }
546  }
547 
548  int tgtTagSize = ( m_eOutputType == DiscretizationType_FV ? 1 : m_nDofsPEl_Dest * m_nDofsPEl_Dest );
549  if( m_remapper->point_cloud_target )
550  {
551  assert( m_nDofsPEl_Dest == 1 );
552  row_gdofmap.resize( m_remapper->m_target_vertices.size(), UINT_MAX );
553  row_dtoc_dofmap.resize( m_remapper->m_target_vertices.size(), UINT_MAX );
554  tgt_soln_gdofs.resize( m_remapper->m_target_vertices.size(), UINT_MAX );
555  rval = m_interface->tag_get_data( m_dofTagDest, m_remapper->m_target_vertices, &tgt_soln_gdofs[0] );MB_CHK_ERR( rval );
556  tgtTagSize = 1;
557  }
558  else
559  {
560  row_gdofmap.resize( m_remapper->m_target_entities.size() * tgtTagSize, UINT_MAX );
561  row_dtoc_dofmap.resize( m_remapper->m_target_entities.size() * tgtTagSize, UINT_MAX );
562  tgt_soln_gdofs.resize( m_remapper->m_target_entities.size() * tgtTagSize, UINT_MAX );
563  rval = m_interface->tag_get_data( m_dofTagDest, m_remapper->m_target_entities, &tgt_soln_gdofs[0] );MB_CHK_ERR( rval );
564  }
565 
566  // Now compute the mapping and store it for the target mesh
567  // To access the GID for each row: row_gdofmap [ row_ldofmap [ 0 : local_ndofs ] ] = GDOF
568  m_nTotDofs_Dest = 0;
569  if( tgtdataGLLNodes == nullptr )
570  { /* we only have a mapping for elements as DoFs */
571  for( unsigned i = 0; i < row_gdofmap.size(); ++i )
572  {
573  assert( tgt_soln_gdofs[i] > 0 );
574  row_gdofmap[i] = tgt_soln_gdofs[i] - 1;
575  row_dtoc_dofmap[i] = i;
576  if( vprint ) std::cout << "Row: " << i << ", " << row_gdofmap[i] << "\n";
577  m_nTotDofs_Dest++;
578  }
579  }
580  else
581  {
582  if( isTgtContinuous ) dgll_cgll_row_ldofmap.resize( m_remapper->m_target_entities.size() * tgtTagSize, false );
583  // Put these remap coefficients into the SparseMatrix map
584  for( unsigned j = 0; j < m_remapper->m_target_entities.size(); j++ )
585  {
586  for( int p = 0; p < m_nDofsPEl_Dest; p++ )
587  {
588  for( int q = 0; q < m_nDofsPEl_Dest; q++ )
589  {
590  const int localDOF = ( *tgtdataGLLNodes )[p][q][j] - 1;
591  const int offsetDOF = j * tgtTagSize + p * m_nDofsPEl_Dest + q;
592  if( isTgtContinuous && !dgll_cgll_row_ldofmap[localDOF] )
593  {
594  m_nTotDofs_Dest++;
595  dgll_cgll_row_ldofmap[localDOF] = true;
596  }
597  if( !isTgtContinuous ) m_nTotDofs_Dest++;
598  assert( tgt_soln_gdofs[offsetDOF] > 0 );
599  row_gdofmap[localDOF] = tgt_soln_gdofs[offsetDOF] - 1;
600  row_dtoc_dofmap[offsetDOF] = localDOF;
601  if( vprint )
602  std::cout << "Row: " << m_remapper->lid_to_gid_tgt[j] << ", " << offsetDOF << ", " << localDOF
603  << ", " << row_gdofmap[offsetDOF] << ", " << m_nTotDofs_Dest << "\n";
604  }
605  }
606  }
607  }
608 
609  // Let us also allocate the local representation of the sparse matrix
610 #if defined( MOAB_HAVE_EIGEN3 ) && defined( VERBOSE )
611  if( vprint )
612  {
613  std::cout << "[" << rank << "]"
614  << "DoFs: row = " << m_nTotDofs_Dest << ", " << row_gdofmap.size() << ", col = " << m_nTotDofs_Src
615  << ", " << m_nTotDofs_SrcCov << ", " << col_gdofmap.size() << "\n";
616  // std::cout << "Max col_dofmap: " << maxcol << ", Min col_dofmap" << mincol << "\n";
617  }
618 #endif
619 
620  // check monotonicity of row_gdofmap and col_gdofmap
621 #ifdef CHECK_INCREASING_DOF
622  for( size_t i = 0; i < row_gdofmap.size() - 1; i++ )
623  {
624  if( row_gdofmap[i] > row_gdofmap[i + 1] )
625  std::cout << " on rank " << rank << " in row_gdofmap[" << i << "]=" << row_gdofmap[i] << " > row_gdofmap["
626  << i + 1 << "]=" << row_gdofmap[i + 1] << " \n";
627  }
628  for( size_t i = 0; i < col_gdofmap.size() - 1; i++ )
629  {
630  if( col_gdofmap[i] > col_gdofmap[i + 1] )
631  std::cout << " on rank " << rank << " in col_gdofmap[" << i << "]=" << col_gdofmap[i] << " > col_gdofmap["
632  << i + 1 << "]=" << col_gdofmap[i + 1] << " \n";
633  }
634 #endif
635 
636  return moab::MB_SUCCESS;
637 }
638 
639 moab::ErrorCode moab::TempestOnlineMap::set_col_dc_dofs( std::vector< int >& values_entities )
640 {
641  // col_gdofmap has global dofs , that should be in the list of values, such that
642  // row_dtoc_dofmap[offsetDOF] = localDOF;
643  // // we need to find col_dtoc_dofmap such that: col_gdofmap[ col_dtoc_dofmap[i] ] == values_entities [i];
644  // we know that col_gdofmap[0..(nbcols-1)] = global_col_dofs -> in values_entities
645  // form first inverse
646 
647  col_dtoc_dofmap.resize( values_entities.size() );
648  for( int j = 0; j < (int)values_entities.size(); j++ )
649  {
650  if( colMap.find( values_entities[j] - 1 ) != colMap.end() )
651  col_dtoc_dofmap[j] = colMap[values_entities[j] - 1];
652  else
653  {
654  col_dtoc_dofmap[j] = -1; // signal that this value should not be used in
655  // std::cout <<"values_entities[j] - 1: " << values_entities[j] - 1 <<" at index j = " << j << " not
656  // found in colMap \n";
657  }
658  }
659  return moab::MB_SUCCESS;
660 }
661 
662 moab::ErrorCode moab::TempestOnlineMap::set_row_dc_dofs( std::vector< int >& values_entities )
663 {
664  // row_dtoc_dofmap = values_entities; // needs to point to local
665  // we need to find row_dtoc_dofmap such that: row_gdofmap[ row_dtoc_dofmap[i] ] == values_entities [i];
666 
667  row_dtoc_dofmap.resize( values_entities.size() );
668  for( int j = 0; j < (int)values_entities.size(); j++ )
669  {
670  if( rowMap.find( values_entities[j] - 1 ) != rowMap.end() )
671  row_dtoc_dofmap[j] = rowMap[values_entities[j] - 1]; // values are 1 based, but rowMap, colMap are not
672  else
673  {
674  row_dtoc_dofmap[j] = -1; // not all values are used
675  // std::cout <<"values_entities[j] - 1: " << values_entities[j] - 1 <<" at index j = " << j << " not
676  // found in rowMap \n";
677  }
678  }
679  return moab::MB_SUCCESS;
680 }
681 ///////////////////////////////////////////////////////////////////////////////
682 
684  std::string strOutputType,
685  const GenerateOfflineMapAlgorithmOptions& mapOptions,
686  const std::string& srcDofTagName,
687  const std::string& tgtDofTagName )
688 {
689  NcError error( NcError::silent_nonfatal );
690 
691  moab::DebugOutput dbgprint( std::cout, rank, 0 );
692  dbgprint.set_prefix( "[TempestOnlineMap]: " );
693  moab::ErrorCode rval;
694 
695  const bool m_bPointCloudSource = ( m_remapper->point_cloud_source );
696  const bool m_bPointCloudTarget = ( m_remapper->point_cloud_target );
697  const bool m_bPointCloud = m_bPointCloudSource || m_bPointCloudTarget;
698 
699  // Build a matrix of source and target discretization so that we know how
700  // to assign the global DoFs in parallel for the mapping weights.
701  // For example,
702  // for FV->FV: the rows represented target DoFs and cols represent source DoFs
703  try
704  {
705  // Check command line parameters (data type arguments)
706  STLStringHelper::ToLower( strInputType );
707  STLStringHelper::ToLower( strOutputType );
708 
709  DiscretizationType eInputType;
710  DiscretizationType eOutputType;
711 
712  if( strInputType == "fv" )
713  {
714  eInputType = DiscretizationType_FV;
715  }
716  else if( strInputType == "cgll" )
717  {
718  eInputType = DiscretizationType_CGLL;
719  }
720  else if( strInputType == "dgll" )
721  {
722  eInputType = DiscretizationType_DGLL;
723  }
724  else if( strInputType == "pcloud" )
725  {
726  eInputType = DiscretizationType_PCLOUD;
727  }
728  else
729  {
730  _EXCEPTION1( "Invalid \"in_type\" value (%s), expected [fv|cgll|dgll]", strInputType.c_str() );
731  }
732 
733  if( strOutputType == "fv" )
734  {
735  eOutputType = DiscretizationType_FV;
736  }
737  else if( strOutputType == "cgll" )
738  {
739  eOutputType = DiscretizationType_CGLL;
740  }
741  else if( strOutputType == "dgll" )
742  {
743  eOutputType = DiscretizationType_DGLL;
744  }
745  else if( strOutputType == "pcloud" )
746  {
747  eOutputType = DiscretizationType_PCLOUD;
748  }
749  else
750  {
751  _EXCEPTION1( "Invalid \"out_type\" value (%s), expected [fv|cgll|dgll]", strOutputType.c_str() );
752  }
753 
754  // set all required input params
755  m_bConserved = !mapOptions.fNoConservation;
756  m_eInputType = eInputType;
757  m_eOutputType = eOutputType;
758 
759  // Method flags
760  std::string strMapAlgorithm( "" );
761  int nMonotoneType = ( mapOptions.fMonotone ) ? ( 1 ) : ( 0 );
762 
763  // Make an index of method arguments
764  std::set< std::string > setMethodStrings;
765  {
766  int iLast = 0;
767  for( size_t i = 0; i <= mapOptions.strMethod.length(); i++ )
768  {
769  if( ( i == mapOptions.strMethod.length() ) || ( mapOptions.strMethod[i] == ';' ) )
770  {
771  std::string strMethodString = mapOptions.strMethod.substr( iLast, i - iLast );
772  STLStringHelper::RemoveWhitespaceInPlace( strMethodString );
773  if( strMethodString.length() > 0 )
774  {
775  setMethodStrings.insert( strMethodString );
776  }
777  iLast = i + 1;
778  }
779  }
780  }
781 
782  for( auto it : setMethodStrings )
783  {
784  // Piecewise constant monotonicity
785  if( it == "mono2" )
786  {
787  if( nMonotoneType != 0 )
788  {
789  _EXCEPTIONT( "Multiple monotonicity specifications found (--mono) or (--method \"mono#\")" );
790  }
791  if( ( m_eInputType == DiscretizationType_FV ) && ( m_eOutputType == DiscretizationType_FV ) )
792  {
793  _EXCEPTIONT( "--method \"mono2\" is only used when remapping to/from CGLL or DGLL grids" );
794  }
795  nMonotoneType = 2;
796 
797  // Piecewise linear monotonicity
798  }
799  else if( it == "mono3" )
800  {
801  if( nMonotoneType != 0 )
802  {
803  _EXCEPTIONT( "Multiple monotonicity specifications found (--mono) or (--method \"mono#\")" );
804  }
805  if( ( m_eInputType == DiscretizationType_FV ) && ( m_eOutputType == DiscretizationType_FV ) )
806  {
807  _EXCEPTIONT( "--method \"mono3\" is only used when remapping to/from CGLL or DGLL grids" );
808  }
809  nMonotoneType = 3;
810 
811  // Volumetric remapping from FV to GLL
812  }
813  else if( it == "volumetric" )
814  {
815  if( ( m_eInputType != DiscretizationType_FV ) || ( m_eOutputType == DiscretizationType_FV ) )
816  {
817  _EXCEPTIONT( "--method \"volumetric\" may only be used for FV->CGLL or FV->DGLL remapping" );
818  }
819  strMapAlgorithm = "volumetric";
820 
821  // Inverse distance mapping
822  }
823  else if( it == "invdist" )
824  {
825  if( ( m_eInputType != DiscretizationType_FV ) || ( m_eOutputType != DiscretizationType_FV ) )
826  {
827  _EXCEPTIONT( "--method \"invdist\" may only be used for FV->FV remapping" );
828  }
829  strMapAlgorithm = "invdist";
830 
831  // Delaunay triangulation mapping
832  }
833  else if( it == "delaunay" )
834  {
835  if( ( m_eInputType != DiscretizationType_FV ) || ( m_eOutputType != DiscretizationType_FV ) )
836  {
837  _EXCEPTIONT( "--method \"delaunay\" may only be used for FV->FV remapping" );
838  }
839  strMapAlgorithm = "delaunay";
840 
841  // Bilinear
842  }
843  else if( it == "bilin" )
844  {
845  if( ( m_eInputType != DiscretizationType_FV ) || ( m_eOutputType != DiscretizationType_FV ) )
846  {
847  _EXCEPTIONT( "--method \"bilin\" may only be used for FV->FV remapping" );
848  }
849  strMapAlgorithm = "fvbilin";
850 
851  // Integrated bilinear (same as mono3 when source grid is CGLL/DGLL)
852  }
853  else if( it == "intbilin" )
854  {
855  if( m_eOutputType != DiscretizationType_FV )
856  {
857  _EXCEPTIONT( "--method \"intbilin\" may only be used when mapping to FV." );
858  }
859  if( m_eInputType == DiscretizationType_FV )
860  {
861  strMapAlgorithm = "fvintbilin";
862  }
863  else
864  {
865  strMapAlgorithm = "mono3";
866  }
867 
868  // Integrated bilinear with generalized Barycentric coordinates
869  }
870  else if( it == "intbilingb" )
871  {
872  if( ( m_eInputType != DiscretizationType_FV ) || ( m_eOutputType != DiscretizationType_FV ) )
873  {
874  _EXCEPTIONT( "--method \"intbilingb\" may only be used for FV->FV remapping" );
875  }
876  strMapAlgorithm = "fvintbilingb";
877  }
878  else
879  {
880  _EXCEPTION1( "Invalid --method argument \"%s\"", it.c_str() );
881  }
882  }
883 
884  m_nDofsPEl_Src =
885  ( m_eInputType == DiscretizationType_FV || m_eInputType == DiscretizationType_PCLOUD ? 1
886  : mapOptions.nPin );
887  m_nDofsPEl_Dest =
888  ( m_eOutputType == DiscretizationType_FV || m_eOutputType == DiscretizationType_PCLOUD ? 1
889  : mapOptions.nPout );
890 
891  rval = SetDOFmapTags( srcDofTagName, tgtDofTagName );MB_CHK_ERR( rval );
892 
893  /// the tag should be created already in the e3sm workflow; if not, create it here
894  Tag areaTag;
895  rval = m_interface->tag_get_handle( "aream", 1, MB_TYPE_DOUBLE, areaTag,
897  if( MB_ALREADY_ALLOCATED == rval )
898  {
899  if( is_root ) dbgprint.printf( 0, "aream tag already defined \n" );
900  }
901 
902  double dTotalAreaInput = 0.0, dTotalAreaOutput = 0.0;
903  if( !m_bPointCloudSource )
904  {
905  // Calculate Input Mesh Face areas
906  if( is_root ) dbgprint.printf( 0, "Calculating input mesh Face areas\n" );
907  double dTotalAreaInput_loc = m_meshInput->CalculateFaceAreas( mapOptions.fSourceConcave );
908  rval = m_interface->tag_set_data( areaTag, m_remapper->m_source_entities, m_meshInput->vecFaceArea );MB_CHK_ERR( rval );
909  dTotalAreaInput = dTotalAreaInput_loc;
910 #ifdef MOAB_HAVE_MPI
911  if( m_pcomm )
912  MPI_Reduce( &dTotalAreaInput_loc, &dTotalAreaInput, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
913 #endif
914  if( is_root ) dbgprint.printf( 0, "Input Mesh Geometric Area: %1.15e\n", dTotalAreaInput );
915 
916  // Input mesh areas
917  m_meshInputCov->CalculateFaceAreas( mapOptions.fSourceConcave );
918  // we do not need to set the area on coverage mesh, only on source and target meshes
919  // rval = m_interface->tag_set_data( areaTag, m_remapper->m_covering_source_entities, m_meshInputCov->vecFaceArea );MB_CHK_ERR( rval );
920  }
921 
922  if( !m_bPointCloudTarget )
923  {
924  // Calculate Output Mesh Face areas
925  if( is_root ) dbgprint.printf( 0, "Calculating output mesh Face areas\n" );
926  double dTotalAreaOutput_loc = m_meshOutput->CalculateFaceAreas( mapOptions.fTargetConcave );
927  dTotalAreaOutput = dTotalAreaOutput_loc;
928 #ifdef MOAB_HAVE_MPI
929  if( m_pcomm )
930  MPI_Reduce( &dTotalAreaOutput_loc, &dTotalAreaOutput, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
931 #endif
932  if( is_root ) dbgprint.printf( 0, "Output Mesh Geometric Area: %1.15e\n", dTotalAreaOutput );
933  rval = m_interface->tag_set_data( areaTag, m_remapper->m_target_entities, m_meshOutput->vecFaceArea );MB_CHK_ERR( rval );
934  }
935 
936  if( !m_bPointCloud )
937  {
938  // Verify that overlap mesh is in the correct order (sanity check)
939  assert( m_meshOverlap->vecSourceFaceIx.size() == m_meshOverlap->vecTargetFaceIx.size() );
940 
941  // Calculate Face areas
942  if( is_root ) dbgprint.printf( 0, "Calculating overlap mesh Face areas\n" );
943  double dTotalAreaOverlap_loc = m_meshOverlap->CalculateFaceAreas( false );
944  double dTotalAreaOverlap = dTotalAreaOverlap_loc;
945 #ifdef MOAB_HAVE_MPI
946  if( m_pcomm )
947  MPI_Reduce( &dTotalAreaOverlap_loc, &dTotalAreaOverlap, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
948 #endif
949  if( is_root ) dbgprint.printf( 0, "Overlap Mesh Area: %1.15e\n", dTotalAreaOverlap );
950 
951  // Correct areas to match the areas calculated in the overlap mesh
952  // if (fCorrectAreas) // In MOAB-TempestRemap, we will always keep this to be true
953  {
954  if( is_root ) dbgprint.printf( 0, "Correcting source/target areas to overlap mesh areas\n" );
955  DataArray1D< double > dSourceArea( m_meshInputCov->faces.size() );
956  DataArray1D< double > dTargetArea( m_meshOutput->faces.size() );
957 
958  assert( m_meshOverlap->vecSourceFaceIx.size() == m_meshOverlap->faces.size() );
959  assert( m_meshOverlap->vecTargetFaceIx.size() == m_meshOverlap->faces.size() );
960  assert( m_meshOverlap->vecFaceArea.GetRows() == m_meshOverlap->faces.size() );
961 
962  assert( m_meshInputCov->vecFaceArea.GetRows() == m_meshInputCov->faces.size() );
963  assert( m_meshOutput->vecFaceArea.GetRows() == m_meshOutput->faces.size() );
964 
965  for( size_t i = 0; i < m_meshOverlap->faces.size(); i++ )
966  {
967  if( m_meshOverlap->vecSourceFaceIx[i] < 0 || m_meshOverlap->vecTargetFaceIx[i] < 0 )
968  continue; // skip this cell since it is ghosted
969 
970  // let us recompute the source/target areas based on overlap mesh areas
971  assert( static_cast< size_t >( m_meshOverlap->vecSourceFaceIx[i] ) < m_meshInputCov->faces.size() );
972  dSourceArea[m_meshOverlap->vecSourceFaceIx[i]] += m_meshOverlap->vecFaceArea[i];
973  assert( static_cast< size_t >( m_meshOverlap->vecTargetFaceIx[i] ) < m_meshOutput->faces.size() );
974  dTargetArea[m_meshOverlap->vecTargetFaceIx[i]] += m_meshOverlap->vecFaceArea[i];
975  }
976 
977  for( size_t i = 0; i < m_meshInputCov->faces.size(); i++ )
978  {
979  if( fabs( dSourceArea[i] - m_meshInputCov->vecFaceArea[i] ) < 1.0e-10 )
980  {
981  m_meshInputCov->vecFaceArea[i] = dSourceArea[i];
982  }
983  }
984  for( size_t i = 0; i < m_meshOutput->faces.size(); i++ )
985  {
986  if( fabs( dTargetArea[i] - m_meshOutput->vecFaceArea[i] ) < 1.0e-10 )
987  {
988  m_meshOutput->vecFaceArea[i] = dTargetArea[i];
989  }
990  }
991  }
992 
993  // Set source mesh areas in map
994  if( !m_bPointCloudSource && eInputType == DiscretizationType_FV )
995  {
996  this->SetSourceAreas( m_meshInputCov->vecFaceArea );
997  if( m_meshInputCov->vecMask.size() )
998  {
999  this->SetSourceMask( m_meshInputCov->vecMask );
1000  }
1001  }
1002 
1003  // Set target mesh areas in map
1004  if( !m_bPointCloudTarget && eOutputType == DiscretizationType_FV )
1005  {
1006  this->SetTargetAreas( m_meshOutput->vecFaceArea );
1007  if( m_meshOutput->vecMask.size() )
1008  {
1009  this->SetTargetMask( m_meshOutput->vecMask );
1010  }
1011  }
1012 
1013  /*
1014  // Recalculate input mesh area from overlap mesh
1015  if (fabs(dTotalAreaOverlap - dTotalAreaInput) > 1.0e-10) {
1016  dbgprint.printf(0, "Overlap mesh only covers a sub-area of the sphere\n");
1017  dbgprint.printf(0, "Recalculating source mesh areas\n");
1018  dTotalAreaInput = m_meshInput->CalculateFaceAreasFromOverlap(m_meshOverlap);
1019  dbgprint.printf(0, "New Input Mesh Geometric Area: %1.15e\n", dTotalAreaInput);
1020  }
1021  */
1022  }
1023 
1024  this->m_pmeshSource = m_meshInputCov;
1025  this->m_pmeshOverlap = m_meshOverlap;
1026 
1027  // Finite volume input / Finite volume output
1028  if( ( eInputType == DiscretizationType_FV ) && ( eOutputType == DiscretizationType_FV ) )
1029  {
1030  // Generate reverse node array and edge map
1031  if( m_meshInputCov->revnodearray.size() == 0 ) m_meshInputCov->ConstructReverseNodeArray();
1032  if( m_meshInputCov->edgemap.size() == 0 ) m_meshInputCov->ConstructEdgeMap( false );
1033 
1034  // Initialize coordinates for map
1035  this->InitializeSourceCoordinatesFromMeshFV( *m_meshInputCov );
1036  this->InitializeTargetCoordinatesFromMeshFV( *m_meshOutput );
1037 
1038  this->m_pdataGLLNodesIn = nullptr;
1039  this->m_pdataGLLNodesOut = nullptr;
1040 
1041  // Finite volume input / Finite element output
1042  rval = this->SetDOFmapAssociation( eInputType, mapOptions.nPin, false, nullptr, nullptr, eOutputType,
1043  mapOptions.nPout, false, nullptr );MB_CHK_ERR( rval );
1044 
1045  // Construct remap for FV-FV
1046  if( is_root ) dbgprint.printf( 0, "Calculating remap weights\n" );
1047 
1048  // Construct OfflineMap
1049  if( strMapAlgorithm == "invdist" )
1050  {
1051  if( is_root ) dbgprint.printf( 0, "Calculating map (invdist)\n" );
1052  if( m_meshInputCov->faces.size() )
1053  LinearRemapFVtoFVInvDist( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, *this );
1054  }
1055  else if( strMapAlgorithm == "delaunay" )
1056  {
1057  if( is_root ) dbgprint.printf( 0, "Calculating map (delaunay)\n" );
1058  if( m_meshInputCov->faces.size() )
1059  LinearRemapTriangulation( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, *this );
1060  }
1061  else if( strMapAlgorithm == "fvintbilin" )
1062  {
1063  if( is_root ) dbgprint.printf( 0, "Calculating map (intbilin)\n" );
1064  if( m_meshInputCov->faces.size() )
1065  LinearRemapIntegratedBilinear( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, *this );
1066  }
1067  else if( strMapAlgorithm == "fvintbilingb" )
1068  {
1069  if( is_root ) dbgprint.printf( 0, "Calculating map (intbilingb)\n" );
1070  if( m_meshInputCov->faces.size() )
1071  LinearRemapIntegratedGeneralizedBarycentric( *m_meshInputCov, *m_meshOutput, *m_meshOverlap,
1072  *this );
1073  }
1074  else if( strMapAlgorithm == "fvbilin" )
1075  {
1076 #ifdef VERBOSE
1077  if( is_root )
1078  {
1079  m_meshInputCov->Write( "SourceMeshMBTR.g" );
1080  m_meshOutput->Write( "TargetMeshMBTR.g" );
1081  }
1082  else
1083  {
1084  m_meshInputCov->Write( "SourceMeshMBTR" + std::to_string( rank ) + ".g" );
1085  m_meshOutput->Write( "TargetMeshMBTR" + std::to_string( rank ) + ".g" );
1086  }
1087 #endif
1088  if( is_root ) dbgprint.printf( 0, "Calculating map (bilin)\n" );
1089  if( m_meshInputCov->faces.size() )
1090  LinearRemapBilinear( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, *this );
1091  }
1092  else
1093  {
1094  if( is_root ) dbgprint.printf( 0, "Calculating conservative FV-FV map\n" );
1095  if( m_meshInputCov->faces.size() )
1096  {
1097 #ifdef USE_NATIVE_TEMPESTREMAP_ROUTINES
1098  LinearRemapFVtoFV( *m_meshInputCov, *m_meshOutput, *m_meshOverlap,
1099  ( mapOptions.fMonotone ) ? ( 1 ) : ( mapOptions.nPin ), *this );
1100 #else
1101  LinearRemapFVtoFV_Tempest_MOAB( ( mapOptions.fMonotone ? 1 : mapOptions.nPin ) );
1102 #endif
1103  }
1104  }
1105  }
1106  else if( eInputType == DiscretizationType_FV )
1107  {
1108  DataArray3D< double > dataGLLJacobian;
1109 
1110  if( is_root ) dbgprint.printf( 0, "Generating output mesh meta data\n" );
1111  double dNumericalArea_loc = GenerateMetaData( *m_meshOutput, mapOptions.nPout, mapOptions.fNoBubble,
1112  dataGLLNodesDest, dataGLLJacobian );
1113 
1114  double dNumericalArea = dNumericalArea_loc;
1115 #ifdef MOAB_HAVE_MPI
1116  if( m_pcomm )
1117  MPI_Reduce( &dNumericalArea_loc, &dNumericalArea, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
1118 #endif
1119  if( is_root ) dbgprint.printf( 0, "Output Mesh Numerical Area: %1.15e\n", dNumericalArea );
1120 
1121  // Initialize coordinates for map
1122  this->InitializeSourceCoordinatesFromMeshFV( *m_meshInputCov );
1123  this->InitializeTargetCoordinatesFromMeshFE( *m_meshOutput, mapOptions.nPout, dataGLLNodesDest );
1124 
1125  this->m_pdataGLLNodesIn = nullptr;
1126  this->m_pdataGLLNodesOut = &dataGLLNodesDest;
1127 
1128  // Generate the continuous Jacobian
1129  bool fContinuous = ( eOutputType == DiscretizationType_CGLL );
1130 
1131  if( eOutputType == DiscretizationType_CGLL )
1132  {
1133  GenerateUniqueJacobian( dataGLLNodesDest, dataGLLJacobian, this->GetTargetAreas() );
1134  }
1135  else
1136  {
1137  GenerateDiscontinuousJacobian( dataGLLJacobian, this->GetTargetAreas() );
1138  }
1139 
1140  // Generate reverse node array and edge map
1141  if( m_meshInputCov->revnodearray.size() == 0 ) m_meshInputCov->ConstructReverseNodeArray();
1142  if( m_meshInputCov->edgemap.size() == 0 ) m_meshInputCov->ConstructEdgeMap( false );
1143 
1144  // Finite volume input / Finite element output
1145  rval = this->SetDOFmapAssociation( eInputType, mapOptions.nPin, false, nullptr, nullptr, eOutputType,
1146  mapOptions.nPout, ( eOutputType == DiscretizationType_CGLL ),
1147  &dataGLLNodesDest );MB_CHK_ERR( rval );
1148 
1149  // Generate remap weights
1150  if( strMapAlgorithm == "volumetric" )
1151  {
1152  if( is_root ) dbgprint.printf( 0, "Calculating remapping weights for FV->GLL (volumetric)\n" );
1153  LinearRemapFVtoGLL_Volumetric( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, dataGLLNodesDest,
1154  dataGLLJacobian, this->GetTargetAreas(), mapOptions.nPin, *this,
1155  nMonotoneType, fContinuous, mapOptions.fNoConservation );
1156  }
1157  else
1158  {
1159  if( is_root ) dbgprint.printf( 0, "Calculating remapping weights for FV->GLL\n" );
1160  LinearRemapFVtoGLL( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, dataGLLNodesDest, dataGLLJacobian,
1161  this->GetTargetAreas(), mapOptions.nPin, *this, nMonotoneType, fContinuous,
1162  mapOptions.fNoConservation );
1163  }
1164  }
1165  else if( ( eInputType == DiscretizationType_PCLOUD ) || ( eOutputType == DiscretizationType_PCLOUD ) )
1166  {
1167  DataArray3D< double > dataGLLJacobian;
1168  if( !m_bPointCloudSource )
1169  {
1170  // Generate reverse node array and edge map
1171  if( m_meshInputCov->revnodearray.size() == 0 ) m_meshInputCov->ConstructReverseNodeArray();
1172  if( m_meshInputCov->edgemap.size() == 0 ) m_meshInputCov->ConstructEdgeMap( false );
1173 
1174  // Initialize coordinates for map
1175  if( eInputType == DiscretizationType_FV )
1176  {
1177  this->InitializeSourceCoordinatesFromMeshFV( *m_meshInputCov );
1178  }
1179  else
1180  {
1181  if( is_root ) dbgprint.printf( 0, "Generating input mesh meta data\n" );
1182  DataArray3D< double > dataGLLJacobianSrc;
1183  GenerateMetaData( *m_meshInputCov, mapOptions.nPin, mapOptions.fNoBubble, dataGLLNodesSrcCov,
1184  dataGLLJacobian );
1185  GenerateMetaData( *m_meshInput, mapOptions.nPin, mapOptions.fNoBubble, dataGLLNodesSrc,
1186  dataGLLJacobianSrc );
1187  }
1188  }
1189  // else { /* Source is a point cloud dataset */ }
1190 
1191  if( !m_bPointCloudTarget )
1192  {
1193  // Generate reverse node array and edge map
1194  if( m_meshOutput->revnodearray.size() == 0 ) m_meshOutput->ConstructReverseNodeArray();
1195  if( m_meshOutput->edgemap.size() == 0 ) m_meshOutput->ConstructEdgeMap( false );
1196 
1197  // Initialize coordinates for map
1198  if( eOutputType == DiscretizationType_FV )
1199  {
1200  this->InitializeSourceCoordinatesFromMeshFV( *m_meshOutput );
1201  }
1202  else
1203  {
1204  if( is_root ) dbgprint.printf( 0, "Generating output mesh meta data\n" );
1205  GenerateMetaData( *m_meshOutput, mapOptions.nPout, mapOptions.fNoBubble, dataGLLNodesDest,
1206  dataGLLJacobian );
1207  }
1208  }
1209  // else { /* Target is a point cloud dataset */ }
1210 
1211  // Finite volume input / Finite element output
1212  rval = this->SetDOFmapAssociation(
1213  eInputType, mapOptions.nPin, ( eInputType == DiscretizationType_CGLL ),
1214  ( m_bPointCloudSource || eInputType == DiscretizationType_FV ? nullptr : &dataGLLNodesSrcCov ),
1215  ( m_bPointCloudSource || eInputType == DiscretizationType_FV ? nullptr : &dataGLLNodesSrc ),
1216  eOutputType, mapOptions.nPout, ( eOutputType == DiscretizationType_CGLL ),
1217  ( m_bPointCloudTarget ? nullptr : &dataGLLNodesDest ) );MB_CHK_ERR( rval );
1218 
1219  // Construct remap
1220  if( is_root ) dbgprint.printf( 0, "Calculating remap weights with Nearest-Neighbor method\n" );
1221  rval = LinearRemapNN_MOAB( true /*use_GID_matching*/, false /*strict_check*/ );MB_CHK_ERR( rval );
1222  }
1223  else if( ( eInputType != DiscretizationType_FV ) && ( eOutputType == DiscretizationType_FV ) )
1224  {
1225  DataArray3D< double > dataGLLJacobianSrc, dataGLLJacobian;
1226 
1227  if( is_root ) dbgprint.printf( 0, "Generating input mesh meta data\n" );
1228  // double dNumericalAreaCov_loc =
1229  GenerateMetaData( *m_meshInputCov, mapOptions.nPin, mapOptions.fNoBubble, dataGLLNodesSrcCov,
1230  dataGLLJacobian );
1231 
1232  double dNumericalArea_loc = GenerateMetaData( *m_meshInput, mapOptions.nPin, mapOptions.fNoBubble,
1233  dataGLLNodesSrc, dataGLLJacobianSrc );
1234 
1235  // if ( is_root ) dbgprint.printf ( 0, "Input Mesh: Coverage Area: %1.15e, Output Area:
1236  // %1.15e\n", dNumericalAreaCov_loc, dTotalAreaOutput_loc );
1237  // assert(dNumericalAreaCov_loc >= dTotalAreaOutput_loc);
1238 
1239  double dNumericalArea = dNumericalArea_loc;
1240 #ifdef MOAB_HAVE_MPI
1241  if( m_pcomm )
1242  MPI_Reduce( &dNumericalArea_loc, &dNumericalArea, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
1243 #endif
1244  if( is_root )
1245  {
1246  dbgprint.printf( 0, "Input Mesh Numerical Area: %1.15e\n", dNumericalArea );
1247 
1248  if( fabs( dNumericalArea - dTotalAreaInput ) > 1.0e-12 )
1249  {
1250  dbgprint.printf( 0, "WARNING: Significant mismatch between input mesh "
1251  "numerical area and geometric area\n" );
1252  }
1253  }
1254 
1255  if( dataGLLNodesSrcCov.GetSubColumns() != m_meshInputCov->faces.size() )
1256  {
1257  _EXCEPTIONT( "Number of element does not match between metadata and "
1258  "input mesh" );
1259  }
1260 
1261  // Initialize coordinates for map
1262  this->InitializeSourceCoordinatesFromMeshFE( *m_meshInputCov, mapOptions.nPin, dataGLLNodesSrcCov );
1263  this->InitializeTargetCoordinatesFromMeshFV( *m_meshOutput );
1264 
1265  // Generate the continuous Jacobian for input mesh
1266  bool fContinuousIn = ( eInputType == DiscretizationType_CGLL );
1267 
1268  if( eInputType == DiscretizationType_CGLL )
1269  {
1270  GenerateUniqueJacobian( dataGLLNodesSrcCov, dataGLLJacobian, this->GetSourceAreas() );
1271  }
1272  else
1273  {
1274  GenerateDiscontinuousJacobian( dataGLLJacobian, this->GetSourceAreas() );
1275  }
1276 
1277  // Finite element input / Finite volume output
1278  rval = this->SetDOFmapAssociation( eInputType, mapOptions.nPin, ( eInputType == DiscretizationType_CGLL ),
1279  &dataGLLNodesSrcCov, &dataGLLNodesSrc, eOutputType, mapOptions.nPout,
1280  false, nullptr );MB_CHK_ERR( rval );
1281 
1282  // Generate remap
1283  if( is_root ) dbgprint.printf( 0, "Calculating remap weights\n" );
1284 
1285  if( strMapAlgorithm == "volumetric" )
1286  {
1287  _EXCEPTIONT( "Unimplemented: Volumetric currently unavailable for"
1288  "GLL input mesh" );
1289  }
1290 
1291  this->m_pdataGLLNodesIn = &dataGLLNodesSrcCov;
1292  this->m_pdataGLLNodesOut = nullptr;
1293 
1294 #ifdef USE_NATIVE_TEMPESTREMAP_ROUTINES
1295  LinearRemapSE4( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, dataGLLNodesSrcCov, dataGLLJacobian,
1296  nMonotoneType, fContinuousIn, mapOptions.fNoConservation, mapOptions.fSparseConstraints,
1297  *this );
1298 #else
1299  LinearRemapSE4_Tempest_MOAB( dataGLLNodesSrcCov, dataGLLJacobian, nMonotoneType, fContinuousIn,
1300  mapOptions.fNoConservation );
1301 #endif
1302  }
1303  else if( ( eInputType != DiscretizationType_FV ) && ( eOutputType != DiscretizationType_FV ) )
1304  {
1305  DataArray3D< double > dataGLLJacobianIn, dataGLLJacobianSrc;
1306  DataArray3D< double > dataGLLJacobianOut;
1307 
1308  // Input metadata
1309  if( is_root ) dbgprint.printf( 0, "Generating input mesh meta data\n" );
1310  GenerateMetaData( *m_meshInputCov, mapOptions.nPin, mapOptions.fNoBubble, dataGLLNodesSrcCov,
1311  dataGLLJacobianIn );
1312 
1313  double dNumericalAreaSrc_loc = GenerateMetaData( *m_meshInput, mapOptions.nPin, mapOptions.fNoBubble,
1314  dataGLLNodesSrc, dataGLLJacobianSrc );
1315  double dNumericalAreaIn = dNumericalAreaSrc_loc;
1316 #ifdef MOAB_HAVE_MPI
1317  if( m_pcomm )
1318  MPI_Reduce( &dNumericalAreaSrc_loc, &dNumericalAreaIn, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
1319 #endif
1320  if( is_root )
1321  {
1322  dbgprint.printf( 0, "Input Mesh Numerical Area: %1.15e\n", dNumericalAreaIn );
1323 
1324  if( fabs( dNumericalAreaIn - dTotalAreaInput ) > 1.0e-12 )
1325  {
1326  dbgprint.printf( 0, "WARNING: Significant mismatch between input mesh "
1327  "numerical area and geometric area\n" );
1328  }
1329  }
1330 
1331  // Output metadata
1332  if( is_root ) dbgprint.printf( 0, "Generating output mesh meta data\n" );
1333  double dNumericalAreaOut_loc = GenerateMetaData( *m_meshOutput, mapOptions.nPout, mapOptions.fNoBubble,
1334  dataGLLNodesDest, dataGLLJacobianOut );
1335  double dNumericalAreaOut = dNumericalAreaOut_loc;
1336 #ifdef MOAB_HAVE_MPI
1337  if( m_pcomm )
1338  MPI_Reduce( &dNumericalAreaOut_loc, &dNumericalAreaOut, 1, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
1339 #endif
1340  if( is_root )
1341  {
1342  dbgprint.printf( 0, "Output Mesh Numerical Area: %1.15e\n", dNumericalAreaOut );
1343 
1344  if( fabs( dNumericalAreaOut - dTotalAreaOutput ) > 1.0e-12 )
1345  {
1346  if( is_root )
1347  dbgprint.printf( 0, "WARNING: Significant mismatch between output mesh "
1348  "numerical area and geometric area\n" );
1349  }
1350  }
1351 
1352  // Initialize coordinates for map
1353  this->InitializeSourceCoordinatesFromMeshFE( *m_meshInputCov, mapOptions.nPin, dataGLLNodesSrcCov );
1354  this->InitializeTargetCoordinatesFromMeshFE( *m_meshOutput, mapOptions.nPout, dataGLLNodesDest );
1355 
1356  // Generate the continuous Jacobian for input mesh
1357  bool fContinuousIn = ( eInputType == DiscretizationType_CGLL );
1358 
1359  if( eInputType == DiscretizationType_CGLL )
1360  {
1361  GenerateUniqueJacobian( dataGLLNodesSrcCov, dataGLLJacobianIn, this->GetSourceAreas() );
1362  }
1363  else
1364  {
1365  GenerateDiscontinuousJacobian( dataGLLJacobianIn, this->GetSourceAreas() );
1366  }
1367 
1368  // Generate the continuous Jacobian for output mesh
1369  bool fContinuousOut = ( eOutputType == DiscretizationType_CGLL );
1370 
1371  if( eOutputType == DiscretizationType_CGLL )
1372  {
1373  GenerateUniqueJacobian( dataGLLNodesDest, dataGLLJacobianOut, this->GetTargetAreas() );
1374  }
1375  else
1376  {
1377  GenerateDiscontinuousJacobian( dataGLLJacobianOut, this->GetTargetAreas() );
1378  }
1379 
1380  // Input Finite Element to Output Finite Element
1381  rval = this->SetDOFmapAssociation( eInputType, mapOptions.nPin, ( eInputType == DiscretizationType_CGLL ),
1382  &dataGLLNodesSrcCov, &dataGLLNodesSrc, eOutputType, mapOptions.nPout,
1383  ( eOutputType == DiscretizationType_CGLL ), &dataGLLNodesDest );MB_CHK_ERR( rval );
1384 
1385  this->m_pdataGLLNodesIn = &dataGLLNodesSrcCov;
1386  this->m_pdataGLLNodesOut = &dataGLLNodesDest;
1387 
1388  // Generate remap
1389  if( is_root ) dbgprint.printf( 0, "Calculating remap weights\n" );
1390 
1391 #ifdef USE_NATIVE_TEMPESTREMAP_ROUTINES
1392  LinearRemapGLLtoGLL_Integrated( *m_meshInputCov, *m_meshOutput, *m_meshOverlap, dataGLLNodesSrcCov,
1393  dataGLLJacobianIn, dataGLLNodesDest, dataGLLJacobianOut,
1394  this->GetTargetAreas(), mapOptions.nPin, mapOptions.nPout, nMonotoneType,
1395  fContinuousIn, fContinuousOut, mapOptions.fSparseConstraints, *this );
1396 #else
1397  LinearRemapGLLtoGLL2_MOAB( dataGLLNodesSrcCov, dataGLLJacobianIn, dataGLLNodesDest, dataGLLJacobianOut,
1398  this->GetTargetAreas(), mapOptions.nPin, mapOptions.nPout, nMonotoneType,
1399  fContinuousIn, fContinuousOut, mapOptions.fNoConservation );
1400 #endif
1401  }
1402  else
1403  {
1404  _EXCEPTIONT( "Not implemented" );
1405  }
1406 
1407 #ifdef MOAB_HAVE_EIGEN3
1408  copy_tempest_sparsemat_to_eigen3();
1409 #endif
1410 
1411 #ifdef MOAB_HAVE_MPI
1412  {
1413  // Remove ghosted entities from overlap set
1414  moab::Range ghostedEnts;
1415  rval = m_remapper->GetOverlapAugmentedEntities( ghostedEnts );MB_CHK_ERR( rval );
1416  moab::EntityHandle m_meshOverlapSet = m_remapper->GetMeshSet( moab::Remapper::OverlapMesh );
1417  rval = m_interface->remove_entities( m_meshOverlapSet, ghostedEnts );MB_CHK_SET_ERR( rval, "Deleting ghosted entities failed" );
1418  }
1419 #endif
1420  // Verify consistency, conservation and monotonicity, globally
1421  if( !mapOptions.fNoCheck )
1422  {
1423  if( is_root ) dbgprint.printf( 0, "Verifying map" );
1424  this->IsConsistent( 1.0e-8 );
1425  if( !mapOptions.fNoConservation ) this->IsConservative( 1.0e-8 );
1426 
1427  if( nMonotoneType != 0 )
1428  {
1429  this->IsMonotone( 1.0e-12 );
1430  }
1431  }
1432  }
1433  catch( Exception& e )
1434  {
1435  dbgprint.printf( 0, "%s", e.ToString().c_str() );
1436  return ( moab::MB_FAILURE );
1437  }
1438  catch( ... )
1439  {
1440  return ( moab::MB_FAILURE );
1441  }
1442  return moab::MB_SUCCESS;
1443 }
1444 
1445 ///////////////////////////////////////////////////////////////////////////////
1446 
1448 {
1449 #ifndef MOAB_HAVE_MPI
1450 
1451  return OfflineMap::IsConsistent( dTolerance );
1452 
1453 #else
1454 
1455  // Get map entries
1456  DataArray1D< int > dataRows;
1457  DataArray1D< int > dataCols;
1458  DataArray1D< double > dataEntries;
1459 
1460  // Calculate row sums
1461  DataArray1D< double > dRowSums;
1462  m_mapRemap.GetEntries( dataRows, dataCols, dataEntries );
1463  dRowSums.Allocate( m_mapRemap.GetRows() );
1464 
1465  for( unsigned i = 0; i < dataRows.GetRows(); i++ )
1466  {
1467  dRowSums[dataRows[i]] += dataEntries[i];
1468  }
1469 
1470  // Verify all row sums are equal to 1
1471  int fConsistent = 0;
1472  for( unsigned i = 0; i < dRowSums.GetRows(); i++ )
1473  {
1474  if( fabs( dRowSums[i] - 1.0 ) > dTolerance )
1475  {
1476  fConsistent++;
1477  int rowGID = row_gdofmap[i];
1478  Announce( "TempestOnlineMap is not consistent in row %i (%1.15e)", rowGID, dRowSums[i] );
1479  }
1480  }
1481 
1482  int ierr;
1483  int fConsistentGlobal = 0;
1484  ierr = MPI_Allreduce( &fConsistent, &fConsistentGlobal, 1, MPI_INT, MPI_SUM, m_pcomm->comm() );
1485  if( ierr != MPI_SUCCESS ) return -1;
1486 
1487  return fConsistentGlobal;
1488 #endif
1489 }
1490 
1491 ///////////////////////////////////////////////////////////////////////////////
1492 
1494 {
1495 #ifndef MOAB_HAVE_MPI
1496 
1497  return OfflineMap::IsConservative( dTolerance );
1498 
1499 #else
1500  // return OfflineMap::IsConservative(dTolerance);
1501 
1502  int ierr;
1503  // Get map entries
1504  DataArray1D< int > dataRows;
1505  DataArray1D< int > dataCols;
1506  DataArray1D< double > dataEntries;
1507  const DataArray1D< double >& dTargetAreas = this->GetTargetAreas();
1508  const DataArray1D< double >& dSourceAreas = this->GetSourceAreas();
1509 
1510  // Calculate column sums
1511  std::vector< int > dColumnsUnique;
1512  std::vector< double > dColumnSums;
1513 
1514  int nColumns = m_mapRemap.GetColumns();
1515  m_mapRemap.GetEntries( dataRows, dataCols, dataEntries );
1516  dColumnSums.resize( m_nTotDofs_SrcCov, 0.0 );
1517  dColumnsUnique.resize( m_nTotDofs_SrcCov, -1 );
1518 
1519  for( unsigned i = 0; i < dataEntries.GetRows(); i++ )
1520  {
1521  dColumnSums[dataCols[i]] += dataEntries[i] * dTargetAreas[dataRows[i]] / dSourceAreas[dataCols[i]];
1522 
1523  assert( dataCols[i] < m_nTotDofs_SrcCov );
1524 
1525  // GID for column DoFs: col_gdofmap[ col_ldofmap [ dataCols[i] ] ]
1526  int colGID = this->GetColGlobalDoF( dataCols[i] ); // col_gdofmap[ col_ldofmap [ dataCols[i] ] ];
1527  // int colGID = col_gdofmap[ col_ldofmap [ dataCols[i] ] ];
1528  dColumnsUnique[dataCols[i]] = colGID;
1529 
1530  // std::cout << "Column dataCols[i]=" << dataCols[i] << " with GID = " << colGID <<
1531  // std::endl;
1532  }
1533 
1534  int rootProc = 0;
1535  std::vector< int > nElementsInProc;
1536  const int nDATA = 3;
1537  nElementsInProc.resize( size * nDATA );
1538  int senddata[nDATA] = { nColumns, m_nTotDofs_SrcCov, m_nTotDofs_Src };
1539  ierr = MPI_Gather( senddata, nDATA, MPI_INT, nElementsInProc.data(), nDATA, MPI_INT, rootProc, m_pcomm->comm() );
1540  if( ierr != MPI_SUCCESS ) return -1;
1541 
1542  int nTotVals = 0, nTotColumns = 0, nTotColumnsUnq = 0;
1543  std::vector< int > dColumnIndices;
1544  std::vector< double > dColumnSourceAreas;
1545  std::vector< double > dColumnSumsTotal;
1546  std::vector< int > displs, rcount;
1547  if( rank == rootProc )
1548  {
1549  displs.resize( size + 1, 0 );
1550  rcount.resize( size, 0 );
1551  int gsum = 0;
1552  for( int ir = 0; ir < size; ++ir )
1553  {
1554  nTotVals += nElementsInProc[ir * nDATA];
1555  nTotColumns += nElementsInProc[ir * nDATA + 1];
1556  nTotColumnsUnq += nElementsInProc[ir * nDATA + 2];
1557 
1558  displs[ir] = gsum;
1559  rcount[ir] = nElementsInProc[ir * nDATA + 1];
1560  gsum += rcount[ir];
1561 
1562  // printf( "%d: nTotColumns: %d, Displs: %d, rcount: %d, gsum = %d\n", ir, nTotColumns, displs[ir], rcount[ir], gsum );
1563  }
1564 
1565  printf( "Total nnz: %d, global source elements = %d\n", nTotVals, gsum );
1566 
1567  dColumnIndices.resize( nTotColumns, -1 );
1568  dColumnSumsTotal.resize( nTotColumns, 0.0 );
1569  // dColumnSourceAreas.resize ( nTotColumns, 0.0 );
1570  }
1571 
1572  // Gather all ColumnSums to root process and accumulate
1573  // We expect that the sums of all columns equate to 1.0 within user specified tolerance
1574  // Need to do a gatherv here since different processes have different number of elements
1575  // MPI_Reduce(&dColumnSums[0], &dColumnSumsTotal[0], m_mapRemap.GetColumns(), MPI_DOUBLE,
1576  // MPI_SUM, 0, m_pcomm->comm());
1577  ierr = MPI_Gatherv( &dColumnsUnique[0], m_nTotDofs_SrcCov, MPI_INT, &dColumnIndices[0], rcount.data(),
1578  displs.data(), MPI_INT, rootProc, m_pcomm->comm() );
1579  if( ierr != MPI_SUCCESS ) return -1;
1580  ierr = MPI_Gatherv( &dColumnSums[0], m_nTotDofs_SrcCov, MPI_DOUBLE, &dColumnSumsTotal[0], rcount.data(),
1581  displs.data(), MPI_DOUBLE, rootProc, m_pcomm->comm() );
1582  if( ierr != MPI_SUCCESS ) return -1;
1583  // ierr = MPI_Gatherv ( &dSourceAreas[0], m_nTotDofs_SrcCov, MPI_DOUBLE, &dColumnSourceAreas[0],
1584  // rcount.data(), displs.data(), MPI_DOUBLE, rootProc, m_pcomm->comm() ); if ( ierr !=
1585  // MPI_SUCCESS ) return -1;
1586 
1587  // Clean out unwanted arrays now
1588  dColumnSums.clear();
1589  dColumnsUnique.clear();
1590 
1591  // Verify all column sums equal the input Jacobian
1592  int fConservative = 0;
1593  if( rank == rootProc )
1594  {
1595  displs[size] = ( nTotColumns );
1596  // std::vector<double> dColumnSumsOnRoot(nTotColumnsUnq, 0.0);
1597  std::map< int, double > dColumnSumsOnRoot;
1598  // std::map<int, double> dColumnSourceAreasOnRoot;
1599  for( int ir = 0; ir < size; ir++ )
1600  {
1601  for( int ips = displs[ir]; ips < displs[ir + 1]; ips++ )
1602  {
1603  if( dColumnIndices[ips] < 0 ) continue;
1604  // printf("%d, %d: dColumnIndices[ips]: %d\n", ir, ips, dColumnIndices[ips]);
1605  assert( dColumnIndices[ips] < nTotColumnsUnq );
1606  dColumnSumsOnRoot[dColumnIndices[ips]] += dColumnSumsTotal[ips]; // / dColumnSourceAreas[ips];
1607  // dColumnSourceAreasOnRoot[ dColumnIndices[ips] ] = dColumnSourceAreas[ips];
1608  // dColumnSourceAreas[ dColumnIndices[ips] ]
1609  }
1610  }
1611 
1612  for( std::map< int, double >::iterator it = dColumnSumsOnRoot.begin(); it != dColumnSumsOnRoot.end(); ++it )
1613  {
1614  // if ( fabs ( it->second - dColumnSourceAreasOnRoot[it->first] ) > dTolerance )
1615  if( fabs( it->second - 1.0 ) > dTolerance )
1616  {
1617  fConservative++;
1618  Announce( "TempestOnlineMap is not conservative in column "
1619  // "%i (%1.15e)", it->first, it->second );
1620  "%i (%1.15e)",
1621  it->first, it->second /* / dColumnSourceAreasOnRoot[it->first] */ );
1622  }
1623  }
1624  }
1625 
1626  // TODO: Just do a broadcast from root instead of a reduction
1627  ierr = MPI_Bcast( &fConservative, 1, MPI_INT, rootProc, m_pcomm->comm() );
1628  if( ierr != MPI_SUCCESS ) return -1;
1629 
1630  return fConservative;
1631 #endif
1632 }
1633 
1634 ///////////////////////////////////////////////////////////////////////////////
1635 
1636 int moab::TempestOnlineMap::IsMonotone( double dTolerance )
1637 {
1638 #ifndef MOAB_HAVE_MPI
1639 
1640  return OfflineMap::IsMonotone( dTolerance );
1641 
1642 #else
1643 
1644  // Get map entries
1645  DataArray1D< int > dataRows;
1646  DataArray1D< int > dataCols;
1647  DataArray1D< double > dataEntries;
1648 
1649  m_mapRemap.GetEntries( dataRows, dataCols, dataEntries );
1650 
1651  // Verify all entries are in the range [0,1]
1652  int fMonotone = 0;
1653  for( unsigned i = 0; i < dataRows.GetRows(); i++ )
1654  {
1655  if( ( dataEntries[i] < -dTolerance ) || ( dataEntries[i] > 1.0 + dTolerance ) )
1656  {
1657  fMonotone++;
1658 
1659  Announce( "TempestOnlineMap is not monotone in entry (%i): %1.15e", i, dataEntries[i] );
1660  }
1661  }
1662 
1663  int ierr;
1664  int fMonotoneGlobal = 0;
1665  ierr = MPI_Allreduce( &fMonotone, &fMonotoneGlobal, 1, MPI_INT, MPI_SUM, m_pcomm->comm() );
1666  if( ierr != MPI_SUCCESS ) return -1;
1667 
1668  return fMonotoneGlobal;
1669 #endif
1670 }
1671 
1672 ///////////////////////////////////////////////////////////////////////////////
1673 
1674 void moab::TempestOnlineMap::ComputeAdjacencyRelations( std::vector< std::unordered_set< int > >& vecAdjFaces,
1675  int nrings,
1676  const Range& entities,
1677  bool useMOABAdjacencies,
1678  Mesh* trMesh )
1679 {
1680  assert( nrings > 0 );
1681  assert( useMOABAdjacencies || trMesh != nullptr );
1682 
1683  const size_t nrows = vecAdjFaces.size();
1684  moab::MeshTopoUtil mtu( m_interface );
1685  for( size_t index = 0; index < nrows; index++ )
1686  {
1687  vecAdjFaces[index].insert( index ); // add self target face first
1688  {
1689  // Compute the adjacent faces to the target face
1690  if( useMOABAdjacencies )
1691  {
1692  moab::Range ents;
1693  // ents.insert( entities.index( entities[index] ) );
1694  ents.insert( entities[index] );
1695  moab::Range adjEnts;
1696  moab::ErrorCode rval = mtu.get_bridge_adjacencies( ents, 0, 2, adjEnts, nrings );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
1697  for( moab::Range::iterator it = adjEnts.begin(); it != adjEnts.end(); ++it )
1698  {
1699  // int adjIndex = m_interface->id_from_handle(*it)-1;
1700  int adjIndex = entities.index( *it );
1701  // printf("rank: %d, Element %lu, entity: %lu, adjIndex %d\n", rank, index, *it, adjIndex);
1702  if( adjIndex >= 0 ) vecAdjFaces[index].insert( adjIndex );
1703  }
1704  }
1705  else
1706  {
1707  /// Vector storing adjacent Faces.
1708  typedef std::pair< int, int > FaceDistancePair;
1709  typedef std::vector< FaceDistancePair > AdjacentFaceVector;
1710  AdjacentFaceVector adjFaces;
1711  Face& face = trMesh->faces[index];
1712  GetAdjacentFaceVectorByEdge( *trMesh, index, nrings * face.edges.size(), adjFaces );
1713 
1714  // Add the adjacent faces to the target face list
1715  for( auto adjFace : adjFaces )
1716  if( adjFace.first >= 0 )
1717  vecAdjFaces[index].insert( adjFace.first ); // map target face to source face
1718  }
1719  }
1720  }
1721 }
1722 
1724  moab::Tag tgtSolutionTag,
1725  bool transpose,
1726  CAASType caasType,
1727  double default_projection )
1728 {
1729  moab::ErrorCode rval;
1730 
1731  std::vector< double > solSTagVals;
1732  std::vector< double > solTTagVals;
1733 
1734  moab::Range sents, tents;
1735  if( m_remapper->point_cloud_source || m_remapper->point_cloud_target )
1736  {
1737  if( m_remapper->point_cloud_source )
1738  {
1739  moab::Range& covSrcEnts = m_remapper->GetMeshVertices( moab::Remapper::CoveringMesh );
1740  solSTagVals.resize( covSrcEnts.size(), default_projection );
1741  sents = covSrcEnts;
1742  }
1743  else
1744  {
1745  moab::Range& covSrcEnts = m_remapper->GetMeshEntities( moab::Remapper::CoveringMesh );
1746  solSTagVals.resize( covSrcEnts.size() * this->GetSourceNDofsPerElement() * this->GetSourceNDofsPerElement(),
1747  default_projection );
1748  sents = covSrcEnts;
1749  }
1750  if( m_remapper->point_cloud_target )
1751  {
1752  moab::Range& tgtEnts = m_remapper->GetMeshVertices( moab::Remapper::TargetMesh );
1753  solTTagVals.resize( tgtEnts.size(), default_projection );
1754  tents = tgtEnts;
1755  }
1756  else
1757  {
1758  moab::Range& tgtEnts = m_remapper->GetMeshEntities( moab::Remapper::TargetMesh );
1759  solTTagVals.resize( tgtEnts.size() * this->GetDestinationNDofsPerElement() *
1760  this->GetDestinationNDofsPerElement(),
1761  default_projection );
1762  tents = tgtEnts;
1763  }
1764  }
1765  else
1766  {
1767  moab::Range& covSrcEnts = m_remapper->GetMeshEntities( moab::Remapper::CoveringMesh );
1768  moab::Range& tgtEnts = m_remapper->GetMeshEntities( moab::Remapper::TargetMesh );
1769  solSTagVals.resize( covSrcEnts.size() * this->GetSourceNDofsPerElement() * this->GetSourceNDofsPerElement(),
1770  default_projection );
1771  solTTagVals.resize( tgtEnts.size() * this->GetDestinationNDofsPerElement() *
1772  this->GetDestinationNDofsPerElement(),
1773  default_projection );
1774 
1775  sents = covSrcEnts;
1776  tents = tgtEnts;
1777  }
1778 
1779  // The tag data is np*np*n_el_src
1780  rval = m_interface->tag_get_data( srcSolutionTag, sents, &solSTagVals[0] );MB_CHK_SET_ERR( rval, "Getting local tag data failed" );
1781 
1782  // Compute the application of weights on the suorce solution data and store it in the
1783  // destination solution vector data Optionally, can also perform the transpose application of
1784  // the weight matrix. Set the 3rd argument to true if this is needed
1785  rval = this->ApplyWeights( solSTagVals, solTTagVals, transpose );MB_CHK_SET_ERR( rval, "Applying remap operator onto source vector data failed" );
1786 
1787  if( caasType != CAAS_NONE )
1788  {
1789  std::string tgtSolutionTagName;
1790  rval = m_interface->tag_get_name( tgtSolutionTag, tgtSolutionTagName );MB_CHK_SET_ERR( rval, "Getting tag name failed" );
1791 
1792  // Perform CAAS iterations iteratively until convergence
1793  constexpr int nmax_caas_iterations = 10;
1794  double mismatch = 1.0;
1795  int caasIteration = 0;
1796  double initialMismatch = 0.0;
1797  while( ( fabs( mismatch / initialMismatch ) > 1e-15 && fabs( mismatch ) > 1e-15 ) &&
1798  caasIteration++ < nmax_caas_iterations ) // iterate until convergence or a maximum of 5 iterations
1799  {
1800  // The tag data is np*np*n_el_dest
1801  rval = m_interface->tag_set_data( tgtSolutionTag, tents, &solTTagVals[0] );MB_CHK_SET_ERR( rval, "Setting local tag data failed" );
1802 
1803  double dMassDiffPostGlobal;
1804  std::pair< double, double > mDefect =
1805  this->ApplyBoundsLimiting( solSTagVals, solTTagVals, caasType, caasIteration, mismatch );
1806 #ifdef MOAB_HAVE_MPI
1807  double dMassDiffPost = mDefect.second;
1808  MPI_Allreduce( &dMassDiffPost, &dMassDiffPostGlobal, 1, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
1809 #else
1810  dMassDiffPostGlobal = mDefect.second;
1811 #endif
1812  if( caasIteration == 1 ) initialMismatch = mDefect.first;
1813  if( m_remapper->verbose && is_root )
1814  {
1815  printf( "Field {%s} -> CAAS iteration: %d, mass defect: %3.4e, post-CAAS: %3.4e\n",
1816  tgtSolutionTagName.c_str(), caasIteration, mDefect.first, dMassDiffPostGlobal );
1817  }
1818  mismatch = dMassDiffPostGlobal;
1819  }
1820  }
1821 
1822  // The tag data is np*np*n_el_dest
1823  rval = m_interface->tag_set_data( tgtSolutionTag, tents, &solTTagVals[0] );MB_CHK_SET_ERR( rval, "Setting local tag data failed" );
1824 
1825  return moab::MB_SUCCESS;
1826 }
1827 
1829  const std::string& solnName,
1831  sample_function testFunction,
1832  moab::Tag* clonedSolnTag,
1833  std::string cloneSolnName )
1834 {
1835  moab::ErrorCode rval;
1836  const bool outputEnabled = ( is_root );
1837  int discOrder;
1838  DiscretizationType discMethod;
1839  // moab::EntityHandle meshset;
1841  Mesh* trmesh;
1842  switch( ctx )
1843  {
1844  case Remapper::SourceMesh:
1845  // meshset = m_remapper->m_covering_source_set;
1846  trmesh = m_remapper->m_covering_source;
1847  entities = ( m_remapper->point_cloud_source ? m_remapper->m_covering_source_vertices
1848  : m_remapper->m_covering_source_entities );
1849  discOrder = m_nDofsPEl_Src;
1850  discMethod = m_eInputType;
1851  break;
1852 
1853  case Remapper::TargetMesh:
1854  // meshset = m_remapper->m_target_set;
1855  trmesh = m_remapper->m_target;
1856  entities =
1857  ( m_remapper->point_cloud_target ? m_remapper->m_target_vertices : m_remapper->m_target_entities );
1858  discOrder = m_nDofsPEl_Dest;
1859  discMethod = m_eOutputType;
1860  break;
1861 
1862  default:
1863  if( outputEnabled )
1864  std::cout << "Invalid context specified for defining an analytical solution tag" << std::endl;
1865  return moab::MB_FAILURE;
1866  }
1867 
1868  // Let us create teh solution tag with appropriate information for name, discretization order
1869  // (DoF space)
1870  rval = m_interface->tag_get_handle( solnName.c_str(), discOrder * discOrder, MB_TYPE_DOUBLE, solnTag,
1871  MB_TAG_DENSE | MB_TAG_CREAT );MB_CHK_ERR( rval );
1872  if( clonedSolnTag != nullptr )
1873  {
1874  if( cloneSolnName.size() == 0 )
1875  {
1876  cloneSolnName = solnName + std::string( "Cloned" );
1877  }
1878  rval = m_interface->tag_get_handle( cloneSolnName.c_str(), discOrder * discOrder, MB_TYPE_DOUBLE,
1879  *clonedSolnTag, MB_TAG_DENSE | MB_TAG_CREAT );MB_CHK_ERR( rval );
1880  }
1881 
1882  // Triangular quadrature rule
1883  const int TriQuadratureOrder = 10;
1884 
1885  if( outputEnabled ) std::cout << "Using triangular quadrature of order " << TriQuadratureOrder << std::endl;
1886 
1887  TriangularQuadratureRule triquadrule( TriQuadratureOrder );
1888 
1889  const int TriQuadraturePoints = triquadrule.GetPoints();
1890 
1891  const DataArray2D< double >& TriQuadratureG = triquadrule.GetG();
1892  const DataArray1D< double >& TriQuadratureW = triquadrule.GetW();
1893 
1894  // Output data
1895  DataArray1D< double > dVar;
1896  DataArray1D< double > dVarMB; // re-arranged local MOAB vector
1897 
1898  // Nodal geometric area
1899  DataArray1D< double > dNodeArea;
1900 
1901  // Calculate element areas
1902  // trmesh->CalculateFaceAreas(fContainsConcaveFaces);
1903 
1904  if( discMethod == DiscretizationType_CGLL || discMethod == DiscretizationType_DGLL )
1905  {
1906  /* Get the spectral points and sample the functionals accordingly */
1907  const bool fGLL = true;
1908  const bool fGLLIntegrate = false;
1909 
1910  // Generate grid metadata
1911  DataArray3D< int > dataGLLNodes;
1912  DataArray3D< double > dataGLLJacobian;
1913 
1914  GenerateMetaData( *trmesh, discOrder, false, dataGLLNodes, dataGLLJacobian );
1915 
1916  // Number of elements
1917  int nElements = trmesh->faces.size();
1918 
1919  // Verify all elements are quadrilaterals
1920  for( int k = 0; k < nElements; k++ )
1921  {
1922  const Face& face = trmesh->faces[k];
1923 
1924  if( face.edges.size() != 4 )
1925  {
1926  _EXCEPTIONT( "Non-quadrilateral face detected; "
1927  "incompatible with --gll" );
1928  }
1929  }
1930 
1931  // Number of unique nodes
1932  int iMaxNode = 0;
1933  for( int i = 0; i < discOrder; i++ )
1934  {
1935  for( int j = 0; j < discOrder; j++ )
1936  {
1937  for( int k = 0; k < nElements; k++ )
1938  {
1939  if( dataGLLNodes[i][j][k] > iMaxNode )
1940  {
1941  iMaxNode = dataGLLNodes[i][j][k];
1942  }
1943  }
1944  }
1945  }
1946 
1947  // Get Gauss-Lobatto quadrature nodes
1948  DataArray1D< double > dG;
1949  DataArray1D< double > dW;
1950 
1951  GaussLobattoQuadrature::GetPoints( discOrder, 0.0, 1.0, dG, dW );
1952 
1953  // Get Gauss quadrature nodes
1954  const int nGaussP = 10;
1955 
1956  DataArray1D< double > dGaussG;
1957  DataArray1D< double > dGaussW;
1958 
1959  GaussQuadrature::GetPoints( nGaussP, 0.0, 1.0, dGaussG, dGaussW );
1960 
1961  // Allocate data
1962  dVar.Allocate( iMaxNode );
1963  dVarMB.Allocate( discOrder * discOrder * nElements );
1964  dNodeArea.Allocate( iMaxNode );
1965 
1966  // Sample data
1967  for( int k = 0; k < nElements; k++ )
1968  {
1969  const Face& face = trmesh->faces[k];
1970 
1971  // Sample data at GLL nodes
1972  if( fGLL )
1973  {
1974  for( int i = 0; i < discOrder; i++ )
1975  {
1976  for( int j = 0; j < discOrder; j++ )
1977  {
1978 
1979  // Apply local map
1980  Node node;
1981  Node dDx1G;
1982  Node dDx2G;
1983 
1984  ApplyLocalMap( face, trmesh->nodes, dG[i], dG[j], node, dDx1G, dDx2G );
1985 
1986  // Sample data at this point
1987  double dNodeLon = atan2( node.y, node.x );
1988  if( dNodeLon < 0.0 )
1989  {
1990  dNodeLon += 2.0 * M_PI;
1991  }
1992  double dNodeLat = asin( node.z );
1993 
1994  double dSample = ( *testFunction )( dNodeLon, dNodeLat );
1995 
1996  dVar[dataGLLNodes[j][i][k] - 1] = dSample;
1997  }
1998  }
1999  // High-order Gaussian integration over basis function
2000  }
2001  else
2002  {
2003  DataArray2D< double > dCoeff( discOrder, discOrder );
2004 
2005  for( int p = 0; p < nGaussP; p++ )
2006  {
2007  for( int q = 0; q < nGaussP; q++ )
2008  {
2009 
2010  // Apply local map
2011  Node node;
2012  Node dDx1G;
2013  Node dDx2G;
2014 
2015  ApplyLocalMap( face, trmesh->nodes, dGaussG[p], dGaussG[q], node, dDx1G, dDx2G );
2016 
2017  // Cross product gives local Jacobian
2018  Node nodeCross = CrossProduct( dDx1G, dDx2G );
2019 
2020  double dJacobian =
2021  sqrt( nodeCross.x * nodeCross.x + nodeCross.y * nodeCross.y + nodeCross.z * nodeCross.z );
2022 
2023  // Find components of quadrature point in basis
2024  // of the first Face
2025  SampleGLLFiniteElement( 0, discOrder, dGaussG[p], dGaussG[q], dCoeff );
2026 
2027  // Sample data at this point
2028  double dNodeLon = atan2( node.y, node.x );
2029  if( dNodeLon < 0.0 )
2030  {
2031  dNodeLon += 2.0 * M_PI;
2032  }
2033  double dNodeLat = asin( node.z );
2034 
2035  double dSample = ( *testFunction )( dNodeLon, dNodeLat );
2036 
2037  // Integrate
2038  for( int i = 0; i < discOrder; i++ )
2039  {
2040  for( int j = 0; j < discOrder; j++ )
2041  {
2042 
2043  double dNodalArea = dCoeff[i][j] * dGaussW[p] * dGaussW[q] * dJacobian;
2044 
2045  dVar[dataGLLNodes[i][j][k] - 1] += dSample * dNodalArea;
2046 
2047  dNodeArea[dataGLLNodes[i][j][k] - 1] += dNodalArea;
2048  }
2049  }
2050  }
2051  }
2052  }
2053  }
2054 
2055  // Divide by area
2056  if( fGLLIntegrate )
2057  {
2058  for( size_t i = 0; i < dVar.GetRows(); i++ )
2059  {
2060  dVar[i] /= dNodeArea[i];
2061  }
2062  }
2063 
2064  // Let us rearrange the data based on DoF ID specification
2065  if( ctx == Remapper::SourceMesh )
2066  {
2067  for( unsigned j = 0; j < entities.size(); j++ )
2068  for( int p = 0; p < discOrder; p++ )
2069  for( int q = 0; q < discOrder; q++ )
2070  {
2071  const int offsetDOF = j * discOrder * discOrder + p * discOrder + q;
2072  dVarMB[offsetDOF] = dVar[col_dtoc_dofmap[offsetDOF]];
2073  }
2074  }
2075  else
2076  {
2077  for( unsigned j = 0; j < entities.size(); j++ )
2078  for( int p = 0; p < discOrder; p++ )
2079  for( int q = 0; q < discOrder; q++ )
2080  {
2081  const int offsetDOF = j * discOrder * discOrder + p * discOrder + q;
2082  dVarMB[offsetDOF] = dVar[row_dtoc_dofmap[offsetDOF]];
2083  }
2084  }
2085 
2086  // Set the tag data
2087  rval = m_interface->tag_set_data( solnTag, entities, &dVarMB[0] );MB_CHK_ERR( rval );
2088  }
2089  else
2090  {
2091  // assert( discOrder == 1 );
2092  if( discMethod == DiscretizationType_FV )
2093  {
2094  /* Compute an element-wise integral to store the sampled solution based on Quadrature
2095  * rules */
2096  // Resize the array
2097  dVar.Allocate( trmesh->faces.size() );
2098 
2099  std::vector< Node >& nodes = trmesh->nodes;
2100 
2101  // Loop through all Faces
2102  for( size_t i = 0; i < trmesh->faces.size(); i++ )
2103  {
2104  const Face& face = trmesh->faces[i];
2105 
2106  // Loop through all sub-triangles
2107  for( size_t j = 0; j < face.edges.size() - 2; j++ )
2108  {
2109 
2110  const Node& node0 = nodes[face[0]];
2111  const Node& node1 = nodes[face[j + 1]];
2112  const Node& node2 = nodes[face[j + 2]];
2113 
2114  // Triangle area
2115  Face faceTri( 3 );
2116  faceTri.SetNode( 0, face[0] );
2117  faceTri.SetNode( 1, face[j + 1] );
2118  faceTri.SetNode( 2, face[j + 2] );
2119 
2120  double dTriangleArea = CalculateFaceArea( faceTri, nodes );
2121 
2122  // Calculate the element average
2123  double dTotalSample = 0.0;
2124 
2125  // Loop through all quadrature points
2126  for( int k = 0; k < TriQuadraturePoints; k++ )
2127  {
2128  Node node( TriQuadratureG[k][0] * node0.x + TriQuadratureG[k][1] * node1.x +
2129  TriQuadratureG[k][2] * node2.x,
2130  TriQuadratureG[k][0] * node0.y + TriQuadratureG[k][1] * node1.y +
2131  TriQuadratureG[k][2] * node2.y,
2132  TriQuadratureG[k][0] * node0.z + TriQuadratureG[k][1] * node1.z +
2133  TriQuadratureG[k][2] * node2.z );
2134 
2135  double dMagnitude = node.Magnitude();
2136  node.x /= dMagnitude;
2137  node.y /= dMagnitude;
2138  node.z /= dMagnitude;
2139 
2140  double dLon = atan2( node.y, node.x );
2141  if( dLon < 0.0 )
2142  {
2143  dLon += 2.0 * M_PI;
2144  }
2145  double dLat = asin( node.z );
2146 
2147  double dSample = ( *testFunction )( dLon, dLat );
2148 
2149  dTotalSample += dSample * TriQuadratureW[k] * dTriangleArea;
2150  }
2151 
2152  dVar[i] += dTotalSample / trmesh->vecFaceArea[i];
2153  }
2154  }
2155  rval = m_interface->tag_set_data( solnTag, entities, &dVar[0] );MB_CHK_ERR( rval );
2156  }
2157  else /* discMethod == DiscretizationType_PCLOUD */
2158  {
2159  /* Get the coordinates of the vertices and sample the functionals accordingly */
2160  std::vector< Node >& nodes = trmesh->nodes;
2161 
2162  // Resize the array
2163  dVar.Allocate( nodes.size() );
2164 
2165  for( size_t j = 0; j < nodes.size(); j++ )
2166  {
2167  Node& node = nodes[j];
2168  double dMagnitude = node.Magnitude();
2169  node.x /= dMagnitude;
2170  node.y /= dMagnitude;
2171  node.z /= dMagnitude;
2172  double dLon = atan2( node.y, node.x );
2173  if( dLon < 0.0 )
2174  {
2175  dLon += 2.0 * M_PI;
2176  }
2177  double dLat = asin( node.z );
2178 
2179  double dSample = ( *testFunction )( dLon, dLat );
2180  dVar[j] = dSample;
2181  }
2182 
2183  rval = m_interface->tag_set_data( solnTag, entities, &dVar[0] );MB_CHK_ERR( rval );
2184  }
2185  }
2186 
2187  return moab::MB_SUCCESS;
2188 }
2189 
2191  moab::Tag& exactTag,
2192  moab::Tag& approxTag,
2193  std::map< std::string, double >& metrics,
2194  bool verbose )
2195 {
2196  moab::ErrorCode rval;
2197  const bool outputEnabled = ( is_root );
2198  int discOrder;
2199  // DiscretizationType discMethod;
2200  // moab::EntityHandle meshset;
2202  // Mesh* trmesh;
2203  switch( ctx )
2204  {
2205  case Remapper::SourceMesh:
2206  // meshset = m_remapper->m_covering_source_set;
2207  // trmesh = m_remapper->m_covering_source;
2208  entities = ( m_remapper->point_cloud_source ? m_remapper->m_covering_source_vertices
2209  : m_remapper->m_covering_source_entities );
2210  discOrder = m_nDofsPEl_Src;
2211  // discMethod = m_eInputType;
2212  break;
2213 
2214  case Remapper::TargetMesh:
2215  // meshset = m_remapper->m_target_set;
2216  // trmesh = m_remapper->m_target;
2217  entities =
2218  ( m_remapper->point_cloud_target ? m_remapper->m_target_vertices : m_remapper->m_target_entities );
2219  discOrder = m_nDofsPEl_Dest;
2220  // discMethod = m_eOutputType;
2221  break;
2222 
2223  default:
2224  if( outputEnabled )
2225  std::cout << "Invalid context specified for defining an analytical solution tag" << std::endl;
2226  return moab::MB_FAILURE;
2227  }
2228 
2229  // Let us create teh solution tag with appropriate information for name, discretization order
2230  // (DoF space)
2231  std::string exactTagName, projTagName;
2232  const int ntotsize = entities.size() * discOrder * discOrder;
2233  std::vector< double > exactSolution( ntotsize, 0.0 ), projSolution( ntotsize, 0.0 );
2234  rval = m_interface->tag_get_name( exactTag, exactTagName );MB_CHK_ERR( rval );
2235  rval = m_interface->tag_get_data( exactTag, entities, &exactSolution[0] );MB_CHK_ERR( rval );
2236  rval = m_interface->tag_get_name( approxTag, projTagName );MB_CHK_ERR( rval );
2237  rval = m_interface->tag_get_data( approxTag, entities, &projSolution[0] );MB_CHK_ERR( rval );
2238 
2239  const auto& ovents = m_remapper->m_overlap_entities;
2240 
2241  std::vector< double > errnorms( 4, 0.0 ), globerrnorms( 4, 0.0 ); // L1Err, L2Err, LinfErr
2242  double sumarea = 0.0;
2243  for( size_t i = 0; i < ovents.size(); ++i )
2244  {
2245  const int srcidx = m_remapper->m_overlap->vecSourceFaceIx[i];
2246  if( srcidx < 0 ) continue; // Skip non-overlapping entities
2247  const int tgtidx = m_remapper->m_overlap->vecTargetFaceIx[i];
2248  if( tgtidx < 0 ) continue; // skip ghost target faces
2249  const double ovarea = m_remapper->m_overlap->vecFaceArea[i];
2250  const double error = fabs( exactSolution[tgtidx] - projSolution[tgtidx] );
2251  errnorms[0] += ovarea * error;
2252  errnorms[1] += ovarea * error * error;
2253  errnorms[3] = ( error > errnorms[3] ? error : errnorms[3] );
2254  sumarea += ovarea;
2255  }
2256  errnorms[2] = sumarea;
2257 #ifdef MOAB_HAVE_MPI
2258  if( m_pcomm )
2259  {
2260  MPI_Reduce( &errnorms[0], &globerrnorms[0], 3, MPI_DOUBLE, MPI_SUM, 0, m_pcomm->comm() );
2261  MPI_Reduce( &errnorms[3], &globerrnorms[3], 1, MPI_DOUBLE, MPI_MAX, 0, m_pcomm->comm() );
2262  }
2263 #else
2264  for( int i = 0; i < 4; ++i )
2265  globerrnorms[i] = errnorms[i];
2266 #endif
2267 
2268  globerrnorms[0] = ( globerrnorms[0] / globerrnorms[2] );
2269  globerrnorms[1] = std::sqrt( globerrnorms[1] / globerrnorms[2] );
2270 
2271  metrics.clear();
2272  metrics["L1Error"] = globerrnorms[0];
2273  metrics["L2Error"] = globerrnorms[1];
2274  metrics["LinfError"] = globerrnorms[3];
2275 
2276  if( verbose && is_root )
2277  {
2278  std::cout << "Error metrics when comparing " << projTagName << " against " << exactTagName << std::endl;
2279  std::cout << "\t Total Intersection area = " << globerrnorms[2] << std::endl;
2280  std::cout << "\t L_1 error = " << globerrnorms[0] << std::endl;
2281  std::cout << "\t L_2 error = " << globerrnorms[1] << std::endl;
2282  std::cout << "\t L_inf error = " << globerrnorms[3] << std::endl;
2283  }
2284 
2285  return moab::MB_SUCCESS;
2286 }