9 #define _USE_MATH_DEFINES
12 #pragma GCC diagnostic push
13 #pragma GCC diagnostic ignored "-Wdeprecated"
14 #pragma GCC diagnostic ignored "-Wsign-compare"
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"
33 #pragma GCC diagnostic pop
41 #include <unordered_set>
44 #define USE_ComputeAdjacencyRelations
64 std::ofstream output_file(
"rowcolindices.txt", std::ios::out );
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";
78 if( use_GID_matching )
80 std::map< unsigned, unsigned > src_gl;
81 for(
unsigned it = 0; it <
col_gdofmap.size(); ++it )
84 std::map< unsigned, unsigned >::iterator iter;
85 for(
unsigned it = 0; it <
row_gdofmap.size(); ++it )
88 iter = src_gl.find( row );
89 if( strict_check && iter == src_gl.end() )
91 std::cout <<
"Searching for global target DOF " << row
92 <<
" but could not find correspondence in source mesh.\n";
95 else if( iter == src_gl.end() )
101 unsigned icol = src_gl[row];
105 m_mapRemap( irow, icol ) = 1.0;
115 return moab::MB_FAILURE;
124 const int TriQuadRuleOrder = 4;
127 if( m_meshInputCov->faces.size() > 0 && m_meshInputCov->revnodearray.size() == 0 )
129 _EXCEPTIONT(
"ReverseNodeArray has not been calculated for m_meshInputCov" );
133 TriangularQuadratureRule triquadrule( TriQuadRuleOrder );
136 #ifdef RECTANGULAR_TRUNCATION
137 int nCoefficients = nOrder * nOrder;
139 #ifdef TRIANGULAR_TRUNCATION
140 int nCoefficients = nOrder * ( nOrder + 1 ) / 2;
144 const int nRequiredFaceSetSize = nCoefficients;
147 const int nFitWeightsExponent = nOrder + 2;
151 dbgprint.set_prefix(
"[LinearRemapFVtoFV_Tempest_MOAB]: " );
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 );
164 const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
166 DataArray2D< double > dIntArray;
167 DataArray1D< double > dConstraint( nCoefficients );
170 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
174 if( ixFirst % outputFrequency == 0 && is_root )
176 dbgprint.printf( 0,
"Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
180 int ixOverlapBegin = ixOverlap;
181 unsigned ixOverlapEnd = ixOverlapBegin;
183 for( ; ixOverlapEnd < m_meshOverlap->faces.size(); ixOverlapEnd++ )
185 if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapEnd] != 0 )
break;
188 unsigned nOverlapFaces = ixOverlapEnd - ixOverlapBegin;
190 if( nOverlapFaces == 0 )
continue;
193 BuildIntegrationArray( *m_meshInputCov, *m_meshOverlap, triquadrule, ixFirst, ixOverlapBegin, ixOverlapEnd,
200 GetAdjacentFaceVectorByEdge( *m_meshInputCov, ixFirst, nRequiredFaceSetSize, vecAdjFaces );
203 int nAdjFaces = vecAdjFaces.size();
206 double dFirstArea = m_meshInputCov->vecFaceArea[ixFirst];
208 for(
int p = 0; p < nCoefficients; p++ )
210 for(
unsigned j = 0; j < nOverlapFaces; j++ )
212 dConstraint[p] += dIntArray[p][j];
214 dConstraint[p] /= dFirstArea;
218 DataArray2D< double > dFitArray;
219 DataArray1D< double > dFitWeights;
220 DataArray2D< double > dFitArrayPlus;
222 BuildFitArray( *m_meshInputCov, triquadrule, ixFirst, vecAdjFaces, nOrder, nFitWeightsExponent, dConstraint,
223 dFitArray, dFitWeights );
226 bool fSuccess = InvertFitArray_Corrected( dConstraint, dFitArray, dFitWeights, dFitArrayPlus );
229 DataArray2D< double > dComposedArray( nAdjFaces, nOverlapFaces );
233 for(
int i = 0; i < nAdjFaces; i++ )
235 for(
size_t j = 0; j < nOverlapFaces; j++ )
237 for(
int k = 0; k < nCoefficients; k++ )
239 dComposedArray( i, j ) += dIntArray( k, j ) * dFitArrayPlus( i, k );
249 dComposedArray.Zero();
250 for(
size_t j = 0; j < nOverlapFaces; j++ )
252 dComposedArray( 0, j ) += dIntArray( 0, j );
257 for(
unsigned i = 0; i < vecAdjFaces.size(); i++ )
259 for(
unsigned j = 0; j < nOverlapFaces; j++ )
261 int& ixFirstFaceLoc = vecAdjFaces[i].first;
262 int& ixSecondFaceLoc = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
268 if( ixSecondFaceLoc < 0 )
continue;
270 m_mapRemap( ixSecondFaceLoc, ixFirstFaceLoc ) +=
271 dComposedArray[i][j] / m_meshOutput->vecFaceArea[ixSecondFaceLoc];
276 ixOverlap += nOverlapFaces;
286 int nrows = m_weightMatrix.rows();
287 int ncols = m_weightMatrix.cols();
288 int NNZ = m_weightMatrix.nonZeros();
292 int arr3[6] = { NNZ, nrows, ncols, -NNZ, -nrows, -ncols };
294 MPI_Reduce( arr3, rarr3, 6, MPI_INT, MPI_MIN, 0, m_pcomm->comm() );
297 MPI_Reduce( arr3, total, 3, MPI_INT, MPI_SUM, 0, m_pcomm->comm() );
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";
303 std::cout <<
"-> Rows: " << nrows <<
", Cols: " << ncols <<
", NNZ: " << NNZ <<
"\n";
307 #ifdef MOAB_HAVE_EIGEN3
308 void moab::TempestOnlineMap::copy_tempest_sparsemat_to_eigen3()
311 #define VERBOSE_ACTIVATED
316 m_weightMatrix.resize( m_nTotDofs_Dest, m_nTotDofs_SrcCov );
317 m_rowVector.resize( m_weightMatrix.rows() );
318 m_colVector.resize( m_weightMatrix.cols() );
321 int locrows = std::max( m_mapRemap.GetRows(), m_nTotDofs_Dest );
322 int loccols = std::max( m_mapRemap.GetColumns(), m_nTotDofs_SrcCov );
324 std::cout << m_weightMatrix.rows() <<
", " << locrows <<
", " << m_weightMatrix.cols() <<
", " << loccols <<
"\n";
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();
335 typedef Eigen::Triplet< double >
Triplet;
336 std::vector< Triplet > tripletList;
337 tripletList.reserve( locvals );
338 for(
size_t iv = 0; iv < locvals; iv++ )
340 tripletList.push_back(
Triplet( lrows[iv], lcols[iv], lvals[iv] ) );
342 m_weightMatrix.setFromTriplets( tripletList.begin(), tripletList.end() );
343 m_weightMatrix.makeCompressed();
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++ )
355 output_file << row_gdofmap[row_ldofmap[lrows[iv]]] <<
" " << col_gdofmap[col_ldofmap[lcols[iv]]] <<
" "
356 << lvals[iv] <<
"\n";
362 #ifdef VERBOSE_ACTIVATED
363 #undef VERBOSE_ACTIVATED
372 template <
typename T >
376 std::vector< size_t > idx( v.size() );
377 std::iota( idx.begin(), idx.end(), 0 );
383 std::stable_sort( idx.begin(), idx.end(), [&v](
size_t i1,
size_t i2 ) { return fabs( v[i1] ) > fabs( v[i2] ); } );
389 std::vector< double >& dataCorrectedField,
390 std::vector< double >& dataLowerBound,
391 std::vector< double >& dataUpperBound,
392 std::vector< double >& dMassDefect )
394 const size_t nrows = dataCorrectedField.size();
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;
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 );
410 ComputeAdjacencyRelations( vecAdjTargetFaces, caasIteration, m_remapper->m_target_entities, useMOABAdjacencies,
411 this->m_remapper->m_target );
417 for(
size_t i = 0; i < nrows; i++ )
421 dataCorrection[
index] = fmax( dataLowerBound[
index], fmin( dataUpperBound[
index], 0.0 ) );
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];
430 #ifndef USE_ComputeAdjacencyRelations
434 if( useMOABAdjacencies )
438 ents.
insert( m_remapper->m_target_entities[
index] );
444 int adjIndex = m_remapper->m_target_entities.index( *it );
446 if( adjIndex >= 0 ) vecAdjTargetFaces[
index].insert( adjIndex );
452 GetAdjacentFaceVectorByEdge( *this->m_remapper->m_target,
index,
453 ( m_output_order + 1 ) * ( m_output_order + 1 ) * ( m_output_order + 1 ),
459 for(
auto adjFace : vecAdjFaces )
460 if( adjFace.first >= 0 )
461 vecAdjTargetFaces[
index].insert( adjFace.first );
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;
475 MPI_Allreduce( localDefects.data(), globalDefects.data(), 4, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
477 dMassL = globalDefects[0];
478 dMassU = globalDefects[1];
479 dMassDiffCum = globalDefects[2];
480 dLMinusU = globalDefects[3];
485 if( fabs( dMassDiffCum ) < 1e-15 || dLMinusU < 1e-15 )
487 for(
size_t i = 0; i < nrows; i++ )
488 dataCorrectedField[i] += dataCorrection[i];
493 if( dMassL > dMassDiffCum )
495 Announce(
"Lower bound mass exceeds target mass by %1.15e: CAAS will need another iteration",
496 dMassL - dMassDiffCum );
497 dMassDiffCum = dMassL;
500 else if( dMassU < dMassDiffCum )
502 Announce(
"Target mass exceeds upper bound mass by %1.15e: CAAS will need another iteration",
503 dMassDiffCum - dMassU );
504 dMassDiffCum = dMassU;
510 for(
size_t i = 0; i < nrows; i++ )
514 const std::unordered_set< int >& neighbors = vecAdjTargetFaces[
index];
515 if( dMassDefect[
index] > 0.0 )
517 double dMassCorrectU = 0.0;
518 for(
auto it : neighbors )
519 dMassCorrectU += dTargetAreas[it] * ( dataUpperBound[it] - dataCorrection[it] );
522 for(
auto it : neighbors )
523 dataCorrection[it] +=
524 dMassDefect[
index] * ( dataUpperBound[it] - dataCorrection[it] ) / dMassCorrectU;
528 double dMassCorrectL = 0.0;
529 for(
auto it : neighbors )
530 dMassCorrectL += dTargetAreas[it] * ( dataCorrection[it] - dataLowerBound[it] );
533 for(
auto it : neighbors )
534 dataCorrection[it] +=
535 dMassDefect[
index] * ( dataCorrection[it] - dataLowerBound[it] ) / dMassCorrectL;
539 for(
size_t i = 0; i < nrows; i++ )
540 dataCorrectedField[i] += dataCorrection[i];
547 std::vector< double >& dataLowerBound,
548 std::vector< double >& dataUpperBound,
551 const size_t nrows = dataCorrectedField.size();
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++ )
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] );
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;
579 MPI_Allreduce( localDefects.data(), globalDefects.data(), 5, MPI_DOUBLE, MPI_SUM, m_pcomm->comm() );
581 dMassL = globalDefects[0];
582 dMassU = globalDefects[1];
583 dMassDiff = globalDefects[2];
584 dMassCorrectL = globalDefects[3];
585 dMassCorrectU = globalDefects[4];
589 if( fabs( dMassDiff ) < 1e-15 || fabs( dLMinusU ) < 1e-15 )
591 for(
size_t i = 0; i < nrows; i++ )
592 dataCorrectedField[i] += dataCorrection[i];
597 if( dMassL > dMassDiff )
599 Announce(
"%d: Lower bound mass exceeds target mass by %1.15e: CAAS will need another iteration", rank,
600 dMassL - dMassDiff );
604 else if( dMassU < dMassDiff )
606 Announce(
"%d: Target mass exceeds upper bound mass by %1.15e: CAAS will need another iteration", rank,
607 dMassDiff - dMassU );
613 DataArray1D< double > dataMassVec( nrows );
614 if( dMassDiff > 0.0 )
616 for(
size_t i = 0; i < nrows; i++ )
618 dataMassVec[i] = ( dataUpperBound[i] - dataCorrection[i] ) / dMassCorrectU;
619 dataCorrection[i] += dMassDiff * dataMassVec[i];
624 for(
size_t i = 0; i < nrows; i++ )
626 dataMassVec[i] = ( dataCorrection[i] - dataLowerBound[i] ) / dMassCorrectL;
627 dataCorrection[i] += dMassDiff * dataMassVec[i];
631 for(
size_t i = 0; i < nrows; i++ )
632 dataCorrectedField[i] += dataCorrection[i];
639 std::vector< double >& dataOutDouble,
646 assert( !dataGLLNodesSrcCov.IsAttached() && !dataGLLNodesDest.IsAttached() );
648 std::pair< double, double > massDefect( 0.0, 0.0 );
651 const size_t nTargetCount = dataOutDouble.size();
652 const DataArray1D< double >& m_dOverlapAreas = this->m_remapper->m_overlap->vecFaceArea;
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 );
662 #undef USE_ComputeAdjacencyRelations
663 constexpr
bool useMOABAdjacencies =
true;
664 #ifdef USE_ComputeAdjacencyRelations
668 if( caasType == CAAS_QLT || caasType == CAAS_LOCAL_ADJACENT )
670 if( useMOABAdjacencies )
673 ComputeAdjacencyRelations( vecSourceOvTarget, caasIteration, m_remapper->m_covering_source_entities,
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" );
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++ )
694 const int ixS = m_meshOverlap->vecSourceFaceIx[i];
695 const int ixT = m_meshOverlap->vecTargetFaceIx[i];
697 if( ixT < 0 )
continue;
699 assert( m_dOverlapAreas[i] > 0.0 );
703 #ifndef USE_ComputeAdjacencyRelations
705 vecSourceOvTarget[ixT].insert( ixS );
706 if( ( caasType == CAAS_QLT || caasType == CAAS_LOCAL_ADJACENT ) )
708 if( useMOABAdjacencies )
711 ents.
insert( m_remapper->m_covering_source_entities[ixS] );
716 int adjIndex = m_remapper->m_covering_source_entities.index( *it );
717 if( adjIndex >= 0 ) vecSourceOvTarget[ixT].insert( adjIndex );
724 GetAdjacentFaceVectorByEdge( *m_meshInputCov, ixS,
725 ( caasIteration ) * ( m_input_order + 1 ) * ( m_input_order + 1 ),
729 for(
size_t iadj = 0; iadj < vecAdjFaces.size(); iadj++ )
730 vecSourceOvTarget[ixT].insert( vecAdjFaces[iadj].
first );
736 dSourceMax = fmax( dSourceMax, dataInDouble[ixS] );
737 dSourceMin = fmin( dSourceMin, dataInDouble[ixS] );
740 dTargetMin = fmin( dTargetMin, dataOutDouble[ixT] );
741 dTargetMax = fmax( dTargetMax, dataOutDouble[ixT] );
743 const double locMassDiff = ( dataInDouble[ixS] * m_dOverlapAreas[i] ) -
744 ( dataOutDouble[ixT] * m_dOverlapAreas[i] );
748 dMassDiff += locMassDiff;
749 massVector[ixT] += locMassDiff;
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;
760 if( caasType == CAAS_GLOBAL )
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,
765 dSourceMin = globalMinMaxDefects[0];
766 dSourceMax = globalMinMaxDefects[2];
767 dTargetMin = globalMinMaxDefects[1];
768 dTargetMax = globalMinMaxDefects[3];
770 if( caasIteration == 1 )
771 MPI_Allreduce( localMinMaxDefects.data() + 4, globalMinMaxDefects.data() + 4, 1, MPI_DOUBLE, MPI_SUM,
774 globalMinMaxDefects[4] = mismatch;
776 dMassDiff = localMinMaxDefects[4];
778 massDefect.first = globalMinMaxDefects[4];
782 massDefect.first = dMassDiff;
787 if( fabs( massDefect.first ) > 1e-20 )
789 if( caasType == CAAS_GLOBAL )
791 for(
size_t i = 0; i < nTargetCount; i++ )
793 dataLowerBound[i] = dSourceMin - dataOutDouble[i];
794 dataUpperBound[i] = dSourceMax - dataOutDouble[i];
800 std::vector< double > vecLocalUpperBound( nTargetCount );
801 std::vector< double > vecLocalLowerBound( nTargetCount );
804 for(
size_t i = 0; i < nTargetCount; i++ )
806 assert( vecSourceOvTarget[i].size() );
809 double dMaxI = -1E10;
812 for(
const auto& srcElem : vecSourceOvTarget[i] )
814 dMinI = fmin( dMinI, dataInDouble[srcElem] );
815 dMaxI = fmax( dMaxI, dataInDouble[srcElem] );
819 vecLocalLowerBound[i] = dMinI;
820 vecLocalUpperBound[i] = dMaxI;
823 for(
size_t i = 0; i < nTargetCount; i++ )
825 dataLowerBound[i] = vecLocalLowerBound[i] - dataOutDouble[i];
826 dataUpperBound[i] = vecLocalUpperBound[i] - dataOutDouble[i];
831 if( fabs( dMassDiff ) > 1e-20 )
833 if( caasType == CAAS_QLT )
834 dMassDiff = QLTLimiter( caasIteration, dataOutDouble, dataLowerBound, dataUpperBound, massVector );
836 CAASLimiter( dataOutDouble, dataLowerBound, dataUpperBound, dMassDiff );
840 double dMassDiffPost = 0.0;
841 for(
size_t i = 0; i < m_meshOverlap->faces.size(); i++ )
843 const int ixS = m_meshOverlap->vecSourceFaceIx[i];
844 const int ixT = m_meshOverlap->vecTargetFaceIx[i];
846 if( ixT < 0 )
continue;
850 dMassDiffPost += ( dataInDouble[ixS] * m_dOverlapAreas[i] ) -
851 ( dataOutDouble[ixT] * m_dOverlapAreas[i] );
854 massDefect.second = dMassDiffPost;
865 const DataArray1D< double >& vecTargetArea,
866 DataArray2D< double >& dCoeff,
868 bool fSparseConstraints =
false );
873 const DataArray1D< double >& vecTargetArea,
874 DataArray2D< double >& dCoeff,
880 const DataArray3D< double >& dataGLLJacobian,
883 bool fNoConservation )
886 int nP = dataGLLNodes.GetRows();
889 const int TriQuadRuleOrder = 4;
892 TriangularQuadratureRule triquadrule( TriQuadRuleOrder );
894 int TriQuadraturePoints = triquadrule.GetPoints();
896 const DataArray2D< double >& TriQuadratureG = triquadrule.GetG();
898 const DataArray1D< double >& TriQuadratureW = triquadrule.GetW();
901 DataArray2D< double > dSampleCoeff( nP, nP );
904 DataArray1D< double > dG;
905 DataArray1D< double > dW;
906 GaussLobattoQuadrature::GetPoints( nP, 0.0, 1.0, dG, dW );
910 dbgprint.set_prefix(
"[LinearRemapSE4_Tempest_MOAB]: " );
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 );
919 SparseMatrix< double >& smatMap = this->GetSparseMatrix();
922 const NodeVector& nodesOverlap = m_meshOverlap->nodes;
923 const NodeVector& nodesFirst = m_meshInputCov->nodes;
926 DataArray1D< double > vecSourceArea( nP * nP );
928 DataArray1D< double > vecTargetArea;
929 DataArray2D< double > dCoeff;
932 std::stringstream sstr;
933 sstr <<
"remapdata_" << rank <<
".txt";
934 std::ofstream output_file( sstr.str() );
940 const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
948 NodeVector nodes( 3 );
949 faceTri.SetNode( 0, 0 );
950 faceTri.SetNode( 1, 1 );
951 faceTri.SetNode( 2, 2 );
954 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
956 const Face& faceFirst = m_meshInputCov->faces[ixFirst];
958 if( faceFirst.edges.size() != 4 )
960 _EXCEPTIONT(
"Only quadrilateral elements allowed for SE remapping" );
964 if( ixFirst % outputFrequency == 0 && is_root )
966 dbgprint.printf( 0,
"Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
974 int nOverlapFaces = 0;
975 size_t ixOverlapTemp = ixOverlap;
976 for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
980 if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 )
break;
986 if( nOverlapFaces == 0 )
continue;
989 DataArray3D< double > dRemapCoeff( nP, nP, nOverlapFaces );
992 for(
int j = 0; j < nOverlapFaces; j++ )
994 const Face& faceOverlap = m_meshOverlap->faces[ixOverlap + j];
995 if( m_meshOverlap->vecFaceArea[ixOverlap + j] < 1.e-16 )
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++ )
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 );
1016 int nbEdges = faceOverlap.edges.size();
1017 int nOverlapTriangles = 1;
1021 nOverlapTriangles = nbEdges;
1022 for(
int k = 0; k < nbEdges; k++ )
1024 const Node& node = nodesOverlap[faceOverlap[k]];
1031 Node node0, node1, node2;
1032 double dTriangleArea;
1035 for(
int k = 0; k < nOverlapTriangles; k++ )
1039 node0 = nodesOverlap[faceOverlap[0]];
1040 node1 = nodesOverlap[faceOverlap[1]];
1041 node2 = nodesOverlap[faceOverlap[2]];
1042 dTriangleArea = CalculateFaceArea( faceOverlap, nodesOverlap );
1047 node1 = nodesOverlap[faceOverlap[k]];
1048 int k1 = ( k + 1 ) % nbEdges;
1049 node2 = nodesOverlap[faceOverlap[k1]];
1053 dTriangleArea = CalculateFaceArea( faceTri, nodes );
1056 for(
int l = 0; l < TriQuadraturePoints; l++ )
1058 Node nodeQuadrature;
1059 nodeQuadrature.x = TriQuadratureG[l][0] * node0.x + TriQuadratureG[l][1] * node1.x +
1060 TriQuadratureG[l][2] * node2.x;
1062 nodeQuadrature.y = TriQuadratureG[l][0] * node0.y + TriQuadratureG[l][1] * node1.y +
1063 TriQuadratureG[l][2] * node2.y;
1065 nodeQuadrature.z = TriQuadratureG[l][0] * node0.z + TriQuadratureG[l][1] * node1.z +
1066 TriQuadratureG[l][2] * node2.z;
1068 nodeQuadrature = nodeQuadrature.Normalized();
1075 ApplyInverseMap( faceFirst, nodesFirst, nodeQuadrature, dAlpha, dBeta );
1078 if( ( dAlpha < -1.0e-13 ) || ( dAlpha > 1.0 + 1.0e-13 ) || ( dBeta < -1.0e-13 ) ||
1079 ( dBeta > 1.0 + 1.0e-13 ) )
1081 _EXCEPTION4(
"Inverse Map for element %d and subtriangle %d out of range "
1083 j, l, dAlpha, dBeta );
1087 SampleGLLFiniteElement( nMonotoneType, nP, dAlpha, dBeta, dSampleCoeff );
1090 for(
int p = 0; p < nP; p++ )
1092 for(
int q = 0; q < nP; q++ )
1094 dRemapCoeff[p][q][j] += TriQuadratureW[l] * dTriangleArea * dSampleCoeff[p][q] /
1095 m_meshOverlap->vecFaceArea[ixOverlap + j];
1103 output_file <<
"[" << m_remapper->lid_to_gid_covsrc[ixFirst] <<
"] \t";
1104 for(
int j = 0; j < nOverlapFaces; j++ )
1106 for(
int p = 0; p < nP; p++ )
1108 for(
int q = 0; q < nP; q++ )
1110 output_file << dRemapCoeff[p][q][j] <<
" ";
1114 output_file << std::endl;
1118 if( !fNoConservation )
1120 double dTargetArea = 0.0;
1121 for(
int j = 0; j < nOverlapFaces; j++ )
1123 dTargetArea += m_meshOverlap->vecFaceArea[ixOverlap + j];
1126 for(
int p = 0; p < nP; p++ )
1128 for(
int q = 0; q < nP; q++ )
1130 vecSourceArea[p * nP + q] = dataGLLJacobian[p][q][ixFirst];
1134 const double areaTolerance = 1e-10;
1136 if( fabs( m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea ) <= areaTolerance )
1138 vecTargetArea.Allocate( nOverlapFaces );
1139 for(
int j = 0; j < nOverlapFaces; j++ )
1141 vecTargetArea[j] = m_meshOverlap->vecFaceArea[ixOverlap + j];
1144 dCoeff.Allocate( nOverlapFaces, nP * nP );
1146 for(
int j = 0; j < nOverlapFaces; j++ )
1148 for(
int p = 0; p < nP; p++ )
1150 for(
int q = 0; q < nP; q++ )
1152 dCoeff[j][p * nP + q] = dRemapCoeff[p][q][j];
1159 else if( m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea > areaTolerance )
1161 double dExtraneousArea = m_meshInputCov->vecFaceArea[ixFirst] - dTargetArea;
1163 vecTargetArea.Allocate( nOverlapFaces + 1 );
1164 for(
int j = 0; j < nOverlapFaces; j++ )
1166 vecTargetArea[j] = m_meshOverlap->vecFaceArea[ixOverlap + j];
1168 vecTargetArea[nOverlapFaces] = dExtraneousArea;
1171 Announce(
"Partial volume: %i (%1.10e / %1.10e)", ixFirst, dTargetArea,
1172 m_meshInputCov->vecFaceArea[ixFirst] );
1174 if( dTargetArea > m_meshInputCov->vecFaceArea[ixFirst] )
1176 _EXCEPTIONT(
"Partial element area exceeds total element area" );
1179 dCoeff.Allocate( nOverlapFaces + 1, nP * nP );
1181 for(
int j = 0; j < nOverlapFaces; j++ )
1183 for(
int p = 0; p < nP; p++ )
1185 for(
int q = 0; q < nP; q++ )
1187 dCoeff[j][p * nP + q] = dRemapCoeff[p][q][j];
1191 for(
int p = 0; p < nP; p++ )
1193 for(
int q = 0; q < nP; q++ )
1195 dCoeff[nOverlapFaces][p * nP + q] = dataGLLJacobian[p][q][ixFirst];
1198 for(
int j = 0; j < nOverlapFaces; j++ )
1200 for(
int p = 0; p < nP; p++ )
1202 for(
int q = 0; q < nP; q++ )
1204 dCoeff[nOverlapFaces][p * nP + q] -=
1205 dRemapCoeff[p][q][j] * m_meshOverlap->vecFaceArea[ixOverlap + j];
1209 for(
int p = 0; p < nP; p++ )
1211 for(
int q = 0; q < nP; q++ )
1213 dCoeff[nOverlapFaces][p * nP + q] /= dExtraneousArea;
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" );
1229 for(
int j = 0; j < nOverlapFaces; j++ )
1231 for(
int p = 0; p < nP; p++ )
1233 for(
int q = 0; q < nP; q++ )
1235 dRemapCoeff[p][q][j] = dCoeff[j][p * nP + q];
1257 for(
int j = 0; j < nOverlapFaces; j++ )
1259 int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
1262 if( ixSecondFace < 0 )
continue;
1264 for(
int p = 0; p < nP; p++ )
1266 for(
int q = 0; q < nP; q++ )
1270 int ixFirstNode = dataGLLNodes[p][q][ixFirst] - 1;
1272 smatMap( ixSecondFace, ixFirstNode ) += dRemapCoeff[p][q][j] *
1273 m_meshOverlap->vecFaceArea[ixOverlap + j] /
1274 m_meshOutput->vecFaceArea[ixSecondFace];
1278 int ixFirstNode = ixFirst * nP * nP + p * nP + q;
1280 smatMap( ixSecondFace, ixFirstNode ) += dRemapCoeff[p][q][j] *
1281 m_meshOverlap->vecFaceArea[ixOverlap + j] /
1282 m_meshOutput->vecFaceArea[ixSecondFace];
1288 ixOverlap += nOverlapFaces;
1291 output_file.flush();
1292 output_file.close();
1301 const DataArray3D< double >& dataGLLJacobianIn,
1302 const DataArray3D< int >& dataGLLNodesOut,
1303 const DataArray3D< double >& dataGLLJacobianOut,
1304 const DataArray1D< double >& dataNodalAreaOut,
1309 bool fContinuousOut,
1310 bool fNoConservation )
1313 TriangularQuadratureRule triquadrule( 8 );
1315 const DataArray2D< double >& dG = triquadrule.GetG();
1316 const DataArray1D< double >& dW = triquadrule.GetW();
1319 SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1322 DataArray2D< double > dSampleCoeffIn( nPin, nPin );
1323 DataArray2D< double > dSampleCoeffOut( nPout, nPout );
1327 dbgprint.set_prefix(
"[LinearRemapGLLtoGLL2_MOAB]: " );
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 );
1336 DataArray3D< double > dGlobalIntArray( nPin * nPin, m_meshOverlap->faces.size(), nPout * nPout );
1339 DataArray1D< int > nAllOverlapFaces( m_meshInputCov->faces.size() );
1342 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1345 int nOverlapFaces = 0;
1346 size_t ixOverlapTemp = ixOverlap;
1347 for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1350 if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 )
1358 nAllOverlapFaces[ixFirst] = nOverlapFaces;
1361 ixOverlap += nAllOverlapFaces[ixFirst];
1365 DataArray2D< double > dGeometricOutputArea( m_meshOutput->faces.size(), nPout * nPout );
1368 DataArray2D< double > dOverlapOutputArea( m_meshOverlap->faces.size(), nPout * nPout );
1373 const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
1375 if( is_root )
dbgprint.printf( 0,
"Building conservative distribution maps\n" );
1383 NodeVector nodes( 3 );
1384 faceTri.SetNode( 0, 0 );
1385 faceTri.SetNode( 1, 1 );
1386 faceTri.SetNode( 2, 2 );
1388 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1392 if( ixFirst % outputFrequency == 0 && is_root )
1394 dbgprint.printf( 0,
"Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1398 const Face& faceFirst = m_meshInputCov->faces[ixFirst];
1400 const NodeVector& nodesFirst = m_meshInputCov->nodes;
1403 int nOverlapFaces = nAllOverlapFaces[ixFirst];
1405 if( !nOverlapFaces )
continue;
1416 for(
int i = 0; i < nOverlapFaces; i++ )
1419 const Face& faceOverlap = m_meshOverlap->faces[ixOverlap + i];
1421 const NodeVector& nodesOverlap = m_meshOverlap->nodes;
1424 int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1427 if( ixSecond < 0 )
continue;
1429 const NodeVector& nodesSecond = m_meshOutput->nodes;
1431 const Face& faceSecond = m_meshOutput->faces[ixSecond];
1433 int nbEdges = faceOverlap.edges.size();
1434 int nOverlapTriangles = 1;
1438 nOverlapTriangles = nbEdges;
1439 for(
int k = 0; k < nbEdges; k++ )
1441 const Node& node = nodesOverlap[faceOverlap[k]];
1448 Node node0, node1, node2;
1452 for(
int j = 0; j < nOverlapTriangles; j++ )
1456 node0 = nodesOverlap[faceOverlap[0]];
1457 node1 = nodesOverlap[faceOverlap[1]];
1458 node2 = nodesOverlap[faceOverlap[2]];
1459 dTriArea = CalculateFaceArea( faceOverlap, nodesOverlap );
1464 node1 = nodesOverlap[faceOverlap[j]];
1465 int j1 = ( j + 1 ) % nbEdges;
1466 node2 = nodesOverlap[faceOverlap[j1]];
1470 dTriArea = CalculateFaceArea( faceTri, nodes );
1473 for(
int k = 0; k < triquadrule.GetPoints(); k++ )
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;
1482 double dMag = sqrt( dX[0] * dX[0] + dX[1] * dX[1] + dX[2] * dX[2] );
1488 Node nodeQuadrature( dX[0], dX[1], dX[2] );
1495 ApplyInverseMap( faceFirst, nodesFirst, nodeQuadrature, dAlphaIn, dBetaIn );
1502 ApplyInverseMap( faceSecond, nodesSecond, nodeQuadrature, dAlphaOut, dBetaOut );
1522 SampleGLLFiniteElement( nMonotoneType, nPin, dAlphaIn, dBetaIn, dSampleCoeffIn );
1525 SampleGLLFiniteElement( nMonotoneType, nPout, dAlphaOut, dBetaOut, dSampleCoeffOut );
1528 for(
int s = 0; s < nPout; s++ )
1530 for(
int t = 0; t < nPout; t++ )
1532 double dNodeArea = dSampleCoeffOut[s][t] * dW[k] * dTriArea;
1534 dOverlapOutputArea[ixOverlap + i][s * nPout + t] += dNodeArea;
1536 dGeometricOutputArea[ixSecond][s * nPout + t] += dNodeArea;
1542 for(
int p = 0; p < nPin; p++ )
1544 for(
int q = 0; q < nPin; q++ )
1547 for(
int s = 0; s < nPout; s++ )
1549 for(
int t = 0; t < nPout; t++ )
1552 dGlobalIntArray[ixp][ixOverlap + i][ixs] +=
1553 dSampleCoeffOut[s][t] * dSampleCoeffIn[p][q] * dW[k] * dTriArea;
1567 DataArray2D< double > dCoeff( nOverlapFaces * nPout * nPout, nPin * nPin );
1569 for(
int i = 0; i < nOverlapFaces; i++ )
1574 for(
int p = 0; p < nPin; p++ )
1576 for(
int q = 0; q < nPin; q++ )
1579 for(
int s = 0; s < nPout; s++ )
1581 for(
int t = 0; t < nPout; t++ )
1583 dCoeff[i * nPout * nPout + ixs][ixp] = dGlobalIntArray[ixp][ixOverlap + i][ixs] /
1584 dOverlapOutputArea[ixOverlap + i][s * nPout + t];
1596 DataArray1D< double > vecSourceArea( nPin * nPin );
1598 for(
int p = 0; p < nPin; p++ )
1600 for(
int q = 0; q < nPin; q++ )
1602 vecSourceArea[p * nPin + q] = dataGLLJacobianIn[p][q][ixFirst];
1607 DataArray1D< double > vecTargetArea( nOverlapFaces * nPout * nPout );
1609 for(
int i = 0; i < nOverlapFaces; i++ )
1613 for(
int s = 0; s < nPout; s++ )
1615 for(
int t = 0; t < nPout; t++ )
1617 vecTargetArea[i * nPout * nPout + ixs] = dOverlapOutputArea[ixOverlap + i][nPout * s + t];
1625 if( !fNoConservation )
1631 for(
int i = 0; i < nOverlapFaces; i++ )
1634 for(
int p = 0; p < nPin; p++ )
1636 for(
int q = 0; q < nPin; q++ )
1639 for(
int s = 0; s < nPout; s++ )
1641 for(
int t = 0; t < nPout; t++ )
1643 dGlobalIntArray[ixp][ixOverlap + i][ixs] =
1644 dCoeff[i * nPout * nPout + ixs][ixp] * dOverlapOutputArea[ixOverlap + i][s * nPout + t];
1657 for(
int i = 0; i < nPin * nPin; i++ )
1659 double dColSum = 0.0;
1660 for(
int j = 0; j < nOverlapFaces * nPout * nPout; j++ )
1662 dColSum += dCoeff[j][i] * vecTargetArea[j];
1664 printf(
"Col %i: %1.15e\n", i, dColSum / vecSourceArea[i] );
1668 for(
int j = 0; j < nOverlapFaces * nPout * nPout; j++ )
1670 double dRowSum = 0.0;
1671 for(
int i = 0; i < nPin * nPin; i++ )
1673 dRowSum += dCoeff[j][i];
1675 printf(
"Row %i: %1.15e\n", j, dRowSum );
1680 ixOverlap += nOverlapFaces;
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() );
1690 for(
size_t ixSecond = 0; ixSecond < m_meshOutput->faces.size(); ixSecond++ )
1692 dRedistributionMaps[ixSecond].Allocate( nPout * nPout, nPout * nPout );
1694 for(
int i = 0; i < nPout * nPout; i++ )
1696 dRedistributionMaps[ixSecond][i][i] = 1.0;
1699 for(
int s = 0; s < nPout * nPout; s++ )
1701 dRedistSourceArea[s] = dGeometricOutputArea[ixSecond][s];
1704 for(
int s = 0; s < nPout * nPout; s++ )
1706 dRedistTargetArea[s] = dataGLLJacobianOut[s / nPout][s % nPout][ixSecond];
1709 if( !fNoConservation )
1712 ( nMonotoneType != 0 ) );
1714 for(
int s = 0; s < nPout * nPout; s++ )
1716 for(
int t = 0; t < nPout * nPout; t++ )
1718 dRedistributionMaps[ixSecond][s][t] *= dRedistTargetArea[s] / dRedistSourceArea[t];
1725 DataArray1D< double > dTotalGeometricArea( dataNodalAreaOut.GetRows() );
1726 for(
size_t ixSecond = 0; ixSecond < m_meshOutput->faces.size(); ixSecond++ )
1728 for(
int s = 0; s < nPout; s++ )
1730 for(
int t = 0; t < nPout; t++ )
1732 dTotalGeometricArea[dataGLLNodesOut[s][t][ixSecond] - 1] +=
1733 dGeometricOutputArea[ixSecond][s * nPout + t];
1741 if( is_root )
dbgprint.printf( 0,
"Assembling map\n" );
1744 DataArray2D< double > dRedistributedOp( nPin * nPin, nPout * nPout );
1746 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1750 if( ixFirst % outputFrequency == 0 && is_root )
1752 dbgprint.printf( 0,
"Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1756 int nOverlapFaces = nAllOverlapFaces[ixFirst];
1758 if( !nOverlapFaces )
continue;
1761 for(
int j = 0; j < nOverlapFaces; j++ )
1763 int ixSecondFace = m_meshOverlap->vecTargetFaceIx[ixOverlap + j];
1766 if( ixSecondFace < 0 )
continue;
1768 dRedistributedOp.Zero();
1769 for(
int p = 0; p < nPin * nPin; p++ )
1771 for(
int s = 0; s < nPout * nPout; s++ )
1773 for(
int t = 0; t < nPout * nPout; t++ )
1775 dRedistributedOp[p][s] +=
1776 dRedistributionMaps[ixSecondFace][s][t] * dGlobalIntArray[p][ixOverlap + j][t];
1782 for(
int p = 0; p < nPin; p++ )
1784 for(
int q = 0; q < nPin; q++ )
1789 ixFirstNode = dataGLLNodesIn[p][q][ixFirst] - 1;
1793 ixFirstNode = ixFirst * nPin * nPin + p * nPin + q;
1797 for(
int s = 0; s < nPout; s++ )
1799 for(
int t = 0; t < nPout; t++ )
1802 if( fContinuousOut )
1804 ixSecondNode = dataGLLNodesOut[s][t][ixSecondFace] - 1;
1806 if( !fNoConservation )
1808 smatMap( ixSecondNode, ixFirstNode ) +=
1809 dRedistributedOp[ixp][ixs] / dataNodalAreaOut[ixSecondNode];
1813 smatMap( ixSecondNode, ixFirstNode ) +=
1814 dRedistributedOp[ixp][ixs] / dTotalGeometricArea[ixSecondNode];
1819 ixSecondNode = ixSecondFace * nPout * nPout + s * nPout + t;
1821 if( !fNoConservation )
1823 smatMap( ixSecondNode, ixFirstNode ) +=
1824 dRedistributedOp[ixp][ixs] / dataGLLJacobianOut[s][t][ixSecondFace];
1828 smatMap( ixSecondNode, ixFirstNode ) +=
1829 dRedistributedOp[ixp][ixs] / dGeometricOutputArea[ixSecondFace][s * nPout + t];
1843 ixOverlap += nOverlapFaces;
1852 const DataArray3D< double >& ,
1853 const DataArray3D< int >& dataGLLNodesOut,
1854 const DataArray3D< double >& ,
1855 const DataArray1D< double >& dataNodalAreaOut,
1860 bool fContinuousOut )
1863 DataArray1D< double > dGL;
1864 DataArray1D< double > dWL;
1866 GaussLobattoQuadrature::GetPoints( nPout, 0.0, 1.0, dGL, dWL );
1869 SparseMatrix< double >& smatMap = this->GetSparseMatrix();
1872 DataArray2D< double > dSampleCoeffIn( nPin, nPin );
1876 dbgprint.set_prefix(
"[LinearRemapGLLtoGLL2_Pointwise_MOAB]: " );
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 );
1885 DataArray1D< int > nAllOverlapFaces( m_meshInputCov->faces.size() );
1889 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1891 size_t ixOverlapTemp = ixOverlap;
1892 for( ; ixOverlapTemp < m_meshOverlap->faces.size(); ixOverlapTemp++ )
1896 if( ixFirst - m_meshOverlap->vecSourceFaceIx[ixOverlapTemp] != 0 )
break;
1898 nAllOverlapFaces[ixFirst]++;
1902 ixOverlap += nAllOverlapFaces[ixFirst];
1906 DataArray1D< bool > fSecondNodeFound( dataNodalAreaOut.GetRows() );
1910 const unsigned outputFrequency = ( m_meshInputCov->faces.size() / 10 ) + 1;
1913 for(
size_t ixFirst = 0; ixFirst < m_meshInputCov->faces.size(); ixFirst++ )
1917 if( ixFirst % outputFrequency == 0 && is_root )
1919 dbgprint.printf( 0,
"Element %zu/%lu\n", ixFirst, m_meshInputCov->faces.size() );
1923 const Face& faceFirst = m_meshInputCov->faces[ixFirst];
1925 const NodeVector& nodesFirst = m_meshInputCov->nodes;
1928 int nOverlapFaces = nAllOverlapFaces[ixFirst];
1931 for(
int i = 0; i < nOverlapFaces; i++ )
1934 int ixSecond = m_meshOverlap->vecTargetFaceIx[ixOverlap + i];
1937 if( ixSecond < 0 )
continue;
1939 const NodeVector& nodesSecond = m_meshOutput->nodes;
1940 const Face& faceSecond = m_meshOutput->faces[ixSecond];
1943 for(
int s = 0; s < nPout; s++ )
1945 for(
int t = 0; t < nPout; t++ )
1947 size_t ixSecondNode;
1948 if( fContinuousOut )
1950 ixSecondNode = dataGLLNodesOut[s][t][ixSecond] - 1;
1954 ixSecondNode = ixSecond * nPout * nPout + s * nPout + t;
1957 if( ixSecondNode >= fSecondNodeFound.GetRows() ) _EXCEPTIONT(
"Logic error" );
1960 if( fSecondNodeFound[ixSecondNode] )
continue;
1967 ApplyLocalMap( faceSecond, nodesSecond, dGL[t], dGL[s], node, dDx1G, dDx2G );
1974 ApplyInverseMap( faceFirst, nodesFirst, node, dAlphaIn, dBetaIn );
1977 if( ( dAlphaIn < -1.0e-10 ) || ( dAlphaIn > 1.0 + 1.0e-10 ) || ( dBetaIn < -1.0e-10 ) ||
1978 ( dBetaIn > 1.0 + 1.0e-10 ) )
1982 fSecondNodeFound[ixSecondNode] =
true;
1985 SampleGLLFiniteElement( nMonotoneType, nPin, dAlphaIn, dBetaIn, dSampleCoeffIn );
1988 for(
int p = 0; p < nPin; p++ )
1990 for(
int q = 0; q < nPin; q++ )
1995 ixFirstNode = dataGLLNodesIn[p][q][ixFirst] - 1;
1999 ixFirstNode = ixFirst * nPin * nPin + p * nPin + q;
2002 smatMap( ixSecondNode, ixFirstNode ) += dSampleCoeffIn[p][q];
2010 ixOverlap += nOverlapFaces;
2014 for(
size_t i = 0; i < fSecondNodeFound.GetRows(); i++ )
2016 if( !fSecondNodeFound[i] )
2018 _EXCEPTION1(
"Can't sample point %i", i );