HYPRE_Int hypre_IJMatrixInitializePETSc(hypre_IJMatrix *matrix) { HYPRE_Int ierr = 0; hypre_ParCSRMatrix *par_matrix = hypre_IJMatrixLocalStorage(matrix); hypre_AuxParCSRMatrix *aux_matrix = hypre_IJMatrixTranslator(matrix); HYPRE_Int local_num_rows = hypre_AuxParCSRMatrixLocalNumRows(aux_matrix); HYPRE_Int local_num_cols = hypre_AuxParCSRMatrixLocalNumCols(aux_matrix); HYPRE_Int *row_space = hypre_AuxParCSRMatrixRowSpace(aux_matrix); HYPRE_Int num_nonzeros = hypre_ParCSRMatrixNumNonzeros(par_matrix); HYPRE_Int local_nnz; HYPRE_Int num_procs, my_id; MPI_Comm comm = hypre_IJMatrixContext(matrix); HYPRE_Int global_num_rows = hypre_IJMatrixM(matrix); hypre_MPI_Comm_size(comm,&num_procs); hypre_MPI_Comm_rank(comm,&my_id); local_nnz = (num_nonzeros/global_num_rows+1)*local_num_rows; if (local_num_rows < 0) hypre_AuxParCSRMatrixLocalNumRows(aux_matrix) = hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(par_matrix)); if (local_num_cols < 0) hypre_AuxParCSRMatrixLocalNumCols(aux_matrix) = hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(par_matrix)); ierr = hypre_AuxParCSRMatrixInitialize(aux_matrix); ierr += hypre_ParCSRMatrixInitialize(par_matrix); return ierr; }
HYPRE_Int hypre_IJMatrixSetLocalSizePETSc(hypre_IJMatrix *matrix, HYPRE_Int local_m, HYPRE_Int local_n) { HYPRE_Int ierr = 0; hypre_AuxParCSRMatrix *aux_data; aux_data = hypre_IJMatrixTranslator(matrix); if (aux_data) { hypre_AuxParCSRMatrixLocalNumRows(aux_data) = local_m; hypre_AuxParCSRMatrixLocalNumCols(aux_data) = local_n; } else { hypre_IJMatrixTranslator(matrix) = hypre_AuxParCSRMatrixCreate(local_m,local_n,NULL); } return ierr; }
/****************************************************************************** * * hypre_IJMatrixSetRowSizesPETSc * *****************************************************************************/ HYPRE_Int hypre_IJMatrixSetRowSizesPETSc(hypre_IJMatrix *matrix, HYPRE_Int *sizes) { HYPRE_Int *row_space; HYPRE_Int local_num_rows; HYPRE_Int i; hypre_AuxParCSRMatrix *aux_matrix; aux_matrix = hypre_IJMatrixTranslator(matrix); if (aux_matrix) local_num_rows = hypre_AuxParCSRMatrixLocalNumRows(aux_matrix); else return -1; row_space = hypre_AuxParCSRMatrixRowSpace(aux_matrix); if (!row_space) row_space = hypre_CTAlloc(HYPRE_Int, local_num_rows); for (i = 0; i < local_num_rows; i++) row_space[i] = sizes[i]; hypre_AuxParCSRMatrixRowSpace(aux_matrix) = row_space; return 0; }
/****************************************************************************** * * hypre_IJMatrixCreatePETSc * * creates AuxParCSRMatrix and ParCSRMatrix if necessary, * generates arrays row_starts and col_starts using either previously * set data local_m and local_n (user defined) or generates them evenly * distributed if not previously defined by user. * *****************************************************************************/ HYPRE_Int hypre_IJMatrixCreatePETSc(hypre_IJMatrix *matrix) { MPI_Comm comm = hypre_IJMatrixContext(matrix); HYPRE_Int global_m = hypre_IJMatrixM(matrix); HYPRE_Int global_n = hypre_IJMatrixN(matrix); hypre_AuxParCSRMatrix *aux_matrix = hypre_IJMatrixTranslator(matrix); HYPRE_Int local_m; HYPRE_Int local_n; HYPRE_Int ierr = 0; HYPRE_Int *row_starts; HYPRE_Int *col_starts; HYPRE_Int num_cols_offd = 0; HYPRE_Int num_nonzeros_diag = 0; HYPRE_Int num_nonzeros_offd = 0; HYPRE_Int num_procs, my_id; HYPRE_Int equal; HYPRE_Int i; hypre_MPI_Comm_size(comm, &num_procs); hypre_MPI_Comm_rank(comm, &my_id); if (aux_matrix) { local_m = hypre_AuxParCSRMatrixLocalNumRows(aux_matrix); local_n = hypre_AuxParCSRMatrixLocalNumCols(aux_matrix); } else { aux_matrix = hypre_AuxParCSRMatrixCreate(-1,-1,NULL); local_m = -1; local_n = -1; hypre_IJMatrixTranslator(matrix) = aux_matrix; } if (local_m < 0) { row_starts = NULL; } else { row_starts = hypre_CTAlloc(HYPRE_Int,num_procs+1); if (my_id == 0 && local_m == global_m) { row_starts[1] = local_m; } else { hypre_MPI_Allgather(&local_m,1,HYPRE_MPI_INT,&row_starts[1],1,HYPRE_MPI_INT,comm); } } if (local_n < 0) { col_starts = NULL; } else { col_starts = hypre_CTAlloc(HYPRE_Int,num_procs+1); if (my_id == 0 && local_n == global_n) { col_starts[1] = local_n; } else { hypre_MPI_Allgather(&local_n,1,HYPRE_MPI_INT,&col_starts[1],1,HYPRE_MPI_INT,comm); } } if (row_starts && col_starts) { equal = 1; for (i=0; i < num_procs; i++) { row_starts[i+1] += row_starts[i]; col_starts[i+1] += col_starts[i]; if (row_starts[i+1] != col_starts[i+1]) equal = 0; } if (equal) { hypre_TFree(col_starts); col_starts = row_starts; } } hypre_IJMatrixLocalStorage(matrix) = hypre_ParCSRMatrixCreate(comm,global_m, global_n,row_starts, col_starts, num_cols_offd, num_nonzeros_diag, num_nonzeros_offd); return ierr; }
/****************************************************************************** * * hypre_IJMatrixInsertRowPETSc * * inserts a row into an IJMatrix, * if diag_i and offd_i are known, those values are inserted directly * into the ParCSRMatrix, * if they are not known, an auxiliary structure, AuxParCSRMatrix is used * *****************************************************************************/ HYPRE_Int hypre_IJMatrixInsertRowPETSc(hypre_IJMatrix *matrix, HYPRE_Int n, HYPRE_Int row, HYPRE_Int *indices, double *coeffs) { HYPRE_Int ierr = 0; hypre_ParCSRMatrix *par_matrix; hypre_AuxParCSRMatrix *aux_matrix; HYPRE_Int *row_starts; HYPRE_Int *col_starts; MPI_Comm comm = hypre_IJMatrixContext(matrix); HYPRE_Int num_procs, my_id; HYPRE_Int row_local; HYPRE_Int col_0, col_n; HYPRE_Int i, temp; HYPRE_Int *indx_diag, *indx_offd; HYPRE_Int **aux_j; HYPRE_Int *local_j; double **aux_data; double *local_data; HYPRE_Int diag_space, offd_space; HYPRE_Int *row_length, *row_space; HYPRE_Int need_aux; HYPRE_Int indx_0; HYPRE_Int diag_indx, offd_indx; hypre_CSRMatrix *diag; HYPRE_Int *diag_i; HYPRE_Int *diag_j; double *diag_data; hypre_CSRMatrix *offd; HYPRE_Int *offd_i; HYPRE_Int *offd_j; double *offd_data; hypre_MPI_Comm_size(comm, &num_procs); hypre_MPI_Comm_rank(comm, &my_id); par_matrix = hypre_IJMatrixLocalStorage( matrix ); aux_matrix = hypre_IJMatrixTranslator(matrix); row_space = hypre_AuxParCSRMatrixRowSpace(aux_matrix); row_length = hypre_AuxParCSRMatrixRowLength(aux_matrix); col_n = hypre_ParCSRMatrixFirstColDiag(par_matrix); row_starts = hypre_ParCSRMatrixRowStarts(par_matrix); col_starts = hypre_ParCSRMatrixColStarts(par_matrix); col_0 = col_starts[my_id]; col_n = col_starts[my_id+1]-1; need_aux = hypre_AuxParCSRMatrixNeedAux(aux_matrix); if (row >= row_starts[my_id] && row < row_starts[my_id+1]) { if (need_aux) { row_local = row - row_starts[my_id]; /* compute local row number */ aux_j = hypre_AuxParCSRMatrixAuxJ(aux_matrix); aux_data = hypre_AuxParCSRMatrixAuxData(aux_matrix); local_j = aux_j[row_local]; local_data = aux_data[row_local]; row_length[row_local] = n; if ( row_space[row_local] < n) { hypre_TFree(local_j); hypre_TFree(local_data); local_j = hypre_CTAlloc(HYPRE_Int,n); local_data = hypre_CTAlloc(double,n); row_space[row_local] = n; } for (i=0; i < n; i++) { local_j[i] = indices[i]; local_data[i] = coeffs[i]; } /* make sure first element is diagonal element, if not, find it and exchange it with first element */ if (local_j[0] != row_local) { for (i=1; i < n; i++) { if (local_j[i] == row_local) { local_j[i] = local_j[0]; local_j[0] = row_local; temp = local_data[0]; local_data[0] = local_data[i]; local_data[i] = temp; break; } } } /* sort data according to column indices, except for first element */ qsort1(local_j,local_data,1,n-1); } else /* insert immediately into data into ParCSRMatrix structure */ {
HYPRE_Int HYPRE_IJMatrixCreate( MPI_Comm comm, HYPRE_Int ilower, HYPRE_Int iupper, HYPRE_Int jlower, HYPRE_Int jupper, HYPRE_IJMatrix *matrix ) { HYPRE_Int *row_partitioning; HYPRE_Int *col_partitioning; HYPRE_Int *info; HYPRE_Int num_procs; HYPRE_Int myid; hypre_IJMatrix *ijmatrix; #ifdef HYPRE_NO_GLOBAL_PARTITION HYPRE_Int row0, col0, rowN, colN; #else HYPRE_Int *recv_buf; HYPRE_Int i, i4; HYPRE_Int square; #endif ijmatrix = hypre_CTAlloc(hypre_IJMatrix, 1); hypre_IJMatrixComm(ijmatrix) = comm; hypre_IJMatrixObject(ijmatrix) = NULL; hypre_IJMatrixTranslator(ijmatrix) = NULL; hypre_IJMatrixObjectType(ijmatrix) = HYPRE_UNITIALIZED; hypre_IJMatrixAssembleFlag(ijmatrix) = 0; hypre_IJMatrixPrintLevel(ijmatrix) = 0; hypre_MPI_Comm_size(comm,&num_procs); hypre_MPI_Comm_rank(comm, &myid); if (ilower > iupper+1 || ilower < 0) { hypre_error_in_arg(2); hypre_TFree(ijmatrix); return hypre_error_flag; } if (iupper < -1) { hypre_error_in_arg(3); hypre_TFree(ijmatrix); return hypre_error_flag; } if (jlower > jupper+1 || jlower < 0) { hypre_error_in_arg(4); hypre_TFree(ijmatrix); return hypre_error_flag; } if (jupper < -1) { hypre_error_in_arg(5); hypre_TFree(ijmatrix); return hypre_error_flag; } #ifdef HYPRE_NO_GLOBAL_PARTITION info = hypre_CTAlloc(HYPRE_Int,2); row_partitioning = hypre_CTAlloc(HYPRE_Int, 2); col_partitioning = hypre_CTAlloc(HYPRE_Int, 2); row_partitioning[0] = ilower; row_partitioning[1] = iupper+1; col_partitioning[0] = jlower; col_partitioning[1] = jupper+1; /* now we need the global number of rows and columns as well as the global first row and column index */ /* proc 0 has the first row and col */ if (myid==0) { info[0] = ilower; info[1] = jlower; } hypre_MPI_Bcast(info, 2, HYPRE_MPI_INT, 0, comm); row0 = info[0]; col0 = info[1]; /* proc (num_procs-1) has the last row and col */ if (myid == (num_procs-1)) { info[0] = iupper; info[1] = jupper; } hypre_MPI_Bcast(info, 2, HYPRE_MPI_INT, num_procs-1, comm); rowN = info[0]; colN = info[1]; hypre_IJMatrixGlobalFirstRow(ijmatrix) = row0; hypre_IJMatrixGlobalFirstCol(ijmatrix) = col0; hypre_IJMatrixGlobalNumRows(ijmatrix) = rowN - row0 + 1; hypre_IJMatrixGlobalNumCols(ijmatrix) = colN - col0 + 1; hypre_TFree(info); #else info = hypre_CTAlloc(HYPRE_Int,4); recv_buf = hypre_CTAlloc(HYPRE_Int,4*num_procs); row_partitioning = hypre_CTAlloc(HYPRE_Int, num_procs+1); info[0] = ilower; info[1] = iupper; info[2] = jlower; info[3] = jupper; /* Generate row- and column-partitioning through information exchange across all processors, check whether the matrix is square, and if the partitionings match. i.e. no overlaps or gaps, if there are overlaps or gaps in the row partitioning or column partitioning , ierr will be set to -9 or -10, respectively */ hypre_MPI_Allgather(info,4,HYPRE_MPI_INT,recv_buf,4,HYPRE_MPI_INT,comm); row_partitioning[0] = recv_buf[0]; square = 1; for (i=0; i < num_procs-1; i++) { i4 = 4*i; if ( recv_buf[i4+1] != (recv_buf[i4+4]-1) ) { hypre_error(HYPRE_ERROR_GENERIC); hypre_TFree(ijmatrix); hypre_TFree(info); hypre_TFree(recv_buf); hypre_TFree(row_partitioning); return hypre_error_flag; } else row_partitioning[i+1] = recv_buf[i4+4]; if ((square && (recv_buf[i4] != recv_buf[i4+2])) || (recv_buf[i4+1] != recv_buf[i4+3]) ) { square = 0; } } i4 = (num_procs-1)*4; row_partitioning[num_procs] = recv_buf[i4+1]+1; if ((recv_buf[i4] != recv_buf[i4+2]) || (recv_buf[i4+1] != recv_buf[i4+3])) square = 0; if (square) col_partitioning = row_partitioning; else { col_partitioning = hypre_CTAlloc(HYPRE_Int,num_procs+1); col_partitioning[0] = recv_buf[2]; for (i=0; i < num_procs-1; i++) { i4 = 4*i; if (recv_buf[i4+3] != recv_buf[i4+6]-1) { hypre_error(HYPRE_ERROR_GENERIC); hypre_TFree(ijmatrix); hypre_TFree(info); hypre_TFree(recv_buf); hypre_TFree(row_partitioning); hypre_TFree(col_partitioning); return hypre_error_flag; } else col_partitioning[i+1] = recv_buf[i4+6]; } col_partitioning[num_procs] = recv_buf[num_procs*4-1]+1; } hypre_IJMatrixGlobalFirstRow(ijmatrix) = row_partitioning[0]; hypre_IJMatrixGlobalFirstCol(ijmatrix) = col_partitioning[0]; hypre_IJMatrixGlobalNumRows(ijmatrix) = row_partitioning[num_procs] - row_partitioning[0]; hypre_IJMatrixGlobalNumCols(ijmatrix) = col_partitioning[num_procs] - col_partitioning[0]; hypre_TFree(info); hypre_TFree(recv_buf); #endif hypre_IJMatrixRowPartitioning(ijmatrix) = row_partitioning; hypre_IJMatrixColPartitioning(ijmatrix) = col_partitioning; *matrix = (HYPRE_IJMatrix) ijmatrix; return hypre_error_flag; }