Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
TempestLinearRemap.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 ///
3 /// \file TempestLinearRemap.cpp
4 /// \author Vijay Mahadevan
5 /// \version Mar 08, 2017
6 ///
7 
8 #ifdef WIN32 /* windows */
9 #define _USE_MATH_DEFINES // For M_PI
10 #endif
11 
12 #pragma GCC diagnostic push
13 #pragma GCC diagnostic ignored "-Wdeprecated-copy-with-user-provided-copy"
14 #pragma GCC diagnostic ignored "-Wsign-compare"
15 
16 #include "Announce.h"
17 #include "DataArray3D.h"
18 #include "FiniteElementTools.h"
19 #include "FiniteVolumeTools.h"
20 #include "GaussLobattoQuadrature.h"
21 #include "TriangularQuadrature.h"
22 #include "MathHelper.h"
23 #include "SparseMatrix.h"
24 #include "OverlapMesh.h"
25 
26 #include "DebugOutput.hpp"
27 #include "moab/AdaptiveKDTree.hpp"
28 
30 #include "moab/TupleList.hpp"
31 #include "moab/MeshTopoUtil.hpp"
32 
33 #pragma GCC diagnostic pop
34 
35 #include <fstream>
36 #include <cmath>
37 #include <cstdlib>
38 #include <sstream>
39 #include <numeric>
40 #include <algorithm>
41 #include <unordered_set>
42 
43 // #define VERBOSE
44 #define USE_ComputeAdjacencyRelations
45 
46 /// <summary>
47 /// Face index and distance metric pair.
48 /// </summary>
49 typedef std::pair< int, int > FaceDistancePair;
50 
51 /// <summary>
52 /// Vector storing adjacent Faces.
53 /// </summary>
54 typedef std::vector< FaceDistancePair > AdjacentFaceVector;
55 
56 ///////////////////////////////////////////////////////////////////////////////
57 
58 moab::ErrorCode moab::TempestOnlineMap::LinearRemapNN_MOAB( bool use_GID_matching, bool strict_check )
59 {
60  /* m_mapRemap size = (m_nTotDofs_Dest X m_nTotDofs_SrcCov) */
61 
62 #ifdef VVERBOSE
63  {
64  std::ofstream output_file( "rowcolindices.txt", std::ios::out );
65  output_file << m_nTotDofs_Dest << " " << m_nTotDofs_SrcCov << " " << row_gdofmap.size() << " "
66  << row_ldofmap.size() << " " << col_gdofmap.size() << " " << col_ldofmap.size() << "\n";
67  output_file << "Rows \n";
68  for( unsigned iv = 0; iv < row_gdofmap.size(); iv++ )
69  output_file << row_gdofmap[iv] << " " << row_dofmap[iv] << "\n";
70  output_file << "Cols \n";
71  for( unsigned iv = 0; iv < col_gdofmap.size(); iv++ )
72  output_file << col_gdofmap[iv] << " " << col_dofmap[iv] << "\n";
73  output_file.flush(); // required here
74  output_file.close();
75  }
76 #endif
77 
78  if( use_GID_matching )
79  {
80  std::map< unsigned, unsigned > src_gl;
81  for( unsigned it = 0; it < col_gdofmap.size(); ++it )
82  src_gl[col_gdofmap[it]] = it;
83 
84  std::map< unsigned, unsigned >::iterator iter;
85  for( unsigned it = 0; it < row_gdofmap.size(); ++it )
86  {
87  unsigned row = row_gdofmap[it];
88  iter = src_gl.find( row );
89  if( strict_check && iter == src_gl.end() )
90  {
91  std::cout << "Searching for global target DOF " << row
92  << " but could not find correspondence in source mesh.\n";
93  assert( false );
94  }
95  else if( iter == src_gl.end() )
96  {
97  continue;
98  }
99  else
100  {
101  unsigned icol = src_gl[row];
102  unsigned irow = it;
103 
104  // Set the permutation matrix in local space
105  m_mapRemap( irow, icol ) = 1.0;
106  }
107  }
108 
109  return moab::MB_SUCCESS;
110  }
111  else
112  {
113  /* Create a Kd-tree to perform local queries to find nearest neighbors */
114 
115  return moab::MB_FAILURE;
116  }
117 }
118 
119 ///////////////////////////////////////////////////////////////////////////////
120 
122 {
123  // Order of triangular quadrature rule
124  const int TriQuadRuleOrder = 4;
125 
126  // Verify ReverseNodeArray has been calculated
127  if( m_meshInputCov->faces.size() > 0 && m_meshInputCov->revnodearray.size() == 0 )
128  {
129  _EXCEPTIONT( "ReverseNodeArray has not been calculated for m_meshInputCov" );
130  }
131 
132  // Triangular quadrature rule
133  TriangularQuadratureRule triquadrule( TriQuadRuleOrder );
134 
135  // Number of coefficients needed at this order
136 #ifdef RECTANGULAR_TRUNCATION
137  int nCoefficients = nOrder * nOrder;
138 #endif
139 #ifdef TRIANGULAR_TRUNCATION
140  int nCoefficients = nOrder * ( nOrder + 1 ) / 2;
141 #endif
142 
143  // Number of faces you need
144  const int nRequiredFaceSetSize = nCoefficients;
145 
146  // Fit weight exponent
147  const int nFitWeightsExponent = nOrder + 2;
148 
149  // Announcements
150  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
151  dbgprint.set_prefix( "[LinearRemapFVtoFV_Tempest_MOAB]: " );
152  if( is_root )
153  {
154  dbgprint.printf( 0, "Finite Volume to Finite Volume Projection\n" );
155  dbgprint.printf( 0, "Triangular quadrature rule order %i\n", TriQuadRuleOrder );
156  dbgprint.printf( 0, "Number of coefficients: %i\n", nCoefficients );
157  dbgprint.printf( 0, "Required adjacency set size: %i\n", nRequiredFaceSetSize );
158  dbgprint.printf( 0, "Fit weights exponent: %i\n", nFitWeightsExponent );
159  }
160 
161  // Current overlap face
162  int ixOverlap = 0;
163 #ifdef VERBOSE
164  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
165 #endif
166  DataArray2D< double > dIntArray;
167  DataArray1D< double > dConstraint( nCoefficients );
168 
169  // Loop through all faces on m_meshInputCov
170  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
171  {
172  // Output every 1000 elements
173 #ifdef VERBOSE
174  if( ixFirst % outputFrequency == 0 && is_root )
175  {
176  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
177  }
178 #endif
179  // Find the set of Faces that overlap faceFirst
180  int ixOverlapBegin = ixOverlap;
181  unsigned ixOverlapEnd = ixOverlapBegin;
182 
183  for( ; ixOverlapEnd < m_meshOverlap->faces.size(); ixOverlapEnd++ )
184  {
185  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapEnd] != 0 ) break;
186  }
187 
188  unsigned nOverlapFaces = ixOverlapEnd - ixOverlapBegin;
189 
190  if( nOverlapFaces == 0 ) continue;
191 
192  // Build integration array
193  BuildIntegrationArray( *m_meshInputCov, *m_meshOverlap, triquadrule, ixFirst, ixOverlapBegin, ixOverlapEnd,
194  nOrder, dIntArray );
195 
196  // Set of Faces to use in building the reconstruction and associated
197  // distance metric.
198  AdjacentFaceVector vecAdjFaces;
199 
200  GetAdjacentFaceVectorByEdge( *m_meshInputCov, ixFirst, nRequiredFaceSetSize, vecAdjFaces );
201 
202  // Number of adjacent Faces
203  int nAdjFaces = vecAdjFaces.size();
204 
205  // Determine the conservative constraint equation
206  double dFirstArea = m_meshInputCov->vecFaceArea[ixFirst];
207  dConstraint.Zero();
208  for( int p = 0; p < nCoefficients; p++ )
209  {
210  for( unsigned j = 0; j < nOverlapFaces; j++ )
211  {
212  dConstraint[p] += dIntArray[p][j];
213  }
214  dConstraint[p] /= dFirstArea;
215  }
216 
217  // Build the fit array from the integration operator
218  DataArray2D< double > dFitArray;
219  DataArray1D< double > dFitWeights;
220  DataArray2D< double > dFitArrayPlus;
221 
222  BuildFitArray( *m_meshInputCov, triquadrule, ixFirst, vecAdjFaces, nOrder, nFitWeightsExponent, dConstraint,
223  dFitArray, dFitWeights );
224 
225  // Compute the inverse fit array
226  bool fSuccess = InvertFitArray_Corrected( dConstraint, dFitArray, dFitWeights, dFitArrayPlus );
227 
228  // Multiply integration array and fit array
229  DataArray2D< double > dComposedArray( nAdjFaces, nOverlapFaces );
230  if( fSuccess )
231  {
232  // Multiply integration array and inverse fit array
233  for( int i = 0; i < nAdjFaces; i++ )
234  {
235  for( size_t j = 0; j < nOverlapFaces; j++ )
236  {
237  for( int k = 0; k < nCoefficients; k++ )
238  {
239  dComposedArray( i, j ) += dIntArray( k, j ) * dFitArrayPlus( i, k );
240  }
241  }
242  }
243 
244  // Unable to invert fit array, drop to 1st order. In this case
245  // dFitArrayPlus(0,0) = 1 and all other entries are zero.
246  }
247  else
248  {
249  dComposedArray.Zero();
250  for( size_t j = 0; j < nOverlapFaces; j++ )
251  {
252  dComposedArray( 0, j ) += dIntArray( 0, j );
253  }
254  }
255 
256  // Put composed array into map
257  for( unsigned i = 0; i < vecAdjFaces.size(); i++ )
258  {
259  for( unsigned j = 0; j < nOverlapFaces; j++ )
260  {
261  int& ixFirstFaceLoc = vecAdjFaces[i].first;
262  int& ixSecondFaceLoc = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
263  // int ixFirstFaceGlob = m_remapper->GetGlobalID(moab::Remapper::SourceMesh,
264  // ixFirstFaceLoc); int ixSecondFaceGlob =
265  // m_remapper->GetGlobalID(moab::Remapper::TargetMesh, ixSecondFaceLoc);
266 
267  // signal to not participate, because it is a ghost target
268  if( ixSecondFaceLoc < 0 ) continue; // do not do anything
269 
270  m_mapRemap( ixSecondFaceLoc, ixFirstFaceLoc ) +=
271  dComposedArray[i][j] / m_meshOutput->vecFaceArea[ixSecondFaceLoc];
272  }
273  }
274 
275  // Increment the current overlap index
276  ixOverlap += nOverlapFaces;
277  }
278 
279  return;
280 }
281 
282 ///////////////////////////////////////////////////////////////////////////////
283 
285 {
286  int nrows = m_weightMatrix.rows(); // Number of rows
287  int ncols = m_weightMatrix.cols(); // Number of columns
288  int NNZ = m_weightMatrix.nonZeros(); // Number of non zero values
289 #ifdef MOAB_HAVE_MPI
290  // find out min/max for NNZ, ncols, nrows
291  // should work on std c++ 11
292  int arr3[6] = { NNZ, nrows, ncols, -NNZ, -nrows, -ncols };
293  int rarr3[6];
294  MPI_Reduce( arr3, rarr3, 6, MPI_INT, MPI_MIN, 0, m_pcomm->comm() );
295 
296  int total[3];
297  MPI_Reduce( arr3, total, 3, MPI_INT, MPI_SUM, 0, m_pcomm->comm() );
298  if( !rank )
299  std::cout << "-> Rows (min/max/sum): (" << rarr3[1] << " / " << -rarr3[4] << " / " << total[1] << "), "
300  << " Cols (min/max/sum): (" << rarr3[2] << " / " << -rarr3[5] << " / " << total[2] << "), "
301  << " NNZ (min/max/sum): (" << rarr3[0] << " / " << -rarr3[3] << " / " << total[0] << ")\n";
302 #else
303  std::cout << "-> Rows: " << nrows << ", Cols: " << ncols << ", NNZ: " << NNZ << "\n";
304 #endif
305 }
306 
307 #ifdef MOAB_HAVE_EIGEN3
308 void moab::TempestOnlineMap::copy_tempest_sparsemat_to_eigen3()
309 {
310 #ifndef VERBOSE
311 #define VERBOSE_ACTIVATED
312 // #define VERBOSE
313 #endif
314 
315  /* Should the columns be the global size of the matrix ? */
316  m_weightMatrix.resize( m_nTotDofs_Dest, m_nTotDofs_SrcCov );
317  m_rowVector.resize( m_weightMatrix.rows() );
318  m_colVector.resize( m_weightMatrix.cols() );
319 
320 #ifdef VERBOSE
321  int locrows = std::max( m_mapRemap.GetRows(), m_nTotDofs_Dest );
322  int loccols = std::max( m_mapRemap.GetColumns(), m_nTotDofs_SrcCov );
323 
324  std::cout << m_weightMatrix.rows() << ", " << locrows << ", " << m_weightMatrix.cols() << ", " << loccols << "\n";
325  // assert(m_weightMatrix.rows() == locrows && m_weightMatrix.cols() == loccols);
326 #endif
327 
328  DataArray1D< int > lrows;
329  DataArray1D< int > lcols;
330  DataArray1D< double > lvals;
331  m_mapRemap.GetEntries( lrows, lcols, lvals );
332  size_t locvals = lvals.GetRows();
333 
334  // first matrix
335  typedef Eigen::Triplet< double > Triplet;
336  std::vector< Triplet > tripletList;
337  tripletList.reserve( locvals );
338  for( size_t iv = 0; iv < locvals; iv++ )
339  {
340  tripletList.push_back( Triplet( lrows[iv], lcols[iv], lvals[iv] ) );
341  }
342  m_weightMatrix.setFromTriplets( tripletList.begin(), tripletList.end() );
343  m_weightMatrix.makeCompressed();
344 
345 #ifdef VERBOSE
346  std::stringstream sstr;
347  sstr << "tempestmatrix.txt.0000" << rank;
348  std::ofstream output_file( sstr.str(), std::ios::out );
349  output_file << "0 " << locrows << " 0 " << loccols << "\n";
350  for( unsigned iv = 0; iv < locvals; iv++ )
351  {
352  // output_file << lrows[iv] << " " << row_ldofmap[lrows[iv]] << " " <<
353  // row_gdofmap[row_ldofmap[lrows[iv]]] << " " << col_gdofmap[col_ldofmap[lcols[iv]]] << " "
354  // << lvals[iv] << "\n";
355  output_file << row_gdofmap[row_ldofmap[lrows[iv]]] << " " << col_gdofmap[col_ldofmap[lcols[iv]]] << " "
356  << lvals[iv] << "\n";
357  }
358  output_file.flush(); // required here
359  output_file.close();
360 #endif
361 
362 #ifdef VERBOSE_ACTIVATED
363 #undef VERBOSE_ACTIVATED
364 #undef VERBOSE
365 #endif
366  return;
367 }
368 
369 ///////////////////////////////////////////////////////////////////////////////
370 
371 template < typename T >
372 static std::vector< size_t > sort_indexes( const std::vector< T >& v )
373 {
374  // initialize original index locations
375  std::vector< size_t > idx( v.size() );
376  std::iota( idx.begin(), idx.end(), 0 );
377 
378  // sort indexes based on comparing values in v
379  // using std::stable_sort instead of std::sort
380  // to avoid unnecessary index re-orderings
381  // when v contains elements of equal values
382  std::stable_sort( idx.begin(), idx.end(), [&v]( size_t i1, size_t i2 ) { return fabs( v[i1] ) > fabs( v[i2] ); } );
383 
384  return idx;
385 }
386 
387 double moab::TempestOnlineMap::QLTLimiter( int caasIteration,
388  std::vector< double >& dataCorrectedField,
389  std::vector< double >& dataLowerBound,
390  std::vector< double >& dataUpperBound,
391  std::vector< double >& dMassDefect )
392 {
393  const size_t nrows = dataCorrectedField.size();
394  double dMassL = 0.0;
395  double dMassU = 0.0;
396  std::vector< double > dataCorrection( nrows );
397  double dMassDiffCum = 0.0;
398  double dLMinusU = fabs( dataUpperBound[0] - dataLowerBound[0] );
399  const DataArray1D< double >& dTargetAreas = this->m_remapper->m_target->vecFaceArea;
400 
401  // std::vector< size_t > sortedIdx = sort_indexes( dMassDefect );
402  std::vector< std::unordered_set< int > > vecAdjTargetFaces( nrows );
403  constexpr bool useMOABAdjacencies = true;
404 #ifdef USE_ComputeAdjacencyRelations
405  if( useMOABAdjacencies )
406  ComputeAdjacencyRelations( vecAdjTargetFaces, caasIteration, m_remapper->m_target_entities,
407  useMOABAdjacencies );
408  else
409  ComputeAdjacencyRelations( vecAdjTargetFaces, caasIteration, m_remapper->m_target_entities, useMOABAdjacencies,
410  this->m_remapper->m_target );
411 #else
412  moab::MeshTopoUtil mtu( m_interface );
413  ;
414 #endif
415 
416  for( size_t i = 0; i < nrows; i++ )
417  {
418  // size_t index = sortedIdx[i];
419  size_t index = i;
420  dataCorrection[index] = fmax( dataLowerBound[index], fmin( dataUpperBound[index], 0.0 ) );
421  // dMassDiff[index] = dMassDefect[index] - dTargetAreas[index] * dataCorrection[index];
422  // dMassDiff[index] = dMassDefect[index];
423 
424  dMassL += dTargetAreas[index] * dataLowerBound[index];
425  dMassU += dTargetAreas[index] * dataUpperBound[index];
426  dLMinusU = fmax( dLMinusU, fabs( dataUpperBound[index] - dataLowerBound[index] ) );
427  dMassDiffCum += dMassDefect[index] - dTargetAreas[index] * dataCorrection[index];
428 
429 #ifndef USE_ComputeAdjacencyRelations
430  vecAdjTargetFaces[index].insert( index ); // add self target face first
431  {
432  // Compute the adjacent faces to the target face
433  if( useMOABAdjacencies )
434  {
435  moab::Range ents;
436  // ents.insert( m_remapper->m_target_entities.index( m_remapper->m_target_entities[index] ) );
437  ents.insert( m_remapper->m_target_entities[index] );
438  moab::Range adjEnts;
439  moab::ErrorCode rval = mtu.get_bridge_adjacencies( ents, 0, 2, adjEnts, caasIteration );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
440  for( moab::Range::iterator it = adjEnts.begin(); it != adjEnts.end(); ++it )
441  {
442  // int adjIndex = m_interface->id_from_handle(*it)-1;
443  int adjIndex = m_remapper->m_target_entities.index( *it );
444  // printf("rank: %d, Element %lu, entity: %lu, adjIndex %d\n", rank, index, *it, adjIndex);
445  if( adjIndex >= 0 ) vecAdjTargetFaces[index].insert( adjIndex );
446  }
447  }
448  else
449  {
450  AdjacentFaceVector vecAdjFaces;
451  GetAdjacentFaceVectorByEdge( *this->m_remapper->m_target, index,
452  ( m_output_order + 1 ) * ( m_output_order + 1 ) * ( m_output_order + 1 ),
453  // ( m_output_order + 1 ) * ( m_output_order + 1 ),
454  // ( 4 ) * ( m_output_order + 1 ) * ( m_output_order + 1 ),
455  vecAdjFaces );
456 
457  // Add the adjacent faces to the target face list
458  for( auto adjFace : vecAdjFaces )
459  if( adjFace.first >= 0 )
460  vecAdjTargetFaces[index].insert( adjFace.first ); // map target face to source face
461  }
462  }
463 #endif
464  }
465 
466 #ifdef MOAB_HAVE_MPI
467  std::vector< double > localDefects( 5, 0.0 ), globalDefects( 5, 0.0 );
468  localDefects[0] = dMassL;
469  localDefects[1] = dMassU;
470  localDefects[2] = dMassDiffCum;
471  localDefects[3] = dLMinusU;
472  // localDefects[4] = dMassCorrectU;
473 
474  MPI_Allreduce( localDefects.data(), globalDefects.data(), 4, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
475 
476  dMassL = globalDefects[0];
477  dMassU = globalDefects[1];
478  dMassDiffCum = globalDefects[2];
479  dLMinusU = globalDefects[3];
480  // dMassCorrectU = globalDefects[4];
481 #endif
482 
483  //If the upper and lower bounds are too close together, just clip
484  if( fabs( dMassDiffCum ) < 1e-15 || dLMinusU < 1e-15 )
485  {
486  for( size_t i = 0; i < nrows; i++ )
487  dataCorrectedField[i] += dataCorrection[i];
488  return dMassDiffCum;
489  }
490  else
491  {
492  if( dMassL > dMassDiffCum )
493  {
494  Announce( "Lower bound mass exceeds target mass by %1.15e: CAAS will need another iteration",
495  dMassL - dMassDiffCum );
496  dMassDiffCum = dMassL;
497  // dMass -= dMassL;
498  }
499  else if( dMassU < dMassDiffCum )
500  {
501  Announce( "Target mass exceeds upper bound mass by %1.15e: CAAS will need another iteration",
502  dMassDiffCum - dMassU );
503  dMassDiffCum = dMassU;
504  // dMass -= dMassU;
505  }
506 
507  // TODO: optimize away dataMassVec by a simple transient double within the loop
508  // DataArray1D< double > dataMassVec( nrows ); //vector of mass redistribution
509  for( size_t i = 0; i < nrows; i++ )
510  {
511  // size_t index = sortedIdx[i];
512  size_t index = i;
513  const std::unordered_set< int >& neighbors = vecAdjTargetFaces[index];
514  if( dMassDefect[index] > 0.0 )
515  {
516  double dMassCorrectU = 0.0;
517  for( auto it : neighbors )
518  dMassCorrectU += dTargetAreas[it] * ( dataUpperBound[it] - dataCorrection[it] );
519 
520  // double dMassDiffCumOld = dMassDefect[index];
521  for( auto it : neighbors )
522  dataCorrection[it] +=
523  dMassDefect[index] * ( dataUpperBound[it] - dataCorrection[it] ) / dMassCorrectU;
524  }
525  else
526  {
527  double dMassCorrectL = 0.0;
528  for( auto it : neighbors )
529  dMassCorrectL += dTargetAreas[it] * ( dataCorrection[it] - dataLowerBound[it] );
530 
531  // double dMassDiffCumOld = dMassDefect[index];
532  for( auto it : neighbors )
533  dataCorrection[it] +=
534  dMassDefect[index] * ( dataCorrection[it] - dataLowerBound[it] ) / dMassCorrectL;
535  }
536  }
537 
538  for( size_t i = 0; i < nrows; i++ )
539  dataCorrectedField[i] += dataCorrection[i];
540  }
541 
542  return dMassDiffCum;
543 }
544 
545 void moab::TempestOnlineMap::CAASLimiter( std::vector< double >& dataCorrectedField,
546  std::vector< double >& dataLowerBound,
547  std::vector< double >& dataUpperBound,
548  double& dMass )
549 {
550  const size_t nrows = dataCorrectedField.size();
551  double dMassL = 0.0;
552  double dMassU = 0.0;
553  std::vector< double > dataCorrection( nrows );
554  const DataArray1D< double >& dTargetAreas = this->m_remapper->m_target->vecFaceArea;
555  double dMassDiff = dMass;
556  double dLMinusU = fabs( dataUpperBound[0] - dataLowerBound[0] );
557  double dMassCorrectU = 0.0;
558  double dMassCorrectL = 0.0;
559  for( size_t i = 0; i < nrows; i++ )
560  {
561  dataCorrection[i] = fmax( dataLowerBound[i], fmin( dataUpperBound[i], 0.0 ) );
562  dMassL += dTargetAreas[i] * dataLowerBound[i];
563  dMassU += dTargetAreas[i] * dataUpperBound[i];
564  dMassDiff -= dTargetAreas[i] * dataCorrection[i];
565  dLMinusU = fmax( dLMinusU, fabs( dataUpperBound[i] - dataLowerBound[i] ) );
566  dMassCorrectL += dTargetAreas[i] * ( dataCorrection[i] - dataLowerBound[i] );
567  dMassCorrectU += dTargetAreas[i] * ( dataUpperBound[i] - dataCorrection[i] );
568  }
569 
570 #ifdef MOAB_HAVE_MPI
571  std::vector< double > localDefects( 5, 0.0 ), globalDefects( 5, 0.0 );
572  localDefects[0] = dMassL;
573  localDefects[1] = dMassU;
574  localDefects[2] = dMassDiff;
575  localDefects[3] = dMassCorrectL;
576  localDefects[4] = dMassCorrectU;
577 
578  MPI_Allreduce( localDefects.data(), globalDefects.data(), 5, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
579 
580  dMassL = globalDefects[0];
581  dMassU = globalDefects[1];
582  dMassDiff = globalDefects[2];
583  dMassCorrectL = globalDefects[3];
584  dMassCorrectU = globalDefects[4];
585 #endif
586 
587  //If the upper and lower bounds are too close together, just clip
588  if( fabs( dMassDiff ) < 1e-15 || fabs( dLMinusU ) < 1e-15 )
589  {
590  for( size_t i = 0; i < nrows; i++ )
591  dataCorrectedField[i] += dataCorrection[i];
592  return;
593  }
594  else
595  {
596  if( dMassL > dMassDiff )
597  {
598  Announce( "%d: Lower bound mass exceeds target mass by %1.15e: CAAS will need another iteration", rank,
599  dMassL - dMassDiff );
600  dMassDiff = dMassL;
601  dMass -= dMassL;
602  }
603  else if( dMassU < dMassDiff )
604  {
605  Announce( "%d: Target mass exceeds upper bound mass by %1.15e: CAAS will need another iteration", rank,
606  dMassDiff - dMassU );
607  dMassDiff = dMassU;
608  dMass -= dMassU;
609  }
610 
611  // TODO: optimize away dataMassVec by a simple transient double within the loop
612  DataArray1D< double > dataMassVec( nrows ); //vector of mass redistribution
613  if( dMassDiff > 0.0 )
614  {
615  for( size_t i = 0; i < nrows; i++ )
616  {
617  dataMassVec[i] = ( dataUpperBound[i] - dataCorrection[i] ) / dMassCorrectU;
618  dataCorrection[i] += dMassDiff * dataMassVec[i];
619  }
620  }
621  else
622  {
623  for( size_t i = 0; i < nrows; i++ )
624  {
625  dataMassVec[i] = ( dataCorrection[i] - dataLowerBound[i] ) / dMassCorrectL;
626  dataCorrection[i] += dMassDiff * dataMassVec[i];
627  }
628  }
629 
630  for( size_t i = 0; i < nrows; i++ )
631  dataCorrectedField[i] += dataCorrection[i];
632  }
633 
634  return;
635 }
636 
637 std::pair< double, double > moab::TempestOnlineMap::ApplyBoundsLimiting( std::vector< double >& dataInDouble,
638  std::vector< double >& dataOutDouble,
639  CAASType caasType,
640  int caasIteration,
641  double mismatch )
642 {
643  // Currently only implemented for FV to FV remapping
644  // We should generalize this to other types of remapping
645  assert( !dataGLLNodesSrcCov.IsAttached() && !dataGLLNodesDest.IsAttached() );
646 
647  std::pair< double, double > massDefect( 0.0, 0.0 );
648 
649  // Check if the source and target data are of the same size
650  const size_t nTargetCount = dataOutDouble.size();
651  const DataArray1D< double >& m_dOverlapAreas = this->m_remapper->m_overlap->vecFaceArea;
652 
653  // Apply the offline map to the data
654  double dMassDiff = 0.0;
655  std::vector< double > x( nTargetCount );
656  std::vector< double > dataLowerBound( nTargetCount );
657  std::vector< double > dataUpperBound( nTargetCount );
658  std::vector< double > massVector( nTargetCount );
659  std::vector< std::unordered_set< int > > vecSourceOvTarget( nTargetCount );
660 
661 #undef USE_ComputeAdjacencyRelations
662  constexpr bool useMOABAdjacencies = true;
663 #ifdef USE_ComputeAdjacencyRelations
664  // Compute the adjacent faces to the source face
665  // However, calling MOAB to do this does not work correctly as we need ixS to be the index
666  // Cannot just iterate over all entities in the source covering mesh
667  if( caasType == CAAS_QLT || caasType == CAAS_LOCAL_ADJACENT )
668  {
669  if( useMOABAdjacencies )
670  {
671  moab::ErrorCode rval =
672  ComputeAdjacencyRelations( vecSourceOvTarget, caasIteration, m_remapper->m_covering_source_entities,
673  useMOABAdjacencies );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
674  }
675  else
676  {
677  moab::ErrorCode rval =
678  ComputeAdjacencyRelations( vecSourceOvTarget, caasIteration, m_remapper->m_covering_source_entities,
679  useMOABAdjacencies, m_meshInputCov );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
680  }
681  }
682 #else
683  moab::MeshTopoUtil mtu( m_interface );
684 #endif
685 
686  // Initialize the bounds on the given source and target data
687  double dSourceMin = dataInDouble[0];
688  double dSourceMax = dataInDouble[0];
689  double dTargetMin = dataOutDouble[0];
690  double dTargetMax = dataOutDouble[0];
691  for( size_t i = 0; i < m_meshOverlap->faces.size(); i++ )
692  {
693  const int ixS = m_meshOverlap->vecSourceFaceIx[i];
694  const int ixT = m_meshOverlap->vecTargetFaceIx[i];
695 
696  if( ixT < 0 ) continue; // skip ghost target faces
697 
698  assert( m_dOverlapAreas[i] > 0.0 );
699  assert( ixS >= 0 );
700  assert( ixT >= 0 );
701 
702 #ifndef USE_ComputeAdjacencyRelations
703  // Compute the adjacent faces to the target face
704  vecSourceOvTarget[ixT].insert( ixS ); // map target face to source face
705  if( ( caasType == CAAS_QLT || caasType == CAAS_LOCAL_ADJACENT ) )
706  {
707  if( useMOABAdjacencies )
708  {
709  moab::Range ents;
710  ents.insert( m_remapper->m_covering_source_entities[ixS] );
711  moab::Range adjEnts;
712  moab::ErrorCode rval = mtu.get_bridge_adjacencies( ents, 0, 2, adjEnts, caasIteration );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
713  for( moab::Range::iterator it = adjEnts.begin(); it != adjEnts.end(); ++it )
714  {
715  int adjIndex = m_remapper->m_covering_source_entities.index( *it );
716  if( adjIndex >= 0 ) vecSourceOvTarget[ixT].insert( adjIndex );
717  }
718  }
719  else
720  {
721  // Compute the adjacent faces to the target face
722  AdjacentFaceVector vecAdjFaces;
723  GetAdjacentFaceVectorByEdge( *m_meshInputCov, ixS,
724  ( caasIteration ) * ( m_input_order + 1 ) * ( m_input_order + 1 ),
725  vecAdjFaces );
726 
727  //Compute min/max over neighboring faces
728  for( size_t iadj = 0; iadj < vecAdjFaces.size(); iadj++ )
729  vecSourceOvTarget[ixT].insert( vecAdjFaces[iadj].first ); // map target face to source face
730  }
731  }
732 #endif
733 
734  // Update the min and max values of the source data
735  dSourceMax = fmax( dSourceMax, dataInDouble[ixS] );
736  dSourceMin = fmin( dSourceMin, dataInDouble[ixS] );
737 
738  // Update the min and max values of the target data
739  dTargetMin = fmin( dTargetMin, dataOutDouble[ixT] );
740  dTargetMax = fmax( dTargetMax, dataOutDouble[ixT] );
741 
742  const double locMassDiff = ( dataInDouble[ixS] * m_dOverlapAreas[i] ) - // source mass
743  ( dataOutDouble[ixT] * m_dOverlapAreas[i] ); // target mass
744 
745  // Update the mass difference between source and target faces
746  // linked to the overlap mesh element
747  dMassDiff += locMassDiff; // target mass
748  massVector[ixT] += locMassDiff;
749  }
750 
751 #ifdef MOAB_HAVE_MPI
752  std::vector< double > localMinMaxDefects( 5, 0.0 ), globalMinMaxDefects( 5, 0.0 );
753  localMinMaxDefects[0] = dSourceMin;
754  localMinMaxDefects[1] = dTargetMin;
755  localMinMaxDefects[2] = dSourceMax;
756  localMinMaxDefects[3] = dTargetMax;
757  localMinMaxDefects[4] = dMassDiff;
758 
759  if( caasType == CAAS_GLOBAL )
760  {
761  MPI_Allreduce( localMinMaxDefects.data(), globalMinMaxDefects.data(), 2, MPI_DOUBLE, MPI_MIN, m_pcomm->comm() );
762  MPI_Allreduce( localMinMaxDefects.data() + 2, globalMinMaxDefects.data() + 2, 2, MPI_DOUBLE, MPI_MAX,
763  m_pcomm->comm() );
764  dSourceMin = globalMinMaxDefects[0];
765  dSourceMax = globalMinMaxDefects[2];
766  dTargetMin = globalMinMaxDefects[1];
767  dTargetMax = globalMinMaxDefects[3];
768  }
769  if( caasIteration == 1 )
770  MPI_Allreduce( localMinMaxDefects.data() + 4, globalMinMaxDefects.data() + 4, 1, MPI_DOUBLE, MPI_SUM,
771  m_pcomm->comm() );
772  else
773  globalMinMaxDefects[4] = mismatch;
774 
775  dMassDiff = localMinMaxDefects[4];
776  // massDefect.first = localMinMaxDefects[4];
777  massDefect.first = globalMinMaxDefects[4];
778 #else
779 
780  // massDefect.first = fabs( dMassDiff / ( dSourceMax - dSourceMin ) );
781  massDefect.first = dMassDiff;
782 #endif
783 
784  // Early exit if the values are monotone already.
785  // if( ( dTargetMax <= dSourceMax && dTargetMin <= dSourceMin ) || fabs( massDefect.first ) < 1e-16 )
786  if( fabs( massDefect.first ) > 1e-20 )
787  {
788  if( caasType == CAAS_GLOBAL )
789  {
790  for( size_t i = 0; i < nTargetCount; i++ )
791  {
792  dataLowerBound[i] = dSourceMin - dataOutDouble[i];
793  dataUpperBound[i] = dSourceMax - dataOutDouble[i];
794  }
795  } // if( caasType == CAAS_GLOBAL )
796  else // caasType == CAAS_LOCAL
797  {
798  // Compute the local min and max values of the target data
799  std::vector< double > vecLocalUpperBound( nTargetCount );
800  std::vector< double > vecLocalLowerBound( nTargetCount );
801  // Loop over the target faces and compute the min and max values
802  // of the source data linked to the target faces
803  for( size_t i = 0; i < nTargetCount; i++ )
804  {
805  assert( vecSourceOvTarget[i].size() );
806 
807  double dMinI = 1E10; // dataInDouble[vecSourceOvTarget[i][0]];
808  double dMaxI = -1E10; // dataInDouble[vecSourceOvTarget[i][0]];
809 
810  // Compute max over intersecting source faces
811  for( const auto& srcElem : vecSourceOvTarget[i] )
812  {
813  dMinI = fmin( dMinI, dataInDouble[srcElem] ); // min over intersecting source faces
814  dMaxI = fmax( dMaxI, dataInDouble[srcElem] ); // max over intersecting source faces
815  }
816 
817  // Update the min and max values of the target data
818  vecLocalLowerBound[i] = dMinI;
819  vecLocalUpperBound[i] = dMaxI;
820  }
821 
822  for( size_t i = 0; i < nTargetCount; i++ )
823  {
824  dataLowerBound[i] = vecLocalLowerBound[i] - dataOutDouble[i];
825  dataUpperBound[i] = vecLocalUpperBound[i] - dataOutDouble[i];
826  }
827  } // caasType == CAAS_LOCAL
828 
829  // Invoke CAAS or QLT application on the map
830  if( fabs( dMassDiff ) > 1e-20 )
831  {
832  if( caasType == CAAS_QLT )
833  dMassDiff = QLTLimiter( caasIteration, dataOutDouble, dataLowerBound, dataUpperBound, massVector );
834  else
835  CAASLimiter( dataOutDouble, dataLowerBound, dataUpperBound, dMassDiff );
836  }
837 
838  // Announce output mass
839  double dMassDiffPost = 0.0;
840  for( size_t i = 0; i < m_meshOverlap->faces.size(); i++ )
841  {
842  const int ixS = m_meshOverlap->vecSourceFaceIx[i];
843  const int ixT = m_meshOverlap->vecTargetFaceIx[i];
844 
845  if( ixT < 0 ) continue; // skip ghost target faces
846 
847  // Update the mass difference between source and target faces
848  // linked to the overlap mesh element
849  dMassDiffPost += ( dataInDouble[ixS] * m_dOverlapAreas[i] ) - // source mass
850  ( dataOutDouble[ixT] * m_dOverlapAreas[i] ); // target mass
851  }
852  // massDefect.second = fabs( dMassDiffPost / ( dSourceMax - dSourceMin ) );
853  massDefect.second = dMassDiffPost;
854  }
855 
856  // Ideally should perform an AllReduce here to get the global mass difference across all processors
857  // But if we satisfy the constraint on every task, essentially, the global mass difference should be zero!
858  return massDefect;
859 }
860 
861 //#define VERBOSE
862 moab::ErrorCode moab::TempestOnlineMap::ApplyWeights( std::vector< double >& srcVals,
863  std::vector< double >& tgtVals,
864  bool transpose )
865 {
866  // Reset the source and target data first
867  m_rowVector.setZero();
868  m_colVector.setZero();
869 
870 #ifdef VERBOSE
871  std::stringstream sstr;
872  static int callId = 0;
873  callId++;
874  sstr << "projection_id_" << callId << "_s_" << size << "_rk_" << rank << ".txt";
875  std::ofstream output_file( sstr.str() );
876 #endif
877  // Perform the actual projection of weights: application of weight matrix onto the source
878  // solution vector
879  if( transpose )
880  {
881  // Permute the source data first
882  for( unsigned i = 0; i < srcVals.size(); ++i )
883  {
884  if( row_dtoc_dofmap[i] >= 0 )
885  m_rowVector( row_dtoc_dofmap[i] ) = srcVals[i]; // permute and set the row (source) vector properly
886  }
887 
888  m_colVector = m_weightMatrix.adjoint() * m_rowVector;
889 
890  // Permute the resulting target data back
891  for( unsigned i = 0; i < tgtVals.size(); ++i )
892  {
893  if( col_dtoc_dofmap[i] >= 0 )
894  tgtVals[i] = m_colVector( col_dtoc_dofmap[i] ); // permute and set the row (source) vector properly
895  }
896  }
897  else
898  {
899  // Permute the source data first
900 #ifdef VERBOSE
901  output_file << "ColVector: " << m_colVector.size() << ", SrcVals: " << srcVals.size()
902  << ", Sizes: " << m_nTotDofs_SrcCov << ", " << col_dtoc_dofmap.size() << "\n";
903 #endif
904  for( unsigned i = 0; i < srcVals.size(); ++i )
905  {
906  if( col_dtoc_dofmap[i] >= 0 )
907  {
908  m_colVector( col_dtoc_dofmap[i] ) = srcVals[i]; // permute and set the row (source) vector properly
909 #ifdef VERBOSE
910  output_file << i << " " << col_gdofmap[col_dtoc_dofmap[i]] + 1 << " " << srcVals[i] << "\n";
911 #endif
912  }
913  }
914 
915  m_rowVector = m_weightMatrix * m_colVector;
916 
917  // Permute the resulting target data back
918 #ifdef VERBOSE
919  output_file << "RowVector: " << m_rowVector.size() << ", TgtVals:" << tgtVals.size()
920  << ", Sizes: " << m_nTotDofs_Dest << ", " << row_gdofmap.size() << "\n";
921 #endif
922  for( unsigned i = 0; i < tgtVals.size(); ++i )
923  {
924  if( row_dtoc_dofmap[i] >= 0 )
925  {
926  tgtVals[i] = m_rowVector( row_dtoc_dofmap[i] ); // permute and set the row (source) vector properly
927 #ifdef VERBOSE
928  output_file << i << " " << row_gdofmap[row_dtoc_dofmap[i]] + 1 << " " << tgtVals[i] << "\n";
929 #endif
930  }
931  }
932  }
933 
934  // if( caasType != CAAS_NONE )
935  // {
936  // constexpr int nmax_caas_iterations = 5;
937  // double mismatch = 1.0;
938  // int caasIteration = 0;
939  // while( mismatch > 1e-15 &&
940  // caasIteration++ < nmax_caas_iterations ) // iterate until convergence or a maximum of 5 iterations
941  // {
942  // std::pair< double, double > mDefect = this->ApplyCAASLimiting( srcVals, tgtVals, caasType );
943  // if( m_remapper->verbose )
944  // printf( "Rank %d: -- Iteration: %d, Net original mass defect: %3.4e, mass defect post-CAAS: %3.4e\n",
945  // m_remapper->rank, caasIteration, mDefect.first, mDefect.second );
946  // mismatch = mDefect.second;
947  // }
948  // }
949 
950 #ifdef VERBOSE
951  output_file.flush(); // required here
952  output_file.close();
953 #endif
954 
955  // All done with matvec application
956  return moab::MB_SUCCESS;
957 }
958 
959 #endif
960 //#undef VERBOSE
961 ///////////////////////////////////////////////////////////////////////////////
962 
963 extern void ForceConsistencyConservation3( const DataArray1D< double >& vecSourceArea,
964  const DataArray1D< double >& vecTargetArea,
965  DataArray2D< double >& dCoeff,
966  bool fMonotone,
967  bool fSparseConstraints = false );
968 
969 ///////////////////////////////////////////////////////////////////////////////
970 
971 extern void ForceIntArrayConsistencyConservation( const DataArray1D< double >& vecSourceArea,
972  const DataArray1D< double >& vecTargetArea,
973  DataArray2D< double >& dCoeff,
974  bool fMonotone );
975 
976 ///////////////////////////////////////////////////////////////////////////////
977 
978 void moab::TempestOnlineMap::LinearRemapSE4_Tempest_MOAB( const DataArray3D< int >& dataGLLNodes,
979  const DataArray3D< double >& dataGLLJacobian,
980  int nMonotoneType,
981  bool fContinuousIn,
982  bool fNoConservation )
983 {
984  // Order of the polynomial interpolant
985  int nP = dataGLLNodes.GetRows();
986 
987  // Order of triangular quadrature rule
988  const int TriQuadRuleOrder = 4;
989 
990  // Triangular quadrature rule
991  TriangularQuadratureRule triquadrule( TriQuadRuleOrder );
992 
993  int TriQuadraturePoints = triquadrule.GetPoints();
994 
995  const DataArray2D< double >& TriQuadratureG = triquadrule.GetG();
996 
997  const DataArray1D< double >& TriQuadratureW = triquadrule.GetW();
998 
999  // Sample coefficients
1000  DataArray2D< double > dSampleCoeff( nP, nP );
1001 
1002  // GLL Quadrature nodes on quadrilateral elements
1003  DataArray1D< double > dG;
1004  DataArray1D< double > dW;
1005  GaussLobattoQuadrature::GetPoints( nP, 0.0, 1.0, dG, dW );
1006 
1007  // Announcements
1008  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
1009  dbgprint.set_prefix( "[LinearRemapSE4_Tempest_MOAB]: " );
1010  if( is_root )
1011  {
1012  dbgprint.printf( 0, "Finite Element to Finite Volume Projection\n" );
1013  dbgprint.printf( 0, "Triangular quadrature rule order %i\n", TriQuadRuleOrder );
1014  dbgprint.printf( 0, "Order of the FE polynomial interpolant: %i\n", nP );
1015  }
1016 
1017  // Get SparseMatrix represntation of the OfflineMap
1018  SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1019 
1020  // NodeVector from m_meshOverlap
1021  const NodeVector& nodesOverlap = m_meshOverlap->nodes;
1022  const NodeVector& nodesFirst = m_meshInputCov->nodes;
1023 
1024  // Vector of source areas
1025  DataArray1D< double > vecSourceArea( nP * nP );
1026 
1027  DataArray1D< double > vecTargetArea;
1028  DataArray2D< double > dCoeff;
1029 
1030 #ifdef VERBOSE
1031  std::stringstream sstr;
1032  sstr << "remapdata_" << rank << ".txt";
1033  std::ofstream output_file( sstr.str() );
1034 #endif
1035 
1036  // Current Overlap Face
1037  int ixOverlap = 0;
1038 #ifdef VERBOSE
1039  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
1040 #endif
1041  // generic triangle used for area computation, for triangles around the center of overlap face;
1042  // used for overlap faces with more than 4 edges;
1043  // nodes array will be set for each triangle;
1044  // these triangles are not part of the mesh structure, they are just temporary during
1045  // aforementioned decomposition.
1046  Face faceTri( 3 );
1047  NodeVector nodes( 3 );
1048  faceTri.SetNode( 0, 0 );
1049  faceTri.SetNode( 1, 1 );
1050  faceTri.SetNode( 2, 2 );
1051 
1052  // Loop over all input Faces
1053  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1054  {
1055  const Face& faceFirst = m_meshInputCov->faces[ixFirst];
1056 
1057  if( faceFirst.edges.size() != 4 )
1058  {
1059  _EXCEPTIONT( "Only quadrilateral elements allowed for SE remapping" );
1060  }
1061 #ifdef VERBOSE
1062  // Announce computation progress
1063  if( ixFirst % outputFrequency == 0 && is_root )
1064  {
1065  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1066  }
1067 #endif
1068  // Need to re-number the overlap elements such that vecSourceFaceIx[a:b] = 0, then 1 and so
1069  // on wrt the input mesh data Then the overlap_end and overlap_begin will be correct.
1070  // However, the relation with MOAB and Tempest will go out of the roof
1071 
1072  // Determine how many overlap Faces and triangles are present
1073  int nOverlapFaces = 0;
1074  size_t ixOverlapTemp = ixOverlap;
1075  for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1076  {
1077  // if( m_meshOverlap->vecTargetFaceIx[ixOverlapTemp] < 0 ) continue; // skip ghost target faces
1078  // const Face & faceOverlap = m_meshOverlap->faces[ixOverlapTemp];
1079  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 ) break;
1080 
1081  nOverlapFaces++;
1082  }
1083 
1084  // No overlaps
1085  if( nOverlapFaces == 0 )
1086  {
1087  continue;
1088  }
1089 
1090  // Allocate remap coefficients array for meshFirst Face
1091  DataArray3D< double > dRemapCoeff( nP, nP, nOverlapFaces );
1092 
1093  // Find the local remap coefficients
1094  for( int j = 0; j < nOverlapFaces; j++ )
1095  {
1096  const Face& faceOverlap = m_meshOverlap->faces[ixOverlap + j];
1097  if( m_meshOverlap->vecFaceArea[ixOverlap + j] < 1.e-16 ) // machine precision
1098  {
1099  Announce( "Very small overlap at index %i area polygon: (%1.10e )", ixOverlap + j,
1100  m_meshOverlap->vecFaceArea[ixOverlap + j] );
1101  int n = faceOverlap.edges.size();
1102  Announce( "Number nodes: %d", n );
1103  for( int k = 0; k < n; k++ )
1104  {
1105  Node nd = nodesOverlap[faceOverlap[k]];
1106  Announce( "Node %d %d : %1.10e %1.10e %1.10e ", k, faceOverlap[k], nd.x, nd.y, nd.z );
1107  }
1108  continue;
1109  }
1110 
1111  // #ifdef VERBOSE
1112  // if ( is_root )
1113  // Announce ( "\tLocal ID: %i/%i = %i, areas = %2.8e", j + ixOverlap, nOverlapFaces,
1114  // m_remapper->lid_to_gid_covsrc[m_meshOverlap->vecSourceFaceIx[ixOverlap + j]],
1115  // m_meshOverlap->vecFaceArea[ixOverlap + j] );
1116  // #endif
1117 
1118  int nbEdges = faceOverlap.edges.size();
1119  int nOverlapTriangles = 1;
1120  Node center; // not used if nbEdges == 3
1121  if( nbEdges > 3 )
1122  { // decompose from center in this case
1123  nOverlapTriangles = nbEdges;
1124  for( int k = 0; k < nbEdges; k++ )
1125  {
1126  const Node& node = nodesOverlap[faceOverlap[k]];
1127  center = center + node;
1128  }
1129  center = center / nbEdges;
1130  center = center.Normalized(); // project back on sphere of radius 1
1131  }
1132 
1133  Node node0, node1, node2;
1134  double dTriangleArea;
1135 
1136  // Loop over all sub-triangles of this Overlap Face
1137  for( int k = 0; k < nOverlapTriangles; k++ )
1138  {
1139  if( nbEdges == 3 ) // will come here only once, nOverlapTriangles == 1 in this case
1140  {
1141  node0 = nodesOverlap[faceOverlap[0]];
1142  node1 = nodesOverlap[faceOverlap[1]];
1143  node2 = nodesOverlap[faceOverlap[2]];
1144  dTriangleArea = CalculateFaceArea( faceOverlap, nodesOverlap );
1145  }
1146  else // decompose polygon in triangles around the center
1147  {
1148  node0 = center;
1149  node1 = nodesOverlap[faceOverlap[k]];
1150  int k1 = ( k + 1 ) % nbEdges;
1151  node2 = nodesOverlap[faceOverlap[k1]];
1152  nodes[0] = center;
1153  nodes[1] = node1;
1154  nodes[2] = node2;
1155  dTriangleArea = CalculateFaceArea( faceTri, nodes );
1156  }
1157  // Coordinates of quadrature Node
1158  for( int l = 0; l < TriQuadraturePoints; l++ )
1159  {
1160  Node nodeQuadrature;
1161  nodeQuadrature.x = TriQuadratureG[l][0] * node0.x + TriQuadratureG[l][1] * node1.x +
1162  TriQuadratureG[l][2] * node2.x;
1163 
1164  nodeQuadrature.y = TriQuadratureG[l][0] * node0.y + TriQuadratureG[l][1] * node1.y +
1165  TriQuadratureG[l][2] * node2.y;
1166 
1167  nodeQuadrature.z = TriQuadratureG[l][0] * node0.z + TriQuadratureG[l][1] * node1.z +
1168  TriQuadratureG[l][2] * node2.z;
1169 
1170  nodeQuadrature = nodeQuadrature.Normalized();
1171 
1172  // Find components of quadrature point in basis
1173  // of the first Face
1174  double dAlpha;
1175  double dBeta;
1176 
1177  ApplyInverseMap( faceFirst, nodesFirst, nodeQuadrature, dAlpha, dBeta );
1178 
1179  // Check inverse map value
1180  if( ( dAlpha < -1.0e-13 ) || ( dAlpha > 1.0 + 1.0e-13 ) || ( dBeta < -1.0e-13 ) ||
1181  ( dBeta > 1.0 + 1.0e-13 ) )
1182  {
1183  _EXCEPTION4( "Inverse Map for element %d and subtriangle %d out of range "
1184  "(%1.5e %1.5e)",
1185  j, l, dAlpha, dBeta );
1186  }
1187 
1188  // Sample the finite element at this point
1189  SampleGLLFiniteElement( nMonotoneType, nP, dAlpha, dBeta, dSampleCoeff );
1190 
1191  // Add sample coefficients to the map if m_meshOverlap->vecFaceArea[ixOverlap + j] > 0
1192  for( int p = 0; p < nP; p++ )
1193  {
1194  for( int q = 0; q < nP; q++ )
1195  {
1196  dRemapCoeff[p][q][j] += TriQuadratureW[l] * dTriangleArea * dSampleCoeff[p][q] /
1197  m_meshOverlap->vecFaceArea[ixOverlap + j];
1198  }
1199  }
1200  }
1201  }
1202  }
1203 
1204 #ifdef VERBOSE
1205  output_file << "[" << m_remapper->lid_to_gid_covsrc[ixFirst] << "] \t";
1206  for( int j = 0; j < nOverlapFaces; j++ )
1207  {
1208  for( int p = 0; p < nP; p++ )
1209  {
1210  for( int q = 0; q < nP; q++ )
1211  {
1212  output_file << dRemapCoeff[p][q][j] << " ";
1213  }
1214  }
1215  }
1216  output_file << std::endl;
1217 #endif
1218 
1219  // Force consistency and conservation
1220  if( !fNoConservation )
1221  {
1222  double dTargetArea = 0.0;
1223  for( int j = 0; j < nOverlapFaces; j++ )
1224  {
1225  dTargetArea += m_meshOverlap->vecFaceArea[ixOverlap + j];
1226  }
1227 
1228  for( int p = 0; p < nP; p++ )
1229  {
1230  for( int q = 0; q < nP; q++ )
1231  {
1232  vecSourceArea[p * nP + q] = dataGLLJacobian[p][q][ixFirst];
1233  }
1234  }
1235 
1236  const double areaTolerance = 1e-10;
1237  // Source elements are completely covered by target volumes
1238  if( fabs( m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea ) <= areaTolerance )
1239  {
1240  vecTargetArea.Allocate( nOverlapFaces );
1241  for( int j = 0; j < nOverlapFaces; j++ )
1242  {
1243  vecTargetArea[j] = m_meshOverlap->vecFaceArea[ixOverlap + j];
1244  }
1245 
1246  dCoeff.Allocate( nOverlapFaces, nP * nP );
1247 
1248  for( int j = 0; j < nOverlapFaces; j++ )
1249  {
1250  for( int p = 0; p < nP; p++ )
1251  {
1252  for( int q = 0; q < nP; q++ )
1253  {
1254  dCoeff[j][p * nP + q] = dRemapCoeff[p][q][j];
1255  }
1256  }
1257  }
1258 
1259  // Target volumes only partially cover source elements
1260  }
1261  else if( m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea > areaTolerance )
1262  {
1263  double dExtraneousArea = m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea;
1264 
1265  vecTargetArea.Allocate( nOverlapFaces + 1 );
1266  for( int j = 0; j < nOverlapFaces; j++ )
1267  {
1268  vecTargetArea[j] = m_meshOverlap->vecFaceArea[ixOverlap + j];
1269  }
1270  vecTargetArea[nOverlapFaces] = dExtraneousArea;
1271 
1272 #ifdef VERBOSE
1273  Announce( "Partial volume: %i (%1.10e / %1.10e)", ixFirst, dTargetArea,
1274  m_meshInputCov->vecFaceArea[ixFirst] );
1275 #endif
1276  if( dTargetArea > m_meshInputCov->vecFaceArea[ixFirst] )
1277  {
1278  _EXCEPTIONT( "Partial element area exceeds total element area" );
1279  }
1280 
1281  dCoeff.Allocate( nOverlapFaces + 1, nP * nP );
1282 
1283  for( int j = 0; j < nOverlapFaces; j++ )
1284  {
1285  for( int p = 0; p < nP; p++ )
1286  {
1287  for( int q = 0; q < nP; q++ )
1288  {
1289  dCoeff[j][p * nP + q] = dRemapCoeff[p][q][j];
1290  }
1291  }
1292  }
1293  for( int p = 0; p < nP; p++ )
1294  {
1295  for( int q = 0; q < nP; q++ )
1296  {
1297  dCoeff[nOverlapFaces][p * nP + q] = dataGLLJacobian[p][q][ixFirst];
1298  }
1299  }
1300  for( int j = 0; j < nOverlapFaces; j++ )
1301  {
1302  for( int p = 0; p < nP; p++ )
1303  {
1304  for( int q = 0; q < nP; q++ )
1305  {
1306  dCoeff[nOverlapFaces][p * nP + q] -=
1307  dRemapCoeff[p][q][j] * m_meshOverlap->vecFaceArea[ixOverlap + j];
1308  }
1309  }
1310  }
1311  for( int p = 0; p < nP; p++ )
1312  {
1313  for( int q = 0; q < nP; q++ )
1314  {
1315  dCoeff[nOverlapFaces][p * nP + q] /= dExtraneousArea;
1316  }
1317  }
1318 
1319  // Source elements only partially cover target volumes
1320  }
1321  else
1322  {
1323  Announce( "Coverage area: %1.10e, and target element area: %1.10e)", ixFirst,
1324  m_meshInputCov->vecFaceArea[ixFirst], dTargetArea );
1325  _EXCEPTIONT( "Target grid must be a subset of source grid" );
1326  }
1327 
1328  ForceConsistencyConservation3( vecSourceArea, vecTargetArea, dCoeff, ( nMonotoneType > 0 )
1329  /*, m_remapper->lid_to_gid_covsrc[ixFirst]*/ );
1330 
1331  for( int j = 0; j < nOverlapFaces; j++ )
1332  {
1333  for( int p = 0; p < nP; p++ )
1334  {
1335  for( int q = 0; q < nP; q++ )
1336  {
1337  dRemapCoeff[p][q][j] = dCoeff[j][p * nP + q];
1338  }
1339  }
1340  }
1341  }
1342 
1343 #ifdef VERBOSE
1344  // output_file << "[" << m_remapper->lid_to_gid_covsrc[ixFirst] << "] \t";
1345  // for ( int j = 0; j < nOverlapFaces; j++ )
1346  // {
1347  // for ( int p = 0; p < nP; p++ )
1348  // {
1349  // for ( int q = 0; q < nP; q++ )
1350  // {
1351  // output_file << dRemapCoeff[p][q][j] << " ";
1352  // }
1353  // }
1354  // }
1355  // output_file << std::endl;
1356 #endif
1357 
1358  // Put these remap coefficients into the SparseMatrix map
1359  for( int j = 0; j < nOverlapFaces; j++ )
1360  {
1361  int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
1362 
1363  // signal to not participate, because it is a ghost target
1364  if( ixSecondFace < 0 ) continue; // do not do anything
1365 
1366  for( int p = 0; p < nP; p++ )
1367  {
1368  for( int q = 0; q < nP; q++ )
1369  {
1370  if( fContinuousIn )
1371  {
1372  int ixFirstNode = dataGLLNodes[p][q][ixFirst] - 1;
1373 
1374  smatMap( ixSecondFace, ixFirstNode ) += dRemapCoeff[p][q][j] *
1375  m_meshOverlap->vecFaceArea[ixOverlap + j] /
1376  m_meshOutput->vecFaceArea[ixSecondFace];
1377  }
1378  else
1379  {
1380  int ixFirstNode = ixFirst * nP * nP + p * nP + q;
1381 
1382  smatMap( ixSecondFace, ixFirstNode ) += dRemapCoeff[p][q][j] *
1383  m_meshOverlap->vecFaceArea[ixOverlap + j] /
1384  m_meshOutput->vecFaceArea[ixSecondFace];
1385  }
1386  }
1387  }
1388  }
1389  // Increment the current overlap index
1390  ixOverlap += nOverlapFaces;
1391  }
1392 #ifdef VERBOSE
1393  output_file.flush(); // required here
1394  output_file.close();
1395 #endif
1396 
1397  return;
1398 }
1399 
1400 ///////////////////////////////////////////////////////////////////////////////
1401 
1402 void moab::TempestOnlineMap::LinearRemapGLLtoGLL2_MOAB( const DataArray3D< int >& dataGLLNodesIn,
1403  const DataArray3D< double >& dataGLLJacobianIn,
1404  const DataArray3D< int >& dataGLLNodesOut,
1405  const DataArray3D< double >& dataGLLJacobianOut,
1406  const DataArray1D< double >& dataNodalAreaOut,
1407  int nPin,
1408  int nPout,
1409  int nMonotoneType,
1410  bool fContinuousIn,
1411  bool fContinuousOut,
1412  bool fNoConservation )
1413 {
1414  // Triangular quadrature rule
1415  TriangularQuadratureRule triquadrule( 8 );
1416 
1417  const DataArray2D< double >& dG = triquadrule.GetG();
1418  const DataArray1D< double >& dW = triquadrule.GetW();
1419 
1420  // Get SparseMatrix represntation of the OfflineMap
1421  SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1422 
1423  // Sample coefficients
1424  DataArray2D< double > dSampleCoeffIn( nPin, nPin );
1425  DataArray2D< double > dSampleCoeffOut( nPout, nPout );
1426 
1427  // Announcemnets
1428  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
1429  dbgprint.set_prefix( "[LinearRemapGLLtoGLL2_MOAB]: " );
1430  if( is_root )
1431  {
1432  dbgprint.printf( 0, "Finite Element to Finite Element Projection\n" );
1433  dbgprint.printf( 0, "Order of the input FE polynomial interpolant: %i\n", nPin );
1434  dbgprint.printf( 0, "Order of the output FE polynomial interpolant: %i\n", nPout );
1435  }
1436 
1437  // Build the integration array for each element on m_meshOverlap
1438  DataArray3D< double > dGlobalIntArray( nPin * nPin, m_meshOverlap->faces.size(), nPout * nPout );
1439 
1440  // Number of overlap Faces per source Face
1441  DataArray1D< int > nAllOverlapFaces( m_meshInputCov->faces.size() );
1442 
1443  int ixOverlap = 0;
1444  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1445  {
1446  // Determine how many overlap Faces and triangles are present
1447  int nOverlapFaces = 0;
1448  size_t ixOverlapTemp = ixOverlap;
1449  for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1450  {
1451  // const Face & faceOverlap = m_meshOverlap->faces[ixOverlapTemp];
1452  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 )
1453  {
1454  break;
1455  }
1456 
1457  nOverlapFaces++;
1458  }
1459 
1460  nAllOverlapFaces[ixFirst] = nOverlapFaces;
1461 
1462  // Increment the current overlap index
1463  ixOverlap += nAllOverlapFaces[ixFirst];
1464  }
1465 
1466  // Geometric area of each output node
1467  DataArray2D< double > dGeometricOutputArea( m_meshOutput->faces.size(), nPout * nPout );
1468 
1469  // Area of each overlap element in the output basis
1470  DataArray2D< double > dOverlapOutputArea( m_meshOverlap->faces.size(), nPout * nPout );
1471 
1472  // Loop through all faces on m_meshInputCov
1473  ixOverlap = 0;
1474 #ifdef VERBOSE
1475  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
1476 #endif
1477  if( is_root ) dbgprint.printf( 0, "Building conservative distribution maps\n" );
1478 
1479  // generic triangle used for area computation, for triangles around the center of overlap face;
1480  // used for overlap faces with more than 4 edges;
1481  // nodes array will be set for each triangle;
1482  // these triangles are not part of the mesh structure, they are just temporary during
1483  // aforementioned decomposition.
1484  Face faceTri( 3 );
1485  NodeVector nodes( 3 );
1486  faceTri.SetNode( 0, 0 );
1487  faceTri.SetNode( 1, 1 );
1488  faceTri.SetNode( 2, 2 );
1489 
1490  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1491  {
1492 #ifdef VERBOSE
1493  // Announce computation progress
1494  if( ixFirst % outputFrequency == 0 && is_root )
1495  {
1496  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1497  }
1498 #endif
1499  // Quantities from the First Mesh
1500  const Face& faceFirst = m_meshInputCov->faces[ixFirst];
1501 
1502  const NodeVector& nodesFirst = m_meshInputCov->nodes;
1503 
1504  // Number of overlapping Faces and triangles
1505  int nOverlapFaces = nAllOverlapFaces[ixFirst];
1506 
1507  if( !nOverlapFaces ) continue;
1508 
1509  // // Calculate total element Jacobian
1510  // double dTotalJacobian = 0.0;
1511  // for (int s = 0; s < nPin; s++) {
1512  // for (int t = 0; t < nPin; t++) {
1513  // dTotalJacobian += dataGLLJacobianIn[s][t][ixFirst];
1514  // }
1515  // }
1516 
1517  // Loop through all Overlap Faces
1518  for( int i = 0; i < nOverlapFaces; i++ )
1519  {
1520  // Quantities from the overlap Mesh
1521  const Face& faceOverlap = m_meshOverlap->faces[ixOverlap + i];
1522 
1523  const NodeVector& nodesOverlap = m_meshOverlap->nodes;
1524 
1525  // Quantities from the Second Mesh
1526  int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1527 
1528  // signal to not participate, because it is a ghost target
1529  if( ixSecond < 0 ) continue; // do not do anything
1530 
1531  const NodeVector& nodesSecond = m_meshOutput->nodes;
1532 
1533  const Face& faceSecond = m_meshOutput->faces[ixSecond];
1534 
1535  int nbEdges = faceOverlap.edges.size();
1536  int nOverlapTriangles = 1;
1537  Node center; // not used if nbEdges == 3
1538  if( nbEdges > 3 )
1539  { // decompose from center in this case
1540  nOverlapTriangles = nbEdges;
1541  for( int k = 0; k < nbEdges; k++ )
1542  {
1543  const Node& node = nodesOverlap[faceOverlap[k]];
1544  center = center + node;
1545  }
1546  center = center / nbEdges;
1547  center = center.Normalized(); // project back on sphere of radius 1
1548  }
1549 
1550  Node node0, node1, node2;
1551  double dTriArea;
1552 
1553  // Loop over all sub-triangles of this Overlap Face
1554  for( int j = 0; j < nOverlapTriangles; j++ )
1555  {
1556  if( nbEdges == 3 ) // will come here only once, nOverlapTriangles == 1 in this case
1557  {
1558  node0 = nodesOverlap[faceOverlap[0]];
1559  node1 = nodesOverlap[faceOverlap[1]];
1560  node2 = nodesOverlap[faceOverlap[2]];
1561  dTriArea = CalculateFaceArea( faceOverlap, nodesOverlap );
1562  }
1563  else // decompose polygon in triangles around the center
1564  {
1565  node0 = center;
1566  node1 = nodesOverlap[faceOverlap[j]];
1567  int j1 = ( j + 1 ) % nbEdges;
1568  node2 = nodesOverlap[faceOverlap[j1]];
1569  nodes[0] = center;
1570  nodes[1] = node1;
1571  nodes[2] = node2;
1572  dTriArea = CalculateFaceArea( faceTri, nodes );
1573  }
1574 
1575  for( int k = 0; k < triquadrule.GetPoints(); k++ )
1576  {
1577  // Get the nodal location of this point
1578  double dX[3];
1579 
1580  dX[0] = dG( k, 0 ) * node0.x + dG( k, 1 ) * node1.x + dG( k, 2 ) * node2.x;
1581  dX[1] = dG( k, 0 ) * node0.y + dG( k, 1 ) * node1.y + dG( k, 2 ) * node2.y;
1582  dX[2] = dG( k, 0 ) * node0.z + dG( k, 1 ) * node1.z + dG( k, 2 ) * node2.z;
1583 
1584  double dMag = sqrt( dX[0] * dX[0] + dX[1] * dX[1] + dX[2] * dX[2] );
1585 
1586  dX[0] /= dMag;
1587  dX[1] /= dMag;
1588  dX[2] /= dMag;
1589 
1590  Node nodeQuadrature( dX[0], dX[1], dX[2] );
1591 
1592  // Find the components of this quadrature point in the basis
1593  // of the first Face.
1594  double dAlphaIn;
1595  double dBetaIn;
1596 
1597  ApplyInverseMap( faceFirst, nodesFirst, nodeQuadrature, dAlphaIn, dBetaIn );
1598 
1599  // Find the components of this quadrature point in the basis
1600  // of the second Face.
1601  double dAlphaOut;
1602  double dBetaOut;
1603 
1604  ApplyInverseMap( faceSecond, nodesSecond, nodeQuadrature, dAlphaOut, dBetaOut );
1605 
1606  /*
1607  // Check inverse map value
1608  if ((dAlphaIn < 0.0) || (dAlphaIn > 1.0) ||
1609  (dBetaIn < 0.0) || (dBetaIn > 1.0)
1610  ) {
1611  _EXCEPTION2("Inverse Map out of range (%1.5e %1.5e)",
1612  dAlphaIn, dBetaIn);
1613  }
1614 
1615  // Check inverse map value
1616  if ((dAlphaOut < 0.0) || (dAlphaOut > 1.0) ||
1617  (dBetaOut < 0.0) || (dBetaOut > 1.0)
1618  ) {
1619  _EXCEPTION2("Inverse Map out of range (%1.5e %1.5e)",
1620  dAlphaOut, dBetaOut);
1621  }
1622  */
1623  // Sample the First finite element at this point
1624  SampleGLLFiniteElement( nMonotoneType, nPin, dAlphaIn, dBetaIn, dSampleCoeffIn );
1625 
1626  // Sample the Second finite element at this point
1627  SampleGLLFiniteElement( nMonotoneType, nPout, dAlphaOut, dBetaOut, dSampleCoeffOut );
1628 
1629  // Overlap output area
1630  for( int s = 0; s < nPout; s++ )
1631  {
1632  for( int t = 0; t < nPout; t++ )
1633  {
1634  double dNodeArea = dSampleCoeffOut[s][t] * dW[k] * dTriArea;
1635 
1636  dOverlapOutputArea[ixOverlap + i][s * nPout + t] += dNodeArea;
1637 
1638  dGeometricOutputArea[ixSecond][s * nPout + t] += dNodeArea;
1639  }
1640  }
1641 
1642  // Compute overlap integral
1643  int ixp = 0;
1644  for( int p = 0; p < nPin; p++ )
1645  {
1646  for( int q = 0; q < nPin; q++ )
1647  {
1648  int ixs = 0;
1649  for( int s = 0; s < nPout; s++ )
1650  {
1651  for( int t = 0; t < nPout; t++ )
1652  {
1653  // Sample the Second finite element at this point
1654  dGlobalIntArray[ixp][ixOverlap + i][ixs] +=
1655  dSampleCoeffOut[s][t] * dSampleCoeffIn[p][q] * dW[k] * dTriArea;
1656 
1657  ixs++;
1658  }
1659  }
1660 
1661  ixp++;
1662  }
1663  }
1664  }
1665  }
1666  }
1667 
1668  // Coefficients
1669  DataArray2D< double > dCoeff( nOverlapFaces * nPout * nPout, nPin * nPin );
1670 
1671  for( int i = 0; i < nOverlapFaces; i++ )
1672  {
1673  // int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1674 
1675  int ixp = 0;
1676  for( int p = 0; p < nPin; p++ )
1677  {
1678  for( int q = 0; q < nPin; q++ )
1679  {
1680  int ixs = 0;
1681  for( int s = 0; s < nPout; s++ )
1682  {
1683  for( int t = 0; t < nPout; t++ )
1684  {
1685  dCoeff[i * nPout * nPout + ixs][ixp] = dGlobalIntArray[ixp][ixOverlap + i][ixs] /
1686  dOverlapOutputArea[ixOverlap + i][s * nPout + t];
1687 
1688  ixs++;
1689  }
1690  }
1691 
1692  ixp++;
1693  }
1694  }
1695  }
1696 
1697  // Source areas
1698  DataArray1D< double > vecSourceArea( nPin * nPin );
1699 
1700  for( int p = 0; p < nPin; p++ )
1701  {
1702  for( int q = 0; q < nPin; q++ )
1703  {
1704  vecSourceArea[p * nPin + q] = dataGLLJacobianIn[p][q][ixFirst];
1705  }
1706  }
1707 
1708  // Target areas
1709  DataArray1D< double > vecTargetArea( nOverlapFaces * nPout * nPout );
1710 
1711  for( int i = 0; i < nOverlapFaces; i++ )
1712  {
1713  // int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1714  int ixs = 0;
1715  for( int s = 0; s < nPout; s++ )
1716  {
1717  for( int t = 0; t < nPout; t++ )
1718  {
1719  vecTargetArea[i * nPout * nPout + ixs] = dOverlapOutputArea[ixOverlap + i][nPout * s + t];
1720 
1721  ixs++;
1722  }
1723  }
1724  }
1725 
1726  // Force consistency and conservation
1727  if( !fNoConservation )
1728  {
1729  ForceIntArrayConsistencyConservation( vecSourceArea, vecTargetArea, dCoeff, ( nMonotoneType != 0 ) );
1730  }
1731 
1732  // Update global coefficients
1733  for( int i = 0; i < nOverlapFaces; i++ )
1734  {
1735  int ixp = 0;
1736  for( int p = 0; p < nPin; p++ )
1737  {
1738  for( int q = 0; q < nPin; q++ )
1739  {
1740  int ixs = 0;
1741  for( int s = 0; s < nPout; s++ )
1742  {
1743  for( int t = 0; t < nPout; t++ )
1744  {
1745  dGlobalIntArray[ixp][ixOverlap + i][ixs] =
1746  dCoeff[i * nPout * nPout + ixs][ixp] * dOverlapOutputArea[ixOverlap + i][s * nPout + t];
1747 
1748  ixs++;
1749  }
1750  }
1751 
1752  ixp++;
1753  }
1754  }
1755  }
1756 
1757 #ifdef VVERBOSE
1758  // Check column sums (conservation)
1759  for( int i = 0; i < nPin * nPin; i++ )
1760  {
1761  double dColSum = 0.0;
1762  for( int j = 0; j < nOverlapFaces * nPout * nPout; j++ )
1763  {
1764  dColSum += dCoeff[j][i] * vecTargetArea[j];
1765  }
1766  printf( "Col %i: %1.15e\n", i, dColSum / vecSourceArea[i] );
1767  }
1768 
1769  // Check row sums (consistency)
1770  for( int j = 0; j < nOverlapFaces * nPout * nPout; j++ )
1771  {
1772  double dRowSum = 0.0;
1773  for( int i = 0; i < nPin * nPin; i++ )
1774  {
1775  dRowSum += dCoeff[j][i];
1776  }
1777  printf( "Row %i: %1.15e\n", j, dRowSum );
1778  }
1779 #endif
1780 
1781  // Increment the current overlap index
1782  ixOverlap += nOverlapFaces;
1783  }
1784 
1785  // Build redistribution map within target element
1786  if( is_root ) dbgprint.printf( 0, "Building redistribution maps on target mesh\n" );
1787  DataArray1D< double > dRedistSourceArea( nPout * nPout );
1788  DataArray1D< double > dRedistTargetArea( nPout * nPout );
1789  std::vector< DataArray2D< double > > dRedistributionMaps;
1790  dRedistributionMaps.resize( m_meshOutput->faces.size() );
1791 
1792  for( size_t ixSecond = 0; ixSecond < m_meshOutput->faces.size(); ixSecond++ )
1793  {
1794  dRedistributionMaps[ixSecond].Allocate( nPout * nPout, nPout * nPout );
1795 
1796  for( int i = 0; i < nPout * nPout; i++ )
1797  {
1798  dRedistributionMaps[ixSecond][i][i] = 1.0;
1799  }
1800 
1801  for( int s = 0; s < nPout * nPout; s++ )
1802  {
1803  dRedistSourceArea[s] = dGeometricOutputArea[ixSecond][s];
1804  }
1805 
1806  for( int s = 0; s < nPout * nPout; s++ )
1807  {
1808  dRedistTargetArea[s] = dataGLLJacobianOut[s / nPout][s % nPout][ixSecond];
1809  }
1810 
1811  if( !fNoConservation )
1812  {
1813  ForceIntArrayConsistencyConservation( dRedistSourceArea, dRedistTargetArea, dRedistributionMaps[ixSecond],
1814  ( nMonotoneType != 0 ) );
1815 
1816  for( int s = 0; s < nPout * nPout; s++ )
1817  {
1818  for( int t = 0; t < nPout * nPout; t++ )
1819  {
1820  dRedistributionMaps[ixSecond][s][t] *= dRedistTargetArea[s] / dRedistSourceArea[t];
1821  }
1822  }
1823  }
1824  }
1825 
1826  // Construct the total geometric area
1827  DataArray1D< double > dTotalGeometricArea( dataNodalAreaOut.GetRows() );
1828  for( size_t ixSecond = 0; ixSecond < m_meshOutput->faces.size(); ixSecond++ )
1829  {
1830  for( int s = 0; s < nPout; s++ )
1831  {
1832  for( int t = 0; t < nPout; t++ )
1833  {
1834  dTotalGeometricArea[dataGLLNodesOut[s][t][ixSecond] - 1] +=
1835  dGeometricOutputArea[ixSecond][s * nPout + t];
1836  }
1837  }
1838  }
1839 
1840  // Compose the integration operator with the output map
1841  ixOverlap = 0;
1842 
1843  if( is_root ) dbgprint.printf( 0, "Assembling map\n" );
1844 
1845  // Map from source DOFs to target DOFs with redistribution applied
1846  DataArray2D< double > dRedistributedOp( nPin * nPin, nPout * nPout );
1847 
1848  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1849  {
1850 #ifdef VERBOSE
1851  // Announce computation progress
1852  if( ixFirst % outputFrequency == 0 && is_root )
1853  {
1854  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1855  }
1856 #endif
1857  // Number of overlapping Faces and triangles
1858  int nOverlapFaces = nAllOverlapFaces[ixFirst];
1859 
1860  if( !nOverlapFaces ) continue;
1861 
1862  // Put composed array into map
1863  for( int j = 0; j < nOverlapFaces; j++ )
1864  {
1865  int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
1866 
1867  // signal to not participate, because it is a ghost target
1868  if( ixSecondFace < 0 ) continue; // do not do anything
1869 
1870  dRedistributedOp.Zero();
1871  for( int p = 0; p < nPin * nPin; p++ )
1872  {
1873  for( int s = 0; s < nPout * nPout; s++ )
1874  {
1875  for( int t = 0; t < nPout * nPout; t++ )
1876  {
1877  dRedistributedOp[p][s] +=
1878  dRedistributionMaps[ixSecondFace][s][t] * dGlobalIntArray[p][ixOverlap + j][t];
1879  }
1880  }
1881  }
1882 
1883  int ixp = 0;
1884  for( int p = 0; p < nPin; p++ )
1885  {
1886  for( int q = 0; q < nPin; q++ )
1887  {
1888  int ixFirstNode;
1889  if( fContinuousIn )
1890  {
1891  ixFirstNode = dataGLLNodesIn[p][q][ixFirst] - 1;
1892  }
1893  else
1894  {
1895  ixFirstNode = ixFirst * nPin * nPin + p * nPin + q;
1896  }
1897 
1898  int ixs = 0;
1899  for( int s = 0; s < nPout; s++ )
1900  {
1901  for( int t = 0; t < nPout; t++ )
1902  {
1903  int ixSecondNode;
1904  if( fContinuousOut )
1905  {
1906  ixSecondNode = dataGLLNodesOut[s][t][ixSecondFace] - 1;
1907 
1908  if( !fNoConservation )
1909  {
1910  smatMap( ixSecondNode, ixFirstNode ) +=
1911  dRedistributedOp[ixp][ixs] / dataNodalAreaOut[ixSecondNode];
1912  }
1913  else
1914  {
1915  smatMap( ixSecondNode, ixFirstNode ) +=
1916  dRedistributedOp[ixp][ixs] / dTotalGeometricArea[ixSecondNode];
1917  }
1918  }
1919  else
1920  {
1921  ixSecondNode = ixSecondFace * nPout * nPout + s * nPout + t;
1922 
1923  if( !fNoConservation )
1924  {
1925  smatMap( ixSecondNode, ixFirstNode ) +=
1926  dRedistributedOp[ixp][ixs] / dataGLLJacobianOut[s][t][ixSecondFace];
1927  }
1928  else
1929  {
1930  smatMap( ixSecondNode, ixFirstNode ) +=
1931  dRedistributedOp[ixp][ixs] / dGeometricOutputArea[ixSecondFace][s * nPout + t];
1932  }
1933  }
1934 
1935  ixs++;
1936  }
1937  }
1938 
1939  ixp++;
1940  }
1941  }
1942  }
1943 
1944  // Increment the current overlap index
1945  ixOverlap += nOverlapFaces;
1946  }
1947 
1948  return;
1949 }
1950 
1951 ///////////////////////////////////////////////////////////////////////////////
1952 
1953 void moab::TempestOnlineMap::LinearRemapGLLtoGLL2_Pointwise_MOAB( const DataArray3D< int >& dataGLLNodesIn,
1954  const DataArray3D< double >& /*dataGLLJacobianIn*/,
1955  const DataArray3D< int >& dataGLLNodesOut,
1956  const DataArray3D< double >& /*dataGLLJacobianOut*/,
1957  const DataArray1D< double >& dataNodalAreaOut,
1958  int nPin,
1959  int nPout,
1960  int nMonotoneType,
1961  bool fContinuousIn,
1962  bool fContinuousOut )
1963 {
1964  // Gauss-Lobatto quadrature within Faces
1965  DataArray1D< double > dGL;
1966  DataArray1D< double > dWL;
1967 
1968  GaussLobattoQuadrature::GetPoints( nPout, 0.0, 1.0, dGL, dWL );
1969 
1970  // Get SparseMatrix represntation of the OfflineMap
1971  SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1972 
1973  // Sample coefficients
1974  DataArray2D< double > dSampleCoeffIn( nPin, nPin );
1975 
1976  // Announcemnets
1977  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
1978  dbgprint.set_prefix( "[LinearRemapGLLtoGLL2_Pointwise_MOAB]: " );
1979  if( is_root )
1980  {
1981  dbgprint.printf( 0, "Finite Element to Finite Element (Pointwise) Projection\n" );
1982  dbgprint.printf( 0, "Order of the input FE polynomial interpolant: %i\n", nPin );
1983  dbgprint.printf( 0, "Order of the output FE polynomial interpolant: %i\n", nPout );
1984  }
1985 
1986  // Number of overlap Faces per source Face
1987  DataArray1D< int > nAllOverlapFaces( m_meshInputCov->faces.size() );
1988 
1989  int ixOverlap = 0;
1990 
1991  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1992  {
1993  size_t ixOverlapTemp = ixOverlap;
1994  for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1995  {
1996  // const Face & faceOverlap = m_meshOverlap->faces[ixOverlapTemp];
1997 
1998  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 ) break;
1999 
2000  nAllOverlapFaces[ixFirst]++;
2001  }
2002 
2003  // Increment the current overlap index
2004  ixOverlap += nAllOverlapFaces[ixFirst];
2005  }
2006 
2007  // Number of times this point was found
2008  DataArray1D< bool > fSecondNodeFound( dataNodalAreaOut.GetRows() );
2009 
2010  ixOverlap = 0;
2011 #ifdef VERBOSE
2012  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
2013 #endif
2014  // Loop through all faces on m_meshInputCov
2015  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
2016  {
2017 #ifdef VERBOSE
2018  // Announce computation progress
2019  if( ixFirst % outputFrequency == 0 && is_root )
2020  {
2021  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
2022  }
2023 #endif
2024  // Quantities from the First Mesh
2025  const Face& faceFirst = m_meshInputCov->faces[ixFirst];
2026 
2027  const NodeVector& nodesFirst = m_meshInputCov->nodes;
2028 
2029  // Number of overlapping Faces and triangles
2030  int nOverlapFaces = nAllOverlapFaces[ixFirst];
2031 
2032  // Loop through all Overlap Faces
2033  for( int i = 0; i < nOverlapFaces; i++ )
2034  {
2035  // Quantities from the Second Mesh
2036  int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
2037 
2038  // signal to not participate, because it is a ghost target
2039  if( ixSecond < 0 ) continue; // do not do anything
2040 
2041  const NodeVector& nodesSecond = m_meshOutput->nodes;
2042  const Face& faceSecond = m_meshOutput->faces[ixSecond];
2043 
2044  // Loop through all nodes on the second face
2045  for( int s = 0; s < nPout; s++ )
2046  {
2047  for( int t = 0; t < nPout; t++ )
2048  {
2049  size_t ixSecondNode;
2050  if( fContinuousOut )
2051  {
2052  ixSecondNode = dataGLLNodesOut[s][t][ixSecond] - 1;
2053  }
2054  else
2055  {
2056  ixSecondNode = ixSecond * nPout * nPout + s * nPout + t;
2057  }
2058 
2059  if( ixSecondNode >= fSecondNodeFound.GetRows() ) _EXCEPTIONT( "Logic error" );
2060 
2061  // Check if this node has been found already
2062  if( fSecondNodeFound[ixSecondNode] ) continue;
2063 
2064  // Check this node
2065  Node node;
2066  Node dDx1G;
2067  Node dDx2G;
2068 
2069  ApplyLocalMap( faceSecond, nodesSecond, dGL[t], dGL[s], node, dDx1G, dDx2G );
2070 
2071  // Find the components of this quadrature point in the basis
2072  // of the first Face.
2073  double dAlphaIn;
2074  double dBetaIn;
2075 
2076  ApplyInverseMap( faceFirst, nodesFirst, node, dAlphaIn, dBetaIn );
2077 
2078  // Check if this node is within the first Face
2079  if( ( dAlphaIn < -1.0e-10 ) || ( dAlphaIn > 1.0 + 1.0e-10 ) || ( dBetaIn < -1.0e-10 ) ||
2080  ( dBetaIn > 1.0 + 1.0e-10 ) )
2081  continue;
2082 
2083  // Node is within the overlap region, mark as found
2084  fSecondNodeFound[ixSecondNode] = true;
2085 
2086  // Sample the First finite element at this point
2087  SampleGLLFiniteElement( nMonotoneType, nPin, dAlphaIn, dBetaIn, dSampleCoeffIn );
2088 
2089  // Add to map
2090  for( int p = 0; p < nPin; p++ )
2091  {
2092  for( int q = 0; q < nPin; q++ )
2093  {
2094  int ixFirstNode;
2095  if( fContinuousIn )
2096  {
2097  ixFirstNode = dataGLLNodesIn[p][q][ixFirst] - 1;
2098  }
2099  else
2100  {
2101  ixFirstNode = ixFirst * nPin * nPin + p * nPin + q;
2102  }
2103 
2104  smatMap( ixSecondNode, ixFirstNode ) += dSampleCoeffIn[p][q];
2105  }
2106  }
2107  }
2108  }
2109  }
2110 
2111  // Increment the current overlap index
2112  ixOverlap += nOverlapFaces;
2113  }
2114 
2115  // Check for missing samples
2116  for( size_t i = 0; i < fSecondNodeFound.GetRows(); i++ )
2117  {
2118  if( !fSecondNodeFound[i] )
2119  {
2120  _EXCEPTION1( "Can't sample point %i", i );
2121  }
2122  }
2123 
2124  return;
2125 }
2126 
2127 ///////////////////////////////////////////////////////////////////////////////