15b89ad90SMark F. Adams /* 25b89ad90SMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 45ef31b24SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> 5f96513f1SMatthew G Knepley 6b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 7b4fbaa2aSMark F. Adams PetscLogEvent gamg_setup_events[NUM_SET]; 8b4fbaa2aSMark F. Adams #endif 9b4fbaa2aSMark F. Adams 10b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 11b4fbaa2aSMark F. Adams 12b8fd24d8SMark F. Adams /*#define GAMG_STAGES*/ 13b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 14b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 15b4fbaa2aSMark F. Adams #endif 16f96513f1SMatthew G Knepley 175b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */ 185b89ad90SMark F. Adams typedef struct gamg_TAG { 195b89ad90SMark F. Adams PetscInt m_dim; 205b89ad90SMark F. Adams PetscInt m_Nlevels; 215b89ad90SMark F. Adams PetscInt m_data_sz; 22d3d6bff4SMark F. Adams PetscInt m_data_rows; 23d3d6bff4SMark F. Adams PetscInt m_data_cols; 24d3d6bff4SMark F. Adams PetscBool m_useSA; 255b89ad90SMark F. Adams PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 265b89ad90SMark F. Adams } PC_GAMG; 275b89ad90SMark F. Adams 285b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 295b89ad90SMark F. Adams /* 305b89ad90SMark F. Adams PCSetCoordinates_GAMG 315b89ad90SMark F. Adams 325b89ad90SMark F. Adams Input Parameter: 335b89ad90SMark F. Adams . pc - the preconditioner context 345b89ad90SMark F. Adams */ 35a92563c5SMark F. Adams EXTERN_C_BEGIN 365b89ad90SMark F. Adams #undef __FUNCT__ 375b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 38eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords ) 395b89ad90SMark F. Adams { 40eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 415b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 426c237d78SBarry Smith PetscErrorCode ierr; 43d3d6bff4SMark F. Adams PetscInt arrsz,bs,my0,kk,ii,jj,nloc,Iend; 44038e3b61SMark F. Adams Mat Amat = a_pc->pmat; 45d3d6bff4SMark F. Adams PetscBool flag; 46785cba28SMark F. Adams char str[64]; 475b89ad90SMark F. Adams 485b89ad90SMark F. Adams PetscFunctionBegin; 49038e3b61SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 50d3d6bff4SMark F. Adams ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 51d3d6bff4SMark F. Adams nloc = (Iend-my0)/bs; 52d3d6bff4SMark F. Adams if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 53038e3b61SMark F. Adams 54785cba28SMark F. Adams ierr = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,64,&flag); CHKERRQ( ierr ); 55d3d6bff4SMark F. Adams pc_gamg->m_useSA = (PetscBool)(flag && strcmp(str,"sa") == 0); 56038e3b61SMark F. Adams 57d3d6bff4SMark F. Adams pc_gamg->m_data_rows = 1; 58d3d6bff4SMark F. Adams if(a_coords == 0) pc_gamg->m_useSA = PETSC_TRUE; /* use SA if no data */ 59d3d6bff4SMark F. Adams if( !pc_gamg->m_useSA ) pc_gamg->m_data_cols = a_ndm; /* coordinates */ 60038e3b61SMark F. Adams else{ /* SA: null space vectors */ 61d3d6bff4SMark F. Adams if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */ 62d3d6bff4SMark F. Adams else if(a_coords != 0) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */ 63d3d6bff4SMark F. Adams else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */ 64d3d6bff4SMark F. Adams pc_gamg->m_data_rows = bs; 65038e3b61SMark F. Adams } 66d3d6bff4SMark F. Adams arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols; 675b89ad90SMark F. Adams 68038e3b61SMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 696c237d78SBarry Smith if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 705b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 71eb07cef2SMark F. Adams ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 725b89ad90SMark F. Adams } 73038e3b61SMark F. Adams for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.; 74eb07cef2SMark F. Adams pc_gamg->m_data[arrsz] = -99.; 75038e3b61SMark F. Adams /* copy data in - column oriented */ 76d3d6bff4SMark F. Adams if( pc_gamg->m_useSA ) { 77d3d6bff4SMark F. Adams const PetscInt M = Iend - my0; 78038e3b61SMark F. Adams for(kk=0;kk<nloc;kk++){ 79038e3b61SMark F. Adams PetscReal *data = &pc_gamg->m_data[kk*bs]; 80d3d6bff4SMark F. Adams if( pc_gamg->m_data_cols==1 ) *data = 1.0; 81038e3b61SMark F. Adams else { 82038e3b61SMark F. Adams for(ii=0;ii<bs;ii++) 83038e3b61SMark F. Adams for(jj=0;jj<bs;jj++) 84038e3b61SMark F. Adams if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */ 85038e3b61SMark F. Adams else data[ii*M + jj] = 0.0; 86eb07cef2SMark F. Adams if( a_coords != 0 ) { 87038e3b61SMark F. Adams if( a_ndm == 2 ){ /* rotational modes */ 88038e3b61SMark F. Adams data += 2*M; 89038e3b61SMark F. Adams data[0] = -a_coords[2*kk+1]; 90038e3b61SMark F. Adams data[1] = a_coords[2*kk]; 915b89ad90SMark F. Adams } 92eb07cef2SMark F. Adams else { 93038e3b61SMark F. Adams data += 3*M; 94038e3b61SMark F. Adams data[0] = 0.0; data[M+0] = a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1]; 95038e3b61SMark F. Adams data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = a_coords[3*kk]; 96038e3b61SMark F. Adams data[2] = a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk]; data[2*M+2] = 0.0; 97038e3b61SMark F. Adams } 98eb07cef2SMark F. Adams } 99eb07cef2SMark F. Adams } 100eb07cef2SMark F. Adams } 101eb07cef2SMark F. Adams } 102eb07cef2SMark F. Adams else { 103038e3b61SMark F. Adams for( kk = 0 ; kk < nloc ; kk++ ){ 104038e3b61SMark F. Adams for( ii = 0 ; ii < a_ndm ; ii++ ) { 105038e3b61SMark F. Adams pc_gamg->m_data[ii*nloc + kk] = a_coords[kk*a_ndm + ii]; 106eb07cef2SMark F. Adams } 107eb07cef2SMark F. Adams } 108038e3b61SMark F. Adams } 109eb07cef2SMark F. Adams assert(pc_gamg->m_data[arrsz] == -99.); 110038e3b61SMark F. Adams 1115b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 112eb07cef2SMark F. Adams pc_gamg->m_dim = a_ndm; 1135b89ad90SMark F. Adams 1145b89ad90SMark F. Adams PetscFunctionReturn(0); 1155b89ad90SMark F. Adams } 116a92563c5SMark F. Adams EXTERN_C_END 1175b89ad90SMark F. Adams 118d3d6bff4SMark F. Adams 119d3d6bff4SMark F. Adams /* -----------------------------------------------------------------------------*/ 120d3d6bff4SMark F. Adams #undef __FUNCT__ 121d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 122d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 123d3d6bff4SMark F. Adams { 124d3d6bff4SMark F. Adams PetscErrorCode ierr; 125d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 126d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 127d3d6bff4SMark F. Adams 128d3d6bff4SMark F. Adams PetscFunctionBegin; 129d3d6bff4SMark F. Adams ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 130d3d6bff4SMark F. Adams pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0; 131d3d6bff4SMark F. Adams PetscFunctionReturn(0); 132d3d6bff4SMark F. Adams } 133d3d6bff4SMark F. Adams 1345b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1355b89ad90SMark F. Adams /* 13611e60469SMark F. Adams partitionLevel 1375b89ad90SMark F. Adams 1385b89ad90SMark F. Adams Input Parameter: 1393530afc2SMark F. Adams . a_Amat_fine - matrix on this fine (k) level 140d3d6bff4SMark F. Adams . a_ndata_rows - size of data to move (coarse grid) 141d3d6bff4SMark F. Adams . a_ndata_cols - size of data to move (coarse grid) 1423530afc2SMark F. Adams In/Output Parameter: 1433530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 144eb07cef2SMark F. Adams . a_coarse_data - data that need to be moved 145afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 14611e60469SMark F. Adams Output Parameter: 1473530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 1485b89ad90SMark F. Adams */ 1495cb416c2SMark F. Adams 150*fd3c6afaSMark F. Adams #define MIN_EQ_PROC 600 15122063be5SMark F. Adams #define TOP_GRID_LIM 1000 15222063be5SMark F. Adams 1535b89ad90SMark F. Adams #undef __FUNCT__ 15411e60469SMark F. Adams #define __FUNCT__ "partitionLevel" 1553530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine, 156d3d6bff4SMark F. Adams PetscInt a_ndata_rows, 157d3d6bff4SMark F. Adams PetscInt a_ndata_cols, 158038e3b61SMark F. Adams PetscInt a_cbs, 1593530afc2SMark F. Adams Mat *a_P_inout, 160eb07cef2SMark F. Adams PetscReal **a_coarse_data, 161afc97cdcSMark F. Adams PetscMPIInt *a_nactive_proc, 1623530afc2SMark F. Adams Mat *a_Amat_crs 16311e60469SMark F. Adams ) 1645b89ad90SMark F. Adams { 1655b89ad90SMark F. Adams PetscErrorCode ierr; 166038e3b61SMark F. Adams Mat Cmat,Pnew,Pold=*a_P_inout; 16711e60469SMark F. Adams IS new_indices,isnum; 1683530afc2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 169afc97cdcSMark F. Adams PetscMPIInt mype,npe; 17022063be5SMark F. Adams PetscInt neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new; 17122063be5SMark F. Adams PetscMPIInt new_npe,nactive,ncrs0; 172737a81a9SMark F. Adams PetscBool flag = PETSC_FALSE; 1735b89ad90SMark F. Adams 1745b89ad90SMark F. Adams PetscFunctionBegin; 17511e60469SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 17611e60469SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 17711e60469SMark F. Adams /* RAP */ 178038e3b61SMark F. Adams ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 179737a81a9SMark F. Adams 180038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 181038e3b61SMark F. Adams ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); 182038e3b61SMark F. Adams ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0); 183038e3b61SMark F. Adams 184737a81a9SMark F. Adams ierr = PetscOptionsHasName(PETSC_NULL,"-pc_gamg_avoid_repartitioning",&flag); 185737a81a9SMark F. Adams CHKERRQ( ierr ); 18622063be5SMark F. Adams if( flag ) { 18722063be5SMark F. Adams *a_Amat_crs = Cmat; /* output */ 18822063be5SMark F. Adams } 18922063be5SMark F. Adams else { 19022063be5SMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 1915ef31b24SMark F. Adams Mat adj; 19222063be5SMark F. Adams const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols; 193b4fbaa2aSMark F. Adams const PetscInt stride0=ncrs0*a_ndata_rows; 19492a756f0SMark F. Adams PetscInt is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx; 195c9a0b8beSMark F. Adams /* create sub communicator */ 196c9a0b8beSMark F. Adams MPI_Comm cm,new_comm; 197afc97cdcSMark F. Adams MPI_Group wg, g2; 198*fd3c6afaSMark F. Adams PetscInt *counts,inpe; 199*fd3c6afaSMark F. Adams PetscMPIInt *ranks; 20022063be5SMark F. Adams IS isscat; 20122063be5SMark F. Adams PetscScalar *array; 20222063be5SMark F. Adams Vec src_crd, dest_crd; 20322063be5SMark F. Adams PetscReal *data = *a_coarse_data; 20422063be5SMark F. Adams VecScatter vecscat; 205b4fbaa2aSMark F. Adams IS isnewproc; 206e33ef3b1SMark F. Adams 20722063be5SMark F. Adams /* get number of PEs to make active, reduce */ 20822063be5SMark F. Adams ierr = MatGetSize( Cmat, &neq, &NN );CHKERRQ(ierr); 20922063be5SMark F. Adams new_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */ 21022063be5SMark F. Adams if( new_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; 21122063be5SMark F. Adams else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */ 21222063be5SMark F. Adams 21392a756f0SMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr); 214*fd3c6afaSMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 21592a756f0SMark F. Adams 216*fd3c6afaSMark F. Adams ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr); 217afc97cdcSMark F. Adams assert(counts[mype]==ncrs0); 218afc97cdcSMark F. Adams /* count real active pes */ 21922063be5SMark F. Adams for( nactive = jj = 0 ; jj < npe ; jj++) { 220afc97cdcSMark F. Adams if( counts[jj] != 0 ) { 22122063be5SMark F. Adams ranks[nactive++] = jj; 222afc97cdcSMark F. Adams } 223afc97cdcSMark F. Adams } 22422063be5SMark F. Adams assert(nactive>=new_npe); 22522063be5SMark F. Adams 22622063be5SMark F. Adams 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); 22722063be5SMark F. Adams *a_nactive_proc = new_npe; /* output */ 2282f03bc48SMark F. Adams 229afc97cdcSMark F. Adams ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr); 23022063be5SMark F. Adams ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr); 231afc97cdcSMark F. Adams ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr); 232b4fbaa2aSMark F. Adams 233afc97cdcSMark F. Adams if( cm != MPI_COMM_NULL ) { 234b4fbaa2aSMark F. Adams assert(ncrs0 != 0); 235c9a0b8beSMark F. Adams ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr); 236c9a0b8beSMark F. Adams ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 237afc97cdcSMark F. Adams } 238b4fbaa2aSMark F. Adams else assert(ncrs0 == 0); 239b4fbaa2aSMark F. Adams 240afc97cdcSMark F. Adams ierr = MPI_Group_free( &wg ); CHKERRQ(ierr); 241afc97cdcSMark F. Adams ierr = MPI_Group_free( &g2 ); CHKERRQ(ierr); 242c9a0b8beSMark F. Adams 2435ef31b24SMark F. Adams /* MatPartitioningApply call MatConvert, which is collective */ 244b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 245b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 246b4fbaa2aSMark F. Adams #endif 247038e3b61SMark F. Adams if( a_cbs == 1) { 248038e3b61SMark F. Adams ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 249eb07cef2SMark F. Adams } 250eb07cef2SMark F. Adams else{ 251038e3b61SMark F. Adams /* make a scalar matrix to partition */ 252eb07cef2SMark F. Adams Mat tMat; 253b4fbaa2aSMark F. Adams PetscInt ncols; 254b4fbaa2aSMark F. Adams const PetscScalar *vals; 255b4fbaa2aSMark F. Adams const PetscInt *idx; 2566876a03eSMark F. Adams MatInfo info; 257b4fbaa2aSMark F. Adams 2586876a03eSMark F. Adams ierr = MatGetInfo(Cmat,MAT_LOCAL,&info); CHKERRQ(ierr); 25955ea7f60SMark F. Adams ncols = (PetscInt)info.nz_used/((ncrs0+1)*a_cbs*a_cbs)+1; 2606876a03eSMark F. Adams 261eb07cef2SMark F. Adams ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 262eb07cef2SMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 2636876a03eSMark F. Adams 2*ncols, PETSC_NULL, ncols, PETSC_NULL, 264eb07cef2SMark F. Adams &tMat ); 2656876a03eSMark F. Adams CHKERRQ(ierr); 266eb07cef2SMark F. Adams 26722063be5SMark F. Adams for ( ii = Istart0; ii < Iend0; ii++ ) { 26822063be5SMark F. Adams PetscInt dest_row = ii/a_cbs; 26922063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 270eb07cef2SMark F. Adams for( jj = 0 ; jj < ncols ; jj++ ){ 271038e3b61SMark F. Adams PetscInt dest_col = idx[jj]/a_cbs; 272eb07cef2SMark F. Adams PetscScalar v = 1.0; 273eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 274eb07cef2SMark F. Adams } 27522063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 276eb07cef2SMark F. Adams } 277eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 278eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 279eb07cef2SMark F. Adams 280b4fbaa2aSMark F. Adams static int llev = 0; 281b4fbaa2aSMark F. Adams if( llev++ == -1 ) { 282b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 283b4fbaa2aSMark F. Adams sprintf(fname,"part_mat_%d.mat",llev); 284b4fbaa2aSMark F. Adams PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 285b4fbaa2aSMark F. Adams ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 286b4fbaa2aSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 287b4fbaa2aSMark F. Adams } 288b4fbaa2aSMark F. Adams 289eb07cef2SMark F. Adams ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 290eb07cef2SMark F. Adams 291eb07cef2SMark F. Adams ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 292eb07cef2SMark F. Adams } 2932f03bc48SMark F. Adams 294afc97cdcSMark F. Adams if( ncrs0 != 0 ){ 295b4fbaa2aSMark F. Adams const PetscInt *is_idx; 296b4fbaa2aSMark F. Adams MatPartitioning mpart; 2975ef31b24SMark F. Adams /* hack to fix global data that pmetis.c uses in 'adj' */ 29822063be5SMark F. Adams for( nactive = jj = 0 ; jj < npe ; jj++ ) { 299afc97cdcSMark F. Adams if( counts[jj] != 0 ) { 30022063be5SMark F. Adams adj->rmap->range[nactive++] = adj->rmap->range[jj]; 301afc97cdcSMark F. Adams } 3025ef31b24SMark F. Adams } 30322063be5SMark F. Adams adj->rmap->range[nactive] = adj->rmap->range[npe]; 3042f03bc48SMark F. Adams 3055ef31b24SMark F. Adams ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr); 3065ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 30711e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 3085ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 30911e60469SMark F. Adams ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 31011e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 31122063be5SMark F. Adams 3125ef31b24SMark F. Adams /* collect IS info */ 3135ef31b24SMark F. Adams ierr = ISGetLocalSize( isnewproc, &is_sz ); CHKERRQ(ierr); 314038e3b61SMark F. Adams ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr); 3155ef31b24SMark F. Adams ierr = ISGetIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 316b4fbaa2aSMark F. Adams /* spread partitioning across machine - best way ??? */ 31722063be5SMark F. Adams NN = npe/new_npe; 318b4fbaa2aSMark F. Adams for( kk = jj = 0 ; kk < is_sz ; kk++ ){ 319afc97cdcSMark F. Adams for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) { 32022063be5SMark F. Adams isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */ 321eb07cef2SMark F. Adams } 3225ef31b24SMark F. Adams } 3235ef31b24SMark F. Adams ierr = ISRestoreIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 3245ef31b24SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 3254624a670SBarry Smith ierr = PetscCommDestroy( &new_comm ); CHKERRQ(ierr); 326b4fbaa2aSMark F. Adams 327b4fbaa2aSMark F. Adams is_sz *= a_cbs; 3285ef31b24SMark F. Adams } 3295ef31b24SMark F. Adams else{ 3305ef31b24SMark F. Adams isnewproc_idx = 0; 3315ef31b24SMark F. Adams is_sz = 0; 3325ef31b24SMark F. Adams } 3335ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 3345ef31b24SMark F. Adams ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc ); 335afc97cdcSMark F. Adams if( isnewproc_idx != 0 ) { 3365ef31b24SMark F. Adams ierr = PetscFree( isnewproc_idx ); CHKERRQ(ierr); 3375ef31b24SMark F. Adams } 338e33ef3b1SMark F. Adams 33911e60469SMark F. Adams /* 34011e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 34111e60469SMark F. Adams */ 34211e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 34311e60469SMark F. Adams /* 34411e60469SMark F. Adams Determine how many elements are assigned to each processor 34511e60469SMark F. Adams */ 346*fd3c6afaSMark F. Adams inpe = npe; 347*fd3c6afaSMark F. Adams ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr); 34811e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 349038e3b61SMark F. Adams ncrs_new = counts[mype]/a_cbs; 35022063be5SMark F. Adams strideNew = ncrs_new*a_ndata_rows; 351b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 352b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 353b4fbaa2aSMark F. Adams #endif 35422063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 35511e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 356d3d6bff4SMark F. Adams ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 35711e60469SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/ 35811e60469SMark F. Adams /* 359d3d6bff4SMark F. Adams There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having 360d3d6bff4SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'a_cbs'. 36111e60469SMark F. Adams */ 36292a756f0SMark F. Adams ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 36311e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 364038e3b61SMark F. Adams for(ii=0,jj=0; ii<ncrs0 ; ii++) { 365d3d6bff4SMark F. Adams PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */ 366d3d6bff4SMark F. Adams for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk; 36711e60469SMark F. Adams } 368038e3b61SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 369d3d6bff4SMark F. Adams ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 37011e60469SMark F. Adams CHKERRQ(ierr); 37192a756f0SMark F. Adams ierr = PetscFree( tidx ); CHKERRQ(ierr); 37211e60469SMark F. Adams /* 37311e60469SMark F. Adams Create a vector to contain the original vertex information for each element 37411e60469SMark F. Adams */ 375d3d6bff4SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 376d3d6bff4SMark F. Adams for( jj=0; jj<a_ndata_cols ; jj++ ) { 377d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs0 ; ii++) { 378d3d6bff4SMark F. Adams for( kk=0; kk<a_ndata_rows ; kk++ ) { 379d3d6bff4SMark F. Adams PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj; 380d3d6bff4SMark F. Adams ierr = VecSetValues( src_crd, 1, &jx, &data[ix], INSERT_VALUES ); CHKERRQ(ierr); 381d3d6bff4SMark F. Adams } 382038e3b61SMark F. Adams } 383eb07cef2SMark F. Adams } 384eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 385eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 38611e60469SMark F. Adams /* 38711e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 38811e60469SMark F. Adams to the correct processor 38911e60469SMark F. Adams */ 39011e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 39111e60469SMark F. Adams CHKERRQ(ierr); 39211e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 39311e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 39411e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 39511e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 39611e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 39711e60469SMark F. Adams /* 39811e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 39911e60469SMark F. Adams */ 400eb07cef2SMark F. Adams ierr = PetscFree( *a_coarse_data ); CHKERRQ(ierr); 401d3d6bff4SMark F. Adams ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data ); CHKERRQ(ierr); 402eb07cef2SMark F. Adams 40311e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 404eb07cef2SMark F. Adams data = *a_coarse_data; 405d3d6bff4SMark F. Adams for( jj=0; jj<a_ndata_cols ; jj++ ) { 406d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs_new ; ii++) { 407d3d6bff4SMark F. Adams for( kk=0; kk<a_ndata_rows ; kk++ ) { 408d3d6bff4SMark F. Adams PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj; 409d3d6bff4SMark F. Adams data[ix] = PetscRealPart(array[jx]); 410d3d6bff4SMark F. Adams array[jx] = 1.e300; 411d3d6bff4SMark F. Adams } 412038e3b61SMark F. Adams } 413038e3b61SMark F. Adams } 41411e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 41511e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 416737a81a9SMark F. Adams 41711e60469SMark F. Adams /* 41811e60469SMark F. Adams Invert for MatGetSubMatrix 41911e60469SMark F. Adams */ 420038e3b61SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr); 42111e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 42211e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 42311e60469SMark F. Adams /* A_crs output */ 424038e3b61SMark F. Adams ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 42511e60469SMark F. Adams CHKERRQ(ierr); 426eb07cef2SMark F. Adams 427038e3b61SMark F. Adams ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 428e33ef3b1SMark F. Adams Cmat = *a_Amat_crs; /* output */ 429038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 430eb07cef2SMark F. Adams 43111e60469SMark F. Adams /* prolongator */ 43211e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 43311e60469SMark F. Adams { 43411e60469SMark F. Adams IS findices; 43511e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 43611e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 43711e60469SMark F. Adams CHKERRQ(ierr); 43811e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 43911e60469SMark F. Adams } 4403530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 4413530afc2SMark F. Adams *a_P_inout = Pnew; /* output */ 44292a756f0SMark F. Adams 44311e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 44492a756f0SMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 44592a756f0SMark F. Adams ierr = PetscFree( ranks ); CHKERRQ(ierr); 446e33ef3b1SMark F. Adams } 4475b89ad90SMark F. Adams 4485b89ad90SMark F. Adams PetscFunctionReturn(0); 4495b89ad90SMark F. Adams } 4505b89ad90SMark F. Adams 4515b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4525b89ad90SMark F. Adams /* 4535b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4545b89ad90SMark F. Adams by setting data structures and options. 4555b89ad90SMark F. Adams 4565b89ad90SMark F. Adams Input Parameter: 4575b89ad90SMark F. Adams . pc - the preconditioner context 4585b89ad90SMark F. Adams 4595b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 4605b89ad90SMark F. Adams 4615b89ad90SMark F. Adams Notes: 4625b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 4635b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 4645b89ad90SMark F. Adams */ 4655b89ad90SMark F. Adams #undef __FUNCT__ 4665b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 467eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc ) 4685b89ad90SMark F. Adams { 4695b89ad90SMark F. Adams PetscErrorCode ierr; 470eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 4715b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 472eb07cef2SMark F. Adams Mat Amat = a_pc->mat, Pmat = a_pc->pmat; 473d3d6bff4SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend; 474eb07cef2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_pc)->comm; 4753530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 476038e3b61SMark F. Adams PetscBool isOK; 477587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 478587fa25dSMark F. Adams PetscReal *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS]; 479737a81a9SMark F. Adams MatInfo info; 4805ef31b24SMark F. Adams 4815b89ad90SMark F. Adams PetscFunctionBegin; 482eb07cef2SMark F. Adams if( a_pc->setupcalled ) { 483eb07cef2SMark F. Adams /* no state data in GAMG to destroy */ 484eb07cef2SMark F. Adams ierr = PCReset_MG( a_pc ); CHKERRQ(ierr); 485eb07cef2SMark F. Adams } 486baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 487baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 4885b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 4895b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 4903530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 491eb07cef2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 492eb07cef2SMark F. Adams nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 493eb07cef2SMark F. Adams 494038e3b61SMark F. Adams /* get data of not around */ 495038e3b61SMark F. Adams if( pc_gamg->m_data == 0 && nloc > 0 ) { 496038e3b61SMark F. Adams ierr = PCSetCoordinates_GAMG( a_pc, -1, 0 ); CHKERRQ( ierr ); 497eb07cef2SMark F. Adams } 498eb07cef2SMark F. Adams data = pc_gamg->m_data; 499038e3b61SMark F. Adams 500eb07cef2SMark F. Adams /* Get A_i and R_i */ 501737a81a9SMark F. Adams ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 502*fd3c6afaSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%ld, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 503*fd3c6afaSMark F. Adams mype,__FUNCT__,0,(long long int)N,(int)pc_gamg->m_data_rows,(int)pc_gamg->m_data_cols, 504*fd3c6afaSMark F. Adams (int)(info.nz_used/(PetscReal)N),(int)npe); 5058f4b7eb5SMark F. Adams for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 506785cba28SMark F. Adams level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM) && (npe==1 || nactivepe>1); 5070205a208SMark F. Adams level++ ){ 5085b89ad90SMark F. Adams level1 = level + 1; 509b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 510b4fbaa2aSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 511b4fbaa2aSMark F. Adams #endif 512b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 513b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 514b4fbaa2aSMark F. Adams #endif 515d3d6bff4SMark F. Adams ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_useSA, 516785cba28SMark F. Adams level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] ); 5175b89ad90SMark F. Adams CHKERRQ(ierr); 518d3d6bff4SMark F. Adams ierr = PetscFree( data ); CHKERRQ( ierr ); 519b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 520b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 521b4fbaa2aSMark F. Adams #endif 522baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 523d3d6bff4SMark F. Adams 524baa4e9faSMark F. Adams if( isOK ) { 525b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 526b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 527b4fbaa2aSMark F. Adams #endif 528d3d6bff4SMark F. Adams ierr = partitionLevel( Aarr[level], pc_gamg->m_useSA ? bs : 1, pc_gamg->m_data_cols, bs, 529eb07cef2SMark F. Adams &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] ); 5303530afc2SMark F. Adams CHKERRQ(ierr); 531b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 532b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 533b4fbaa2aSMark F. Adams #endif 5343530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 535737a81a9SMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 536*fd3c6afaSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%ld, bs=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 537*fd3c6afaSMark F. Adams mype,__FUNCT__,(int)level1,(long long int)N,(int)bs,(int)pc_gamg->m_data_cols, 538*fd3c6afaSMark F. Adams (int)(info.nz_used/(PetscReal)N),(int)nactivepe); 539e33ef3b1SMark F. Adams /* coarse grids with SA can have zero row/cols from singleton aggregates */ 540e33ef3b1SMark F. Adams /* aggregation method can probably gaurrentee this does not happen! - be safe for now */ 541737a81a9SMark F. Adams 542e33ef3b1SMark F. Adams if( PETSC_TRUE ){ 543785cba28SMark F. Adams Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; 544785cba28SMark F. Adams v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */ 545e33ef3b1SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); 546785cba28SMark F. Adams nloceq = Iend-Istart; 547e33ef3b1SMark F. Adams ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); 548e33ef3b1SMark F. Adams ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); 549e33ef3b1SMark F. Adams ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); 550785cba28SMark F. Adams for(kk=0;kk<nloceq;kk++){ 551e33ef3b1SMark F. Adams if(data_arr[kk]==0.0) { 552785cba28SMark F. Adams id = kk + Istart; 553e33ef3b1SMark F. Adams ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); 554e33ef3b1SMark F. Adams CHKERRQ(ierr); 555785cba28SMark F. Adams PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added diag to zero (%d) on level %d \n",mype,__FUNCT__,id,level); 556e33ef3b1SMark F. Adams } 557e33ef3b1SMark F. Adams } 558e33ef3b1SMark F. Adams ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); 559e33ef3b1SMark F. Adams ierr = VecDestroy( &diag ); CHKERRQ(ierr); 560e33ef3b1SMark F. Adams ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 561e33ef3b1SMark F. Adams ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 562e33ef3b1SMark F. Adams } 563baa4e9faSMark F. Adams } 564baa4e9faSMark F. Adams else{ 565be544d3cSMark F. Adams coarse_data = 0; 566baa4e9faSMark F. Adams break; 567baa4e9faSMark F. Adams } 568eb07cef2SMark F. Adams data = coarse_data; 569b4fbaa2aSMark F. Adams 570b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 571b4fbaa2aSMark F. Adams ierr = PetscLogStagePop(); CHKERRQ( ierr ); 572b4fbaa2aSMark F. Adams #endif 5735b89ad90SMark F. Adams } 574be544d3cSMark F. Adams if( coarse_data ) { 575eb07cef2SMark F. Adams ierr = PetscFree( coarse_data ); CHKERRQ( ierr ); 576be544d3cSMark F. Adams } 577baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 5785b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 5795b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 5805b89ad90SMark F. Adams fine_level = level; 581eb07cef2SMark F. Adams ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 5825b89ad90SMark F. Adams 583fc4362bfSMark F. Adams /* set default smoothers */ 584587fa25dSMark F. Adams for ( lidx=1, level = pc_gamg->m_Nlevels-2; 585587fa25dSMark F. Adams lidx <= fine_level; 586587fa25dSMark F. Adams lidx++, level--) { 5875745e0f5SMark F. Adams PetscReal emax, emin; 5885b89ad90SMark F. Adams KSP smoother; PC subpc; 589587fa25dSMark F. Adams ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 5905b89ad90SMark F. Adams ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 591038e3b61SMark F. Adams if( emaxs[level] > 0.0 ) emax=emaxs[level]; 592587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 593587fa25dSMark F. Adams KSP eksp; Mat Lmat = Aarr[level]; 594fc4362bfSMark F. Adams Vec bb, xx; PC pc; 595038e3b61SMark F. Adams 5965745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 5975745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 598fc4362bfSMark F. Adams { 599fc4362bfSMark F. Adams PetscRandom rctx; 600fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 601fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 602fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 603fc4362bfSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 6045b89ad90SMark F. Adams } 605fc4362bfSMark F. Adams ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 606e33ef3b1SMark F. Adams ierr = KSPSetType( eksp, KSPCG ); CHKERRQ(ierr); 607fc4362bfSMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 6085745e0f5SMark F. Adams ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr ); 609fc4362bfSMark F. Adams ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 6101fddbf69SMark F. Adams ierr = PCSetType( pc, PCPBJACOBI ); CHKERRQ(ierr); /* should be same as above */ 611038e3b61SMark F. Adams ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 612fc4362bfSMark F. Adams CHKERRQ(ierr); 613e33ef3b1SMark F. Adams //ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr); 614fc4362bfSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 615fc4362bfSMark F. Adams 616fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 617fc4362bfSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 618fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 619fc4362bfSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 620fc4362bfSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 621fc4362bfSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 6221fddbf69SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s max eigen=%e min=%e PC=%s\n",__FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER); 623fc4362bfSMark F. Adams } 624038e3b61SMark F. Adams { 625038e3b61SMark F. Adams PetscInt N1, N0, tt; 626038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 627038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 628785cba28SMark F. Adams emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */ 629038e3b61SMark F. Adams emax *= 1.05; 6301fddbf69SMark F. Adams 631038e3b61SMark F. Adams } 632038e3b61SMark F. Adams 633587fa25dSMark F. Adams ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN ); 634fc4362bfSMark F. Adams ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 6350e1b4bd6SMark F. Adams /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */ 6365745e0f5SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 637e33ef3b1SMark F. Adams ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); 6385745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 6395745e0f5SMark F. Adams } 6405745e0f5SMark F. Adams { 6415745e0f5SMark F. Adams KSP smoother; /* coarse grid */ 6425745e0f5SMark F. Adams Mat Lmat = Aarr[pc_gamg->m_Nlevels-1]; 643eb07cef2SMark F. Adams ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 6445745e0f5SMark F. Adams ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); 6455745e0f5SMark F. Adams CHKERRQ(ierr); 6465745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 647fc4362bfSMark F. Adams } 648737a81a9SMark F. Adams 649fc4362bfSMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 650eb07cef2SMark F. Adams ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr); 6515b89ad90SMark F. Adams { 6525b89ad90SMark F. Adams PetscBool galerkin; 653eb07cef2SMark F. Adams ierr = PCMGGetGalerkin( a_pc, &galerkin); CHKERRQ(ierr); 6545b89ad90SMark F. Adams if(galerkin){ 6555b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 6565b89ad90SMark F. Adams } 6575b89ad90SMark F. Adams } 6585745e0f5SMark F. Adams 659b4fbaa2aSMark F. Adams /* set interpolation between the levels, clean up */ 660587fa25dSMark F. Adams for (lidx=0,level=pc_gamg->m_Nlevels-1; 661587fa25dSMark F. Adams lidx<fine_level; 662587fa25dSMark F. Adams lidx++, level--){ 663587fa25dSMark F. Adams ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr); 6646c237d78SBarry Smith if( !PETSC_TRUE ) { 66511e60469SMark F. Adams PetscViewer viewer; char fname[32]; 666*fd3c6afaSMark F. Adams sprintf(fname,"Pmat_%d.m",(int)level); 66711e60469SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 6685b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 669038e3b61SMark F. Adams ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr); 6705b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 671*fd3c6afaSMark F. Adams sprintf(fname,"Amat_%d.m",(int)level); 672e33ef3b1SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 673e33ef3b1SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 674e33ef3b1SMark F. Adams ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr); 675e33ef3b1SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 6765b89ad90SMark F. Adams } 677587fa25dSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 678587fa25dSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 6795b89ad90SMark F. Adams } 6805745e0f5SMark F. Adams 6815b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 682eb07cef2SMark F. Adams a_pc->setupcalled = 0; 683eb07cef2SMark F. Adams ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr); 684e33ef3b1SMark F. Adams 6855b89ad90SMark F. Adams PetscFunctionReturn(0); 6865b89ad90SMark F. Adams } 6875b89ad90SMark F. Adams 688eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 6895b89ad90SMark F. Adams /* 6905b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 6915b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 6925b89ad90SMark F. Adams 6935b89ad90SMark F. Adams Input Parameter: 6945b89ad90SMark F. Adams . pc - the preconditioner context 6955b89ad90SMark F. Adams 6965b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 6975b89ad90SMark F. Adams */ 6985b89ad90SMark F. Adams #undef __FUNCT__ 6995b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 7005b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 7015b89ad90SMark F. Adams { 7025b89ad90SMark F. Adams PetscErrorCode ierr; 7035b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7045b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 7055b89ad90SMark F. Adams 7065b89ad90SMark F. Adams PetscFunctionBegin; 7075b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 7085b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 7095b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 7105b89ad90SMark F. Adams PetscFunctionReturn(0); 7115b89ad90SMark F. Adams } 7125b89ad90SMark F. Adams 7135b89ad90SMark F. Adams #undef __FUNCT__ 7145b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 7155b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 7165b89ad90SMark F. Adams { 7175b89ad90SMark F. Adams /* PetscErrorCode ierr; */ 7185b89ad90SMark F. Adams /* PC_MG *mg = (PC_MG*)pc->data; */ 7195b89ad90SMark F. Adams /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 7205b89ad90SMark F. Adams /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 7215b89ad90SMark F. Adams 7225b89ad90SMark F. Adams PetscFunctionBegin; 7235b89ad90SMark F. Adams PetscFunctionReturn(0); 7245b89ad90SMark F. Adams } 7255b89ad90SMark F. Adams 7265b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 7275b89ad90SMark F. Adams /* 7285b89ad90SMark F. Adams PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 7295b89ad90SMark F. Adams 7305b89ad90SMark F. Adams Input Parameter: 7315b89ad90SMark F. Adams . pc - the preconditioner context 7325b89ad90SMark F. Adams 7335b89ad90SMark F. Adams Application Interface Routine: PCCreate() 7345b89ad90SMark F. Adams 7355b89ad90SMark F. Adams */ 7365b89ad90SMark F. Adams /* MC 7375b89ad90SMark F. Adams PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 7385b89ad90SMark F. Adams fine grid discretization matrix and coordinates on the fine grid. 7395b89ad90SMark F. Adams 7405b89ad90SMark F. Adams Options Database Key: 7415b89ad90SMark F. Adams Multigrid options(inherited) 7425b89ad90SMark F. Adams + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 7435b89ad90SMark F. Adams . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 7445b89ad90SMark F. Adams . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 7455b89ad90SMark F. Adams -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 7465b89ad90SMark F. Adams GAMG options: 7475b89ad90SMark F. Adams 7485b89ad90SMark F. Adams Level: intermediate 7495b89ad90SMark F. Adams Concepts: multigrid 7505b89ad90SMark F. Adams 7515b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 7525b89ad90SMark F. Adams PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 7535b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 7545b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 7555b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 7565b89ad90SMark F. Adams M */ 7575b89ad90SMark F. Adams 7585b89ad90SMark F. Adams EXTERN_C_BEGIN 7595b89ad90SMark F. Adams #undef __FUNCT__ 7605b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 7615b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 7625b89ad90SMark F. Adams { 7635b89ad90SMark F. Adams PetscErrorCode ierr; 7645b89ad90SMark F. Adams PC_GAMG *pc_gamg; 7655b89ad90SMark F. Adams PC_MG *mg; 7665ef31b24SMark F. Adams PetscClassId cookie; 7675b89ad90SMark F. Adams 7685b89ad90SMark F. Adams PetscFunctionBegin; 7695b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 7705b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 7715b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 7725b89ad90SMark F. Adams 7735b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 7745b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 7755b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 7765b89ad90SMark F. Adams mg->innerctx = pc_gamg; 7775b89ad90SMark F. Adams 7785b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 7795b89ad90SMark F. Adams 7805b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 7815b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 7825b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 7835b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 7845b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 7855b89ad90SMark F. Adams 7865b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 7875b89ad90SMark F. Adams "PCSetCoordinates_C", 7885b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 7895b89ad90SMark F. Adams PCSetCoordinates_GAMG);CHKERRQ(ierr); 790b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 791785cba28SMark F. Adams static int count = 0; 792785cba28SMark F. Adams if( count++ == 0 ) { 7935ef31b24SMark F. Adams PetscClassIdRegister("GAMG Setup",&cookie); 794b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); 795b4fbaa2aSMark F. Adams PetscLogEventRegister(" make graph", cookie, &gamg_setup_events[SET3]); 796b4fbaa2aSMark F. Adams PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); 797b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); 798b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); 799b4fbaa2aSMark F. Adams PetscLogEventRegister(" search & set", cookie, &gamg_setup_events[FIND_V]); 800b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: init", cookie, &gamg_setup_events[SET7]); 801b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); */ 802b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); 803b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); 804b4fbaa2aSMark F. Adams PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_events[SET12]); 805b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ 806b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ 807b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ 8085ef31b24SMark F. Adams 809b4fbaa2aSMark F. Adams /* create timer stages */ 810b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 811b4fbaa2aSMark F. Adams { 812b4fbaa2aSMark F. Adams char str[32]; 813b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d (finest)",0); 814b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[0]); 815b4fbaa2aSMark F. Adams PetscInt lidx; 816b4fbaa2aSMark F. Adams for (lidx=1;lidx<9;lidx++){ 817b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 818b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[lidx]); 819b4fbaa2aSMark F. Adams } 820b4fbaa2aSMark F. Adams } 821b4fbaa2aSMark F. Adams #endif 822b4fbaa2aSMark F. Adams } 823b4fbaa2aSMark F. Adams #endif 8245b89ad90SMark F. Adams PetscFunctionReturn(0); 8255b89ad90SMark F. Adams } 8265b89ad90SMark F. Adams EXTERN_C_END 827