HYPRE_Int hypre_BoomerAMGCoarsenCGC (hypre_ParCSRMatrix *S,HYPRE_Int numberofgrids,HYPRE_Int coarsen_type,HYPRE_Int *CF_marker) /* CGC algorithm * ==================================================================================================== * coupling : the strong couplings * numberofgrids : the number of grids * coarsen_type : the coarsening type * gridpartition : the grid partition * =====================================================================================================*/ { HYPRE_Int j,/*p,*/mpisize,mpirank,/*rstart,rend,*/choice,*coarse,ierr=0; HYPRE_Int *vertexrange = NULL; HYPRE_Int *vertexrange_all = NULL; HYPRE_Int *CF_marker_offd = NULL; HYPRE_Int num_variables = hypre_CSRMatrixNumRows (hypre_ParCSRMatrixDiag(S)); /* HYPRE_Int num_cols_offd = hypre_CSRMatrixNumCols (hypre_ParCSRMatrixOffd (S)); */ /* HYPRE_Int *col_map_offd = hypre_ParCSRMatrixColMapOffd (S); */ /* HYPRE_Real wall_time; */ HYPRE_IJMatrix ijG; hypre_ParCSRMatrix *G; hypre_CSRMatrix *Gseq; MPI_Comm comm = hypre_ParCSRMatrixComm(S); hypre_MPI_Comm_size (comm,&mpisize); hypre_MPI_Comm_rank (comm,&mpirank); #if 0 if (!mpirank) { wall_time = time_getWallclockSeconds(); hypre_printf ("Starting CGC preparation\n"); } #endif AmgCGCPrepare (S,numberofgrids,CF_marker,&CF_marker_offd,coarsen_type,&vertexrange); #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds() - wall_time; hypre_printf ("Finished CGC preparation, wall_time = %f s\n",wall_time); wall_time = time_getWallclockSeconds(); hypre_printf ("Starting CGC matrix assembly\n"); } #endif AmgCGCGraphAssemble (S,vertexrange,CF_marker,CF_marker_offd,coarsen_type,&ijG); #if 0 HYPRE_IJMatrixPrint (ijG,"graph.txt"); #endif HYPRE_IJMatrixGetObject (ijG,(void**)&G); #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds() - wall_time; hypre_printf ("Finished CGC matrix assembly, wall_time = %f s\n",wall_time); wall_time = time_getWallclockSeconds(); hypre_printf ("Starting CGC matrix communication\n"); } #endif #ifdef HYPRE_NO_GLOBAL_PARTITION { /* classical CGC does not really make sense in combination with HYPRE_NO_GLOBAL_PARTITION, but anyway, here it is: */ HYPRE_Int nlocal = vertexrange[1]-vertexrange[0]; vertexrange_all = hypre_CTAlloc (HYPRE_Int,mpisize+1); hypre_MPI_Allgather (&nlocal,1,HYPRE_MPI_INT,vertexrange_all+1,1,HYPRE_MPI_INT,comm); vertexrange_all[0]=0; for (j=2;j<=mpisize;j++) vertexrange_all[j]+=vertexrange_all[j-1]; } #else vertexrange_all = vertexrange; #endif Gseq = hypre_ParCSRMatrixToCSRMatrixAll (G); #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds() - wall_time; hypre_printf ("Finished CGC matrix communication, wall_time = %f s\n",wall_time); } #endif if (Gseq) { /* BM Aug 31, 2006: Gseq==NULL if G has no local rows */ #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds(); hypre_printf ("Starting CGC election\n"); } #endif AmgCGCChoose (Gseq,vertexrange_all,mpisize,&coarse); #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds() - wall_time; hypre_printf ("Finished CGC election, wall_time = %f s\n",wall_time); } #endif #if 0 /* debugging */ if (!mpirank) { for (j=0;j<mpisize;j++) hypre_printf ("Processor %d, choice = %d of range %d - %d\n",j,coarse[j],vertexrange_all[j]+1,vertexrange_all[j+1]); } fflush(stdout); #endif #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds(); hypre_printf ("Starting CGC CF assignment\n"); } #endif choice = coarse[mpirank]; for (j=0;j<num_variables;j++) { if (CF_marker[j]==choice) CF_marker[j] = C_PT; else CF_marker[j] = F_PT; } hypre_CSRMatrixDestroy (Gseq); hypre_TFree (coarse); } else for (j=0;j<num_variables;j++) CF_marker[j] = F_PT; #if 0 if (!mpirank) { wall_time = time_getWallclockSeconds() - wall_time; hypre_printf ("Finished CGC CF assignment, wall_time = %f s\n",wall_time); } #endif #if 0 /* debugging */ if (!mpirank) { wall_time = time_getWallclockSeconds(); hypre_printf ("Starting CGC cleanup\n"); } #endif HYPRE_IJMatrixDestroy (ijG); hypre_TFree (vertexrange); #ifdef HYPRE_NO_GLOBAL_PARTITION hypre_TFree (vertexrange_all); #endif hypre_TFree (CF_marker_offd); #if 0 if (!mpirank) { wall_time = time_getWallclockSeconds() - wall_time; hypre_printf ("Finished CGC cleanup, wall_time = %f s\n",wall_time); } #endif return(ierr); }
HYPRE_Int hypre_BoomerAMGRelaxT( hypre_ParCSRMatrix *A, hypre_ParVector *f, HYPRE_Int *cf_marker, HYPRE_Int relax_type, HYPRE_Int relax_points, HYPRE_Real relax_weight, hypre_ParVector *u, hypre_ParVector *Vtemp ) { hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A); HYPRE_Real *A_diag_data = hypre_CSRMatrixData(A_diag); HYPRE_Int *A_diag_i = hypre_CSRMatrixI(A_diag); HYPRE_Int n_global= hypre_ParCSRMatrixGlobalNumRows(A); HYPRE_Int n = hypre_CSRMatrixNumRows(A_diag); HYPRE_Int first_index = hypre_ParVectorFirstIndex(u); hypre_Vector *u_local = hypre_ParVectorLocalVector(u); HYPRE_Real *u_data = hypre_VectorData(u_local); hypre_Vector *Vtemp_local = hypre_ParVectorLocalVector(Vtemp); HYPRE_Real *Vtemp_data = hypre_VectorData(Vtemp_local); hypre_CSRMatrix *A_CSR; HYPRE_Int *A_CSR_i; HYPRE_Int *A_CSR_j; HYPRE_Real *A_CSR_data; hypre_Vector *f_vector; HYPRE_Real *f_vector_data; HYPRE_Int i; HYPRE_Int jj; HYPRE_Int column; HYPRE_Int relax_error = 0; HYPRE_Real *A_mat; HYPRE_Real *b_vec; HYPRE_Real zero = 0.0; /*----------------------------------------------------------------------- * Switch statement to direct control based on relax_type: * relax_type = 7 -> Jacobi (uses ParMatvec) * relax_type = 9 -> Direct Solve *-----------------------------------------------------------------------*/ switch (relax_type) { case 7: /* Jacobi (uses ParMatvec) */ { /*----------------------------------------------------------------- * Copy f into temporary vector. *-----------------------------------------------------------------*/ hypre_ParVectorCopy(f,Vtemp); /*----------------------------------------------------------------- * Perform MatvecT Vtemp=f-A^Tu *-----------------------------------------------------------------*/ hypre_ParCSRMatrixMatvecT(-1.0,A, u, 1.0, Vtemp); for (i = 0; i < n; i++) { /*----------------------------------------------------------- * If diagonal is nonzero, relax point i; otherwise, skip it. *-----------------------------------------------------------*/ if (A_diag_data[A_diag_i[i]] != zero) { u_data[i] += relax_weight * Vtemp_data[i] / A_diag_data[A_diag_i[i]]; } } } break; case 9: /* Direct solve: use gaussian elimination */ { /*----------------------------------------------------------------- * Generate CSR matrix from ParCSRMatrix A *-----------------------------------------------------------------*/ if (n) { A_CSR = hypre_ParCSRMatrixToCSRMatrixAll(A); f_vector = hypre_ParVectorToVectorAll(f); A_CSR_i = hypre_CSRMatrixI(A_CSR); A_CSR_j = hypre_CSRMatrixJ(A_CSR); A_CSR_data = hypre_CSRMatrixData(A_CSR); f_vector_data = hypre_VectorData(f_vector); A_mat = hypre_CTAlloc(HYPRE_Real, n_global*n_global); b_vec = hypre_CTAlloc(HYPRE_Real, n_global); /*--------------------------------------------------------------- * Load transpose of CSR matrix into A_mat. *---------------------------------------------------------------*/ for (i = 0; i < n_global; i++) { for (jj = A_CSR_i[i]; jj < A_CSR_i[i+1]; jj++) { column = A_CSR_j[jj]; A_mat[column*n_global+i] = A_CSR_data[jj]; } b_vec[i] = f_vector_data[i]; } relax_error = gselim(A_mat,b_vec,n_global); for (i = 0; i < n; i++) { u_data[i] = b_vec[first_index+i]; } hypre_TFree(A_mat); hypre_TFree(b_vec); hypre_CSRMatrixDestroy(A_CSR); A_CSR = NULL; hypre_SeqVectorDestroy(f_vector); f_vector = NULL; } } break; } return(relax_error); }
HYPRE_Int main( HYPRE_Int argc, char *argv[] ) { hypre_CSRMatrix *matrix; hypre_CSRMatrix *matrix1; hypre_ParCSRMatrix *par_matrix; hypre_Vector *x_local; hypre_Vector *y_local; hypre_Vector *y2_local; hypre_ParVector *x; hypre_ParVector *x2; hypre_ParVector *y; hypre_ParVector *y2; HYPRE_Int vecstride_x, idxstride_x, vecstride_y, idxstride_y; HYPRE_Int num_procs, my_id; HYPRE_Int local_size; HYPRE_Int num_vectors; HYPRE_Int global_num_rows, global_num_cols; HYPRE_Int first_index; HYPRE_Int i, j, ierr=0; double *data, *data2; HYPRE_Int *row_starts, *col_starts; char file_name[80]; /* Initialize MPI */ hypre_MPI_Init(&argc, &argv); hypre_MPI_Comm_size(hypre_MPI_COMM_WORLD, &num_procs); hypre_MPI_Comm_rank(hypre_MPI_COMM_WORLD, &my_id); hypre_printf(" my_id: %d num_procs: %d\n", my_id, num_procs); if (my_id == 0) { matrix = hypre_CSRMatrixRead("input"); hypre_printf(" read input\n"); } row_starts = NULL; col_starts = NULL; par_matrix = hypre_CSRMatrixToParCSRMatrix(hypre_MPI_COMM_WORLD, matrix, row_starts, col_starts); hypre_printf(" converted\n"); matrix1 = hypre_ParCSRMatrixToCSRMatrixAll(par_matrix); hypre_sprintf(file_name,"matrix1.%d",my_id); if (matrix1) hypre_CSRMatrixPrint(matrix1, file_name); hypre_ParCSRMatrixPrint(par_matrix,"matrix"); hypre_ParCSRMatrixPrintIJ(par_matrix,0,0,"matrixIJ"); par_matrix = hypre_ParCSRMatrixRead(hypre_MPI_COMM_WORLD,"matrix"); global_num_cols = hypre_ParCSRMatrixGlobalNumCols(par_matrix); hypre_printf(" global_num_cols %d\n", global_num_cols); global_num_rows = hypre_ParCSRMatrixGlobalNumRows(par_matrix); col_starts = hypre_ParCSRMatrixColStarts(par_matrix); first_index = col_starts[my_id]; local_size = col_starts[my_id+1] - first_index; num_vectors = 3; x = hypre_ParMultiVectorCreate( hypre_MPI_COMM_WORLD, global_num_cols, col_starts, num_vectors ); hypre_ParVectorSetPartitioningOwner(x,0); hypre_ParVectorInitialize(x); x_local = hypre_ParVectorLocalVector(x); data = hypre_VectorData(x_local); vecstride_x = hypre_VectorVectorStride(x_local); idxstride_x = hypre_VectorIndexStride(x_local); for ( j=0; j<num_vectors; ++j ) for (i=0; i < local_size; i++) data[i*idxstride_x + j*vecstride_x] = first_index+i+1 + 100*j; x2 = hypre_ParMultiVectorCreate( hypre_MPI_COMM_WORLD, global_num_cols, col_starts, num_vectors ); hypre_ParVectorSetPartitioningOwner(x2,0); hypre_ParVectorInitialize(x2); hypre_ParVectorSetConstantValues(x2,2.0); row_starts = hypre_ParCSRMatrixRowStarts(par_matrix); first_index = row_starts[my_id]; local_size = row_starts[my_id+1] - first_index; y = hypre_ParMultiVectorCreate( hypre_MPI_COMM_WORLD, global_num_rows, row_starts, num_vectors ); hypre_ParVectorSetPartitioningOwner(y,0); hypre_ParVectorInitialize(y); y_local = hypre_ParVectorLocalVector(y); y2 = hypre_ParMultiVectorCreate( hypre_MPI_COMM_WORLD, global_num_rows, row_starts, num_vectors ); hypre_ParVectorSetPartitioningOwner(y2,0); hypre_ParVectorInitialize(y2); y2_local = hypre_ParVectorLocalVector(y2); data2 = hypre_VectorData(y2_local); vecstride_y = hypre_VectorVectorStride(y2_local); idxstride_y = hypre_VectorIndexStride(y2_local); for ( j=0; j<num_vectors; ++j ) for (i=0; i < local_size; i++) data2[i*idxstride_y+j*vecstride_y] = first_index+i+1 + 100*j; hypre_ParVectorSetConstantValues(y,1.0); hypre_printf(" initialized vectors, first_index=%i\n", first_index); hypre_ParVectorPrint(x, "vectorx"); hypre_ParVectorPrint(y, "vectory"); hypre_MatvecCommPkgCreate(par_matrix); hypre_ParCSRMatrixMatvec ( 1.0, par_matrix, x, 1.0, y); hypre_printf(" did matvec\n"); hypre_ParVectorPrint(y, "result"); ierr = hypre_ParCSRMatrixMatvecT ( 1.0, par_matrix, y2, 1.0, x2); hypre_printf(" did matvecT %d\n", ierr); hypre_ParVectorPrint(x2, "transp"); hypre_ParCSRMatrixDestroy(par_matrix); hypre_ParVectorDestroy(x); hypre_ParVectorDestroy(x2); hypre_ParVectorDestroy(y); hypre_ParVectorDestroy(y2); if (my_id == 0) hypre_CSRMatrixDestroy(matrix); if (matrix1) hypre_CSRMatrixDestroy(matrix1); /* Finalize MPI */ hypre_MPI_Finalize(); return 0; }