int HYPRE_StructMatrixCreate( MPI_Comm comm, HYPRE_StructGrid grid, HYPRE_StructStencil stencil, HYPRE_StructMatrix *matrix ) { *matrix = hypre_StructMatrixCreate(comm, grid, stencil); return 0; }
hypre_StructMatrix * hypre_PFMGCreateCoarseOp5( hypre_StructMatrix *R, hypre_StructMatrix *A, hypre_StructMatrix *P, hypre_StructGrid *coarse_grid, HYPRE_Int cdir ) { hypre_StructMatrix *RAP; hypre_Index *RAP_stencil_shape; hypre_StructStencil *RAP_stencil; HYPRE_Int RAP_stencil_size; HYPRE_Int RAP_stencil_dim; HYPRE_Int RAP_num_ghost[] = {1, 1, 1, 1, 1, 1}; hypre_Index index_temp; HYPRE_Int j, i; HYPRE_Int stencil_rank; RAP_stencil_dim = 2; /*----------------------------------------------------------------------- * Define RAP_stencil *-----------------------------------------------------------------------*/ stencil_rank = 0; /*----------------------------------------------------------------------- * non-symmetric case *-----------------------------------------------------------------------*/ if (!hypre_StructMatrixSymmetric(A)) { /*-------------------------------------------------------------------- * 5 point coarse grid stencil *--------------------------------------------------------------------*/ RAP_stencil_size = 5; RAP_stencil_shape = hypre_CTAlloc(hypre_Index, RAP_stencil_size); for (j = -1; j < 2; j++) { for (i = -1; i < 2; i++) { /*-------------------------------------------------------------- * Storage for 5 elements (c,w,e,n,s) *--------------------------------------------------------------*/ if (i*j == 0) { hypre_SetIndex(index_temp,i,j,0); MapIndex(index_temp, cdir, RAP_stencil_shape[stencil_rank]); stencil_rank++; } } } } /*----------------------------------------------------------------------- * symmetric case *-----------------------------------------------------------------------*/ else { /*-------------------------------------------------------------------- * 5 point coarse grid stencil * Only store the lower triangular part + diagonal = 3 entries, * lower triangular means the lower triangular part on the matrix * in the standard lexicographic ordering. *--------------------------------------------------------------------*/ RAP_stencil_size = 3; RAP_stencil_shape = hypre_CTAlloc(hypre_Index, RAP_stencil_size); for (j = -1; j < 1; j++) { for (i = -1; i < 1; i++) { /*-------------------------------------------------------------- * Store 3 elements in (c,w,s) *--------------------------------------------------------------*/ if( i*j == 0 ) { hypre_SetIndex(index_temp,i,j,0); MapIndex(index_temp, cdir, RAP_stencil_shape[stencil_rank]); stencil_rank++; } } } } RAP_stencil = hypre_StructStencilCreate(RAP_stencil_dim, RAP_stencil_size, RAP_stencil_shape); RAP = hypre_StructMatrixCreate(hypre_StructMatrixComm(A), coarse_grid, RAP_stencil); hypre_StructStencilDestroy(RAP_stencil); /*----------------------------------------------------------------------- * Coarse operator in symmetric iff fine operator is *-----------------------------------------------------------------------*/ hypre_StructMatrixSymmetric(RAP) = hypre_StructMatrixSymmetric(A); /*----------------------------------------------------------------------- * Set number of ghost points - one one each boundary *-----------------------------------------------------------------------*/ hypre_StructMatrixSetNumGhost(RAP, RAP_num_ghost); return RAP; }
hypre_StructMatrix * hypre_SMG2CreateRAPOp( hypre_StructMatrix *R, hypre_StructMatrix *A, hypre_StructMatrix *PT, hypre_StructGrid *coarse_grid ) { hypre_StructMatrix *RAP; hypre_Index *RAP_stencil_shape; hypre_StructStencil *RAP_stencil; int RAP_stencil_size; int RAP_stencil_dim; int RAP_num_ghost[] = {1, 1, 1, 1, 0, 0}; int j, i; int stencil_rank; RAP_stencil_dim = 2; /*----------------------------------------------------------------------- * Define RAP_stencil *-----------------------------------------------------------------------*/ stencil_rank = 0; /*----------------------------------------------------------------------- * non-symmetric case *-----------------------------------------------------------------------*/ if (!hypre_StructMatrixSymmetric(A)) { /*-------------------------------------------------------------------- * 5 or 9 point fine grid stencil produces 9 point RAP *--------------------------------------------------------------------*/ RAP_stencil_size = 9; RAP_stencil_shape = hypre_CTAlloc(hypre_Index, RAP_stencil_size); for (j = -1; j < 2; j++) { for (i = -1; i < 2; i++) { /*-------------------------------------------------------------- * Storage for 9 elements (c,w,e,n,s,sw,se,nw,se) *--------------------------------------------------------------*/ hypre_SetIndex(RAP_stencil_shape[stencil_rank],i,j,0); stencil_rank++; } } } /*----------------------------------------------------------------------- * symmetric case *-----------------------------------------------------------------------*/ else { /*-------------------------------------------------------------------- * 5 or 9 point fine grid stencil produces 9 point RAP * Only store the lower triangular part + diagonal = 5 entries, * lower triangular means the lower triangular part on the matrix * in the standard lexicalgraphic ordering. *--------------------------------------------------------------------*/ RAP_stencil_size = 5; RAP_stencil_shape = hypre_CTAlloc(hypre_Index, RAP_stencil_size); for (j = -1; j < 1; j++) { for (i = -1; i < 2; i++) { /*-------------------------------------------------------------- * Store 5 elements in (c,w,s,sw,se) *--------------------------------------------------------------*/ if( i+j <=0 ) { hypre_SetIndex(RAP_stencil_shape[stencil_rank],i,j,0); stencil_rank++; } } } } RAP_stencil = hypre_StructStencilCreate(RAP_stencil_dim, RAP_stencil_size, RAP_stencil_shape); RAP = hypre_StructMatrixCreate(hypre_StructMatrixComm(A), coarse_grid, RAP_stencil); hypre_StructStencilDestroy(RAP_stencil); /*----------------------------------------------------------------------- * Coarse operator in symmetric iff fine operator is *-----------------------------------------------------------------------*/ hypre_StructMatrixSymmetric(RAP) = hypre_StructMatrixSymmetric(A); /*----------------------------------------------------------------------- * Set number of ghost points *-----------------------------------------------------------------------*/ if (hypre_StructMatrixSymmetric(A)) { RAP_num_ghost[1] = 0; RAP_num_ghost[3] = 0; } hypre_StructMatrixSetNumGhost(RAP, RAP_num_ghost); return RAP; }
int hypre_SStructPMatrixCreate( MPI_Comm comm, hypre_SStructPGrid *pgrid, hypre_SStructStencil **stencils, hypre_SStructPMatrix **pmatrix_ptr ) { hypre_SStructPMatrix *pmatrix; int nvars; int **smaps; hypre_StructStencil ***sstencils; hypre_StructMatrix ***smatrices; int **symmetric; hypre_StructStencil *sstencil; int *vars; hypre_Index *sstencil_shape; int sstencil_size; int new_dim; int *new_sizes; hypre_Index **new_shapes; int size; hypre_StructGrid *sgrid; int vi, vj; int i, j, k; pmatrix = hypre_TAlloc(hypre_SStructPMatrix, 1); hypre_SStructPMatrixComm(pmatrix) = comm; hypre_SStructPMatrixPGrid(pmatrix) = pgrid; hypre_SStructPMatrixStencils(pmatrix) = stencils; nvars = hypre_SStructPGridNVars(pgrid); hypre_SStructPMatrixNVars(pmatrix) = nvars; /* create sstencils */ smaps = hypre_TAlloc(int *, nvars); sstencils = hypre_TAlloc(hypre_StructStencil **, nvars); new_sizes = hypre_TAlloc(int, nvars); new_shapes = hypre_TAlloc(hypre_Index *, nvars); size = 0; for (vi = 0; vi < nvars; vi++) { sstencils[vi] = hypre_TAlloc(hypre_StructStencil *, nvars); for (vj = 0; vj < nvars; vj++) { sstencils[vi][vj] = NULL; new_sizes[vj] = 0; } sstencil = hypre_SStructStencilSStencil(stencils[vi]); vars = hypre_SStructStencilVars(stencils[vi]); sstencil_shape = hypre_StructStencilShape(sstencil); sstencil_size = hypre_StructStencilSize(sstencil); smaps[vi] = hypre_TAlloc(int, sstencil_size); for (i = 0; i < sstencil_size; i++) { j = vars[i]; new_sizes[j]++; } for (vj = 0; vj < nvars; vj++) { if (new_sizes[vj]) { new_shapes[vj] = hypre_TAlloc(hypre_Index, new_sizes[vj]); new_sizes[vj] = 0; } } for (i = 0; i < sstencil_size; i++) { j = vars[i]; k = new_sizes[j]; hypre_CopyIndex(sstencil_shape[i], new_shapes[j][k]); smaps[vi][i] = k; new_sizes[j]++; } new_dim = hypre_StructStencilDim(sstencil); for (vj = 0; vj < nvars; vj++) { if (new_sizes[vj]) { sstencils[vi][vj] = hypre_StructStencilCreate(new_dim, new_sizes[vj], new_shapes[vj]); } size = hypre_max(size, new_sizes[vj]); } } hypre_SStructPMatrixSMaps(pmatrix) = smaps; hypre_SStructPMatrixSStencils(pmatrix) = sstencils; hypre_TFree(new_sizes); hypre_TFree(new_shapes); /* create smatrices */ smatrices = hypre_TAlloc(hypre_StructMatrix **, nvars); for (vi = 0; vi < nvars; vi++) { smatrices[vi] = hypre_TAlloc(hypre_StructMatrix *, nvars); for (vj = 0; vj < nvars; vj++) { smatrices[vi][vj] = NULL; if (sstencils[vi][vj] != NULL) { sgrid = hypre_SStructPGridSGrid(pgrid, vi); smatrices[vi][vj] = hypre_StructMatrixCreate(comm, sgrid, sstencils[vi][vj]); } } } hypre_SStructPMatrixSMatrices(pmatrix) = smatrices; /* create symmetric */ symmetric = hypre_TAlloc(int *, nvars); for (vi = 0; vi < nvars; vi++) { symmetric[vi] = hypre_TAlloc(int, nvars); for (vj = 0; vj < nvars; vj++) { symmetric[vi][vj] = 0; } } hypre_SStructPMatrixSymmetric(pmatrix) = symmetric; hypre_SStructPMatrixSEntriesSize(pmatrix) = size; hypre_SStructPMatrixSEntries(pmatrix) = hypre_TAlloc(int, size); hypre_SStructPMatrixRefCount(pmatrix) = 1; *pmatrix_ptr = pmatrix; return hypre_error_flag; }