Mesh Oriented datABase  (version 5.6.0)
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"
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 #endif
369 
370 ///////////////////////////////////////////////////////////////////////////////
371 
372 template < typename T >
373 static std::vector< size_t > sort_indexes( const std::vector< T >& v )
374 {
375  // initialize original index locations
376  std::vector< size_t > idx( v.size() );
377  std::iota( idx.begin(), idx.end(), 0 );
378 
379  // sort indexes based on comparing values in v
380  // using std::stable_sort instead of std::sort
381  // to avoid unnecessary index re-orderings
382  // when v contains elements of equal values
383  std::stable_sort( idx.begin(), idx.end(), [&v]( size_t i1, size_t i2 ) { return fabs( v[i1] ) > fabs( v[i2] ); } );
384 
385  return idx;
386 }
387 
388 double moab::TempestOnlineMap::QLTLimiter( int caasIteration,
389  std::vector< double >& dataCorrectedField,
390  std::vector< double >& dataLowerBound,
391  std::vector< double >& dataUpperBound,
392  std::vector< double >& dMassDefect )
393 {
394  const size_t nrows = dataCorrectedField.size();
395  double dMassL = 0.0;
396  double dMassU = 0.0;
397  std::vector< double > dataCorrection( nrows );
398  double dMassDiffCum = 0.0;
399  double dLMinusU = fabs( dataUpperBound[0] - dataLowerBound[0] );
400  const DataArray1D< double >& dTargetAreas = this->m_remapper->m_target->vecFaceArea;
401 
402  // std::vector< size_t > sortedIdx = sort_indexes( dMassDefect );
403  std::vector< std::unordered_set< int > > vecAdjTargetFaces( nrows );
404  constexpr bool useMOABAdjacencies = true;
405 #ifdef USE_ComputeAdjacencyRelations
406  if( useMOABAdjacencies )
407  ComputeAdjacencyRelations( vecAdjTargetFaces, caasIteration, m_remapper->m_target_entities,
408  useMOABAdjacencies );
409  else
410  ComputeAdjacencyRelations( vecAdjTargetFaces, caasIteration, m_remapper->m_target_entities, useMOABAdjacencies,
411  this->m_remapper->m_target );
412 #else
413  moab::MeshTopoUtil mtu( m_interface );
414  ;
415 #endif
416 
417  for( size_t i = 0; i < nrows; i++ )
418  {
419  // size_t index = sortedIdx[i];
420  size_t index = i;
421  dataCorrection[index] = fmax( dataLowerBound[index], fmin( dataUpperBound[index], 0.0 ) );
422  // dMassDiff[index] = dMassDefect[index] - dTargetAreas[index] * dataCorrection[index];
423  // dMassDiff[index] = dMassDefect[index];
424 
425  dMassL += dTargetAreas[index] * dataLowerBound[index];
426  dMassU += dTargetAreas[index] * dataUpperBound[index];
427  dLMinusU = fmax( dLMinusU, fabs( dataUpperBound[index] - dataLowerBound[index] ) );
428  dMassDiffCum += dMassDefect[index] - dTargetAreas[index] * dataCorrection[index];
429 
430 #ifndef USE_ComputeAdjacencyRelations
431  vecAdjTargetFaces[index].insert( index ); // add self target face first
432  {
433  // Compute the adjacent faces to the target face
434  if( useMOABAdjacencies )
435  {
436  moab::Range ents;
437  // ents.insert( m_remapper->m_target_entities.index( m_remapper->m_target_entities[index] ) );
438  ents.insert( m_remapper->m_target_entities[index] );
439  moab::Range adjEnts;
440  moab::ErrorCode rval = mtu.get_bridge_adjacencies( ents, 0, 2, adjEnts, caasIteration );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
441  for( moab::Range::iterator it = adjEnts.begin(); it != adjEnts.end(); ++it )
442  {
443  // int adjIndex = m_interface->id_from_handle(*it)-1;
444  int adjIndex = m_remapper->m_target_entities.index( *it );
445  // printf("rank: %d, Element %lu, entity: %lu, adjIndex %d\n", rank, index, *it, adjIndex);
446  if( adjIndex >= 0 ) vecAdjTargetFaces[index].insert( adjIndex );
447  }
448  }
449  else
450  {
451  AdjacentFaceVector vecAdjFaces;
452  GetAdjacentFaceVectorByEdge( *this->m_remapper->m_target, index,
453  ( m_output_order + 1 ) * ( m_output_order + 1 ) * ( m_output_order + 1 ),
454  // ( m_output_order + 1 ) * ( m_output_order + 1 ),
455  // ( 4 ) * ( m_output_order + 1 ) * ( m_output_order + 1 ),
456  vecAdjFaces );
457 
458  // Add the adjacent faces to the target face list
459  for( auto adjFace : vecAdjFaces )
460  if( adjFace.first >= 0 )
461  vecAdjTargetFaces[index].insert( adjFace.first ); // map target face to source face
462  }
463  }
464 #endif
465  }
466 
467 #ifdef MOAB_HAVE_MPI
468  std::vector< double > localDefects( 5, 0.0 ), globalDefects( 5, 0.0 );
469  localDefects[0] = dMassL;
470  localDefects[1] = dMassU;
471  localDefects[2] = dMassDiffCum;
472  localDefects[3] = dLMinusU;
473  // localDefects[4] = dMassCorrectU;
474 
475  MPI_Allreduce( localDefects.data(), globalDefects.data(), 4, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
476 
477  dMassL = globalDefects[0];
478  dMassU = globalDefects[1];
479  dMassDiffCum = globalDefects[2];
480  dLMinusU = globalDefects[3];
481  // dMassCorrectU = globalDefects[4];
482 #endif
483 
484  //If the upper and lower bounds are too close together, just clip
485  if( fabs( dMassDiffCum ) < 1e-15 || dLMinusU < 1e-15 )
486  {
487  for( size_t i = 0; i < nrows; i++ )
488  dataCorrectedField[i] += dataCorrection[i];
489  return dMassDiffCum;
490  }
491  else
492  {
493  if( dMassL > dMassDiffCum )
494  {
495  Announce( "Lower bound mass exceeds target mass by %1.15e: CAAS will need another iteration",
496  dMassL - dMassDiffCum );
497  dMassDiffCum = dMassL;
498  // dMass -= dMassL;
499  }
500  else if( dMassU < dMassDiffCum )
501  {
502  Announce( "Target mass exceeds upper bound mass by %1.15e: CAAS will need another iteration",
503  dMassDiffCum - dMassU );
504  dMassDiffCum = dMassU;
505  // dMass -= dMassU;
506  }
507 
508  // TODO: optimize away dataMassVec by a simple transient double within the loop
509  // DataArray1D< double > dataMassVec( nrows ); //vector of mass redistribution
510  for( size_t i = 0; i < nrows; i++ )
511  {
512  // size_t index = sortedIdx[i];
513  size_t index = i;
514  const std::unordered_set< int >& neighbors = vecAdjTargetFaces[index];
515  if( dMassDefect[index] > 0.0 )
516  {
517  double dMassCorrectU = 0.0;
518  for( auto it : neighbors )
519  dMassCorrectU += dTargetAreas[it] * ( dataUpperBound[it] - dataCorrection[it] );
520 
521  // double dMassDiffCumOld = dMassDefect[index];
522  for( auto it : neighbors )
523  dataCorrection[it] +=
524  dMassDefect[index] * ( dataUpperBound[it] - dataCorrection[it] ) / dMassCorrectU;
525  }
526  else
527  {
528  double dMassCorrectL = 0.0;
529  for( auto it : neighbors )
530  dMassCorrectL += dTargetAreas[it] * ( dataCorrection[it] - dataLowerBound[it] );
531 
532  // double dMassDiffCumOld = dMassDefect[index];
533  for( auto it : neighbors )
534  dataCorrection[it] +=
535  dMassDefect[index] * ( dataCorrection[it] - dataLowerBound[it] ) / dMassCorrectL;
536  }
537  }
538 
539  for( size_t i = 0; i < nrows; i++ )
540  dataCorrectedField[i] += dataCorrection[i];
541  }
542 
543  return dMassDiffCum;
544 }
545 
546 void moab::TempestOnlineMap::CAASLimiter( std::vector< double >& dataCorrectedField,
547  std::vector< double >& dataLowerBound,
548  std::vector< double >& dataUpperBound,
549  double& dMass )
550 {
551  const size_t nrows = dataCorrectedField.size();
552  double dMassL = 0.0;
553  double dMassU = 0.0;
554  std::vector< double > dataCorrection( nrows );
555  const DataArray1D< double >& dTargetAreas = this->m_remapper->m_target->vecFaceArea;
556  double dMassDiff = dMass;
557  double dLMinusU = fabs( dataUpperBound[0] - dataLowerBound[0] );
558  double dMassCorrectU = 0.0;
559  double dMassCorrectL = 0.0;
560  for( size_t i = 0; i < nrows; i++ )
561  {
562  dataCorrection[i] = fmax( dataLowerBound[i], fmin( dataUpperBound[i], 0.0 ) );
563  dMassL += dTargetAreas[i] * dataLowerBound[i];
564  dMassU += dTargetAreas[i] * dataUpperBound[i];
565  dMassDiff -= dTargetAreas[i] * dataCorrection[i];
566  dLMinusU = fmax( dLMinusU, fabs( dataUpperBound[i] - dataLowerBound[i] ) );
567  dMassCorrectL += dTargetAreas[i] * ( dataCorrection[i] - dataLowerBound[i] );
568  dMassCorrectU += dTargetAreas[i] * ( dataUpperBound[i] - dataCorrection[i] );
569  }
570 
571 #ifdef MOAB_HAVE_MPI
572  std::vector< double > localDefects( 5, 0.0 ), globalDefects( 5, 0.0 );
573  localDefects[0] = dMassL;
574  localDefects[1] = dMassU;
575  localDefects[2] = dMassDiff;
576  localDefects[3] = dMassCorrectL;
577  localDefects[4] = dMassCorrectU;
578 
579  MPI_Allreduce( localDefects.data(), globalDefects.data(), 5, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
580 
581  dMassL = globalDefects[0];
582  dMassU = globalDefects[1];
583  dMassDiff = globalDefects[2];
584  dMassCorrectL = globalDefects[3];
585  dMassCorrectU = globalDefects[4];
586 #endif
587 
588  //If the upper and lower bounds are too close together, just clip
589  if( fabs( dMassDiff ) < 1e-15 || fabs( dLMinusU ) < 1e-15 )
590  {
591  for( size_t i = 0; i < nrows; i++ )
592  dataCorrectedField[i] += dataCorrection[i];
593  return;
594  }
595  else
596  {
597  if( dMassL > dMassDiff )
598  {
599  Announce( "%d: Lower bound mass exceeds target mass by %1.15e: CAAS will need another iteration", rank,
600  dMassL - dMassDiff );
601  dMassDiff = dMassL;
602  dMass -= dMassL;
603  }
604  else if( dMassU < dMassDiff )
605  {
606  Announce( "%d: Target mass exceeds upper bound mass by %1.15e: CAAS will need another iteration", rank,
607  dMassDiff - dMassU );
608  dMassDiff = dMassU;
609  dMass -= dMassU;
610  }
611 
612  // TODO: optimize away dataMassVec by a simple transient double within the loop
613  DataArray1D< double > dataMassVec( nrows ); //vector of mass redistribution
614  if( dMassDiff > 0.0 )
615  {
616  for( size_t i = 0; i < nrows; i++ )
617  {
618  dataMassVec[i] = ( dataUpperBound[i] - dataCorrection[i] ) / dMassCorrectU;
619  dataCorrection[i] += dMassDiff * dataMassVec[i];
620  }
621  }
622  else
623  {
624  for( size_t i = 0; i < nrows; i++ )
625  {
626  dataMassVec[i] = ( dataCorrection[i] - dataLowerBound[i] ) / dMassCorrectL;
627  dataCorrection[i] += dMassDiff * dataMassVec[i];
628  }
629  }
630 
631  for( size_t i = 0; i < nrows; i++ )
632  dataCorrectedField[i] += dataCorrection[i];
633  }
634 
635  return;
636 }
637 
638 std::pair< double, double > moab::TempestOnlineMap::ApplyBoundsLimiting( std::vector< double >& dataInDouble,
639  std::vector< double >& dataOutDouble,
640  CAASType caasType,
641  int caasIteration,
642  double mismatch )
643 {
644  // Currently only implemented for FV to FV remapping
645  // We should generalize this to other types of remapping
646  assert( !dataGLLNodesSrcCov.IsAttached() && !dataGLLNodesDest.IsAttached() );
647 
648  std::pair< double, double > massDefect( 0.0, 0.0 );
649 
650  // Check if the source and target data are of the same size
651  const size_t nTargetCount = dataOutDouble.size();
652  const DataArray1D< double >& m_dOverlapAreas = this->m_remapper->m_overlap->vecFaceArea;
653 
654  // Apply the offline map to the data
655  double dMassDiff = 0.0;
656  std::vector< double > x( nTargetCount );
657  std::vector< double > dataLowerBound( nTargetCount );
658  std::vector< double > dataUpperBound( nTargetCount );
659  std::vector< double > massVector( nTargetCount );
660  std::vector< std::unordered_set< int > > vecSourceOvTarget( nTargetCount );
661 
662 #undef USE_ComputeAdjacencyRelations
663  constexpr bool useMOABAdjacencies = true;
664 #ifdef USE_ComputeAdjacencyRelations
665  // Compute the adjacent faces to the source face
666  // However, calling MOAB to do this does not work correctly as we need ixS to be the index
667  // Cannot just iterate over all entities in the source covering mesh
668  if( caasType == CAAS_QLT || caasType == CAAS_LOCAL_ADJACENT )
669  {
670  if( useMOABAdjacencies )
671  {
672  moab::ErrorCode rval =
673  ComputeAdjacencyRelations( vecSourceOvTarget, caasIteration, m_remapper->m_covering_source_entities,
674  useMOABAdjacencies );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
675  }
676  else
677  {
678  moab::ErrorCode rval =
679  ComputeAdjacencyRelations( vecSourceOvTarget, caasIteration, m_remapper->m_covering_source_entities,
680  useMOABAdjacencies, m_meshInputCov );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
681  }
682  }
683 #else
684  moab::MeshTopoUtil mtu( m_interface );
685 #endif
686 
687  // Initialize the bounds on the given source and target data
688  double dSourceMin = dataInDouble[0];
689  double dSourceMax = dataInDouble[0];
690  double dTargetMin = dataOutDouble[0];
691  double dTargetMax = dataOutDouble[0];
692  for( size_t i = 0; i < m_meshOverlap->faces.size(); i++ )
693  {
694  const int ixS = m_meshOverlap->vecSourceFaceIx[i];
695  const int ixT = m_meshOverlap->vecTargetFaceIx[i];
696 
697  if( ixT < 0 ) continue; // skip ghost target faces
698 
699  assert( m_dOverlapAreas[i] > 0.0 );
700  assert( ixS >= 0 );
701  assert( ixT >= 0 );
702 
703 #ifndef USE_ComputeAdjacencyRelations
704  // Compute the adjacent faces to the target face
705  vecSourceOvTarget[ixT].insert( ixS ); // map target face to source face
706  if( ( caasType == CAAS_QLT || caasType == CAAS_LOCAL_ADJACENT ) )
707  {
708  if( useMOABAdjacencies )
709  {
710  moab::Range ents;
711  ents.insert( m_remapper->m_covering_source_entities[ixS] );
712  moab::Range adjEnts;
713  moab::ErrorCode rval = mtu.get_bridge_adjacencies( ents, 0, 2, adjEnts, caasIteration );MB_CHK_SET_ERR_CONT( rval, "Failed to get adjacent faces" );
714  for( moab::Range::iterator it = adjEnts.begin(); it != adjEnts.end(); ++it )
715  {
716  int adjIndex = m_remapper->m_covering_source_entities.index( *it );
717  if( adjIndex >= 0 ) vecSourceOvTarget[ixT].insert( adjIndex );
718  }
719  }
720  else
721  {
722  // Compute the adjacent faces to the target face
723  AdjacentFaceVector vecAdjFaces;
724  GetAdjacentFaceVectorByEdge( *m_meshInputCov, ixS,
725  ( caasIteration ) * ( m_input_order + 1 ) * ( m_input_order + 1 ),
726  vecAdjFaces );
727 
728  //Compute min/max over neighboring faces
729  for( size_t iadj = 0; iadj < vecAdjFaces.size(); iadj++ )
730  vecSourceOvTarget[ixT].insert( vecAdjFaces[iadj].first ); // map target face to source face
731  }
732  }
733 #endif
734 
735  // Update the min and max values of the source data
736  dSourceMax = fmax( dSourceMax, dataInDouble[ixS] );
737  dSourceMin = fmin( dSourceMin, dataInDouble[ixS] );
738 
739  // Update the min and max values of the target data
740  dTargetMin = fmin( dTargetMin, dataOutDouble[ixT] );
741  dTargetMax = fmax( dTargetMax, dataOutDouble[ixT] );
742 
743  const double locMassDiff = ( dataInDouble[ixS] * m_dOverlapAreas[i] ) - // source mass
744  ( dataOutDouble[ixT] * m_dOverlapAreas[i] ); // target mass
745 
746  // Update the mass difference between source and target faces
747  // linked to the overlap mesh element
748  dMassDiff += locMassDiff; // target mass
749  massVector[ixT] += locMassDiff;
750  }
751 
752 #ifdef MOAB_HAVE_MPI
753  std::vector< double > localMinMaxDefects( 5, 0.0 ), globalMinMaxDefects( 5, 0.0 );
754  localMinMaxDefects[0] = dSourceMin;
755  localMinMaxDefects[1] = dTargetMin;
756  localMinMaxDefects[2] = dSourceMax;
757  localMinMaxDefects[3] = dTargetMax;
758  localMinMaxDefects[4] = dMassDiff;
759 
760  if( caasType == CAAS_GLOBAL )
761  {
762  MPI_Allreduce( localMinMaxDefects.data(), globalMinMaxDefects.data(), 2, MPI_DOUBLE, MPI_MIN, m_pcomm->comm() );
763  MPI_Allreduce( localMinMaxDefects.data() + 2, globalMinMaxDefects.data() + 2, 2, MPI_DOUBLE, MPI_MAX,
764  m_pcomm->comm() );
765  dSourceMin = globalMinMaxDefects[0];
766  dSourceMax = globalMinMaxDefects[2];
767  dTargetMin = globalMinMaxDefects[1];
768  dTargetMax = globalMinMaxDefects[3];
769  }
770  if( caasIteration == 1 )
771  MPI_Allreduce( localMinMaxDefects.data() + 4, globalMinMaxDefects.data() + 4, 1, MPI_DOUBLE, MPI_SUM,
772  m_pcomm->comm() );
773  else
774  globalMinMaxDefects[4] = mismatch;
775 
776  dMassDiff = localMinMaxDefects[4];
777  // massDefect.first = localMinMaxDefects[4];
778  massDefect.first = globalMinMaxDefects[4];
779 #else
780 
781  // massDefect.first = fabs( dMassDiff / ( dSourceMax - dSourceMin ) );
782  massDefect.first = dMassDiff;
783 #endif
784 
785  // Early exit if the values are monotone already.
786  // if( ( dTargetMax <= dSourceMax && dTargetMin <= dSourceMin ) || fabs( massDefect.first ) < 1e-16 )
787  if( fabs( massDefect.first ) > 1e-20 )
788  {
789  if( caasType == CAAS_GLOBAL )
790  {
791  for( size_t i = 0; i < nTargetCount; i++ )
792  {
793  dataLowerBound[i] = dSourceMin - dataOutDouble[i];
794  dataUpperBound[i] = dSourceMax - dataOutDouble[i];
795  }
796  } // if( caasType == CAAS_GLOBAL )
797  else // caasType == CAAS_LOCAL
798  {
799  // Compute the local min and max values of the target data
800  std::vector< double > vecLocalUpperBound( nTargetCount );
801  std::vector< double > vecLocalLowerBound( nTargetCount );
802  // Loop over the target faces and compute the min and max values
803  // of the source data linked to the target faces
804  for( size_t i = 0; i < nTargetCount; i++ )
805  {
806  assert( vecSourceOvTarget[i].size() );
807 
808  double dMinI = 1E10; // dataInDouble[vecSourceOvTarget[i][0]];
809  double dMaxI = -1E10; // dataInDouble[vecSourceOvTarget[i][0]];
810 
811  // Compute max over intersecting source faces
812  for( const auto& srcElem : vecSourceOvTarget[i] )
813  {
814  dMinI = fmin( dMinI, dataInDouble[srcElem] ); // min over intersecting source faces
815  dMaxI = fmax( dMaxI, dataInDouble[srcElem] ); // max over intersecting source faces
816  }
817 
818  // Update the min and max values of the target data
819  vecLocalLowerBound[i] = dMinI;
820  vecLocalUpperBound[i] = dMaxI;
821  }
822 
823  for( size_t i = 0; i < nTargetCount; i++ )
824  {
825  dataLowerBound[i] = vecLocalLowerBound[i] - dataOutDouble[i];
826  dataUpperBound[i] = vecLocalUpperBound[i] - dataOutDouble[i];
827  }
828  } // caasType == CAAS_LOCAL
829 
830  // Invoke CAAS or QLT application on the map
831  if( fabs( dMassDiff ) > 1e-20 )
832  {
833  if( caasType == CAAS_QLT )
834  dMassDiff = QLTLimiter( caasIteration, dataOutDouble, dataLowerBound, dataUpperBound, massVector );
835  else
836  CAASLimiter( dataOutDouble, dataLowerBound, dataUpperBound, dMassDiff );
837  }
838 
839  // Announce output mass
840  double dMassDiffPost = 0.0;
841  for( size_t i = 0; i < m_meshOverlap->faces.size(); i++ )
842  {
843  const int ixS = m_meshOverlap->vecSourceFaceIx[i];
844  const int ixT = m_meshOverlap->vecTargetFaceIx[i];
845 
846  if( ixT < 0 ) continue; // skip ghost target faces
847 
848  // Update the mass difference between source and target faces
849  // linked to the overlap mesh element
850  dMassDiffPost += ( dataInDouble[ixS] * m_dOverlapAreas[i] ) - // source mass
851  ( dataOutDouble[ixT] * m_dOverlapAreas[i] ); // target mass
852  }
853  // massDefect.second = fabs( dMassDiffPost / ( dSourceMax - dSourceMin ) );
854  massDefect.second = dMassDiffPost;
855  }
856 
857  // Ideally should perform an AllReduce here to get the global mass difference across all processors
858  // But if we satisfy the constraint on every task, essentially, the global mass difference should be zero!
859  return massDefect;
860 }
861 
862 ///////////////////////////////////////////////////////////////////////////////
863 
864 extern void ForceConsistencyConservation3( const DataArray1D< double >& vecSourceArea,
865  const DataArray1D< double >& vecTargetArea,
866  DataArray2D< double >& dCoeff,
867  bool fMonotone,
868  bool fSparseConstraints = false );
869 
870 ///////////////////////////////////////////////////////////////////////////////
871 
872 extern void ForceIntArrayConsistencyConservation( const DataArray1D< double >& vecSourceArea,
873  const DataArray1D< double >& vecTargetArea,
874  DataArray2D< double >& dCoeff,
875  bool fMonotone );
876 
877 ///////////////////////////////////////////////////////////////////////////////
878 
879 void moab::TempestOnlineMap::LinearRemapSE4_Tempest_MOAB( const DataArray3D< int >& dataGLLNodes,
880  const DataArray3D< double >& dataGLLJacobian,
881  int nMonotoneType,
882  bool fContinuousIn,
883  bool fNoConservation )
884 {
885  // Order of the polynomial interpolant
886  int nP = dataGLLNodes.GetRows();
887 
888  // Order of triangular quadrature rule
889  const int TriQuadRuleOrder = 4;
890 
891  // Triangular quadrature rule
892  TriangularQuadratureRule triquadrule( TriQuadRuleOrder );
893 
894  int TriQuadraturePoints = triquadrule.GetPoints();
895 
896  const DataArray2D< double >& TriQuadratureG = triquadrule.GetG();
897 
898  const DataArray1D< double >& TriQuadratureW = triquadrule.GetW();
899 
900  // Sample coefficients
901  DataArray2D< double > dSampleCoeff( nP, nP );
902 
903  // GLL Quadrature nodes on quadrilateral elements
904  DataArray1D< double > dG;
905  DataArray1D< double > dW;
906  GaussLobattoQuadrature::GetPoints( nP, 0.0, 1.0, dG, dW );
907 
908  // Announcements
909  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
910  dbgprint.set_prefix( "[LinearRemapSE4_Tempest_MOAB]: " );
911  if( is_root )
912  {
913  dbgprint.printf( 0, "Finite Element to Finite Volume Projection\n" );
914  dbgprint.printf( 0, "Triangular quadrature rule order %i\n", TriQuadRuleOrder );
915  dbgprint.printf( 0, "Order of the FE polynomial interpolant: %i\n", nP );
916  }
917 
918  // Get SparseMatrix represntation of the OfflineMap
919  SparseMatrix< double >& smatMap = this->GetSparseMatrix();
920 
921  // NodeVector from m_meshOverlap
922  const NodeVector& nodesOverlap = m_meshOverlap->nodes;
923  const NodeVector& nodesFirst = m_meshInputCov->nodes;
924 
925  // Vector of source areas
926  DataArray1D< double > vecSourceArea( nP * nP );
927 
928  DataArray1D< double > vecTargetArea;
929  DataArray2D< double > dCoeff;
930 
931 #ifdef VERBOSE
932  std::stringstream sstr;
933  sstr << "remapdata_" << rank << ".txt";
934  std::ofstream output_file( sstr.str() );
935 #endif
936 
937  // Current Overlap Face
938  int ixOverlap = 0;
939 #ifdef VERBOSE
940  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
941 #endif
942  // generic triangle used for area computation, for triangles around the center of overlap face;
943  // used for overlap faces with more than 4 edges;
944  // nodes array will be set for each triangle;
945  // these triangles are not part of the mesh structure, they are just temporary during
946  // aforementioned decomposition.
947  Face faceTri( 3 );
948  NodeVector nodes( 3 );
949  faceTri.SetNode( 0, 0 );
950  faceTri.SetNode( 1, 1 );
951  faceTri.SetNode( 2, 2 );
952 
953  // Loop over all input Faces
954  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
955  {
956  const Face& faceFirst = m_meshInputCov->faces[ixFirst];
957 
958  if( faceFirst.edges.size() != 4 )
959  {
960  _EXCEPTIONT( "Only quadrilateral elements allowed for SE remapping" );
961  }
962 #ifdef VERBOSE
963  // Announce computation progress
964  if( ixFirst % outputFrequency == 0 && is_root )
965  {
966  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
967  }
968 #endif
969  // Need to re-number the overlap elements such that vecSourceFaceIx[a:b] = 0, then 1 and so
970  // on wrt the input mesh data Then the overlap_end and overlap_begin will be correct.
971  // However, the relation with MOAB and Tempest will go out of the roof
972 
973  // Determine how many overlap Faces and triangles are present
974  int nOverlapFaces = 0;
975  size_t ixOverlapTemp = ixOverlap;
976  for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
977  {
978  // if( m_meshOverlap->vecTargetFaceIx[ixOverlapTemp] < 0 ) continue; // skip ghost target faces
979  // const Face & faceOverlap = m_meshOverlap->faces[ixOverlapTemp];
980  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 ) break;
981 
982  nOverlapFaces++;
983  }
984 
985  // No overlaps
986  if( nOverlapFaces == 0 ) continue;
987 
988  // Allocate remap coefficients array for meshFirst Face
989  DataArray3D< double > dRemapCoeff( nP, nP, nOverlapFaces );
990 
991  // Find the local remap coefficients
992  for( int j = 0; j < nOverlapFaces; j++ )
993  {
994  const Face& faceOverlap = m_meshOverlap->faces[ixOverlap + j];
995  if( m_meshOverlap->vecFaceArea[ixOverlap + j] < 1.e-16 ) // machine precision
996  {
997  Announce( "Very small overlap at index %i area polygon: (%1.10e )", ixOverlap + j,
998  m_meshOverlap->vecFaceArea[ixOverlap + j] );
999  int n = faceOverlap.edges.size();
1000  Announce( "Number nodes: %d", n );
1001  for( int k = 0; k < n; k++ )
1002  {
1003  Node nd = nodesOverlap[faceOverlap[k]];
1004  Announce( "Node %d %d : %1.10e %1.10e %1.10e ", k, faceOverlap[k], nd.x, nd.y, nd.z );
1005  }
1006  continue;
1007  }
1008 
1009  // #ifdef VERBOSE
1010  // if ( is_root )
1011  // Announce ( "\tLocal ID: %i/%i = %i, areas = %2.8e", j + ixOverlap, nOverlapFaces,
1012  // m_remapper->lid_to_gid_covsrc[m_meshOverlap->vecSourceFaceIx[ixOverlap + j]],
1013  // m_meshOverlap->vecFaceArea[ixOverlap + j] );
1014  // #endif
1015 
1016  int nbEdges = faceOverlap.edges.size();
1017  int nOverlapTriangles = 1;
1018  Node center; // not used if nbEdges == 3
1019  if( nbEdges > 3 )
1020  { // decompose from center in this case
1021  nOverlapTriangles = nbEdges;
1022  for( int k = 0; k < nbEdges; k++ )
1023  {
1024  const Node& node = nodesOverlap[faceOverlap[k]];
1025  center = center + node;
1026  }
1027  center = center / nbEdges;
1028  center = center.Normalized(); // project back on sphere of radius 1
1029  }
1030 
1031  Node node0, node1, node2;
1032  double dTriangleArea;
1033 
1034  // Loop over all sub-triangles of this Overlap Face
1035  for( int k = 0; k < nOverlapTriangles; k++ )
1036  {
1037  if( nbEdges == 3 ) // will come here only once, nOverlapTriangles == 1 in this case
1038  {
1039  node0 = nodesOverlap[faceOverlap[0]];
1040  node1 = nodesOverlap[faceOverlap[1]];
1041  node2 = nodesOverlap[faceOverlap[2]];
1042  dTriangleArea = CalculateFaceArea( faceOverlap, nodesOverlap );
1043  }
1044  else // decompose polygon in triangles around the center
1045  {
1046  node0 = center;
1047  node1 = nodesOverlap[faceOverlap[k]];
1048  int k1 = ( k + 1 ) % nbEdges;
1049  node2 = nodesOverlap[faceOverlap[k1]];
1050  nodes[0] = center;
1051  nodes[1] = node1;
1052  nodes[2] = node2;
1053  dTriangleArea = CalculateFaceArea( faceTri, nodes );
1054  }
1055  // Coordinates of quadrature Node
1056  for( int l = 0; l < TriQuadraturePoints; l++ )
1057  {
1058  Node nodeQuadrature;
1059  nodeQuadrature.x = TriQuadratureG[l][0] * node0.x + TriQuadratureG[l][1] * node1.x +
1060  TriQuadratureG[l][2] * node2.x;
1061 
1062  nodeQuadrature.y = TriQuadratureG[l][0] * node0.y + TriQuadratureG[l][1] * node1.y +
1063  TriQuadratureG[l][2] * node2.y;
1064 
1065  nodeQuadrature.z = TriQuadratureG[l][0] * node0.z + TriQuadratureG[l][1] * node1.z +
1066  TriQuadratureG[l][2] * node2.z;
1067 
1068  nodeQuadrature = nodeQuadrature.Normalized();
1069 
1070  // Find components of quadrature point in basis
1071  // of the first Face
1072  double dAlpha;
1073  double dBeta;
1074 
1075  ApplyInverseMap( faceFirst, nodesFirst, nodeQuadrature, dAlpha, dBeta );
1076 
1077  // Check inverse map value
1078  if( ( dAlpha < -1.0e-13 ) || ( dAlpha > 1.0 + 1.0e-13 ) || ( dBeta < -1.0e-13 ) ||
1079  ( dBeta > 1.0 + 1.0e-13 ) )
1080  {
1081  _EXCEPTION4( "Inverse Map for element %d and subtriangle %d out of range "
1082  "(%1.5e %1.5e)",
1083  j, l, dAlpha, dBeta );
1084  }
1085 
1086  // Sample the finite element at this point
1087  SampleGLLFiniteElement( nMonotoneType, nP, dAlpha, dBeta, dSampleCoeff );
1088 
1089  // Add sample coefficients to the map if m_meshOverlap->vecFaceArea[ixOverlap + j] > 0
1090  for( int p = 0; p < nP; p++ )
1091  {
1092  for( int q = 0; q < nP; q++ )
1093  {
1094  dRemapCoeff[p][q][j] += TriQuadratureW[l] * dTriangleArea * dSampleCoeff[p][q] /
1095  m_meshOverlap->vecFaceArea[ixOverlap + j];
1096  }
1097  }
1098  }
1099  }
1100  }
1101 
1102 #ifdef VERBOSE
1103  output_file << "[" << m_remapper->lid_to_gid_covsrc[ixFirst] << "] \t";
1104  for( int j = 0; j < nOverlapFaces; j++ )
1105  {
1106  for( int p = 0; p < nP; p++ )
1107  {
1108  for( int q = 0; q < nP; q++ )
1109  {
1110  output_file << dRemapCoeff[p][q][j] << " ";
1111  }
1112  }
1113  }
1114  output_file << std::endl;
1115 #endif
1116 
1117  // Force consistency and conservation
1118  if( !fNoConservation )
1119  {
1120  double dTargetArea = 0.0;
1121  for( int j = 0; j < nOverlapFaces; j++ )
1122  {
1123  dTargetArea += m_meshOverlap->vecFaceArea[ixOverlap + j];
1124  }
1125 
1126  for( int p = 0; p < nP; p++ )
1127  {
1128  for( int q = 0; q < nP; q++ )
1129  {
1130  vecSourceArea[p * nP + q] = dataGLLJacobian[p][q][ixFirst];
1131  }
1132  }
1133 
1134  const double areaTolerance = 1e-10;
1135  // Source elements are completely covered by target volumes
1136  if( fabs( m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea ) <= areaTolerance )
1137  {
1138  vecTargetArea.Allocate( nOverlapFaces );
1139  for( int j = 0; j < nOverlapFaces; j++ )
1140  {
1141  vecTargetArea[j] = m_meshOverlap->vecFaceArea[ixOverlap + j];
1142  }
1143 
1144  dCoeff.Allocate( nOverlapFaces, nP * nP );
1145 
1146  for( int j = 0; j < nOverlapFaces; j++ )
1147  {
1148  for( int p = 0; p < nP; p++ )
1149  {
1150  for( int q = 0; q < nP; q++ )
1151  {
1152  dCoeff[j][p * nP + q] = dRemapCoeff[p][q][j];
1153  }
1154  }
1155  }
1156 
1157  // Target volumes only partially cover source elements
1158  }
1159  else if( m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea > areaTolerance )
1160  {
1161  double dExtraneousArea = m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea;
1162 
1163  vecTargetArea.Allocate( nOverlapFaces + 1 );
1164  for( int j = 0; j < nOverlapFaces; j++ )
1165  {
1166  vecTargetArea[j] = m_meshOverlap->vecFaceArea[ixOverlap + j];
1167  }
1168  vecTargetArea[nOverlapFaces] = dExtraneousArea;
1169 
1170 #ifdef VERBOSE
1171  Announce( "Partial volume: %i (%1.10e / %1.10e)", ixFirst, dTargetArea,
1172  m_meshInputCov->vecFaceArea[ixFirst] );
1173 #endif
1174  if( dTargetArea > m_meshInputCov->vecFaceArea[ixFirst] )
1175  {
1176  _EXCEPTIONT( "Partial element area exceeds total element area" );
1177  }
1178 
1179  dCoeff.Allocate( nOverlapFaces + 1, nP * nP );
1180 
1181  for( int j = 0; j < nOverlapFaces; j++ )
1182  {
1183  for( int p = 0; p < nP; p++ )
1184  {
1185  for( int q = 0; q < nP; q++ )
1186  {
1187  dCoeff[j][p * nP + q] = dRemapCoeff[p][q][j];
1188  }
1189  }
1190  }
1191  for( int p = 0; p < nP; p++ )
1192  {
1193  for( int q = 0; q < nP; q++ )
1194  {
1195  dCoeff[nOverlapFaces][p * nP + q] = dataGLLJacobian[p][q][ixFirst];
1196  }
1197  }
1198  for( int j = 0; j < nOverlapFaces; j++ )
1199  {
1200  for( int p = 0; p < nP; p++ )
1201  {
1202  for( int q = 0; q < nP; q++ )
1203  {
1204  dCoeff[nOverlapFaces][p * nP + q] -=
1205  dRemapCoeff[p][q][j] * m_meshOverlap->vecFaceArea[ixOverlap + j];
1206  }
1207  }
1208  }
1209  for( int p = 0; p < nP; p++ )
1210  {
1211  for( int q = 0; q < nP; q++ )
1212  {
1213  dCoeff[nOverlapFaces][p * nP + q] /= dExtraneousArea;
1214  }
1215  }
1216 
1217  // Source elements only partially cover target volumes
1218  }
1219  else
1220  {
1221  Announce( "Coverage area: %1.10e, and target element area: %1.10e)", ixFirst,
1222  m_meshInputCov->vecFaceArea[ixFirst], dTargetArea );
1223  _EXCEPTIONT( "Target grid must be a subset of source grid" );
1224  }
1225 
1226  ForceConsistencyConservation3( vecSourceArea, vecTargetArea, dCoeff, ( nMonotoneType > 0 )
1227  /*, m_remapper->lid_to_gid_covsrc[ixFirst]*/ );
1228 
1229  for( int j = 0; j < nOverlapFaces; j++ )
1230  {
1231  for( int p = 0; p < nP; p++ )
1232  {
1233  for( int q = 0; q < nP; q++ )
1234  {
1235  dRemapCoeff[p][q][j] = dCoeff[j][p * nP + q];
1236  }
1237  }
1238  }
1239  }
1240 
1241 #ifdef VERBOSE
1242  // output_file << "[" << m_remapper->lid_to_gid_covsrc[ixFirst] << "] \t";
1243  // for ( int j = 0; j < nOverlapFaces; j++ )
1244  // {
1245  // for ( int p = 0; p < nP; p++ )
1246  // {
1247  // for ( int q = 0; q < nP; q++ )
1248  // {
1249  // output_file << dRemapCoeff[p][q][j] << " ";
1250  // }
1251  // }
1252  // }
1253  // output_file << std::endl;
1254 #endif
1255 
1256  // Put these remap coefficients into the SparseMatrix map
1257  for( int j = 0; j < nOverlapFaces; j++ )
1258  {
1259  int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
1260 
1261  // signal to not participate, because it is a ghost target
1262  if( ixSecondFace < 0 ) continue; // do not do anything
1263 
1264  for( int p = 0; p < nP; p++ )
1265  {
1266  for( int q = 0; q < nP; q++ )
1267  {
1268  if( fContinuousIn )
1269  {
1270  int ixFirstNode = dataGLLNodes[p][q][ixFirst] - 1;
1271 
1272  smatMap( ixSecondFace, ixFirstNode ) += dRemapCoeff[p][q][j] *
1273  m_meshOverlap->vecFaceArea[ixOverlap + j] /
1274  m_meshOutput->vecFaceArea[ixSecondFace];
1275  }
1276  else
1277  {
1278  int ixFirstNode = ixFirst * nP * nP + p * nP + q;
1279 
1280  smatMap( ixSecondFace, ixFirstNode ) += dRemapCoeff[p][q][j] *
1281  m_meshOverlap->vecFaceArea[ixOverlap + j] /
1282  m_meshOutput->vecFaceArea[ixSecondFace];
1283  }
1284  }
1285  }
1286  }
1287  // Increment the current overlap index
1288  ixOverlap += nOverlapFaces;
1289  }
1290 #ifdef VERBOSE
1291  output_file.flush(); // required here
1292  output_file.close();
1293 #endif
1294 
1295  return;
1296 }
1297 
1298 ///////////////////////////////////////////////////////////////////////////////
1299 
1300 void moab::TempestOnlineMap::LinearRemapGLLtoGLL2_MOAB( const DataArray3D< int >& dataGLLNodesIn,
1301  const DataArray3D< double >& dataGLLJacobianIn,
1302  const DataArray3D< int >& dataGLLNodesOut,
1303  const DataArray3D< double >& dataGLLJacobianOut,
1304  const DataArray1D< double >& dataNodalAreaOut,
1305  int nPin,
1306  int nPout,
1307  int nMonotoneType,
1308  bool fContinuousIn,
1309  bool fContinuousOut,
1310  bool fNoConservation )
1311 {
1312  // Triangular quadrature rule
1313  TriangularQuadratureRule triquadrule( 8 );
1314 
1315  const DataArray2D< double >& dG = triquadrule.GetG();
1316  const DataArray1D< double >& dW = triquadrule.GetW();
1317 
1318  // Get SparseMatrix represntation of the OfflineMap
1319  SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1320 
1321  // Sample coefficients
1322  DataArray2D< double > dSampleCoeffIn( nPin, nPin );
1323  DataArray2D< double > dSampleCoeffOut( nPout, nPout );
1324 
1325  // Announcemnets
1326  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
1327  dbgprint.set_prefix( "[LinearRemapGLLtoGLL2_MOAB]: " );
1328  if( is_root )
1329  {
1330  dbgprint.printf( 0, "Finite Element to Finite Element Projection\n" );
1331  dbgprint.printf( 0, "Order of the input FE polynomial interpolant: %i\n", nPin );
1332  dbgprint.printf( 0, "Order of the output FE polynomial interpolant: %i\n", nPout );
1333  }
1334 
1335  // Build the integration array for each element on m_meshOverlap
1336  DataArray3D< double > dGlobalIntArray( nPin * nPin, m_meshOverlap->faces.size(), nPout * nPout );
1337 
1338  // Number of overlap Faces per source Face
1339  DataArray1D< int > nAllOverlapFaces( m_meshInputCov->faces.size() );
1340 
1341  int ixOverlap = 0;
1342  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1343  {
1344  // Determine how many overlap Faces and triangles are present
1345  int nOverlapFaces = 0;
1346  size_t ixOverlapTemp = ixOverlap;
1347  for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1348  {
1349  // const Face & faceOverlap = m_meshOverlap->faces[ixOverlapTemp];
1350  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 )
1351  {
1352  break;
1353  }
1354 
1355  nOverlapFaces++;
1356  }
1357 
1358  nAllOverlapFaces[ixFirst] = nOverlapFaces;
1359 
1360  // Increment the current overlap index
1361  ixOverlap += nAllOverlapFaces[ixFirst];
1362  }
1363 
1364  // Geometric area of each output node
1365  DataArray2D< double > dGeometricOutputArea( m_meshOutput->faces.size(), nPout * nPout );
1366 
1367  // Area of each overlap element in the output basis
1368  DataArray2D< double > dOverlapOutputArea( m_meshOverlap->faces.size(), nPout * nPout );
1369 
1370  // Loop through all faces on m_meshInputCov
1371  ixOverlap = 0;
1372 #ifdef VERBOSE
1373  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
1374 #endif
1375  if( is_root ) dbgprint.printf( 0, "Building conservative distribution maps\n" );
1376 
1377  // generic triangle used for area computation, for triangles around the center of overlap face;
1378  // used for overlap faces with more than 4 edges;
1379  // nodes array will be set for each triangle;
1380  // these triangles are not part of the mesh structure, they are just temporary during
1381  // aforementioned decomposition.
1382  Face faceTri( 3 );
1383  NodeVector nodes( 3 );
1384  faceTri.SetNode( 0, 0 );
1385  faceTri.SetNode( 1, 1 );
1386  faceTri.SetNode( 2, 2 );
1387 
1388  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1389  {
1390 #ifdef VERBOSE
1391  // Announce computation progress
1392  if( ixFirst % outputFrequency == 0 && is_root )
1393  {
1394  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1395  }
1396 #endif
1397  // Quantities from the First Mesh
1398  const Face& faceFirst = m_meshInputCov->faces[ixFirst];
1399 
1400  const NodeVector& nodesFirst = m_meshInputCov->nodes;
1401 
1402  // Number of overlapping Faces and triangles
1403  int nOverlapFaces = nAllOverlapFaces[ixFirst];
1404 
1405  if( !nOverlapFaces ) continue;
1406 
1407  // // Calculate total element Jacobian
1408  // double dTotalJacobian = 0.0;
1409  // for (int s = 0; s < nPin; s++) {
1410  // for (int t = 0; t < nPin; t++) {
1411  // dTotalJacobian += dataGLLJacobianIn[s][t][ixFirst];
1412  // }
1413  // }
1414 
1415  // Loop through all Overlap Faces
1416  for( int i = 0; i < nOverlapFaces; i++ )
1417  {
1418  // Quantities from the overlap Mesh
1419  const Face& faceOverlap = m_meshOverlap->faces[ixOverlap + i];
1420 
1421  const NodeVector& nodesOverlap = m_meshOverlap->nodes;
1422 
1423  // Quantities from the Second Mesh
1424  int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1425 
1426  // signal to not participate, because it is a ghost target
1427  if( ixSecond < 0 ) continue; // do not do anything
1428 
1429  const NodeVector& nodesSecond = m_meshOutput->nodes;
1430 
1431  const Face& faceSecond = m_meshOutput->faces[ixSecond];
1432 
1433  int nbEdges = faceOverlap.edges.size();
1434  int nOverlapTriangles = 1;
1435  Node center; // not used if nbEdges == 3
1436  if( nbEdges > 3 )
1437  { // decompose from center in this case
1438  nOverlapTriangles = nbEdges;
1439  for( int k = 0; k < nbEdges; k++ )
1440  {
1441  const Node& node = nodesOverlap[faceOverlap[k]];
1442  center = center + node;
1443  }
1444  center = center / nbEdges;
1445  center = center.Normalized(); // project back on sphere of radius 1
1446  }
1447 
1448  Node node0, node1, node2;
1449  double dTriArea;
1450 
1451  // Loop over all sub-triangles of this Overlap Face
1452  for( int j = 0; j < nOverlapTriangles; j++ )
1453  {
1454  if( nbEdges == 3 ) // will come here only once, nOverlapTriangles == 1 in this case
1455  {
1456  node0 = nodesOverlap[faceOverlap[0]];
1457  node1 = nodesOverlap[faceOverlap[1]];
1458  node2 = nodesOverlap[faceOverlap[2]];
1459  dTriArea = CalculateFaceArea( faceOverlap, nodesOverlap );
1460  }
1461  else // decompose polygon in triangles around the center
1462  {
1463  node0 = center;
1464  node1 = nodesOverlap[faceOverlap[j]];
1465  int j1 = ( j + 1 ) % nbEdges;
1466  node2 = nodesOverlap[faceOverlap[j1]];
1467  nodes[0] = center;
1468  nodes[1] = node1;
1469  nodes[2] = node2;
1470  dTriArea = CalculateFaceArea( faceTri, nodes );
1471  }
1472 
1473  for( int k = 0; k < triquadrule.GetPoints(); k++ )
1474  {
1475  // Get the nodal location of this point
1476  double dX[3];
1477 
1478  dX[0] = dG( k, 0 ) * node0.x + dG( k, 1 ) * node1.x + dG( k, 2 ) * node2.x;
1479  dX[1] = dG( k, 0 ) * node0.y + dG( k, 1 ) * node1.y + dG( k, 2 ) * node2.y;
1480  dX[2] = dG( k, 0 ) * node0.z + dG( k, 1 ) * node1.z + dG( k, 2 ) * node2.z;
1481 
1482  double dMag = sqrt( dX[0] * dX[0] + dX[1] * dX[1] + dX[2] * dX[2] );
1483 
1484  dX[0] /= dMag;
1485  dX[1] /= dMag;
1486  dX[2] /= dMag;
1487 
1488  Node nodeQuadrature( dX[0], dX[1], dX[2] );
1489 
1490  // Find the components of this quadrature point in the basis
1491  // of the first Face.
1492  double dAlphaIn;
1493  double dBetaIn;
1494 
1495  ApplyInverseMap( faceFirst, nodesFirst, nodeQuadrature, dAlphaIn, dBetaIn );
1496 
1497  // Find the components of this quadrature point in the basis
1498  // of the second Face.
1499  double dAlphaOut;
1500  double dBetaOut;
1501 
1502  ApplyInverseMap( faceSecond, nodesSecond, nodeQuadrature, dAlphaOut, dBetaOut );
1503 
1504  /*
1505  // Check inverse map value
1506  if ((dAlphaIn < 0.0) || (dAlphaIn > 1.0) ||
1507  (dBetaIn < 0.0) || (dBetaIn > 1.0)
1508  ) {
1509  _EXCEPTION2("Inverse Map out of range (%1.5e %1.5e)",
1510  dAlphaIn, dBetaIn);
1511  }
1512 
1513  // Check inverse map value
1514  if ((dAlphaOut < 0.0) || (dAlphaOut > 1.0) ||
1515  (dBetaOut < 0.0) || (dBetaOut > 1.0)
1516  ) {
1517  _EXCEPTION2("Inverse Map out of range (%1.5e %1.5e)",
1518  dAlphaOut, dBetaOut);
1519  }
1520  */
1521  // Sample the First finite element at this point
1522  SampleGLLFiniteElement( nMonotoneType, nPin, dAlphaIn, dBetaIn, dSampleCoeffIn );
1523 
1524  // Sample the Second finite element at this point
1525  SampleGLLFiniteElement( nMonotoneType, nPout, dAlphaOut, dBetaOut, dSampleCoeffOut );
1526 
1527  // Overlap output area
1528  for( int s = 0; s < nPout; s++ )
1529  {
1530  for( int t = 0; t < nPout; t++ )
1531  {
1532  double dNodeArea = dSampleCoeffOut[s][t] * dW[k] * dTriArea;
1533 
1534  dOverlapOutputArea[ixOverlap + i][s * nPout + t] += dNodeArea;
1535 
1536  dGeometricOutputArea[ixSecond][s * nPout + t] += dNodeArea;
1537  }
1538  }
1539 
1540  // Compute overlap integral
1541  int ixp = 0;
1542  for( int p = 0; p < nPin; p++ )
1543  {
1544  for( int q = 0; q < nPin; q++ )
1545  {
1546  int ixs = 0;
1547  for( int s = 0; s < nPout; s++ )
1548  {
1549  for( int t = 0; t < nPout; t++ )
1550  {
1551  // Sample the Second finite element at this point
1552  dGlobalIntArray[ixp][ixOverlap + i][ixs] +=
1553  dSampleCoeffOut[s][t] * dSampleCoeffIn[p][q] * dW[k] * dTriArea;
1554 
1555  ixs++;
1556  }
1557  }
1558 
1559  ixp++;
1560  }
1561  }
1562  }
1563  }
1564  }
1565 
1566  // Coefficients
1567  DataArray2D< double > dCoeff( nOverlapFaces * nPout * nPout, nPin * nPin );
1568 
1569  for( int i = 0; i < nOverlapFaces; i++ )
1570  {
1571  // int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1572 
1573  int ixp = 0;
1574  for( int p = 0; p < nPin; p++ )
1575  {
1576  for( int q = 0; q < nPin; q++ )
1577  {
1578  int ixs = 0;
1579  for( int s = 0; s < nPout; s++ )
1580  {
1581  for( int t = 0; t < nPout; t++ )
1582  {
1583  dCoeff[i * nPout * nPout + ixs][ixp] = dGlobalIntArray[ixp][ixOverlap + i][ixs] /
1584  dOverlapOutputArea[ixOverlap + i][s * nPout + t];
1585 
1586  ixs++;
1587  }
1588  }
1589 
1590  ixp++;
1591  }
1592  }
1593  }
1594 
1595  // Source areas
1596  DataArray1D< double > vecSourceArea( nPin * nPin );
1597 
1598  for( int p = 0; p < nPin; p++ )
1599  {
1600  for( int q = 0; q < nPin; q++ )
1601  {
1602  vecSourceArea[p * nPin + q] = dataGLLJacobianIn[p][q][ixFirst];
1603  }
1604  }
1605 
1606  // Target areas
1607  DataArray1D< double > vecTargetArea( nOverlapFaces * nPout * nPout );
1608 
1609  for( int i = 0; i < nOverlapFaces; i++ )
1610  {
1611  // int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1612  int ixs = 0;
1613  for( int s = 0; s < nPout; s++ )
1614  {
1615  for( int t = 0; t < nPout; t++ )
1616  {
1617  vecTargetArea[i * nPout * nPout + ixs] = dOverlapOutputArea[ixOverlap + i][nPout * s + t];
1618 
1619  ixs++;
1620  }
1621  }
1622  }
1623 
1624  // Force consistency and conservation
1625  if( !fNoConservation )
1626  {
1627  ForceIntArrayConsistencyConservation( vecSourceArea, vecTargetArea, dCoeff, ( nMonotoneType != 0 ) );
1628  }
1629 
1630  // Update global coefficients
1631  for( int i = 0; i < nOverlapFaces; i++ )
1632  {
1633  int ixp = 0;
1634  for( int p = 0; p < nPin; p++ )
1635  {
1636  for( int q = 0; q < nPin; q++ )
1637  {
1638  int ixs = 0;
1639  for( int s = 0; s < nPout; s++ )
1640  {
1641  for( int t = 0; t < nPout; t++ )
1642  {
1643  dGlobalIntArray[ixp][ixOverlap + i][ixs] =
1644  dCoeff[i * nPout * nPout + ixs][ixp] * dOverlapOutputArea[ixOverlap + i][s * nPout + t];
1645 
1646  ixs++;
1647  }
1648  }
1649 
1650  ixp++;
1651  }
1652  }
1653  }
1654 
1655 #ifdef VVERBOSE
1656  // Check column sums (conservation)
1657  for( int i = 0; i < nPin * nPin; i++ )
1658  {
1659  double dColSum = 0.0;
1660  for( int j = 0; j < nOverlapFaces * nPout * nPout; j++ )
1661  {
1662  dColSum += dCoeff[j][i] * vecTargetArea[j];
1663  }
1664  printf( "Col %i: %1.15e\n", i, dColSum / vecSourceArea[i] );
1665  }
1666 
1667  // Check row sums (consistency)
1668  for( int j = 0; j < nOverlapFaces * nPout * nPout; j++ )
1669  {
1670  double dRowSum = 0.0;
1671  for( int i = 0; i < nPin * nPin; i++ )
1672  {
1673  dRowSum += dCoeff[j][i];
1674  }
1675  printf( "Row %i: %1.15e\n", j, dRowSum );
1676  }
1677 #endif
1678 
1679  // Increment the current overlap index
1680  ixOverlap += nOverlapFaces;
1681  }
1682 
1683  // Build redistribution map within target element
1684  if( is_root ) dbgprint.printf( 0, "Building redistribution maps on target mesh\n" );
1685  DataArray1D< double > dRedistSourceArea( nPout * nPout );
1686  DataArray1D< double > dRedistTargetArea( nPout * nPout );
1687  std::vector< DataArray2D< double > > dRedistributionMaps;
1688  dRedistributionMaps.resize( m_meshOutput->faces.size() );
1689 
1690  for( size_t ixSecond = 0; ixSecond < m_meshOutput->faces.size(); ixSecond++ )
1691  {
1692  dRedistributionMaps[ixSecond].Allocate( nPout * nPout, nPout * nPout );
1693 
1694  for( int i = 0; i < nPout * nPout; i++ )
1695  {
1696  dRedistributionMaps[ixSecond][i][i] = 1.0;
1697  }
1698 
1699  for( int s = 0; s < nPout * nPout; s++ )
1700  {
1701  dRedistSourceArea[s] = dGeometricOutputArea[ixSecond][s];
1702  }
1703 
1704  for( int s = 0; s < nPout * nPout; s++ )
1705  {
1706  dRedistTargetArea[s] = dataGLLJacobianOut[s / nPout][s % nPout][ixSecond];
1707  }
1708 
1709  if( !fNoConservation )
1710  {
1711  ForceIntArrayConsistencyConservation( dRedistSourceArea, dRedistTargetArea, dRedistributionMaps[ixSecond],
1712  ( nMonotoneType != 0 ) );
1713 
1714  for( int s = 0; s < nPout * nPout; s++ )
1715  {
1716  for( int t = 0; t < nPout * nPout; t++ )
1717  {
1718  dRedistributionMaps[ixSecond][s][t] *= dRedistTargetArea[s] / dRedistSourceArea[t];
1719  }
1720  }
1721  }
1722  }
1723 
1724  // Construct the total geometric area
1725  DataArray1D< double > dTotalGeometricArea( dataNodalAreaOut.GetRows() );
1726  for( size_t ixSecond = 0; ixSecond < m_meshOutput->faces.size(); ixSecond++ )
1727  {
1728  for( int s = 0; s < nPout; s++ )
1729  {
1730  for( int t = 0; t < nPout; t++ )
1731  {
1732  dTotalGeometricArea[dataGLLNodesOut[s][t][ixSecond] - 1] +=
1733  dGeometricOutputArea[ixSecond][s * nPout + t];
1734  }
1735  }
1736  }
1737 
1738  // Compose the integration operator with the output map
1739  ixOverlap = 0;
1740 
1741  if( is_root ) dbgprint.printf( 0, "Assembling map\n" );
1742 
1743  // Map from source DOFs to target DOFs with redistribution applied
1744  DataArray2D< double > dRedistributedOp( nPin * nPin, nPout * nPout );
1745 
1746  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1747  {
1748 #ifdef VERBOSE
1749  // Announce computation progress
1750  if( ixFirst % outputFrequency == 0 && is_root )
1751  {
1752  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1753  }
1754 #endif
1755  // Number of overlapping Faces and triangles
1756  int nOverlapFaces = nAllOverlapFaces[ixFirst];
1757 
1758  if( !nOverlapFaces ) continue;
1759 
1760  // Put composed array into map
1761  for( int j = 0; j < nOverlapFaces; j++ )
1762  {
1763  int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
1764 
1765  // signal to not participate, because it is a ghost target
1766  if( ixSecondFace < 0 ) continue; // do not do anything
1767 
1768  dRedistributedOp.Zero();
1769  for( int p = 0; p < nPin * nPin; p++ )
1770  {
1771  for( int s = 0; s < nPout * nPout; s++ )
1772  {
1773  for( int t = 0; t < nPout * nPout; t++ )
1774  {
1775  dRedistributedOp[p][s] +=
1776  dRedistributionMaps[ixSecondFace][s][t] * dGlobalIntArray[p][ixOverlap + j][t];
1777  }
1778  }
1779  }
1780 
1781  int ixp = 0;
1782  for( int p = 0; p < nPin; p++ )
1783  {
1784  for( int q = 0; q < nPin; q++ )
1785  {
1786  int ixFirstNode;
1787  if( fContinuousIn )
1788  {
1789  ixFirstNode = dataGLLNodesIn[p][q][ixFirst] - 1;
1790  }
1791  else
1792  {
1793  ixFirstNode = ixFirst * nPin * nPin + p * nPin + q;
1794  }
1795 
1796  int ixs = 0;
1797  for( int s = 0; s < nPout; s++ )
1798  {
1799  for( int t = 0; t < nPout; t++ )
1800  {
1801  int ixSecondNode;
1802  if( fContinuousOut )
1803  {
1804  ixSecondNode = dataGLLNodesOut[s][t][ixSecondFace] - 1;
1805 
1806  if( !fNoConservation )
1807  {
1808  smatMap( ixSecondNode, ixFirstNode ) +=
1809  dRedistributedOp[ixp][ixs] / dataNodalAreaOut[ixSecondNode];
1810  }
1811  else
1812  {
1813  smatMap( ixSecondNode, ixFirstNode ) +=
1814  dRedistributedOp[ixp][ixs] / dTotalGeometricArea[ixSecondNode];
1815  }
1816  }
1817  else
1818  {
1819  ixSecondNode = ixSecondFace * nPout * nPout + s * nPout + t;
1820 
1821  if( !fNoConservation )
1822  {
1823  smatMap( ixSecondNode, ixFirstNode ) +=
1824  dRedistributedOp[ixp][ixs] / dataGLLJacobianOut[s][t][ixSecondFace];
1825  }
1826  else
1827  {
1828  smatMap( ixSecondNode, ixFirstNode ) +=
1829  dRedistributedOp[ixp][ixs] / dGeometricOutputArea[ixSecondFace][s * nPout + t];
1830  }
1831  }
1832 
1833  ixs++;
1834  }
1835  }
1836 
1837  ixp++;
1838  }
1839  }
1840  }
1841 
1842  // Increment the current overlap index
1843  ixOverlap += nOverlapFaces;
1844  }
1845 
1846  return;
1847 }
1848 
1849 ///////////////////////////////////////////////////////////////////////////////
1850 
1851 void moab::TempestOnlineMap::LinearRemapGLLtoGLL2_Pointwise_MOAB( const DataArray3D< int >& dataGLLNodesIn,
1852  const DataArray3D< double >& /*dataGLLJacobianIn*/,
1853  const DataArray3D< int >& dataGLLNodesOut,
1854  const DataArray3D< double >& /*dataGLLJacobianOut*/,
1855  const DataArray1D< double >& dataNodalAreaOut,
1856  int nPin,
1857  int nPout,
1858  int nMonotoneType,
1859  bool fContinuousIn,
1860  bool fContinuousOut )
1861 {
1862  // Gauss-Lobatto quadrature within Faces
1863  DataArray1D< double > dGL;
1864  DataArray1D< double > dWL;
1865 
1866  GaussLobattoQuadrature::GetPoints( nPout, 0.0, 1.0, dGL, dWL );
1867 
1868  // Get SparseMatrix represntation of the OfflineMap
1869  SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1870 
1871  // Sample coefficients
1872  DataArray2D< double > dSampleCoeffIn( nPin, nPin );
1873 
1874  // Announcemnets
1875  moab::DebugOutput dbgprint( std::cout, this->rank, 0 );
1876  dbgprint.set_prefix( "[LinearRemapGLLtoGLL2_Pointwise_MOAB]: " );
1877  if( is_root )
1878  {
1879  dbgprint.printf( 0, "Finite Element to Finite Element (Pointwise) Projection\n" );
1880  dbgprint.printf( 0, "Order of the input FE polynomial interpolant: %i\n", nPin );
1881  dbgprint.printf( 0, "Order of the output FE polynomial interpolant: %i\n", nPout );
1882  }
1883 
1884  // Number of overlap Faces per source Face
1885  DataArray1D< int > nAllOverlapFaces( m_meshInputCov->faces.size() );
1886 
1887  int ixOverlap = 0;
1888 
1889  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1890  {
1891  size_t ixOverlapTemp = ixOverlap;
1892  for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1893  {
1894  // const Face & faceOverlap = m_meshOverlap->faces[ixOverlapTemp];
1895 
1896  if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 ) break;
1897 
1898  nAllOverlapFaces[ixFirst]++;
1899  }
1900 
1901  // Increment the current overlap index
1902  ixOverlap += nAllOverlapFaces[ixFirst];
1903  }
1904 
1905  // Number of times this point was found
1906  DataArray1D< bool > fSecondNodeFound( dataNodalAreaOut.GetRows() );
1907 
1908  ixOverlap = 0;
1909 #ifdef VERBOSE
1910  const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
1911 #endif
1912  // Loop through all faces on m_meshInputCov
1913  for( size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1914  {
1915 #ifdef VERBOSE
1916  // Announce computation progress
1917  if( ixFirst % outputFrequency == 0 && is_root )
1918  {
1919  dbgprint.printf( 0, "Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1920  }
1921 #endif
1922  // Quantities from the First Mesh
1923  const Face& faceFirst = m_meshInputCov->faces[ixFirst];
1924 
1925  const NodeVector& nodesFirst = m_meshInputCov->nodes;
1926 
1927  // Number of overlapping Faces and triangles
1928  int nOverlapFaces = nAllOverlapFaces[ixFirst];
1929 
1930  // Loop through all Overlap Faces
1931  for( int i = 0; i < nOverlapFaces; i++ )
1932  {
1933  // Quantities from the Second Mesh
1934  int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1935 
1936  // signal to not participate, because it is a ghost target
1937  if( ixSecond < 0 ) continue; // do not do anything
1938 
1939  const NodeVector& nodesSecond = m_meshOutput->nodes;
1940  const Face& faceSecond = m_meshOutput->faces[ixSecond];
1941 
1942  // Loop through all nodes on the second face
1943  for( int s = 0; s < nPout; s++ )
1944  {
1945  for( int t = 0; t < nPout; t++ )
1946  {
1947  size_t ixSecondNode;
1948  if( fContinuousOut )
1949  {
1950  ixSecondNode = dataGLLNodesOut[s][t][ixSecond] - 1;
1951  }
1952  else
1953  {
1954  ixSecondNode = ixSecond * nPout * nPout + s * nPout + t;
1955  }
1956 
1957  if( ixSecondNode >= fSecondNodeFound.GetRows() ) _EXCEPTIONT( "Logic error" );
1958 
1959  // Check if this node has been found already
1960  if( fSecondNodeFound[ixSecondNode] ) continue;
1961 
1962  // Check this node
1963  Node node;
1964  Node dDx1G;
1965  Node dDx2G;
1966 
1967  ApplyLocalMap( faceSecond, nodesSecond, dGL[t], dGL[s], node, dDx1G, dDx2G );
1968 
1969  // Find the components of this quadrature point in the basis
1970  // of the first Face.
1971  double dAlphaIn;
1972  double dBetaIn;
1973 
1974  ApplyInverseMap( faceFirst, nodesFirst, node, dAlphaIn, dBetaIn );
1975 
1976  // Check if this node is within the first Face
1977  if( ( dAlphaIn < -1.0e-10 ) || ( dAlphaIn > 1.0 + 1.0e-10 ) || ( dBetaIn < -1.0e-10 ) ||
1978  ( dBetaIn > 1.0 + 1.0e-10 ) )
1979  continue;
1980 
1981  // Node is within the overlap region, mark as found
1982  fSecondNodeFound[ixSecondNode] = true;
1983 
1984  // Sample the First finite element at this point
1985  SampleGLLFiniteElement( nMonotoneType, nPin, dAlphaIn, dBetaIn, dSampleCoeffIn );
1986 
1987  // Add to map
1988  for( int p = 0; p < nPin; p++ )
1989  {
1990  for( int q = 0; q < nPin; q++ )
1991  {
1992  int ixFirstNode;
1993  if( fContinuousIn )
1994  {
1995  ixFirstNode = dataGLLNodesIn[p][q][ixFirst] - 1;
1996  }
1997  else
1998  {
1999  ixFirstNode = ixFirst * nPin * nPin + p * nPin + q;
2000  }
2001 
2002  smatMap( ixSecondNode, ixFirstNode ) += dSampleCoeffIn[p][q];
2003  }
2004  }
2005  }
2006  }
2007  }
2008 
2009  // Increment the current overlap index
2010  ixOverlap += nOverlapFaces;
2011  }
2012 
2013  // Check for missing samples
2014  for( size_t i = 0; i < fSecondNodeFound.GetRows(); i++ )
2015  {
2016  if( !fSecondNodeFound[i] )
2017  {
2018  _EXCEPTION1( "Can't sample point %i", i );
2019  }
2020  }
2021 
2022  return;
2023 }
2024 
2025 ///////////////////////////////////////////////////////////////////////////////