31 #define moab_hypre_assert( a, b ) \
34 #define moab_hypre_assert_t( a, b ) \
38 std::cout << "HYPRE Error: " << b << std::endl; \
43 template <
typename TargetT,
typename SourceT >
44 static TargetT* DuplicateAs(
const SourceT* array,
int size,
bool cplusplus =
true )
46 TargetT* target_array = cplusplus ?
new TargetT[
size]
48 : hypre_TAlloc( TargetT,
size );
50 for(
int i = 0; i <
size; i++ )
52 target_array[i] = array[i];
58 void HypreParMatrix::Init()
74 HYPRE_Int* row_starts,
79 resize( glob_size, row_starts, nnz_pr );
82 void HypreParMatrix::resize( HYPRE_Int glob_size,
83 HYPRE_Int* row_starts,
84 HYPRE_Int nnz_pr_diag,
85 HYPRE_Int nnz_pr_offdiag )
90 HYPRE_IJMatrixCreate( pcomm->comm(), row_starts[0], row_starts[1], row_starts[0], row_starts[1], &A );
92 HYPRE_IJMatrixSetObjectType( A, HYPRE_PARCSR );
94 int locsize = row_starts[1] - row_starts[0] + 1;
96 MPI_Comm_size( pcomm->comm(), &num_procs );
97 MPI_Comm_rank( pcomm->comm(), &rank );
99 if( nnz_pr_diag != 0 || nnz_pr_offdiag != 0 )
103 std::vector< HYPRE_Int > m_nnz_pr_diag( locsize, nnz_pr_diag );
104 std::vector< HYPRE_Int > m_onz_pr_diag( locsize, nnz_pr_offdiag );
105 HYPRE_IJMatrixSetDiagOffdSizes( A, &m_nnz_pr_diag[0], &m_onz_pr_diag[0] );
109 std::vector< HYPRE_Int > m_nnz_pr_diag( locsize, nnz_pr_diag );
110 HYPRE_IJMatrixSetRowSizes( A, &m_nnz_pr_diag[0] );
115 HYPRE_IJMatrixInitialize( A );
117 HYPRE_IJMatrixGetObject( A, (
void**)&A_parcsr );
120 HYPRE_MatvecFunctions matvec_fn;
121 interpreter = hypre_CTAlloc( mv_InterfaceInterpreter, 1 );
122 HYPRE_ParCSRSetupInterpreter( interpreter );
123 HYPRE_ParCSRSetupMatvec( &matvec_fn );
126 height = GetNumRows();
127 width = GetNumCols();
130 void HypreParMatrix::resize( HYPRE_Int global_num_rows,
131 HYPRE_Int global_num_cols,
132 HYPRE_Int* row_starts,
133 HYPRE_Int* col_starts,
134 HYPRE_Int* nnz_pr_diag,
135 HYPRE_Int* onz_pr_diag,
136 HYPRE_Int nnz_pr_offdiag )
141 HYPRE_IJMatrixCreate( pcomm->comm(), row_starts[0], row_starts[1], col_starts[0], col_starts[1], &A );
143 HYPRE_IJMatrixSetObjectType( A, HYPRE_PARCSR );
146 MPI_Comm_size( pcomm->comm(), &num_procs );
151 if( nnz_pr_diag == NULL && onz_pr_diag == NULL )
153 std::cout <<
"Parameter nnz_pr_diag and onz_pr_diag cannot be NULL.\n";
154 moab_hypre_assert_t( nnz_pr_diag,
true );
157 HYPRE_IJMatrixSetDiagOffdSizes( A, nnz_pr_diag, onz_pr_diag );
161 HYPRE_IJMatrixSetRowSizes( A, nnz_pr_diag );
166 HYPRE_IJMatrixSetMaxOffProcElmts( A, nnz_pr_offdiag );
170 HYPRE_IJMatrixInitialize( A );
172 HYPRE_IJMatrixGetObject( A, (
void**)&A_parcsr );
175 HYPRE_MatvecFunctions matvec_fn;
176 interpreter = hypre_CTAlloc( mv_InterfaceInterpreter, 1 );
177 HYPRE_ParCSRSetupInterpreter( interpreter );
178 HYPRE_ParCSRSetupMatvec( &matvec_fn );
179 gnrows = global_num_rows;
180 gncols = global_num_cols;
181 height = GetNumRows();
182 width = GetNumCols();
187 HYPRE_Int global_num_rows,
188 HYPRE_Int global_num_cols,
189 HYPRE_Int* row_starts,
190 HYPRE_Int* col_starts,
191 HYPRE_Int nnz_pr_diag,
192 HYPRE_Int onz_pr_diag,
193 HYPRE_Int nnz_pr_offdiag )
197 std::vector< HYPRE_Int > m_nnz_pr_diag( row_starts[1] - row_starts[0], nnz_pr_diag );
198 std::vector< HYPRE_Int > m_onz_pr_diag( row_starts[1] - row_starts[0], onz_pr_diag );
199 resize( global_num_rows, global_num_cols, row_starts, col_starts, &m_nnz_pr_diag[0], &m_onz_pr_diag[0],
203 void HypreParMatrix::MakeRef(
const HypreParMatrix& master )
208 A_parcsr = master.A_parcsr;
210 height = master.GetNumRows();
211 width = master.GetNumCols();
214 HYPRE_IJMatrix HypreParMatrix::StealData()
219 moab_hypre_assert( diagOwner == -1 && offdOwner == -1 && colMapOwner == -1,
"" );
220 moab_hypre_assert( ParCSROwner,
"" );
221 HYPRE_IJMatrix
R = A;
229 void HypreParMatrix::StealData( HypreParMatrix& X )
234 moab_hypre_assert( diagOwner == -1 && offdOwner == -1 && colMapOwner == -1,
"" );
235 moab_hypre_assert( ParCSROwner,
"" );
236 moab_hypre_assert( X.diagOwner == -1 && X.offdOwner == -1 && X.colMapOwner == -1,
"" );
237 moab_hypre_assert( X.ParCSROwner,
"" );
239 A_parcsr = X.A_parcsr;
248 void HypreParMatrix::GetDiag( Eigen::VectorXd& diag )
const
250 int size = this->height;
253 for(
int j = 0; j <
size; j++ )
255 diag( j ) = A_parcsr->diag->data[A_parcsr->diag->i[j]];
256 moab_hypre_assert( A_parcsr->diag->j[A_parcsr->diag->i[j]] == j,
257 "the first entry in each row must be the diagonal one" );
261 static void MakeWrapper(
const hypre_CSRMatrix* mat, Eigen::SparseMatrix< double, Eigen::RowMajor >& wrapper )
263 HYPRE_Int
nr = hypre_CSRMatrixNumRows( mat );
264 HYPRE_Int nc = hypre_CSRMatrixNumCols( mat );
265 Eigen::SparseMatrix< double, Eigen::RowMajor > tmp(
nr, nc );
266 typedef Eigen::Triplet< double > T;
267 HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros( mat );
271 rindx = hypre_CSRMatrixI( mat );
272 cindx = hypre_CSRMatrixJ( mat );
273 dindx = hypre_CSRMatrixData( mat );
275 rindx = DuplicateAs< int >( hypre_CSRMatrixI( mat ),
nr + 1 );
276 cindx = DuplicateAs< int >( hypre_CSRMatrixJ( mat ), nnz );
277 dindx = hypre_CSRMatrixData( mat );
279 std::vector< T > tripletList( nnz );
281 for(
int i = 0; i < nnz; ++i )
283 tripletList.push_back( T( rindx[i], cindx[i], dindx[i] ) );
286 tmp.setFromTriplets( tripletList.begin(), tripletList.end() );
290 void HypreParMatrix::GetDiag( Eigen::SparseMatrix< double, Eigen::RowMajor >& diag )
const
292 MakeWrapper( A_parcsr->diag, diag );
295 void HypreParMatrix::GetOffd( Eigen::SparseMatrix< double, Eigen::RowMajor >& offd, HYPRE_Int*& cmap )
const
297 MakeWrapper( A_parcsr->offd, offd );
298 cmap = A_parcsr->col_map_offd;
335 HYPRE_Int HypreParMatrix::Mult( HypreParVector& x, HypreParVector& y,
double a,
double b )
337 return hypre_ParCSRMatrixMatvec( a, A_parcsr, x.x_par, b, y.x_par );
340 void HypreParMatrix::Mult(
double a,
const HypreParVector& x,
double b, HypreParVector& y )
const
342 hypre_ParCSRMatrixMatvec( a, A_parcsr, x.x_par, b, y.x_par );
345 void HypreParMatrix::MultTranspose(
double a,
const HypreParVector& x,
double b, HypreParVector& y )
const
347 hypre_ParCSRMatrixMatvecT( a, A_parcsr, y.x_par, b, x.x_par );
350 HYPRE_Int HypreParMatrix::Mult( HYPRE_ParVector x, HYPRE_ParVector y,
double a,
double b )
352 return hypre_ParCSRMatrixMatvec( a, A_parcsr, (hypre_ParVector*)x, b, (hypre_ParVector*)y );
355 HYPRE_Int HypreParMatrix::MultTranspose( HypreParVector& x, HypreParVector& y,
double a,
double b )
357 return hypre_ParCSRMatrixMatvecT( a, A_parcsr, x.x_par, b, y.x_par );
438 void HypreParMatrix::ScaleRows(
const Eigen::VectorXd& diag )
440 if( hypre_CSRMatrixNumRows( A_parcsr->diag ) != hypre_CSRMatrixNumRows( A_parcsr->offd ) )
445 if( hypre_CSRMatrixNumRows( A_parcsr->diag ) != diag.size() )
447 MB_SET_ERR_RET(
"Note the Eigen::VectorXd diag is not of compatible dimensions with A\n" );
450 int size = this->height;
451 double* Adiag_data = hypre_CSRMatrixData( A_parcsr->diag );
452 HYPRE_Int* Adiag_i = hypre_CSRMatrixI( A_parcsr->diag );
453 double* Aoffd_data = hypre_CSRMatrixData( A_parcsr->offd );
454 HYPRE_Int* Aoffd_i = hypre_CSRMatrixI( A_parcsr->offd );
458 for(
int i( 0 ); i <
size; ++i )
462 for( jj = Adiag_i[i]; jj < Adiag_i[i + 1]; ++jj )
464 Adiag_data[jj] *= val;
467 for( jj = Aoffd_i[i]; jj < Aoffd_i[i + 1]; ++jj )
469 Aoffd_data[jj] *= val;
474 void HypreParMatrix::InvScaleRows(
const Eigen::VectorXd& diag )
476 if( hypre_CSRMatrixNumRows( A_parcsr->diag ) != hypre_CSRMatrixNumRows( A_parcsr->offd ) )
481 if( hypre_CSRMatrixNumRows( A_parcsr->diag ) != diag.size() )
483 MB_SET_ERR_RET(
"Note the Eigen::VectorXd diag is not of compatible dimensions with A_parcsr\n" );
486 int size = this->height;
487 double* Adiag_data = hypre_CSRMatrixData( A_parcsr->diag );
488 HYPRE_Int* Adiag_i = hypre_CSRMatrixI( A_parcsr->diag );
489 double* Aoffd_data = hypre_CSRMatrixData( A_parcsr->offd );
490 HYPRE_Int* Aoffd_i = hypre_CSRMatrixI( A_parcsr->offd );
494 for(
int i( 0 ); i <
size; ++i )
498 if( 0.0 == diag( i ) )
504 val = 1. / diag( i );
506 for( jj = Adiag_i[i]; jj < Adiag_i[i + 1]; ++jj )
508 Adiag_data[jj] *= val;
511 for( jj = Aoffd_i[i]; jj < Aoffd_i[i + 1]; ++jj )
513 Aoffd_data[jj] *= val;
518 void HypreParMatrix::operator*=(
double s )
520 if( hypre_CSRMatrixNumRows( A_parcsr->diag ) != hypre_CSRMatrixNumRows( A_parcsr->offd ) )
525 HYPRE_Int
size = hypre_CSRMatrixNumRows( A_parcsr->diag );
527 double* Adiag_data = hypre_CSRMatrixData( A_parcsr->diag );
528 HYPRE_Int* Adiag_i = hypre_CSRMatrixI( A_parcsr->diag );
530 for( jj = 0; jj < Adiag_i[
size]; ++jj )
535 double* Aoffd_data = hypre_CSRMatrixData( A_parcsr->offd );
536 HYPRE_Int* Aoffd_i = hypre_CSRMatrixI( A_parcsr->offd );
538 for( jj = 0; jj < Aoffd_i[
size]; ++jj )
544 HYPRE_Int HypreParMatrix::GetValues(
const HYPRE_Int nrows,
548 HYPRE_Complex* values )
550 return HYPRE_IJMatrixGetValues( A, nrows, ncols, rows, cols, values );
553 HYPRE_Int HypreParMatrix::SetValues(
const HYPRE_Int nrows,
555 const HYPRE_Int* rows,
556 const HYPRE_Int* cols,
557 const HYPRE_Complex* values )
559 return HYPRE_IJMatrixSetValues( A, nrows, ncols, rows, cols, values );
562 HYPRE_Int HypreParMatrix::AddToValues(
const HYPRE_Int nrows,
564 const HYPRE_Int* rows,
565 const HYPRE_Int* cols,
566 const HYPRE_Complex* values )
568 return HYPRE_IJMatrixAddToValues( A, nrows, ncols, rows, cols, values );
571 HYPRE_Int HypreParMatrix::verbosity(
const HYPRE_Int level )
573 return HYPRE_IJMatrixSetPrintLevel( A, level );
576 HYPRE_Int HypreParMatrix::FinalizeAssembly()
578 return HYPRE_IJMatrixAssemble( A );
581 static void get_sorted_rows_cols(
const std::vector< HYPRE_Int >& rows_cols, std::vector< HYPRE_Int >& hypre_sorted )
583 hypre_sorted.resize( rows_cols.size() );
584 std::copy( rows_cols.begin(), rows_cols.end(), hypre_sorted.begin() );
585 std::sort( hypre_sorted.begin(), hypre_sorted.end() );
588 void HypreParMatrix::Threshold(
double threshold )
592 hypre_CSRMatrix* csr_A;
593 hypre_CSRMatrix* csr_A_wo_z;
594 hypre_ParCSRMatrix* parcsr_A_ptr;
595 int* row_starts = NULL;
596 int* col_starts = NULL;
601 MPI_Comm_size( pcomm->comm(), &num_procs );
602 ierr += hypre_ParCSRMatrixGetLocalRange( A_parcsr, &row_start, &row_end, &col_start, &col_end );
603 row_starts = hypre_ParCSRMatrixRowStarts( A_parcsr );
604 col_starts = hypre_ParCSRMatrixColStarts( A_parcsr );
605 parcsr_A_ptr = hypre_ParCSRMatrixCreate( pcomm->comm(), row_starts[num_procs], col_starts[num_procs], row_starts,
606 col_starts, 0, 0, 0 );
607 csr_A = hypre_MergeDiagAndOffd( A_parcsr );
608 csr_A_wo_z = hypre_CSRMatrixDeleteZeros( csr_A, threshold );
612 if( csr_A_wo_z == NULL )
618 ierr += hypre_CSRMatrixDestroy( csr_A );
621 ierr += GenerateDiagAndOffd( csr_A_wo_z, parcsr_A_ptr, col_start, col_end );
622 ierr += hypre_CSRMatrixDestroy( csr_A_wo_z );
623 ierr += hypre_ParCSRMatrixDestroy( A_parcsr );
624 hypre_IJMatrixObject( A ) = parcsr_A_ptr;
626 HYPRE_IJMatrixGetObject( A, (
void**)A_parcsr );
629 void HypreParMatrix::EliminateRowsCols(
const std::vector< HYPRE_Int >& rows_cols,
630 const HypreParVector& X,
633 std::vector< HYPRE_Int > rc_sorted;
634 get_sorted_rows_cols( rows_cols, rc_sorted );
635 internal::hypre_ParCSRMatrixEliminateAXB( A_parcsr, rc_sorted.size(), rc_sorted.data(), X.x_par, B.x_par );
638 HypreParMatrix* HypreParMatrix::EliminateRowsCols(
const std::vector< HYPRE_Int >& rows_cols )
640 std::vector< HYPRE_Int > rc_sorted;
641 get_sorted_rows_cols( rows_cols, rc_sorted );
642 hypre_ParCSRMatrix* Ae;
643 internal::hypre_ParCSRMatrixEliminateAAe( A_parcsr, &Ae, rc_sorted.size(), rc_sorted.data() );
644 HypreParMatrix* tmpMat =
new HypreParMatrix( pcomm, M(), N(), RowPart(), ColPart(), 0, 0 );
645 hypre_IJMatrixObject( tmpMat->A ) = Ae;
647 HYPRE_IJMatrixGetObject( tmpMat->A, (
void**)tmpMat->A_parcsr );
651 void HypreParMatrix::Print(
const char* fname, HYPRE_Int offi, HYPRE_Int offj )
653 hypre_ParCSRMatrixPrintIJ( A_parcsr, offi, offj, fname );
656 void HypreParMatrix::Read(
const char* fname )
660 HYPRE_Int base_i, base_j;
661 hypre_ParCSRMatrixReadIJ( pcomm->comm(), fname, &base_i, &base_j, &A_parcsr );
662 hypre_ParCSRMatrixSetNumNonzeros( A_parcsr );
663 hypre_MatvecCommPkgCreate( A_parcsr );
664 height = GetNumRows();
665 width = GetNumCols();
668 void HypreParMatrix::Destroy()
672 hypre_TFree( interpreter );
682 HYPRE_IJMatrixDestroy( A );
686 HypreParMatrix* ParMult( HypreParMatrix* A, HypreParMatrix* B )
688 hypre_ParCSRMatrix* ab;
689 ab = hypre_ParMatmul( A->A_parcsr, B->A_parcsr );
690 hypre_MatvecCommPkgCreate( ab );
691 HypreParMatrix* tmpMat =
new HypreParMatrix( A->pcomm, A->M(), A->N(), A->RowPart(), A->ColPart(), 0, 0 );
692 hypre_IJMatrixObject( tmpMat->A ) = ab;
694 HYPRE_IJMatrixGetObject( tmpMat->A, (
void**)tmpMat->A_parcsr );
698 HypreParMatrix* RAP( HypreParMatrix* A, HypreParMatrix* P )
700 HYPRE_Int P_owns_its_col_starts = hypre_ParCSRMatrixOwnsColStarts( (hypre_ParCSRMatrix*)( P->A_parcsr ) );
701 hypre_ParCSRMatrix* rap;
702 hypre_BoomerAMGBuildCoarseOperator( P->A_parcsr, A->A_parcsr, P->A_parcsr, &rap );
703 hypre_ParCSRMatrixSetNumNonzeros( rap );
707 hypre_ParCSRMatrixSetRowStartsOwner( rap, 0 );
708 hypre_ParCSRMatrixSetColStartsOwner( rap, 0 );
710 if( P_owns_its_col_starts )
712 hypre_ParCSRMatrixSetColStartsOwner( P->A_parcsr, 1 );
715 HypreParMatrix* tmpMat =
new HypreParMatrix( A->pcomm, A->M(), A->N(), A->RowPart(), A->ColPart(), 0, 0 );
716 hypre_IJMatrixObject( tmpMat->A ) = rap;
718 HYPRE_IJMatrixGetObject( tmpMat->A, (
void**)tmpMat->A_parcsr );
722 HypreParMatrix* RAP( HypreParMatrix* Rt, HypreParMatrix* A, HypreParMatrix* P )
724 HYPRE_Int P_owns_its_col_starts = hypre_ParCSRMatrixOwnsColStarts( (hypre_ParCSRMatrix*)( P->A_parcsr ) );
725 HYPRE_Int Rt_owns_its_col_starts = hypre_ParCSRMatrixOwnsColStarts( (hypre_ParCSRMatrix*)( Rt->A_parcsr ) );
726 hypre_ParCSRMatrix* rap;
727 hypre_BoomerAMGBuildCoarseOperator( Rt->A_parcsr, A->A_parcsr, P->A_parcsr, &rap );
728 hypre_ParCSRMatrixSetNumNonzeros( rap );
731 if( !P_owns_its_col_starts )
735 hypre_ParCSRMatrixSetColStartsOwner( rap, 0 );
738 if( !Rt_owns_its_col_starts )
742 hypre_ParCSRMatrixSetRowStartsOwner( rap, 0 );
745 HypreParMatrix* tmpMat =
new HypreParMatrix( A->pcomm, A->M(), A->N(), A->RowPart(), A->ColPart(), 0, 0 );
746 hypre_IJMatrixObject( tmpMat->A ) = rap;
748 HYPRE_IJMatrixGetObject( tmpMat->A, (
void**)tmpMat->A_parcsr );
752 void EliminateBC( HypreParMatrix& A,
754 const std::vector< int >& ess_dof_list,
755 const HypreParVector& X,
759 Ae.Mult( -1.0, X, 1.0, B );
760 hypre_CSRMatrix* A_diag = hypre_ParCSRMatrixDiag( (hypre_ParCSRMatrix*)A.A_parcsr );
761 double* data = hypre_CSRMatrixData( A_diag );
762 HYPRE_Int* I = hypre_CSRMatrixI( A_diag );
764 HYPRE_Int* J = hypre_CSRMatrixJ( A_diag );
765 hypre_CSRMatrix* A_offd = hypre_ParCSRMatrixOffd( (hypre_ParCSRMatrix*)A.A_parcsr );
766 HYPRE_Int* I_offd = hypre_CSRMatrixI( A_offd );
767 double* data_offd = hypre_CSRMatrixData( A_offd );
769 std::vector< HYPRE_Complex > bdata( ess_dof_list.size() ), xdata( ess_dof_list.size() );
770 B.GetValues( ess_dof_list.size(), ess_dof_list.data(), bdata.data() );
771 X.GetValues( ess_dof_list.size(), ess_dof_list.data(), xdata.data() );
773 for(
size_t i = 0; i < ess_dof_list.size(); i++ )
775 int r = ess_dof_list[i];
776 bdata[r] = data[I[r]] * xdata[r];
784 MB_SET_ERR_RET(
"the diagonal entry must be the first entry in the row!" );
787 for(
int j = I[r] + 1; j < I[r + 1]; j++ )
795 for(
int j = I_offd[r]; j < I_offd[r + 1]; j++ )
797 if( data_offd[j] != 0.0 )
810 int ParCSRRelax_Taubin( hypre_ParCSRMatrix* A,
820 hypre_CSRMatrix* A_diag = hypre_ParCSRMatrixDiag( A );
821 HYPRE_Int num_rows = hypre_CSRMatrixNumRows( A_diag );
822 double* u_data = hypre_VectorData( hypre_ParVectorLocalVector( u ) );
823 double* r_data = hypre_VectorData( hypre_ParVectorLocalVector( r ) );
825 for(
int i = 0; i < N; i++ )
828 hypre_ParVectorCopy( f, r );
829 hypre_ParCSRMatrixMatvec( -1.0, A, u, 1.0, r );
831 ( 0 == ( i % 2 ) ) ? coef = lambda : coef = mu;
833 for( HYPRE_Int j = 0; j < num_rows; j++ )
835 u_data[j] += coef * r_data[j] / max_eig;
845 int ParCSRRelax_FIR( hypre_ParCSRMatrix* A,
854 hypre_ParVector* x3 )
857 hypre_CSRMatrix* A_diag = hypre_ParCSRMatrixDiag( A );
858 HYPRE_Int num_rows = hypre_CSRMatrixNumRows( A_diag );
859 double* u_data = hypre_VectorData( hypre_ParVectorLocalVector( u ) );
860 double* x0_data = hypre_VectorData( hypre_ParVectorLocalVector( x0 ) );
861 double* x1_data = hypre_VectorData( hypre_ParVectorLocalVector( x1 ) );
862 double* x2_data = hypre_VectorData( hypre_ParVectorLocalVector( x2 ) );
863 double* x3_data = hypre_VectorData( hypre_ParVectorLocalVector( x3 ) );
864 hypre_ParVectorCopy( u, x0 );
866 hypre_ParVectorCopy( f, x1 );
867 hypre_ParCSRMatrixMatvec( -1.0, A, x0, 1.0, x1 );
869 for( HYPRE_Int i = 0; i < num_rows; i++ )
871 x1_data[i] /= -max_eig;
875 for( HYPRE_Int i = 0; i < num_rows; i++ )
877 x1_data[i] = x0_data[i] - x1_data[i];
881 for( HYPRE_Int i = 0; i < num_rows; i++ )
883 x3_data[i] = fir_coeffs[0] * x0_data[i] + fir_coeffs[1] * x1_data[i];
886 for(
int n = 2; n <= poly_order; n++ )
889 hypre_ParVectorCopy( f, x2 );
890 hypre_ParCSRMatrixMatvec( -1.0, A, x1, 1.0, x2 );
892 for( HYPRE_Int i = 0; i < num_rows; i++ )
894 x2_data[i] /= -max_eig;
902 for( HYPRE_Int i = 0; i < num_rows; i++ )
904 x2_data[i] = ( x1_data[i] - x0_data[i] ) + ( x1_data[i] - 2 * x2_data[i] );
905 x3_data[i] += fir_coeffs[n] * x2_data[i];
906 x0_data[i] = x1_data[i];
907 x1_data[i] = x2_data[i];
911 for( HYPRE_Int i = 0; i < num_rows; i++ )
913 u_data[i] = x3_data[i];