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 6f96513f1SMatthew G Knepley PetscLogEvent gamg_setup_stages[NUM_SET]; 7f96513f1SMatthew G Knepley 85b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */ 95b89ad90SMark F. Adams typedef struct gamg_TAG { 105b89ad90SMark F. Adams PetscInt m_dim; 115b89ad90SMark F. Adams PetscInt m_Nlevels; 125b89ad90SMark F. Adams PetscInt m_data_sz; 135b89ad90SMark F. Adams PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 145b89ad90SMark F. Adams } PC_GAMG; 155b89ad90SMark F. Adams 16fc4362bfSMark F. Adams #define TOP_GRID_LIM 100 175ef31b24SMark F. Adams 185b89ad90SMark F. Adams /* -----------------------------------------------------------------------------*/ 195b89ad90SMark F. Adams #undef __FUNCT__ 205b89ad90SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 215b89ad90SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 225b89ad90SMark F. Adams { 235b89ad90SMark F. Adams PetscErrorCode ierr; 245b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 255b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 265b89ad90SMark F. Adams 275b89ad90SMark F. Adams PetscFunctionBegin; 285b89ad90SMark F. Adams ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 29eb07cef2SMark F. Adams pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0; 305b89ad90SMark F. Adams PetscFunctionReturn(0); 315b89ad90SMark F. Adams } 325b89ad90SMark F. Adams 335b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 345b89ad90SMark F. Adams /* 355b89ad90SMark F. Adams PCSetCoordinates_GAMG 365b89ad90SMark F. Adams 375b89ad90SMark F. Adams Input Parameter: 385b89ad90SMark F. Adams . pc - the preconditioner context 395b89ad90SMark F. Adams */ 40a92563c5SMark F. Adams EXTERN_C_BEGIN 415b89ad90SMark F. Adams #undef __FUNCT__ 425b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 43eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords ) 445b89ad90SMark F. Adams { 45eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 465b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 476c237d78SBarry Smith PetscErrorCode ierr; 48038e3b61SMark F. Adams PetscInt arrsz,bs,my0,kk,ii,jj,nloc,sz,M; 49038e3b61SMark F. Adams Mat Amat = a_pc->pmat; 50eb07cef2SMark F. Adams PetscBool useSA = PETSC_FALSE, flag; 51eb07cef2SMark F. Adams char str[16]; 525b89ad90SMark F. Adams 535b89ad90SMark F. Adams PetscFunctionBegin; 54038e3b61SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 55038e3b61SMark F. Adams ierr = MatGetOwnershipRange( Amat, &my0, &kk ); CHKERRQ(ierr); 56038e3b61SMark F. Adams nloc = (kk-my0)/bs; M = kk - my0; 57038e3b61SMark F. Adams if((kk-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 58038e3b61SMark F. Adams 59eb07cef2SMark F. Adams ierr = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,16,&flag); CHKERRQ( ierr ); 60eb07cef2SMark F. Adams useSA = (PetscBool)(flag && strcmp(str,"sa") == 0); 61038e3b61SMark F. Adams 62038e3b61SMark F. Adams if(a_coords == 0) useSA = PETSC_TRUE; /* use SA if no data */ 63038e3b61SMark F. Adams if( !useSA ) sz = a_ndm; /* coordinates */ 64038e3b61SMark F. Adams else{ /* SA: null space vectors */ 65038e3b61SMark F. Adams if(a_coords != 0) sz = (a_ndm==2 ? 3*bs : 6*bs); /* elasticity */ 66038e3b61SMark F. Adams else sz = bs*bs; /* no data, force SA with constant null space vectors */ 67038e3b61SMark F. Adams if(a_coords != 0 && bs==1 ) sz = 1; /* scalar problem with coords and SA (not needed) */ 68038e3b61SMark F. Adams } 69eb07cef2SMark F. Adams arrsz = nloc*sz; 705b89ad90SMark F. Adams 71038e3b61SMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 726c237d78SBarry Smith if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 735b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 74eb07cef2SMark F. Adams ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 755b89ad90SMark F. Adams } 76038e3b61SMark F. Adams for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.; 77eb07cef2SMark F. Adams pc_gamg->m_data[arrsz] = -99.; 78038e3b61SMark F. Adams /* copy data in - column oriented */ 79eb07cef2SMark F. Adams if( useSA ) { 80038e3b61SMark F. Adams for(kk=0;kk<nloc;kk++){ 81038e3b61SMark F. Adams PetscReal *data = &pc_gamg->m_data[kk*bs]; 82038e3b61SMark F. Adams if( sz==1 ) *data = 1.0; 83038e3b61SMark F. Adams else { 84038e3b61SMark F. Adams for(ii=0;ii<bs;ii++) 85038e3b61SMark F. Adams for(jj=0;jj<bs;jj++) 86038e3b61SMark F. Adams if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */ 87038e3b61SMark F. Adams else data[ii*M + jj] = 0.0; 88eb07cef2SMark F. Adams if( a_coords != 0 ) { 89038e3b61SMark F. Adams if( a_ndm == 2 ){ /* rotational modes */ 90038e3b61SMark F. Adams data += 2*M; 91038e3b61SMark F. Adams data[0] = -a_coords[2*kk+1]; 92038e3b61SMark F. Adams data[1] = a_coords[2*kk]; 935b89ad90SMark F. Adams } 94eb07cef2SMark F. Adams else { 95038e3b61SMark F. Adams data += 3*M; 96038e3b61SMark F. Adams data[0] = 0.0; data[M+0] = a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1]; 97038e3b61SMark F. Adams data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = a_coords[3*kk]; 98038e3b61SMark F. Adams data[2] = a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk]; data[2*M+2] = 0.0; 99038e3b61SMark F. Adams } 100eb07cef2SMark F. Adams } 101eb07cef2SMark F. Adams } 102eb07cef2SMark F. Adams } 103eb07cef2SMark F. Adams } 104eb07cef2SMark F. Adams else { 105038e3b61SMark F. Adams for( kk = 0 ; kk < nloc ; kk++ ){ 106038e3b61SMark F. Adams for( ii = 0 ; ii < a_ndm ; ii++ ) { 107038e3b61SMark F. Adams pc_gamg->m_data[ii*nloc + kk] = a_coords[kk*a_ndm + ii]; 108eb07cef2SMark F. Adams } 109eb07cef2SMark F. Adams } 110038e3b61SMark F. Adams } 111038e3b61SMark F. Adams 112eb07cef2SMark F. Adams assert(pc_gamg->m_data[arrsz] == -99.); 113038e3b61SMark F. Adams for(kk=0;kk<arrsz;kk++){ 114038e3b61SMark F. Adams assert(pc_gamg->m_data[kk] != -999.); // debug 115038e3b61SMark F. Adams } 116038e3b61SMark F. Adams 1175b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 118eb07cef2SMark F. Adams pc_gamg->m_dim = a_ndm; 1195b89ad90SMark F. Adams 1205b89ad90SMark F. Adams PetscFunctionReturn(0); 1215b89ad90SMark F. Adams } 122a92563c5SMark F. Adams EXTERN_C_END 1235b89ad90SMark F. Adams 1245b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1255b89ad90SMark F. Adams /* 12611e60469SMark F. Adams partitionLevel 1275b89ad90SMark F. Adams 1285b89ad90SMark F. Adams Input Parameter: 1293530afc2SMark F. Adams . a_Amat_fine - matrix on this fine (k) level 130038e3b61SMark F. Adams . a_data_sz - size of data to move (coarse grid) 1313530afc2SMark F. Adams In/Output Parameter: 1323530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 133eb07cef2SMark F. Adams . a_coarse_data - data that need to be moved 1343530afc2SMark F. Adams . a_active_proc - number of active procs 13511e60469SMark F. Adams Output Parameter: 1363530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 1375b89ad90SMark F. Adams */ 1385cb416c2SMark F. Adams 1395b89ad90SMark F. Adams #undef __FUNCT__ 14011e60469SMark F. Adams #define __FUNCT__ "partitionLevel" 1413530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine, 142eb07cef2SMark F. Adams PetscInt a_data_sz, 143038e3b61SMark F. Adams PetscInt a_cbs, 1443530afc2SMark F. Adams Mat *a_P_inout, 145eb07cef2SMark F. Adams PetscReal **a_coarse_data, 1463530afc2SMark F. Adams PetscMPIInt *a_active_proc, 1473530afc2SMark F. Adams Mat *a_Amat_crs 14811e60469SMark F. Adams ) 1495b89ad90SMark F. Adams { 1505b89ad90SMark F. Adams PetscErrorCode ierr; 151038e3b61SMark F. Adams Mat Cmat,Pnew,Pold=*a_P_inout; 15211e60469SMark F. Adams IS new_indices,isnum; 1533530afc2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 1545ef31b24SMark F. Adams PetscMPIInt nactive_procs,mype,npe; 155038e3b61SMark F. Adams PetscInt Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,fbs; 156*a6828334SMark F. Adams PetscInt neq,NN; 157e33ef3b1SMark F. Adams PetscMPIInt new_npe,targ_npe; 1585b89ad90SMark F. Adams 1595b89ad90SMark F. Adams PetscFunctionBegin; 16011e60469SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 16111e60469SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 162038e3b61SMark F. Adams ierr = MatGetBlockSize( a_Amat_fine, &fbs ); CHKERRQ(ierr); 16311e60469SMark F. Adams /* RAP */ 164038e3b61SMark F. Adams ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 165038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 166acadaa72SMark F. Adams 167038e3b61SMark F. Adams ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); 168038e3b61SMark F. Adams ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0); 169038e3b61SMark F. Adams 170038e3b61SMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 171038e3b61SMark F. Adams ierr = MatGetSize( Cmat, &neq, &NN );CHKERRQ(ierr); 1723530afc2SMark F. Adams #define MIN_EQ_PROC 100 1735ef31b24SMark F. Adams nactive_procs = *a_active_proc; 1743530afc2SMark F. Adams targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */ 1755ef31b24SMark F. Adams if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */ 1765ef31b24SMark F. Adams else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */ 1773530afc2SMark F. Adams else { 1785ef31b24SMark F. Adams PetscMPIInt factstart,fact; 1793530afc2SMark F. Adams new_npe = -9999; 1805ef31b24SMark F. Adams factstart = nactive_procs; 1813530afc2SMark F. Adams for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */ 1825ef31b24SMark F. Adams if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) { 1835ef31b24SMark F. Adams new_npe = nactive_procs/fact; 1843530afc2SMark F. Adams } 1853530afc2SMark F. Adams } 186038e3b61SMark F. Adams if(new_npe == -9999) new_npe = nactive_procs; 1873530afc2SMark F. Adams } 188e33ef3b1SMark F. Adams 189e33ef3b1SMark F. Adams if( PETSC_TRUE ) { /* partition: get 'isnewproc' */ 19011e60469SMark F. Adams MatPartitioning mpart; 1915ef31b24SMark F. Adams Mat adj; 1925ef31b24SMark F. Adams const PetscInt *is_idx; 193*a6828334SMark F. Adams PetscInt is_sz,kk,jj,ii,old_fact=(npe/nactive_procs),new_fact=(npe/new_npe),*isnewproc_idx,counts[npe]; 194c9a0b8beSMark F. Adams /* create sub communicator */ 195c9a0b8beSMark F. Adams MPI_Comm cm,new_comm; 196*a6828334SMark F. Adams PetscInt membershipKey = mype % old_fact; 197e33ef3b1SMark F. Adams IS isnewproc; 198e33ef3b1SMark F. Adams 199e33ef3b1SMark F. Adams *a_active_proc = new_npe; /* output for next time */ 200c9a0b8beSMark F. Adams ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr); 201c9a0b8beSMark F. Adams ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr); 202c9a0b8beSMark F. Adams ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 203c9a0b8beSMark F. Adams 2045ef31b24SMark F. Adams /* MatPartitioningApply call MatConvert, which is collective */ 205038e3b61SMark F. Adams if( a_cbs==1) { 206038e3b61SMark F. Adams ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 207eb07cef2SMark F. Adams } 208eb07cef2SMark F. Adams else { 209038e3b61SMark F. Adams /* make a scalar matrix to partition */ 210eb07cef2SMark F. Adams Mat tMat; 211eb07cef2SMark F. Adams PetscInt Ii,ncols; const PetscScalar *vals; const PetscInt *idx; 212eb07cef2SMark F. Adams ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 213eb07cef2SMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 214eb07cef2SMark F. Adams 25, PETSC_NULL, 10, PETSC_NULL, 215eb07cef2SMark F. Adams &tMat ); 216eb07cef2SMark F. Adams 217eb07cef2SMark F. Adams for ( Ii = Istart0; Ii < Iend0; Ii++ ) { 218038e3b61SMark F. Adams PetscInt dest_row = Ii/a_cbs; 219038e3b61SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr); 220eb07cef2SMark F. Adams for( jj = 0 ; jj < ncols ; jj++ ){ 221038e3b61SMark F. Adams PetscInt dest_col = idx[jj]/a_cbs; 222eb07cef2SMark F. Adams PetscScalar v = 1.0; 223eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 224eb07cef2SMark F. Adams } 225038e3b61SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr); 226eb07cef2SMark F. Adams } 227eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 228eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 229eb07cef2SMark F. Adams 230eb07cef2SMark F. Adams ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 231eb07cef2SMark F. Adams 232eb07cef2SMark F. Adams ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 233eb07cef2SMark F. Adams } 2345ef31b24SMark F. Adams if( membershipKey == 0 ) { 235038e3b61SMark F. Adams if(ncrs0==0)SETERRQ(wcomm,PETSC_ERR_LIB,"zero local nodes -- increase min"); 2365ef31b24SMark F. Adams /* hack to fix global data that pmetis.c uses in 'adj' */ 2375ef31b24SMark F. Adams for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) { 2385ef31b24SMark F. Adams adj->rmap->range[jj] = adj->rmap->range[kk]; 2395ef31b24SMark F. Adams } 2405ef31b24SMark F. Adams ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr); 2415ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 24211e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 2435ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr); 24411e60469SMark F. Adams ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 24511e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 2465ef31b24SMark F. Adams /* collect IS info */ 2475ef31b24SMark F. Adams ierr = ISGetLocalSize( isnewproc, &is_sz ); CHKERRQ(ierr); 248038e3b61SMark F. Adams ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr); 2495ef31b24SMark F. Adams ierr = ISGetIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 250eb07cef2SMark F. Adams /* spread partitioning across machine - probably the right thing to do but machine spec. */ 251eb07cef2SMark F. Adams for(kk=0,jj=0;kk<is_sz;kk++){ 252038e3b61SMark F. Adams for(ii=0 ; ii<a_cbs ; ii++, jj++ ) { /* expand for equation level by 'bs' */ 253eb07cef2SMark F. Adams isnewproc_idx[jj] = is_idx[kk] * new_fact; 254eb07cef2SMark F. Adams } 2555ef31b24SMark F. Adams } 2565ef31b24SMark F. Adams ierr = ISRestoreIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 2575ef31b24SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 258038e3b61SMark F. Adams is_sz *= a_cbs; 2595ef31b24SMark F. Adams } 2605ef31b24SMark F. Adams else { 2615ef31b24SMark F. Adams isnewproc_idx = 0; 2625ef31b24SMark F. Adams is_sz = 0; 2635ef31b24SMark F. Adams } 2645ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 2655ef31b24SMark F. Adams ierr = MPI_Comm_free( &new_comm ); CHKERRQ(ierr); 2665ef31b24SMark F. Adams ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc ); 2675ef31b24SMark F. Adams if( membershipKey == 0 ) { 2685ef31b24SMark F. Adams ierr = PetscFree( isnewproc_idx ); CHKERRQ(ierr); 2695ef31b24SMark F. Adams } 270e33ef3b1SMark F. Adams 27111e60469SMark F. Adams /* 27211e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 27311e60469SMark F. Adams */ 27411e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 27511e60469SMark F. Adams /* 27611e60469SMark F. Adams Determine how many elements are assigned to each processor 27711e60469SMark F. Adams */ 27811e60469SMark F. Adams ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr); 27911e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 280038e3b61SMark F. Adams ncrs_new = counts[mype]/a_cbs; 2815ef31b24SMark F. Adams 28211e60469SMark F. Adams { /* Create a vector to contain the newly ordered element information */ 28311e60469SMark F. Adams const PetscInt *idx; 284038e3b61SMark F. Adams PetscInt ii,jj,kk; 28511e60469SMark F. Adams IS isscat; 28611e60469SMark F. Adams PetscScalar *array; 28711e60469SMark F. Adams Vec src_crd, dest_crd; 288eb07cef2SMark F. Adams PetscReal *data = *a_coarse_data; 28911e60469SMark F. Adams VecScatter vecscat; 290eb07cef2SMark F. Adams PetscInt tidx[ncrs0*a_data_sz]; 291acadaa72SMark F. Adams 29211e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 293eb07cef2SMark F. Adams ierr = VecSetSizes( dest_crd, a_data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 29411e60469SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/ 29511e60469SMark F. Adams /* 296eb07cef2SMark F. Adams There are 'a_data_sz' data items per node, (one can think of the vectors of having a 297038e3b61SMark F. Adams block size of 'a_data_sz'). Note, ISs are expanded into equation space by 'a_cbs'. 29811e60469SMark F. Adams */ 29911e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 300038e3b61SMark F. Adams for(ii=0,jj=0; ii<ncrs0 ; ii++) { 301038e3b61SMark F. Adams PetscInt id = idx[ii*a_cbs]/a_cbs; 302038e3b61SMark F. Adams for( kk=0; kk<a_data_sz ; kk++, jj++) tidx[jj] = id*a_data_sz + kk; 30311e60469SMark F. Adams } 304038e3b61SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 305eb07cef2SMark F. Adams ierr = ISCreateGeneral( wcomm, a_data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 30611e60469SMark F. Adams CHKERRQ(ierr); 30711e60469SMark F. Adams /* 30811e60469SMark F. Adams Create a vector to contain the original vertex information for each element 30911e60469SMark F. Adams */ 310eb07cef2SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, a_data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 311038e3b61SMark F. Adams for(ii=0, jj=0; ii<ncrs0 ; ii++) { 312038e3b61SMark F. Adams for( kk=0; kk < a_data_sz ; kk++, jj++) { 313038e3b61SMark F. Adams ierr = VecSetValues( src_crd, 1, &jj, &data[kk*ncrs0 + ii], INSERT_VALUES ); CHKERRQ(ierr); 314038e3b61SMark F. Adams } 315eb07cef2SMark F. Adams } 316eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 317eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 31811e60469SMark F. Adams /* 31911e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 32011e60469SMark F. Adams to the correct processor 32111e60469SMark F. Adams */ 32211e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 32311e60469SMark F. Adams CHKERRQ(ierr); 32411e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 32511e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 32611e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 32711e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 32811e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 32911e60469SMark F. Adams /* 33011e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 33111e60469SMark F. Adams */ 332eb07cef2SMark F. Adams ierr = PetscFree( *a_coarse_data ); CHKERRQ(ierr); 333eb07cef2SMark F. Adams ierr = PetscMalloc( a_data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data ); CHKERRQ(ierr); 334038e3b61SMark F. Adams VecGetLocalSize( dest_crd, &ii ); assert(ii==a_data_sz*ncrs_new); 335eb07cef2SMark F. Adams 33611e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 337eb07cef2SMark F. Adams data = *a_coarse_data; 338038e3b61SMark F. Adams for(ii=0, jj=0; ii<ncrs_new ; ii++) { 339038e3b61SMark F. Adams for( kk=0; kk < a_data_sz ; kk++, jj++) { 340038e3b61SMark F. Adams data[kk*ncrs_new + ii] = PetscRealPart(array[jj]); 341038e3b61SMark F. Adams } 342038e3b61SMark F. Adams } 34311e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 34411e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 34511e60469SMark F. Adams } 34611e60469SMark F. Adams /* 34711e60469SMark F. Adams Invert for MatGetSubMatrix 34811e60469SMark F. Adams */ 349038e3b61SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr); 35011e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 35111e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 35211e60469SMark F. Adams /* A_crs output */ 353038e3b61SMark F. Adams ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 35411e60469SMark F. Adams CHKERRQ(ierr); 355eb07cef2SMark F. Adams 356038e3b61SMark F. Adams ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 357e33ef3b1SMark F. Adams Cmat = *a_Amat_crs; /* output */ 358038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 359eb07cef2SMark F. Adams 36011e60469SMark F. Adams /* prolongator */ 36111e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 36211e60469SMark F. Adams { 36311e60469SMark F. Adams IS findices; 36411e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 36511e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 36611e60469SMark F. Adams CHKERRQ(ierr); 36711e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 36811e60469SMark F. Adams } 3693530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 3703530afc2SMark F. Adams *a_P_inout = Pnew; /* output */ 37111e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 372e33ef3b1SMark F. Adams } 373e33ef3b1SMark F. Adams else { 374e33ef3b1SMark F. Adams *a_Amat_crs = Cmat; /* output */ 375e33ef3b1SMark F. Adams } 3765b89ad90SMark F. Adams 3775b89ad90SMark F. Adams PetscFunctionReturn(0); 3785b89ad90SMark F. Adams } 3795b89ad90SMark F. Adams 380e33ef3b1SMark F. Adams #define GAMG_MAXLEVELS 30 3815ef31b24SMark F. Adams #if defined(PETSC_USE_LOG) 3825ef31b24SMark F. Adams PetscLogStage gamg_stages[20]; 3835ef31b24SMark F. Adams #endif 3845b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 3855b89ad90SMark F. Adams /* 3865b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 3875b89ad90SMark F. Adams by setting data structures and options. 3885b89ad90SMark F. Adams 3895b89ad90SMark F. Adams Input Parameter: 3905b89ad90SMark F. Adams . pc - the preconditioner context 3915b89ad90SMark F. Adams 3925b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 3935b89ad90SMark F. Adams 3945b89ad90SMark F. Adams Notes: 3955b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 3965b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 3975b89ad90SMark F. Adams */ 3985b89ad90SMark F. Adams #undef __FUNCT__ 3995b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 400eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc ) 4015b89ad90SMark F. Adams { 4025b89ad90SMark F. Adams PetscErrorCode ierr; 403eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 4045b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 405eb07cef2SMark F. Adams Mat Amat = a_pc->mat, Pmat = a_pc->pmat; 406eb07cef2SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, data_sz, Istart, Iend; 407eb07cef2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_pc)->comm; 4083530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 409038e3b61SMark F. Adams PetscBool isOK; 410587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 411587fa25dSMark F. Adams PetscReal *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS]; 4125ef31b24SMark F. Adams 4135b89ad90SMark F. Adams PetscFunctionBegin; 414eb07cef2SMark F. Adams if( a_pc->setupcalled ) { 415eb07cef2SMark F. Adams /* no state data in GAMG to destroy */ 416eb07cef2SMark F. Adams ierr = PCReset_MG( a_pc ); CHKERRQ(ierr); 417eb07cef2SMark F. Adams } 418baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 419baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 4205b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 4215b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 4223530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 423eb07cef2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 424eb07cef2SMark F. Adams nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 425eb07cef2SMark F. Adams 426038e3b61SMark F. Adams /* get data of not around */ 427038e3b61SMark F. Adams if( pc_gamg->m_data == 0 && nloc > 0 ) { 428038e3b61SMark F. Adams ierr = PCSetCoordinates_GAMG( a_pc, -1, 0 ); CHKERRQ( ierr ); 429eb07cef2SMark F. Adams } 430eb07cef2SMark F. Adams data = pc_gamg->m_data; 431038e3b61SMark F. Adams data_sz = pc_gamg->m_data_sz/nloc; assert(data!=0 || data_sz==0); 432038e3b61SMark F. Adams 433eb07cef2SMark F. Adams /* Get A_i and R_i */ 434e33ef3b1SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, data size=%d m_data_sz=%d nloc=%d, np=%d\n",0,__FUNCT__,0,N,data_sz,pc_gamg->m_data_sz,nloc,npe); 4358f4b7eb5SMark F. Adams for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 436eb07cef2SMark F. Adams level < GAMG_MAXLEVELS-1 && (level==0 || M/bs>TOP_GRID_LIM) && (npe==1 || nactivepe>1); 4370205a208SMark F. Adams level++ ) { 4385b89ad90SMark F. Adams level1 = level + 1; 4395ef31b24SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr); 440038e3b61SMark F. Adams ierr = createProlongation( Aarr[level], data, pc_gamg->m_dim, &bs, &data_sz, 441587fa25dSMark F. Adams &Parr[level1], &coarse_data, &isOK, &emaxs[level] ); 4425b89ad90SMark F. Adams CHKERRQ(ierr); 4435ef31b24SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr); 4445ef31b24SMark F. Adams 445eb07cef2SMark F. Adams ierr = PetscFree( data ); CHKERRQ( ierr ); 446baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 447baa4e9faSMark F. Adams if( isOK ) { 4485ef31b24SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr); 449038e3b61SMark F. Adams ierr = partitionLevel( Aarr[level], data_sz, bs, 450eb07cef2SMark F. Adams &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] ); 4513530afc2SMark F. Adams CHKERRQ(ierr); 4525ef31b24SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr); 4533530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 454eb07cef2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, %d data size, %d active pes\n",0,__FUNCT__,level1,N,data_sz,nactivepe); 455e33ef3b1SMark F. Adams /* coarse grids with SA can have zero row/cols from singleton aggregates */ 456e33ef3b1SMark F. Adams /* aggregation method can probably gaurrentee this does not happen! - be safe for now */ 457e33ef3b1SMark F. Adams if( PETSC_TRUE ){ 458e33ef3b1SMark F. Adams Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloc,id; 459e33ef3b1SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); 460e33ef3b1SMark F. Adams nloc = Iend-Istart; 461e33ef3b1SMark F. Adams ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); 462e33ef3b1SMark F. Adams ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); 463e33ef3b1SMark F. Adams ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); 464e33ef3b1SMark F. Adams for(kk=0;kk<nloc;kk++){ 465e33ef3b1SMark F. Adams if(data_arr[kk]==0.0) { 466e33ef3b1SMark F. Adams id = kk + Istart; v = 1.0; 467e33ef3b1SMark F. Adams ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); 468e33ef3b1SMark F. Adams CHKERRQ(ierr); 469e33ef3b1SMark F. Adams PetscPrintf(PETSC_COMM_SELF,"[%d]%s warning: added diag to zero (%d)\n",mype,__FUNCT__,id); 470e33ef3b1SMark F. Adams } 471e33ef3b1SMark F. Adams } 472e33ef3b1SMark F. Adams ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); 473e33ef3b1SMark F. Adams ierr = VecDestroy( &diag ); CHKERRQ(ierr); 474e33ef3b1SMark F. Adams ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 475e33ef3b1SMark F. Adams ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 476e33ef3b1SMark F. Adams } 477baa4e9faSMark F. Adams } 478baa4e9faSMark F. Adams else{ 479baa4e9faSMark F. Adams break; 480baa4e9faSMark F. Adams } 481eb07cef2SMark F. Adams data = coarse_data; 4825b89ad90SMark F. Adams } 483eb07cef2SMark F. Adams ierr = PetscFree( coarse_data ); CHKERRQ( ierr ); 484baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 4855b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 4865b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 4875b89ad90SMark F. Adams fine_level = level; 488eb07cef2SMark F. Adams ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 4895b89ad90SMark F. Adams 490fc4362bfSMark F. Adams /* set default smoothers */ 491587fa25dSMark F. Adams for ( lidx=1, level = pc_gamg->m_Nlevels-2; 492587fa25dSMark F. Adams lidx <= fine_level; 493587fa25dSMark F. Adams lidx++, level--) { 4945745e0f5SMark F. Adams PetscReal emax, emin; 4955b89ad90SMark F. Adams KSP smoother; PC subpc; 496587fa25dSMark F. Adams ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 4975b89ad90SMark F. Adams ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 498038e3b61SMark F. Adams if( emaxs[level] > 0.0 ) emax = emaxs[level]; 499587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 500587fa25dSMark F. Adams KSP eksp; Mat Lmat = Aarr[level]; 501fc4362bfSMark F. Adams Vec bb, xx; PC pc; 502038e3b61SMark F. Adams 5035745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 5045745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 505fc4362bfSMark F. Adams { 506fc4362bfSMark F. Adams PetscRandom rctx; 507fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 508fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 509fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 510fc4362bfSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 5115b89ad90SMark F. Adams } 512fc4362bfSMark F. Adams ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 513e33ef3b1SMark F. Adams ierr = KSPSetType( eksp, KSPCG ); CHKERRQ(ierr); 514fc4362bfSMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 5155745e0f5SMark F. Adams ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr ); 516fc4362bfSMark F. Adams ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 517e33ef3b1SMark F. Adams ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as above */ 518038e3b61SMark F. Adams ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 519fc4362bfSMark F. Adams CHKERRQ(ierr); 520e33ef3b1SMark F. Adams //ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr); 521fc4362bfSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 522fc4362bfSMark F. Adams 523fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 524fc4362bfSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 525fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 526fc4362bfSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 527fc4362bfSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 528fc4362bfSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 529e33ef3b1SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"%s max eigen = %e min = %e\n",__FUNCT__,emax,emin); 530fc4362bfSMark F. Adams } 531038e3b61SMark F. Adams { 532038e3b61SMark F. Adams PetscInt N1, N0, tt; 533038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 534038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 535e33ef3b1SMark F. Adams emin = .5*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */ 536038e3b61SMark F. Adams emax *= 1.05; 537e33ef3b1SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t%s eigen reduction factor = %g\n",__FUNCT__,emax/emin); 538038e3b61SMark F. Adams } 539038e3b61SMark F. Adams 540587fa25dSMark F. Adams ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN ); 541fc4362bfSMark F. Adams ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 542e33ef3b1SMark F. Adams ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); 5435745e0f5SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 544e33ef3b1SMark F. Adams ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); 5455745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 5465745e0f5SMark F. Adams } 5475745e0f5SMark F. Adams { 5485745e0f5SMark F. Adams KSP smoother; /* coarse grid */ 5495745e0f5SMark F. Adams Mat Lmat = Aarr[pc_gamg->m_Nlevels-1]; 550eb07cef2SMark F. Adams ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 5515745e0f5SMark F. Adams ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); 5525745e0f5SMark F. Adams CHKERRQ(ierr); 5535745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 554fc4362bfSMark F. Adams } 555fc4362bfSMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 556eb07cef2SMark F. Adams ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr); 5575b89ad90SMark F. Adams { 5585b89ad90SMark F. Adams PetscBool galerkin; 559eb07cef2SMark F. Adams ierr = PCMGGetGalerkin( a_pc, &galerkin); CHKERRQ(ierr); 5605b89ad90SMark F. Adams if(galerkin){ 5615b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 5625b89ad90SMark F. Adams } 5635b89ad90SMark F. Adams } 5645745e0f5SMark F. Adams 5655745e0f5SMark F. Adams /* set interpolation between the levels, create timer stages, clean up */ 5668f4b7eb5SMark F. Adams if( PETSC_FALSE ) { 5675ef31b24SMark F. Adams char str[32]; 5685ef31b24SMark F. Adams sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1); 5695ef31b24SMark F. Adams PetscLogStageRegister(str, &gamg_stages[fine_level]); 5705ef31b24SMark F. Adams } 571587fa25dSMark F. Adams for (lidx=0,level=pc_gamg->m_Nlevels-1; 572587fa25dSMark F. Adams lidx<fine_level; 573587fa25dSMark F. Adams lidx++, level--){ 574587fa25dSMark F. Adams ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr); 5756c237d78SBarry Smith if( !PETSC_TRUE ) { 57611e60469SMark F. Adams PetscViewer viewer; char fname[32]; 577038e3b61SMark F. Adams sprintf(fname,"Pmat_%d.m",level); 57811e60469SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 5795b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 580038e3b61SMark F. Adams ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr); 5815b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 582e33ef3b1SMark F. Adams sprintf(fname,"Amat_%d.m",level); 583e33ef3b1SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 584e33ef3b1SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 585e33ef3b1SMark F. Adams ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr); 586e33ef3b1SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 5875b89ad90SMark F. Adams } 588587fa25dSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 589587fa25dSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 5908f4b7eb5SMark F. Adams if( PETSC_FALSE ) { 5915ef31b24SMark F. Adams char str[32]; 592587fa25dSMark F. Adams sprintf(str,"MG Level %d (%d)",lidx+1,level-1); 593587fa25dSMark F. Adams PetscLogStageRegister(str, &gamg_stages[level-1]); 594a92563c5SMark F. Adams } 5955b89ad90SMark F. Adams } 5965745e0f5SMark F. Adams 5975b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 598eb07cef2SMark F. Adams a_pc->setupcalled = 0; 599eb07cef2SMark F. Adams ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr); 600e33ef3b1SMark F. Adams 6015b89ad90SMark F. Adams PetscFunctionReturn(0); 6025b89ad90SMark F. Adams } 6035b89ad90SMark F. Adams 604eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 6055b89ad90SMark F. Adams /* 6065b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 6075b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 6085b89ad90SMark F. Adams 6095b89ad90SMark F. Adams Input Parameter: 6105b89ad90SMark F. Adams . pc - the preconditioner context 6115b89ad90SMark F. Adams 6125b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 6135b89ad90SMark F. Adams */ 6145b89ad90SMark F. Adams #undef __FUNCT__ 6155b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 6165b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 6175b89ad90SMark F. Adams { 6185b89ad90SMark F. Adams PetscErrorCode ierr; 6195b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 6205b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 6215b89ad90SMark F. Adams 6225b89ad90SMark F. Adams PetscFunctionBegin; 6235b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 6245b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 6255b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 6265b89ad90SMark F. Adams PetscFunctionReturn(0); 6275b89ad90SMark F. Adams } 6285b89ad90SMark F. Adams 6295b89ad90SMark F. Adams #undef __FUNCT__ 6305b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 6315b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 6325b89ad90SMark F. Adams { 6335b89ad90SMark F. Adams /* PetscErrorCode ierr; */ 6345b89ad90SMark F. Adams /* PC_MG *mg = (PC_MG*)pc->data; */ 6355b89ad90SMark F. Adams /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 6365b89ad90SMark F. Adams /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 6375b89ad90SMark F. Adams 6385b89ad90SMark F. Adams PetscFunctionBegin; 6395b89ad90SMark F. Adams PetscFunctionReturn(0); 6405b89ad90SMark F. Adams } 6415b89ad90SMark F. Adams 6425b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 6435b89ad90SMark F. Adams /* 6445b89ad90SMark F. Adams PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 6455b89ad90SMark F. Adams 6465b89ad90SMark F. Adams Input Parameter: 6475b89ad90SMark F. Adams . pc - the preconditioner context 6485b89ad90SMark F. Adams 6495b89ad90SMark F. Adams Application Interface Routine: PCCreate() 6505b89ad90SMark F. Adams 6515b89ad90SMark F. Adams */ 6525b89ad90SMark F. Adams /* MC 6535b89ad90SMark F. Adams PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 6545b89ad90SMark F. Adams fine grid discretization matrix and coordinates on the fine grid. 6555b89ad90SMark F. Adams 6565b89ad90SMark F. Adams Options Database Key: 6575b89ad90SMark F. Adams Multigrid options(inherited) 6585b89ad90SMark F. Adams + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 6595b89ad90SMark F. Adams . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 6605b89ad90SMark F. Adams . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 6615b89ad90SMark F. Adams -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 6625b89ad90SMark F. Adams GAMG options: 6635b89ad90SMark F. Adams 6645b89ad90SMark F. Adams Level: intermediate 6655b89ad90SMark F. Adams Concepts: multigrid 6665b89ad90SMark F. Adams 6675b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 6685b89ad90SMark F. Adams PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 6695b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 6705b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 6715b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 6725b89ad90SMark F. Adams M */ 6735b89ad90SMark F. Adams 6745b89ad90SMark F. Adams EXTERN_C_BEGIN 6755b89ad90SMark F. Adams #undef __FUNCT__ 6765b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 6775b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 6785b89ad90SMark F. Adams { 6795b89ad90SMark F. Adams PetscErrorCode ierr; 6805b89ad90SMark F. Adams PC_GAMG *pc_gamg; 6815b89ad90SMark F. Adams PC_MG *mg; 6825ef31b24SMark F. Adams PetscClassId cookie; 6835b89ad90SMark F. Adams 6845b89ad90SMark F. Adams PetscFunctionBegin; 6855b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 6865b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 6875b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 6885b89ad90SMark F. Adams 6895b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 6905b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 6915b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 6925b89ad90SMark F. Adams mg->innerctx = pc_gamg; 6935b89ad90SMark F. Adams 6945b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 6955b89ad90SMark F. Adams 6965b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 6975b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 6985b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 6995b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 7005b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 7015b89ad90SMark F. Adams 7025b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 7035b89ad90SMark F. Adams "PCSetCoordinates_C", 7045b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 7055b89ad90SMark F. Adams PCSetCoordinates_GAMG);CHKERRQ(ierr); 7065ef31b24SMark F. Adams 7075ef31b24SMark F. Adams PetscClassIdRegister("GAMG Setup",&cookie); 7085ef31b24SMark F. Adams PetscLogEventRegister("GAMG-createProl", cookie, &gamg_setup_stages[SET1]); 7095ef31b24SMark F. Adams PetscLogEventRegister("GAMG-partLevel", cookie, &gamg_setup_stages[SET2]); 7105ef31b24SMark F. Adams PetscLogEventRegister("GAMG-MIS Graph", cookie, &gamg_setup_stages[SET3]); 7115ef31b24SMark F. Adams PetscLogEventRegister("GAMG-MIS-Agg", cookie, &gamg_setup_stages[SET4]); 7125ef31b24SMark F. Adams PetscLogEventRegister("GAMG-growSupp", cookie, &gamg_setup_stages[SET5]); 7135ef31b24SMark F. Adams PetscLogEventRegister("GAMG-tri-Prol", cookie, &gamg_setup_stages[SET6]); 714e33ef3b1SMark F. Adams PetscLogEventRegister("GAMG-sa-init", cookie, &gamg_setup_stages[SET7]); 715e33ef3b1SMark F. Adams PetscLogEventRegister("GAMG-sa-frmProl0", cookie, &gamg_setup_stages[SET8]); 716e33ef3b1SMark F. Adams PetscLogEventRegister("GAMG-sa-smooth", cookie, &gamg_setup_stages[SET9]); 7175ef31b24SMark F. Adams 7185b89ad90SMark F. Adams PetscFunctionReturn(0); 7195b89ad90SMark F. Adams } 7205b89ad90SMark F. Adams EXTERN_C_END 721