/* GAMG geometric-algebric multiogrid PC - Mark Adams 2011 */ #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ #include "private/matimpl.h" #if defined PETSC_USE_LOG PetscLogEvent gamg_setup_events[NUM_SET]; #endif #define GAMG_MAXLEVELS 30 /*#define GAMG_STAGES*/ #if (defined PETSC_USE_LOG && defined GAMG_STAGES) static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; #endif /* Private context for the GAMG preconditioner */ static PetscBool s_avoid_repart = PETSC_FALSE; static PetscInt s_min_eq_proc = 600; static PetscReal s_threshold = 0.001; typedef struct gamg_TAG { PetscInt m_dim; PetscInt m_Nlevels; PetscInt m_data_sz; PetscInt m_data_rows; PetscInt m_data_cols; PetscInt m_count; PetscInt m_method; /* 0: geomg; 1: plain agg 'pa'; 2: smoothed agg 'sa' */ PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ char m_type[64]; } PC_GAMG; /* -------------------------------------------------------------------------- */ /* PCSetCoordinates_GAMG Input Parameter: . pc - the preconditioner context */ EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCSetCoordinates_GAMG" PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords ) { PC_MG *mg = (PC_MG*)a_pc->data; PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; PetscErrorCode ierr; PetscInt arrsz,bs,my0,kk,ii,jj,nloc,Iend; Mat Amat = a_pc->pmat; PetscFunctionBegin; PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 ); ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); nloc = (Iend-my0)/bs; if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); pc_gamg->m_data_rows = 1; if(a_coords==0 && pc_gamg->m_method==0) pc_gamg->m_method = 2; /* use SA if no coords */ if( pc_gamg->m_method==0 ) pc_gamg->m_data_cols = a_ndm; /* coordinates */ else{ /* SA: null space vectors */ if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */ else if(a_coords != 0 ) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */ else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */ pc_gamg->m_data_rows = bs; } arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols; /* create data - syntactic sugar that should be refactored at some point */ if (pc_gamg->m_data==0 || (pc_gamg->m_data_sz != arrsz)) { ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); } for(kk=0;kkm_data[kk] = -999.; pc_gamg->m_data[arrsz] = -99.; /* copy data in - column oriented */ if( pc_gamg->m_method != 0 ) { const PetscInt M = Iend - my0; for(kk=0;kkm_data[kk*bs]; if( pc_gamg->m_data_cols==1 ) *data = 1.0; else { for(ii=0;iim_data[ii*nloc + kk] = a_coords[kk*a_ndm + ii]; } } } assert(pc_gamg->m_data[arrsz] == -99.); pc_gamg->m_data_sz = arrsz; pc_gamg->m_dim = a_ndm; PetscFunctionReturn(0); } EXTERN_C_END /* ----------------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "PCReset_GAMG" PetscErrorCode PCReset_GAMG(PC pc) { PetscErrorCode ierr; PC_MG *mg = (PC_MG*)pc->data; PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; PetscFunctionBegin; if( pc_gamg->m_data != 0 ) { /* this should not happen, cleaned up in SetUp */ ierr = PetscFree(pc_gamg->m_data); CHKERRQ(ierr); } pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0; PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ /* partitionLevel Input Parameter: . a_Amat_fine - matrix on this fine (k) level . a_ndata_rows - size of data to move (coarse grid) . a_ndata_cols - size of data to move (coarse grid) In/Output Parameter: . a_P_inout - prolongation operator to the next level (k-1) . a_coarse_data - data that need to be moved . a_nactive_proc - number of active procs Output Parameter: . a_Amat_crs - coarse matrix that is created (k-1) */ #define TOP_GRID_LIM 2*s_min_eq_proc /* this will happen anyway */ #undef __FUNCT__ #define __FUNCT__ "partitionLevel" PetscErrorCode partitionLevel( Mat a_Amat_fine, PetscInt a_ndata_rows, PetscInt a_ndata_cols, PetscInt a_cbs, Mat *a_P_inout, PetscReal **a_coarse_data, PetscMPIInt *a_nactive_proc, Mat *a_Amat_crs ) { PetscErrorCode ierr; Mat Cmat,Pnew,Pold=*a_P_inout; IS new_indices,isnum; MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; PetscMPIInt mype,npe,new_npe,nactive; PetscInt neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0; PetscFunctionBegin; ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); /* RAP */ ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0); if( s_avoid_repart ) { *a_Amat_crs = Cmat; /* output */ } else { /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ Mat adj; const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols; const PetscInt stride0=ncrs0*a_ndata_rows; PetscInt is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx; /* create sub communicator */ MPI_Comm cm,new_comm; MPI_Group wg, g2; PetscInt *counts,inpe; PetscMPIInt *ranks; IS isscat; PetscScalar *array; Vec src_crd, dest_crd; PetscReal *data = *a_coarse_data; VecScatter vecscat; IS isnewproc; /* get number of PEs to make active, reduce */ ierr = MatGetSize( Cmat, &neq, &NN ); CHKERRQ(ierr); new_npe = neq/s_min_eq_proc; /* hardwire min. number of eq/proc */ if( new_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */ ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr); ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr); assert(counts[mype]==ncrs0); /* count real active pes */ for( nactive = jj = 0 ; jj < npe ; jj++) { if( counts[jj] != 0 ) { ranks[nactive++] = jj; } } if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */ #ifdef VERBOSE PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s npe (active): %d --> %d. new npe = %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,nactive,new_npe,neq); #endif *a_nactive_proc = new_npe; /* output */ ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr); ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr); ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr); if( cm != MPI_COMM_NULL ) { assert(ncrs0 != 0); ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr); ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); } else assert(ncrs0 == 0); ierr = MPI_Group_free( &wg ); CHKERRQ(ierr); ierr = MPI_Group_free( &g2 ); CHKERRQ(ierr); /* MatPartitioningApply call MatConvert, which is collective */ #if defined PETSC_USE_LOG ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); #endif if( a_cbs == 1) { ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); } else{ /* make a scalar matrix to partition */ Mat tMat; PetscInt ncols,jj,Ii; const PetscScalar *vals; const PetscInt *idx; PetscInt *d_nnz; ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) { ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); d_nnz[jj] = ncols/a_cbs; if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0; ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); } ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, PETSC_DETERMINE, PETSC_DETERMINE, 0, d_nnz, 0, d_nnz, &tMat ); CHKERRQ(ierr); ierr = PetscFree( d_nnz ); CHKERRQ(ierr); for ( ii = Istart0; ii < Iend0; ii++ ) { PetscInt dest_row = ii/a_cbs; ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); for( jj = 0 ; jj < ncols ; jj++ ){ PetscInt dest_col = idx[jj]/a_cbs; PetscScalar v = 1.0; ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); } ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); } ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); static int llev = 0; if( llev++ == -1 ) { PetscViewer viewer; char fname[32]; sprintf(fname,"part_mat_%d.mat",llev); PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); ierr = MatView( tMat, viewer ); CHKERRQ(ierr); ierr = PetscViewerDestroy( &viewer ); } ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); ierr = MatDestroy( &tMat ); CHKERRQ(ierr); } if( ncrs0 != 0 ){ const PetscInt *is_idx; MatPartitioning mpart; /* hack to fix global data that pmetis.c uses in 'adj' */ for( nactive = jj = 0 ; jj < npe ; jj++ ) { if( counts[jj] != 0 ) { adj->rmap->range[nactive++] = adj->rmap->range[jj]; } } adj->rmap->range[nactive] = adj->rmap->range[npe]; ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr); ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); /* collect IS info */ ierr = ISGetLocalSize( isnewproc, &is_sz ); CHKERRQ(ierr); ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr); ierr = ISGetIndices( isnewproc, &is_idx ); CHKERRQ(ierr); /* spread partitioning across machine - best way ??? */ NN = 1; /*npe/new_npe;*/ for( kk = jj = 0 ; kk < is_sz ; kk++ ){ for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) { isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */ } } ierr = ISRestoreIndices( isnewproc, &is_idx ); CHKERRQ(ierr); ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); ierr = PetscCommDestroy( &new_comm ); CHKERRQ(ierr); is_sz *= a_cbs; } else{ isnewproc_idx = 0; is_sz = 0; } ierr = MatDestroy( &adj ); CHKERRQ(ierr); ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc ); if( isnewproc_idx != 0 ) { ierr = PetscFree( isnewproc_idx ); CHKERRQ(ierr); } /* Create an index set from the isnewproc index set to indicate the mapping TO */ ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); /* Determine how many elements are assigned to each processor */ inpe = npe; ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr); ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); ncrs_new = counts[mype]/a_cbs; strideNew = ncrs_new*a_ndata_rows; #if defined PETSC_USE_LOG ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); #endif /* Create a vector to contain the newly ordered element information */ ierr = VecCreate( wcomm, &dest_crd ); ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ /* There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having a block size of ...). Note, ISs are expanded into equation space by 'a_cbs'. */ ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); for(ii=0,jj=0; iiele */ ierr = PetscFree( *a_coarse_data ); CHKERRQ(ierr); ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data ); CHKERRQ(ierr); ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); data = *a_coarse_data; for( jj=0; jjdata; PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; PC_MG_Levels **mglevels = mg->levels; Mat Amat = a_pc->mat, Pmat = a_pc->pmat; PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend; MPI_Comm wcomm = ((PetscObject)a_pc)->comm; PetscMPIInt mype,npe,nactivepe; PetscBool isOK; Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; PetscReal *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS]; MatInfo info; PetscFunctionBegin; pc_gamg->m_count++; if( a_pc->setupcalled > 0 ) { /* just do Galerkin grids */ Mat B,dA,dB; /* PCSetUp_MG seems to insists on setting this to GMRES */ ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); /* currently only handle case where mat and pmat are the same on coarser levels */ ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); /* (re)set to get dirty flag */ ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr); for (level=pc_gamg->m_Nlevels-2; level>-1; level--) { ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); /* the first time through the matrix structure has changed from repartitioning */ if( pc_gamg->m_count == 2 ) { ierr = MatDestroy( &B ); CHKERRQ(ierr); ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); mglevels[level]->A = B; } else { ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); } ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); dB = B; /* setup KSP/PC */ ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr); } #define PRINT_MATS !PETSC_TRUE /* plot levels - A */ if( PRINT_MATS ) { for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){ PetscViewer viewer; char fname[32]; KSP smoother; Mat Tmat, TTm; ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr); sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level); ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); ierr = MatView( Tmat, viewer ); CHKERRQ(ierr); ierr = PetscViewerDestroy( &viewer ); } } a_pc->setupcalled = 2; PetscFunctionReturn(0); } ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); /* GAMG requires input of fine-grid matrix. It determines nlevels. */ ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); /* get data of not around */ if( pc_gamg->m_data == 0 && nloc > 0 ) { ierr = PCSetCoordinates_GAMG( a_pc, -1, 0 ); CHKERRQ( ierr ); } data = pc_gamg->m_data; /* Get A_i and R_i */ ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); #ifdef VERBOSE PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols, (int)(info.nz_used/(PetscReal)N),npe); #endif for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM) && (npe==1 || nactivepe>1); level++ ){ level1 = level + 1; #if (defined PETSC_USE_LOG && defined GAMG_STAGES) ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); #endif #if defined PETSC_USE_LOG ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); #endif ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_method, level, s_threshold, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] ); CHKERRQ(ierr); ierr = PetscFree( data ); CHKERRQ( ierr ); #if defined PETSC_USE_LOG ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); #endif if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ if( isOK ) { #if defined PETSC_USE_LOG ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); #endif ierr = partitionLevel( Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs, &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] ); CHKERRQ(ierr); #if defined PETSC_USE_LOG ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); #endif ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); #ifdef VERBOSE PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols, (int)(info.nz_used/(PetscReal)N),nactivepe); #endif /* coarse grids with SA can have zero row/cols from singleton aggregates */ /* aggregation method should gaurrentee this does not happen! */ #ifdef VERBOSE if( PETSC_TRUE ){ Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */ ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); nloceq = Iend-Istart; ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); for(kk=0;kkm_data = 0; /* destroyed coordinate data */ pc_gamg->m_Nlevels = level + 1; fine_level = level; ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); /* set default smoothers */ for ( lidx=1, level = pc_gamg->m_Nlevels-2; lidx <= fine_level; lidx++, level--) { PetscReal emax, emin; KSP smoother; PC subpc; ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); if( emaxs[level] > 0.0 ) emax=emaxs[level]; else{ /* eigen estimate 'emax' */ KSP eksp; Mat Lmat = Aarr[level]; Vec bb, xx; PC pc; ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); { PetscRandom rctx; ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); } ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); ierr = KSPSetType( eksp, KSPCG ); CHKERRQ(ierr); ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as eigen estimates op. */ ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); CHKERRQ(ierr); ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); ierr = VecDestroy( &xx ); CHKERRQ(ierr); ierr = VecDestroy( &bb ); CHKERRQ(ierr); ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); #ifdef VERBOSE PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s max eigen=%e min=%e PC=%s\n",__FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER); #endif } { PetscInt N1, N0, tt; ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */ emax *= 1.05; } ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */ ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); } { /* coarse grid */ KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; Mat Lmat = Aarr[pc_gamg->m_Nlevels-1]; ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); ierr = PCSetUp( subpc ); CHKERRQ(ierr); ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); assert(ii==1); ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); } /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr); { PetscBool galerkin; ierr = PCMGGetGalerkin( a_pc, &galerkin); CHKERRQ(ierr); if(galerkin){ SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); } } /* plot levels - R/P */ if( PRINT_MATS ) { for (level=pc_gamg->m_Nlevels-1;level>0;level--){ PetscViewer viewer; char fname[32]; sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level); ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr); ierr = PetscViewerDestroy( &viewer ); sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level); ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr); ierr = PetscViewerDestroy( &viewer ); } } /* set interpolation between the levels, clean up */ for (lidx=0,level=pc_gamg->m_Nlevels-1; lidxsetupcalled = 0; ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr ); a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ { KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); ierr = KSPSetUp( smoother ); CHKERRQ(ierr); } PetscFunctionReturn(0); } /* ------------------------------------------------------------------------- */ /* PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner that was created with PCCreate_GAMG(). Input Parameter: . pc - the preconditioner context Application Interface Routine: PCDestroy() */ #undef __FUNCT__ #define __FUNCT__ "PCDestroy_GAMG" PetscErrorCode PCDestroy_GAMG(PC pc) { PetscErrorCode ierr; PC_MG *mg = (PC_MG*)pc->data; PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; PetscFunctionBegin; ierr = PCReset_GAMG(pc);CHKERRQ(ierr); ierr = PetscFree(pc_gamg);CHKERRQ(ierr); ierr = PCDestroy_MG(pc);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PCGAMGSetProcEqLim" /*@ PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction. Not Collective on PC Input Parameters: . pc - the preconditioner context Options Database Key: . -pc_gamg_process_eq_limit Level: intermediate Concepts: Unstructured multrigrid preconditioner .seealso: () @*/ PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(pc,PC_CLASSID,1); ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) { /* PC_MG *mg = (PC_MG*)pc->data; */ /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ PetscFunctionBegin; if(n>0) s_min_eq_proc = n; PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "PCGAMGAvoidRepartitioning" /*@ PCGAMGAvoidRepartitioning - Do not repartition the coarse grids Collective on PC Input Parameters: . pc - the preconditioner context Options Database Key: . -pc_gamg_avoid_repartitioning Level: intermediate Concepts: Unstructured multrigrid preconditioner .seealso: () @*/ PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(pc,PC_CLASSID,1); ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG" PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n) { /* PC_MG *mg = (PC_MG*)pc->data; */ /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ PetscFunctionBegin; s_avoid_repart = n; PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "PCGAMGSetThreshold" /*@ PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph Not collective on PC Input Parameters: . pc - the preconditioner context Options Database Key: . -pc_gamg_threshold Level: intermediate Concepts: Unstructured multrigrid preconditioner .seealso: () @*/ PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(pc,PC_CLASSID,1); ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCGAMGSetThreshold_GAMG" PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) { /* PC_MG *mg = (PC_MG*)pc->data; */ /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ PetscFunctionBegin; s_threshold = n; PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "PCGAMGSetSolverType" /*@ PCGAMGSetSolverType - Set solution method. Collective on PC Input Parameters: . pc - the preconditioner context Options Database Key: . -pc_gamg_type Level: intermediate Concepts: Unstructured multrigrid preconditioner .seealso: () @*/ PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz ) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(pc,PC_CLASSID,1); ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz)); CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCGAMGSetSolverType_GAMG" PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz ) { PC_MG *mg = (PC_MG*)pc->data; PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; PetscFunctionBegin; if(sz < 64) strcpy(pc_gamg->m_type,str); PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "PCSetFromOptions_GAMG" PetscErrorCode PCSetFromOptions_GAMG(PC pc) { PetscErrorCode ierr; PC_MG *mg = (PC_MG*)pc->data; PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; PetscBool flag; PetscFunctionBegin; ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); { ierr = PetscOptionsString("-pc_gamg_type", "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)", "PCGAMGSetSolverType", pc_gamg->m_type, pc_gamg->m_type, 64, &flag ); CHKERRQ(ierr); if (flag && strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2; else if (flag && strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1; else pc_gamg->m_method = 0; /* common (static) variable */ ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning", "Do not repartion coarse grids (false)", "PCGAMGAvoidRepartitioning", s_avoid_repart, &s_avoid_repart, &flag); CHKERRQ(ierr); /* common (static) variable */ ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", "Limit (goal) on number of equations per process on coarse grids", "PCGAMGSetProcEqLim", s_min_eq_proc, &s_min_eq_proc, &flag ); CHKERRQ(ierr); /* common (static) variable */ ierr = PetscOptionsReal("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", s_threshold, &s_threshold, &flag ); CHKERRQ(ierr); } ierr = PetscOptionsTail();CHKERRQ(ierr); PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ /* PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG Input Parameter: . pc - the preconditioner context Application Interface Routine: PCCreate() */ /* MC PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide fine grid discretization matrix and coordinates on the fine grid. Options Database Key: Multigrid options(inherited) + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) -pc_mg_type : (one of) additive multiplicative full cascade kascade GAMG options: Level: intermediate Concepts: multigrid .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() M */ EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCCreate_GAMG" PetscErrorCode PCCreate_GAMG(PC pc) { PetscErrorCode ierr; PC_GAMG *pc_gamg; PC_MG *mg; PetscClassId cookie; PetscFunctionBegin; /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); /* create a supporting struct and attach it to pc */ ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0; mg = (PC_MG*)pc->data; mg->innerctx = pc_gamg; pc_gamg->m_Nlevels = -1; /* overwrite the pointers of PCMG by the functions of PCGAMG */ pc->ops->setfromoptions = PCSetFromOptions_GAMG; pc->ops->setup = PCSetUp_GAMG; pc->ops->reset = PCReset_GAMG; pc->ops->destroy = PCDestroy_GAMG; ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, "PCSetCoordinates_C", "PCSetCoordinates_GAMG", PCSetCoordinates_GAMG); CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, "PCGAMGSetProcEqLim_C", "PCGAMGSetProcEqLim_GAMG", PCGAMGSetProcEqLim_GAMG); CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, "PCGAMGAvoidRepartitioning_C", "PCGAMGAvoidRepartitioning_GAMG", PCGAMGAvoidRepartitioning_GAMG); CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, "PCGAMGSetThreshold_C", "PCGAMGSetThreshold_GAMG", PCGAMGSetThreshold_GAMG); CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, "PCGAMGSetSolverType_C", "PCGAMGSetSolverType_GAMG", PCGAMGSetSolverType_GAMG); CHKERRQ(ierr); #if defined PETSC_USE_LOG static int count = 0; if( count++ == 0 ) { PetscClassIdRegister("GAMG Setup",&cookie); PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); PetscLogEventRegister(" Graph", cookie, &gamg_setup_events[GRAPH]); PetscLogEventRegister(" G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); PetscLogEventRegister(" G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); PetscLogEventRegister(" G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); PetscLogEventRegister(" search&set", cookie, &gamg_setup_events[FIND_V]); PetscLogEventRegister(" SA: init", cookie, &gamg_setup_events[SET7]); /* PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); */ PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_events[SET12]); /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ /* create timer stages */ #if defined GAMG_STAGES { char str[32]; sprintf(str,"MG Level %d (finest)",0); PetscLogStageRegister(str, &gamg_stages[0]); PetscInt lidx; for (lidx=1;lidx<9;lidx++){ sprintf(str,"MG Level %d",lidx); PetscLogStageRegister(str, &gamg_stages[lidx]); } } #endif } #endif PetscFunctionReturn(0); } EXTERN_C_END