/*-------------------------------------------------------------------------- * hypre_SStructPMatrixInitialize *--------------------------------------------------------------------------*/ int hypre_SStructPMatrixInitialize( hypre_SStructPMatrix *pmatrix ) { int nvars = hypre_SStructPMatrixNVars(pmatrix); int **symmetric = hypre_SStructPMatrixSymmetric(pmatrix); hypre_SStructPGrid *pgrid = hypre_SStructPMatrixPGrid(pmatrix); HYPRE_SStructVariable *vartypes = hypre_SStructPGridVarTypes(pgrid); int ndim = hypre_SStructPGridNDim(pgrid); int num_ghost[6]= {1, 1, 1, 1, 1, 1}; hypre_StructMatrix *smatrix; hypre_StructGrid *sgrid; hypre_Index varoffset; int vi, vj, d; for (vi = 0; vi < nvars; vi++) { /* use variable vi add_numghost */ sgrid= hypre_SStructPGridSGrid(pgrid, vi); hypre_SStructVariableGetOffset(vartypes[vi], ndim, varoffset); for (vj = 0; vj < nvars; vj++) { smatrix = hypre_SStructPMatrixSMatrix(pmatrix, vi, vj); if (smatrix != NULL) { HYPRE_StructMatrixSetSymmetric(smatrix, symmetric[vi][vj]); hypre_StructMatrixSetNumGhost(smatrix, num_ghost); for (d = 0; d < 3; d++) { hypre_StructMatrixAddNumGhost(smatrix)[2*d]= hypre_IndexD(varoffset, d); hypre_StructMatrixAddNumGhost(smatrix)[2*d+1]= hypre_IndexD(varoffset, d); } hypre_StructMatrixInitialize(smatrix); } } } return hypre_error_flag; }
HYPRE_Int hypre_AMR_CFCoarsen( hypre_SStructMatrix * A, hypre_SStructMatrix * fac_A, hypre_Index refine_factors, HYPRE_Int level ) { MPI_Comm comm = hypre_SStructMatrixComm(A); hypre_SStructGraph *graph = hypre_SStructMatrixGraph(A); HYPRE_Int graph_type = hypre_SStructGraphObjectType(graph); hypre_SStructGrid *grid = hypre_SStructGraphGrid(graph); HYPRE_Int nUventries = hypre_SStructGraphNUVEntries(graph); HYPRE_IJMatrix ij_A = hypre_SStructMatrixIJMatrix(A); HYPRE_Int matrix_type= hypre_SStructMatrixObjectType(A); HYPRE_Int ndim = hypre_SStructMatrixNDim(A); hypre_SStructPMatrix *A_pmatrix; hypre_StructMatrix *smatrix_var; hypre_StructStencil *stencils; HYPRE_Int stencil_size; hypre_Index stencil_shape_i; hypre_Index loop_size; hypre_Box refined_box; double **a_ptrs; hypre_Box *A_dbox; HYPRE_Int part_crse= level-1; HYPRE_Int part_fine= level; hypre_BoxManager *fboxman; hypre_BoxManEntry **boxman_entries, *boxman_entry; HYPRE_Int nboxman_entries; hypre_Box boxman_entry_box; hypre_BoxArrayArray ***fgrid_cinterface_extents; hypre_StructGrid *cgrid; hypre_BoxArray *cgrid_boxes; hypre_Box *cgrid_box; hypre_Index node_extents; hypre_Index stridec, stridef; hypre_BoxArrayArray *cinterface_arrays; hypre_BoxArray *cinterface_array; hypre_Box *fgrid_cinterface; HYPRE_Int centre; HYPRE_Int ci, fi, boxi; HYPRE_Int max_stencil_size= 27; HYPRE_Int false= 0; HYPRE_Int true = 1; HYPRE_Int found; HYPRE_Int *stencil_ranks, *rank_stencils; HYPRE_Int rank, startrank; double *vals; HYPRE_Int i, j, iA; HYPRE_Int nvars, var1; hypre_Index lindex, zero_index; hypre_Index index1, index2; hypre_Index index_temp; hypre_SStructUVEntry *Uventry; HYPRE_Int nUentries, cnt1; HYPRE_Int box_array_size; HYPRE_Int *ncols, *rows, *cols; HYPRE_Int *temp1, *temp2; HYPRE_Int myid; hypre_MPI_Comm_rank(comm, &myid); hypre_SetIndex(zero_index, 0, 0, 0); /*-------------------------------------------------------------------------- * Task: Coarsen the CF interface connections of A into fac_A so that * fac_A will have the stencil coefficients extending into a coarsened * fbox. The centre coefficient is constructed to preserve the row sum. *--------------------------------------------------------------------------*/ if (graph_type == HYPRE_SSTRUCT) { startrank = hypre_SStructGridGhstartRank(grid); } if (graph_type == HYPRE_PARCSR) { startrank = hypre_SStructGridStartRank(grid); } /*-------------------------------------------------------------------------- * Fine grid strides by the refinement factors. *--------------------------------------------------------------------------*/ hypre_SetIndex(stridec, 1, 1, 1); for (i= 0; i< ndim; i++) { stridef[i]= refine_factors[i]; } for (i= ndim; i< 3; i++) { stridef[i]= 1; } /*-------------------------------------------------------------------------- * Determine the c/f interface index boxes: fgrid_cinterface_extents. * These are between fpart= level and cpart= (level-1). The * fgrid_cinterface_extents are indexed by cboxes, but fboxes that * abutt a given cbox must be considered. Moreover, for each fbox, * we can have a c/f interface from a number of different stencil * directions- i.e., we have a boxarrayarray for each cbox, each * fbox leading to a boxarray. * * Algo.: For each cbox: * 1) refine & stretch by a unit in each dimension. * 2) boxman_intersect with the fgrid boxman to get all fboxes contained * or abutting this cbox. * 3) get the fgrid_cinterface_extents for each of these fboxes. * * fgrid_cinterface_extents[var1][ci] *--------------------------------------------------------------------------*/ A_pmatrix= hypre_SStructMatrixPMatrix(fac_A, part_crse); nvars = hypre_SStructPMatrixNVars(A_pmatrix); fgrid_cinterface_extents= hypre_TAlloc(hypre_BoxArrayArray **, nvars); for (var1= 0; var1< nvars; var1++) { fboxman= hypre_SStructGridBoxManager(grid, part_fine, var1); stencils= hypre_SStructPMatrixSStencil(A_pmatrix, var1, var1); cgrid= hypre_SStructPGridSGrid(hypre_SStructPMatrixPGrid(A_pmatrix), var1); cgrid_boxes= hypre_StructGridBoxes(cgrid); fgrid_cinterface_extents[var1]= hypre_TAlloc(hypre_BoxArrayArray *, hypre_BoxArraySize(cgrid_boxes)); hypre_ForBoxI(ci, cgrid_boxes) { cgrid_box= hypre_BoxArrayBox(cgrid_boxes, ci); hypre_StructMapCoarseToFine(hypre_BoxIMin(cgrid_box), zero_index, refine_factors, hypre_BoxIMin(&refined_box)); hypre_SetIndex(index1, refine_factors[0]-1, refine_factors[1]-1, refine_factors[2]-1); hypre_StructMapCoarseToFine(hypre_BoxIMax(cgrid_box), index1, refine_factors, hypre_BoxIMax(&refined_box)); /*------------------------------------------------------------------------ * Stretch the refined_box so that a BoxManIntersect will get abutting * fboxes. *------------------------------------------------------------------------*/ for (i= 0; i< ndim; i++) { hypre_BoxIMin(&refined_box)[i]-= 1; hypre_BoxIMax(&refined_box)[i]+= 1; } hypre_BoxManIntersect(fboxman, hypre_BoxIMin(&refined_box), hypre_BoxIMax(&refined_box), &boxman_entries, &nboxman_entries); fgrid_cinterface_extents[var1][ci]= hypre_BoxArrayArrayCreate(nboxman_entries); /*------------------------------------------------------------------------ * Get the fgrid_cinterface_extents using var1-var1 stencil (only like- * variables couple). *------------------------------------------------------------------------*/ if (stencils != NULL) { for (i= 0; i< nboxman_entries; i++) { hypre_BoxManEntryGetExtents(boxman_entries[i], hypre_BoxIMin(&boxman_entry_box), hypre_BoxIMax(&boxman_entry_box)); hypre_CFInterfaceExtents2(&boxman_entry_box, cgrid_box, stencils, refine_factors, hypre_BoxArrayArrayBoxArray(fgrid_cinterface_extents[var1][ci], i) ); } } hypre_TFree(boxman_entries); } /* hypre_ForBoxI(ci, cgrid_boxes) */ } /* for (var1= 0; var1< nvars; var1++) */
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; }
/*-------------------------------------------------------------------------- * hypre_FacZeroCFSten: Zeroes the coarse stencil coefficients that reach * into an underlying coarsened refinement box. * Algo: For each cbox * { * 1) refine cbox and expand by one in each direction * 2) boxman_intersect with the fboxman * 3) loop over intersection boxes to see if stencil * reaches over. * } *--------------------------------------------------------------------------*/ HYPRE_Int hypre_FacZeroCFSten( hypre_SStructPMatrix *Af, hypre_SStructPMatrix *Ac, hypre_SStructGrid *grid, HYPRE_Int fine_part, hypre_Index rfactors ) { hypre_BoxManager *fboxman; hypre_BoxManEntry **boxman_entries; HYPRE_Int nboxman_entries; hypre_SStructPGrid *p_cgrid; hypre_Box fgrid_box; hypre_StructGrid *cgrid; hypre_BoxArray *cgrid_boxes; hypre_Box *cgrid_box; hypre_Box scaled_box; hypre_Box *shift_ibox; hypre_StructMatrix *smatrix; hypre_StructStencil *stencils; HYPRE_Int stencil_size; hypre_Index refine_factors, upper_shift; hypre_Index stride; hypre_Index stencil_shape; hypre_Index zero_index, ilower, iupper; HYPRE_Int nvars, var1, var2; HYPRE_Int ndim; hypre_Box *ac_dbox; HYPRE_Real *ac_ptr; hypre_Index loop_size; HYPRE_Int iac; HYPRE_Int ci, i, j; HYPRE_Int abs_shape; HYPRE_Int ierr = 0; p_cgrid = hypre_SStructPMatrixPGrid(Ac); nvars = hypre_SStructPMatrixNVars(Ac); ndim = hypre_SStructPGridNDim(p_cgrid); hypre_BoxInit(&fgrid_box, ndim); hypre_BoxInit(&scaled_box, ndim); hypre_ClearIndex(zero_index); hypre_ClearIndex(stride); hypre_ClearIndex(upper_shift); for (i= 0; i< ndim; i++) { stride[i]= 1; upper_shift[i]= rfactors[i]-1; } hypre_CopyIndex(rfactors, refine_factors); if (ndim < 3) { for (i= ndim; i< 3; i++) { refine_factors[i]= 1; } } for (var1= 0; var1< nvars; var1++) { cgrid= hypre_SStructPGridSGrid(hypre_SStructPMatrixPGrid(Ac), var1); cgrid_boxes= hypre_StructGridBoxes(cgrid); fboxman= hypre_SStructGridBoxManager(grid, fine_part, var1); /*------------------------------------------------------------------ * For each parent coarse box find all fboxes that may be connected * through a stencil entry- refine this box, expand it by one * in each direction, and boxman_intersect with fboxman *------------------------------------------------------------------*/ hypre_ForBoxI(ci, cgrid_boxes) { cgrid_box= hypre_BoxArrayBox(cgrid_boxes, ci); hypre_StructMapCoarseToFine(hypre_BoxIMin(cgrid_box), zero_index, refine_factors, hypre_BoxIMin(&scaled_box)); hypre_StructMapCoarseToFine(hypre_BoxIMax(cgrid_box), upper_shift, refine_factors, hypre_BoxIMax(&scaled_box)); hypre_SubtractIndexes(hypre_BoxIMin(&scaled_box), stride, 3, hypre_BoxIMin(&scaled_box)); hypre_AddIndexes(hypre_BoxIMax(&scaled_box), stride, 3, hypre_BoxIMax(&scaled_box)); hypre_BoxManIntersect(fboxman, hypre_BoxIMin(&scaled_box), hypre_BoxIMax(&scaled_box), &boxman_entries, &nboxman_entries); for (var2= 0; var2< nvars; var2++) { stencils= hypre_SStructPMatrixSStencil(Ac, var1, var2); if (stencils != NULL) { stencil_size= hypre_StructStencilSize(stencils); smatrix = hypre_SStructPMatrixSMatrix(Ac, var1, var2); ac_dbox = hypre_BoxArrayBox(hypre_StructMatrixDataSpace(smatrix), ci); /*--------------------------------------------------------- * Find the stencil coefficients that must be zeroed off. * Loop over all possible boxes. *---------------------------------------------------------*/ for (i= 0; i< stencil_size; i++) { hypre_CopyIndex(hypre_StructStencilElement(stencils, i), stencil_shape); AbsStencilShape(stencil_shape, abs_shape); if (abs_shape) /* non-centre stencils are zeroed */ { /* look for connecting fboxes that must be zeroed. */ for (j= 0; j< nboxman_entries; j++) { hypre_BoxManEntryGetExtents(boxman_entries[j], ilower, iupper); hypre_BoxSetExtents(&fgrid_box, ilower, iupper); shift_ibox= hypre_CF_StenBox(&fgrid_box, cgrid_box, stencil_shape, refine_factors, ndim); if ( hypre_BoxVolume(shift_ibox) ) { ac_ptr= hypre_StructMatrixExtractPointerByIndex(smatrix, ci, stencil_shape); hypre_BoxGetSize(shift_ibox, loop_size); hypre_BoxLoop1Begin(ndim, loop_size, ac_dbox, hypre_BoxIMin(shift_ibox), stride, iac); #ifdef HYPRE_USING_OPENMP #pragma omp parallel for private(HYPRE_BOX_PRIVATE,iac) HYPRE_SMP_SCHEDULE #endif hypre_BoxLoop1For(iac) { ac_ptr[iac] = 0.0; } hypre_BoxLoop1End(iac); } /* if ( hypre_BoxVolume(shift_ibox) ) */ hypre_BoxDestroy(shift_ibox); } /* for (j= 0; j< nboxman_entries; j++) */ } /* if (abs_shape) */ } /* for (i= 0; i< stencil_size; i++) */ } /* if (stencils != NULL) */ } /* for (var2= 0; var2< nvars; var2++) */ hypre_TFree(boxman_entries); } /* hypre_ForBoxI ci */