15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4313a3e24SSatish Balay #include "private/matimpl.h" 5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 65a9b9e01SMark F. Adams #include <private/kspimpl.h> 7f96513f1SMatthew G Knepley 8b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 9b4fbaa2aSMark F. Adams PetscLogEvent gamg_setup_events[NUM_SET]; 10b4fbaa2aSMark F. Adams #endif 11b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 12b4fbaa2aSMark F. Adams 13b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 14b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 15b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 16b4fbaa2aSMark F. Adams #endif 17f96513f1SMatthew G Knepley 189d5b6da9SMark F. Adams static PetscFList GAMGList = 0; 199d5b6da9SMark F. Adams 20d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 21d3d6bff4SMark F. Adams #undef __FUNCT__ 22d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 23d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 24d3d6bff4SMark F. Adams { 25d3d6bff4SMark F. Adams PetscErrorCode ierr; 26d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 27d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 28d3d6bff4SMark F. Adams 29d3d6bff4SMark F. Adams PetscFunctionBegin; 309d5b6da9SMark F. Adams if( pc_gamg->data != 0 ) { /* this should not happen, cleaned up in SetUp */ 319d5b6da9SMark F. Adams ierr = PetscFree(pc_gamg->data); CHKERRQ(ierr); 3258471d46SMark F. Adams } 339d5b6da9SMark F. Adams pc_gamg->data = 0; pc_gamg->data_sz = 0; 34d3d6bff4SMark F. Adams PetscFunctionReturn(0); 35d3d6bff4SMark F. Adams } 36d3d6bff4SMark F. Adams 375b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 385b89ad90SMark F. Adams /* 39a147abb0SMark F. Adams createLevel: create coarse op with RAP. repartition and/or reduce number 40a147abb0SMark F. Adams of active processors. 415b89ad90SMark F. Adams 425b89ad90SMark F. Adams Input Parameter: 439d5b6da9SMark F. Adams . pc - parameters 449d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 459d5b6da9SMark F. Adams . cbs - coarse block size 469d5b6da9SMark F. Adams . isLast - 473530afc2SMark F. Adams In/Output Parameter: 483530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 49afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 5011e60469SMark F. Adams Output Parameter: 513530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 525b89ad90SMark F. Adams */ 535cb416c2SMark F. Adams 545b89ad90SMark F. Adams #undef __FUNCT__ 558263b398SMark F. Adams #define __FUNCT__ "createLevel" 569d5b6da9SMark F. Adams PetscErrorCode createLevel( const PC pc, 579d5b6da9SMark F. Adams const Mat Amat_fine, 589d5b6da9SMark F. Adams const PetscInt cbs, 599d5b6da9SMark F. Adams const PetscBool isLast, 603530afc2SMark F. Adams Mat *a_P_inout, 61afc97cdcSMark F. Adams PetscMPIInt *a_nactive_proc, 62389730f3SMark F. Adams Mat *a_Amat_crs 6311e60469SMark F. Adams ) 645b89ad90SMark F. Adams { 659d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 66486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 679d5b6da9SMark F. Adams const PetscBool repart = pc_gamg->repart; 689d5b6da9SMark F. Adams const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 695b89ad90SMark F. Adams PetscErrorCode ierr; 70038e3b61SMark F. Adams Mat Cmat,Pnew,Pold=*a_P_inout; 7111e60469SMark F. Adams IS new_indices,isnum; 729d5b6da9SMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat_fine)->comm; 735a9b9e01SMark F. Adams PetscMPIInt mype,npe,new_npe,nactive = *a_nactive_proc; 745a9b9e01SMark F. Adams PetscInt neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0,rfactor; 755b89ad90SMark F. Adams 765b89ad90SMark F. Adams PetscFunctionBegin; 7711e60469SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 7811e60469SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 79f6536408SMark F. Adams 8011e60469SMark F. Adams /* RAP */ 8196434bc1SMark F. Adams #ifdef USE_R 8296434bc1SMark F. Adams /* make R wih brute force for now */ 8396434bc1SMark F. Adams ierr = MatTranspose( Pold, Pnew ); 8496434bc1SMark F. Adams ierr = MatDestroy( &Pold ); CHKERRQ(ierr); 859d5b6da9SMark F. Adams ierr = MatRARt( Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 8696434bc1SMark F. Adams Pold = Pnew; 8796434bc1SMark F. Adams #else 889d5b6da9SMark F. Adams ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 8996434bc1SMark F. Adams #endif 909d5b6da9SMark F. Adams ierr = MatSetBlockSize( Cmat, cbs ); CHKERRQ(ierr); 91038e3b61SMark F. Adams ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); 929d5b6da9SMark F. Adams ncrs0 = (Iend0-Istart0)/cbs; assert((Iend0-Istart0)%cbs == 0); 93038e3b61SMark F. Adams 94f852f58cSMark F. Adams /* get number of PEs to make active, reduce */ 95f852f58cSMark F. Adams ierr = MatGetSize( Cmat, &neq, &NN ); CHKERRQ(ierr); 9681708dccSMark F. Adams new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 97f852f58cSMark F. Adams if( new_npe == 0 || neq < coarse_max ) new_npe = 1; 985a9b9e01SMark F. Adams else if (new_npe >= nactive ) new_npe = nactive; /* no change, rare */ 999d5b6da9SMark F. Adams if( isLast ) new_npe = 1; 100f852f58cSMark F. Adams 1015a9b9e01SMark F. Adams if( !repart && new_npe==nactive ) { 102a147abb0SMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */ 10322063be5SMark F. Adams } 10422063be5SMark F. Adams else { 10522063be5SMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 106*c8b0795cSMark F. Adams const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,data_sz=ndata_rows*ndata_cols; 107*c8b0795cSMark F. Adams const PetscInt stride0=ncrs0*pc_gamg->data_cell_rows; 1085a9b9e01SMark F. Adams PetscInt *counts,is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx,inpe,targetPE; 1095a9b9e01SMark F. Adams IS isnewproc; 1105a9b9e01SMark F. Adams VecScatter vecscat; 11122063be5SMark F. Adams PetscScalar *array; 11222063be5SMark F. Adams Vec src_crd, dest_crd; 1135a9b9e01SMark F. Adams IS isscat; 114e33ef3b1SMark F. Adams 115fd3c6afaSMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 116b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 117b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 118b4fbaa2aSMark F. Adams #endif 1195a9b9e01SMark F. Adams if( repart ) { 1205a9b9e01SMark F. Adams /* create sub communicator */ 1215a9b9e01SMark F. Adams Mat adj; 1225a9b9e01SMark F. Adams 1239d5b6da9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq); 1245a9b9e01SMark F. Adams 1255a9b9e01SMark F. Adams /* MatPartitioningApply call MatConvert, which is collective */ 1269d5b6da9SMark F. Adams if( cbs == 1 ) { 127038e3b61SMark F. Adams ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 128eb07cef2SMark F. Adams } 129eb07cef2SMark F. Adams else{ 130038e3b61SMark F. Adams /* make a scalar matrix to partition */ 131eb07cef2SMark F. Adams Mat tMat; 13258471d46SMark F. Adams PetscInt ncols,jj,Ii; 133b4fbaa2aSMark F. Adams const PetscScalar *vals; 134b4fbaa2aSMark F. Adams const PetscInt *idx; 1355f8cf99dSMark F. Adams PetscInt *d_nnz, *o_nnz; 1365ee9c036SSatish Balay static int llev = 0; 137b4fbaa2aSMark F. Adams 13858471d46SMark F. Adams ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); 1395f8cf99dSMark F. Adams ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr); 1409d5b6da9SMark F. Adams for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += cbs, jj++ ) { 14158471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 1429d5b6da9SMark F. Adams d_nnz[jj] = ncols/cbs; 1439d5b6da9SMark F. Adams o_nnz[jj] = ncols/cbs; 14458471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 1455f8cf99dSMark F. Adams if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0; 1469d5b6da9SMark F. Adams if( o_nnz[jj] > (neq/cbs-ncrs0) ) o_nnz[jj] = neq/cbs-ncrs0; 14758471d46SMark F. Adams } 1486876a03eSMark F. Adams 149eb07cef2SMark F. Adams ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 150eb07cef2SMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 1515f8cf99dSMark F. Adams 0, d_nnz, 0, o_nnz, 152eb07cef2SMark F. Adams &tMat ); 1536876a03eSMark F. Adams CHKERRQ(ierr); 15458471d46SMark F. Adams ierr = PetscFree( d_nnz ); CHKERRQ(ierr); 1555f8cf99dSMark F. Adams ierr = PetscFree( o_nnz ); CHKERRQ(ierr); 156eb07cef2SMark F. Adams 15722063be5SMark F. Adams for ( ii = Istart0; ii < Iend0; ii++ ) { 1589d5b6da9SMark F. Adams PetscInt dest_row = ii/cbs; 15922063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 160eb07cef2SMark F. Adams for( jj = 0 ; jj < ncols ; jj++ ){ 1619d5b6da9SMark F. Adams PetscInt dest_col = idx[jj]/cbs; 162eb07cef2SMark F. Adams PetscScalar v = 1.0; 163eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 164eb07cef2SMark F. Adams } 16522063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 166eb07cef2SMark F. Adams } 167eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169eb07cef2SMark F. Adams 170b4fbaa2aSMark F. Adams if( llev++ == -1 ) { 171b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 172b4fbaa2aSMark F. Adams sprintf(fname,"part_mat_%d.mat",llev); 173b4fbaa2aSMark F. Adams PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 174b4fbaa2aSMark F. Adams ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 175b4fbaa2aSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 176b4fbaa2aSMark F. Adams } 177b4fbaa2aSMark F. Adams 178eb07cef2SMark F. Adams ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 179eb07cef2SMark F. Adams 180eb07cef2SMark F. Adams ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 181eb07cef2SMark F. Adams } 182f150b916SMark F. Adams 1835a9b9e01SMark F. Adams { /* partition: get newproc_idx, set is_sz */ 1845a9b9e01SMark F. Adams char prefix[256]; 1855a9b9e01SMark F. Adams const char *pcpre; 186b4fbaa2aSMark F. Adams const PetscInt *is_idx; 187b4fbaa2aSMark F. Adams MatPartitioning mpart; 188a4b7d37bSMark F. Adams IS proc_is; 1892f03bc48SMark F. Adams 1905a9b9e01SMark F. Adams ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr); 1915ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 1929d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc,&pcpre);CHKERRQ(ierr); 19359a0be82SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr); 19459a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 19511e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 1965ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 197a4b7d37bSMark F. Adams ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr); 19811e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 1995a9b9e01SMark F. Adams 2005ef31b24SMark F. Adams /* collect IS info */ 201a4b7d37bSMark F. Adams ierr = ISGetLocalSize( proc_is, &is_sz ); CHKERRQ(ierr); 2029d5b6da9SMark F. Adams ierr = PetscMalloc( cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr); 203a4b7d37bSMark F. Adams ierr = ISGetIndices( proc_is, &is_idx ); CHKERRQ(ierr); 204a147abb0SMark F. Adams NN = 1; /* bring to "front" of machine */ 205a147abb0SMark F. Adams /*NN = npe/new_npe;*/ /* spread partitioning across machine */ 206b4fbaa2aSMark F. Adams for( kk = jj = 0 ; kk < is_sz ; kk++ ){ 2079d5b6da9SMark F. Adams for( ii = 0 ; ii < cbs ; ii++, jj++ ){ 2088263b398SMark F. Adams newproc_idx[jj] = is_idx[kk] * NN; /* distribution */ 209eb07cef2SMark F. Adams } 2105ef31b24SMark F. Adams } 211a4b7d37bSMark F. Adams ierr = ISRestoreIndices( proc_is, &is_idx ); CHKERRQ(ierr); 212a4b7d37bSMark F. Adams ierr = ISDestroy( &proc_is ); CHKERRQ(ierr); 213b4fbaa2aSMark F. Adams 2149d5b6da9SMark F. Adams is_sz *= cbs; 2155ef31b24SMark F. Adams } 2165ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 2175a9b9e01SMark F. Adams 2188263b398SMark F. Adams ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc ); 2195a9b9e01SMark F. Adams CHKERRQ(ierr); 2208263b398SMark F. Adams if( newproc_idx != 0 ) { 2218263b398SMark F. Adams ierr = PetscFree( newproc_idx ); CHKERRQ(ierr); 2225ef31b24SMark F. Adams } 2235a9b9e01SMark F. Adams } 2245a9b9e01SMark F. Adams else { /* simple aggreagtion of parts */ 2255a9b9e01SMark F. Adams /* find factor */ 2265a9b9e01SMark F. Adams if( new_npe == 1 ) rfactor = npe; /* easy */ 2275a9b9e01SMark F. Adams else { 2285a9b9e01SMark F. Adams PetscReal best_fact = 0.; 2295a9b9e01SMark F. Adams jj = -1; 2305a9b9e01SMark F. Adams for( kk = 1 ; kk <= npe ; kk++ ){ 2315a9b9e01SMark F. Adams if( npe%kk==0 ) { /* a candidate */ 2325a9b9e01SMark F. Adams PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe; 2335a9b9e01SMark F. Adams if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 2345a9b9e01SMark F. Adams if( fact > best_fact ) { 2355a9b9e01SMark F. Adams best_fact = fact; jj = kk; 2365a9b9e01SMark F. Adams } 2375a9b9e01SMark F. Adams } 2385a9b9e01SMark F. Adams } 2395a9b9e01SMark F. Adams if(jj!=-1) rfactor = jj; 2405a9b9e01SMark F. Adams else rfactor = 1; /* prime? */ 2415a9b9e01SMark F. Adams } 2425a9b9e01SMark F. Adams new_npe = npe/rfactor; 2435a9b9e01SMark F. Adams 2445a9b9e01SMark F. Adams if( new_npe==nactive ) { 2455a9b9e01SMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */ 2465a9b9e01SMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 2475a9b9e01SMark F. Adams *a_nactive_proc = new_npe; /* output */ 2489d5b6da9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s aggregate processors noop: new_npe=%d, neq=%d\n",mype,__FUNCT__,new_npe,neq); 2495a9b9e01SMark F. Adams PetscFunctionReturn(0); 2505a9b9e01SMark F. Adams } 2515a9b9e01SMark F. Adams 2525a9b9e01SMark F. Adams /* reduce - isnewproc */ 2535a9b9e01SMark F. Adams targetPE = mype/rfactor; 2545a9b9e01SMark F. Adams 2559d5b6da9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s aggregate processors: npe: %d --> %d, neq=%d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq); 2569d5b6da9SMark F. Adams is_sz = ncrs0*cbs; 2575a9b9e01SMark F. Adams ierr = ISCreateStride( wcomm, is_sz, targetPE, 0, &isnewproc ); 2585a9b9e01SMark F. Adams CHKERRQ(ierr); 2595a9b9e01SMark F. Adams } 260e33ef3b1SMark F. Adams 26111e60469SMark F. Adams /* 26211e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 26311e60469SMark F. Adams */ 26411e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 26511e60469SMark F. Adams /* 26611e60469SMark F. Adams Determine how many elements are assigned to each processor 26711e60469SMark F. Adams */ 268fd3c6afaSMark F. Adams inpe = npe; 269fd3c6afaSMark F. Adams ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr); 27011e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 2719d5b6da9SMark F. Adams ncrs_new = counts[mype]/cbs; 2729d5b6da9SMark F. Adams strideNew = ncrs_new*ndata_rows; 273b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 274b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 275b4fbaa2aSMark F. Adams #endif 27622063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 27711e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 278d3d6bff4SMark F. Adams ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 279470e26f8SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ 28011e60469SMark F. Adams /* 2819d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 2829d5b6da9SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cbs'. 28311e60469SMark F. Adams */ 28492a756f0SMark F. Adams ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 28511e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 286038e3b61SMark F. Adams for(ii=0,jj=0; ii<ncrs0 ; ii++) { 2879d5b6da9SMark F. Adams PetscInt id = idx[ii*cbs]/cbs; /* get node back */ 288d3d6bff4SMark F. Adams for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk; 28911e60469SMark F. Adams } 290038e3b61SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 291d3d6bff4SMark F. Adams ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 29211e60469SMark F. Adams CHKERRQ(ierr); 29392a756f0SMark F. Adams ierr = PetscFree( tidx ); CHKERRQ(ierr); 29411e60469SMark F. Adams /* 29511e60469SMark F. Adams Create a vector to contain the original vertex information for each element 29611e60469SMark F. Adams */ 297d3d6bff4SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 2989d5b6da9SMark F. Adams for( jj=0; jj<ndata_cols ; jj++ ) { 299d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs0 ; ii++) { 3009d5b6da9SMark F. Adams for( kk=0; kk<ndata_rows ; kk++ ) { 3019d5b6da9SMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*ndata_cols + jj; 302*c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 303676e1743SMark F. Adams ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES ); CHKERRQ(ierr); 304d3d6bff4SMark F. Adams } 305038e3b61SMark F. Adams } 306eb07cef2SMark F. Adams } 307eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 308eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 30911e60469SMark F. Adams /* 31011e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 31111e60469SMark F. Adams to the correct processor 31211e60469SMark F. Adams */ 31311e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 31411e60469SMark F. Adams CHKERRQ(ierr); 31511e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 31611e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 31711e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 31811e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 31911e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 32011e60469SMark F. Adams /* 32111e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 32211e60469SMark F. Adams */ 323*c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 324*c8b0795cSMark F. Adams ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 325*c8b0795cSMark F. Adams pc_gamg->data_sz = data_sz*ncrs_new; 32611e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 3279d5b6da9SMark F. Adams for( jj=0; jj<ndata_cols ; jj++ ) { 328d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs_new ; ii++) { 3299d5b6da9SMark F. Adams for( kk=0; kk<ndata_rows ; kk++ ) { 3309d5b6da9SMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*ndata_cols + jj; 331*c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 332d3d6bff4SMark F. Adams array[jx] = 1.e300; 333d3d6bff4SMark F. Adams } 334038e3b61SMark F. Adams } 335038e3b61SMark F. Adams } 33611e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 33711e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 338ed3f9983SMark F. Adams #if defined PETSC_USE_LOG 339ed3f9983SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 340ed3f9983SMark F. Adams #endif 34111e60469SMark F. Adams /* 34211e60469SMark F. Adams Invert for MatGetSubMatrix 34311e60469SMark F. Adams */ 3449d5b6da9SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new*cbs, &new_indices ); CHKERRQ(ierr); 34511e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 34611e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 347ed3f9983SMark F. Adams #if defined PETSC_USE_LOG 348ed3f9983SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 349ed3f9983SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 350ed3f9983SMark F. Adams #endif 35111e60469SMark F. Adams /* A_crs output */ 352038e3b61SMark F. Adams ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 35311e60469SMark F. Adams CHKERRQ(ierr); 354eb07cef2SMark F. Adams 355038e3b61SMark F. Adams ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 356e33ef3b1SMark F. Adams Cmat = *a_Amat_crs; /* output */ 3579d5b6da9SMark F. Adams ierr = MatSetBlockSize( Cmat, cbs ); CHKERRQ(ierr); 358ed3f9983SMark F. Adams #if defined PETSC_USE_LOG 359ed3f9983SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 360ed3f9983SMark F. Adams #endif 36111e60469SMark F. Adams /* prolongator */ 36211e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 36311e60469SMark F. Adams { 36411e60469SMark F. Adams IS findices; 365ed3f9983SMark F. Adams #if defined PETSC_USE_LOG 366ed3f9983SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 367ed3f9983SMark F. Adams #endif 36811e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 36996434bc1SMark F. Adams #ifdef USE_R 37096434bc1SMark F. Adams ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew ); 37196434bc1SMark F. Adams #else 37211e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 37396434bc1SMark F. Adams #endif 37411e60469SMark F. Adams CHKERRQ(ierr); 37511e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 376ed3f9983SMark F. Adams #if defined PETSC_USE_LOG 377ed3f9983SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 378ed3f9983SMark F. Adams #endif 37911e60469SMark F. Adams } 3803530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 381a147abb0SMark F. Adams *a_P_inout = Pnew; /* output - repartitioned */ 38211e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 38392a756f0SMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 384e33ef3b1SMark F. Adams } 3855b89ad90SMark F. Adams 3865a9b9e01SMark F. Adams *a_nactive_proc = new_npe; /* output */ 3875a9b9e01SMark F. Adams 388*c8b0795cSMark F. Adams if( !PETSC_TRUE ) { 389*c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 390*c8b0795cSMark F. Adams if(llev==0) { 391*c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 392*c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 393*c8b0795cSMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 394*c8b0795cSMark F. Adams ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr); 395*c8b0795cSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 396*c8b0795cSMark F. Adams } 397*c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 398*c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 399*c8b0795cSMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 400*c8b0795cSMark F. Adams ierr = MatView(Cmat, viewer ); CHKERRQ(ierr); 401*c8b0795cSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 402*c8b0795cSMark F. Adams } 403*c8b0795cSMark F. Adams 4045b89ad90SMark F. Adams PetscFunctionReturn(0); 4055b89ad90SMark F. Adams } 4065b89ad90SMark F. Adams 4075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4085b89ad90SMark F. Adams /* 4095b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4105b89ad90SMark F. Adams by setting data structures and options. 4115b89ad90SMark F. Adams 4125b89ad90SMark F. Adams Input Parameter: 4135b89ad90SMark F. Adams . pc - the preconditioner context 4145b89ad90SMark F. Adams 4155b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 4165b89ad90SMark F. Adams 4175b89ad90SMark F. Adams Notes: 4185b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 4195b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 4205b89ad90SMark F. Adams */ 4215b89ad90SMark F. Adams #undef __FUNCT__ 4225b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 4239d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc ) 4245b89ad90SMark F. Adams { 4255b89ad90SMark F. Adams PetscErrorCode ierr; 4269d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 4275b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 42858471d46SMark F. Adams PC_MG_Levels **mglevels = mg->levels; 4299d5b6da9SMark F. Adams Mat Amat = pc->mat, Pmat = pc->pmat; 430d3d6bff4SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend; 4319d5b6da9SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 4323530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 433587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 434*c8b0795cSMark F. Adams PetscReal emaxs[GAMG_MAXLEVELS]; 435737a81a9SMark F. Adams MatInfo info; 4365ef31b24SMark F. Adams 4375b89ad90SMark F. Adams PetscFunctionBegin; 4389d5b6da9SMark F. Adams pc_gamg->setup_count++; 439fca9ed99SMark F. Adams 4409d5b6da9SMark F. Adams if( pc->setupcalled > 0 ) { 44103a628feSMark F. Adams /* just do Galerkin grids */ 44258471d46SMark F. Adams Mat B,dA,dB; 44358471d46SMark F. Adams 44458471d46SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 44558471d46SMark F. Adams ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 44658471d46SMark F. Adams 4479d5b6da9SMark F. Adams if( pc_gamg->Nlevels > 1 ) { 44858471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 4499d5b6da9SMark F. Adams ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 45058471d46SMark F. Adams /* (re)set to get dirty flag */ 4519d5b6da9SMark F. Adams ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 4529d5b6da9SMark F. Adams ierr = KSPSetUp( mglevels[pc_gamg->Nlevels-1]->smoothd ); CHKERRQ(ierr); 45358471d46SMark F. Adams 4549d5b6da9SMark F. Adams for (level=pc_gamg->Nlevels-2; level>-1; level--) { 45558471d46SMark F. Adams ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 45603a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 4579d5b6da9SMark F. Adams if( pc_gamg->setup_count == 2 ) { 45803a628feSMark F. Adams ierr = MatDestroy( &B ); CHKERRQ(ierr); 45903a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 46003a628feSMark F. Adams mglevels[level]->A = B; 46103a628feSMark F. Adams } 46203a628feSMark F. Adams else { 46358471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 46403a628feSMark F. Adams } 46558471d46SMark F. Adams ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 46658471d46SMark F. Adams dB = B; 46758471d46SMark F. Adams /* setup KSP/PC */ 46858471d46SMark F. Adams ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr); 46958471d46SMark F. Adams } 4705f8cf99dSMark F. Adams } 47158471d46SMark F. Adams 4729d5b6da9SMark F. Adams pc->setupcalled = 2; 47303a628feSMark F. Adams 47458471d46SMark F. Adams PetscFunctionReturn(0); 475eb07cef2SMark F. Adams } 476f6536408SMark F. Adams 477baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 478baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 4795b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 4805b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 4813530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 482eb07cef2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 483eb07cef2SMark F. Adams nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 484eb07cef2SMark F. Adams 4859d5b6da9SMark F. Adams if( pc_gamg->data == 0 && nloc > 0 ) { 4869d5b6da9SMark F. Adams if(!pc_gamg->createdefaultdata){ 487*c8b0795cSMark F. Adams SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!"); 488eb07cef2SMark F. Adams } 4899d5b6da9SMark F. Adams ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr); 4909d5b6da9SMark F. Adams } 491038e3b61SMark F. Adams 492eb07cef2SMark F. Adams /* Get A_i and R_i */ 493737a81a9SMark F. Adams ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 494*c8b0795cSMark F. Adams if (pc_gamg->verbose) { 495*c8b0795cSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 496*c8b0795cSMark F. Adams mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 4972c047e2dSMark F. Adams (int)(info.nz_used/(PetscReal)N),npe); 498*c8b0795cSMark F. Adams } 4998f4b7eb5SMark F. Adams for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 5009d5b6da9SMark F. Adams level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 5010205a208SMark F. Adams level++ ){ 5025b89ad90SMark F. Adams level1 = level + 1; 503b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 504b4fbaa2aSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 505b4fbaa2aSMark F. Adams #endif 506b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 507b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr); 508b4fbaa2aSMark F. Adams #endif 509*c8b0795cSMark F. Adams { /* construct prolongator */ 510*c8b0795cSMark F. Adams Mat Gmat; assert(pc_gamg->graph); 511*c8b0795cSMark F. Adams ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr); 512*c8b0795cSMark F. Adams 513*c8b0795cSMark F. Adams IS selected, llist; assert(pc_gamg->coarsen); 514*c8b0795cSMark F. Adams ierr = pc_gamg->coarsen( pc, Gmat, &selected, &llist ); CHKERRQ(ierr); 515*c8b0795cSMark F. Adams 516*c8b0795cSMark F. Adams ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, selected, llist, &Parr[level1] ); 5175b89ad90SMark F. Adams CHKERRQ(ierr); 518*c8b0795cSMark F. Adams 519*c8b0795cSMark F. Adams if( Parr[level1] ){ 5209d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 521*c8b0795cSMark F. Adams if( pc_gamg->col_bs_id != -1 ){ 5229d5b6da9SMark F. Adams PetscBool flag; 5239d5b6da9SMark F. Adams ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag ); 524*c8b0795cSMark F. Adams CHKERRQ( ierr ); assert(flag); 5259d5b6da9SMark F. Adams } 526*c8b0795cSMark F. Adams } 527*c8b0795cSMark F. Adams 528*c8b0795cSMark F. Adams if( pc_gamg->optprol && Parr[level1] ){ 529*c8b0795cSMark F. Adams /* smooth */ 530*c8b0795cSMark F. Adams ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr); 531*c8b0795cSMark F. Adams } 532*c8b0795cSMark F. Adams 533*c8b0795cSMark F. Adams ierr = MatDestroy( &Gmat ); CHKERRQ(ierr); 534*c8b0795cSMark F. Adams ierr = ISDestroy( &selected ); CHKERRQ(ierr); 535*c8b0795cSMark F. Adams ierr = ISDestroy( &llist ); CHKERRQ(ierr); 536*c8b0795cSMark F. Adams } 537*c8b0795cSMark F. Adams #if defined PETSC_USE_LOG 538*c8b0795cSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 539*c8b0795cSMark F. Adams #endif 5409d5b6da9SMark F. Adams /* cache eigen estimate */ 5419d5b6da9SMark F. Adams if( pc_gamg->emax_id != -1 ){ 5429d5b6da9SMark F. Adams PetscBool flag; 543f73f8d2cSSatish Balay ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag ); 5449d5b6da9SMark F. Adams CHKERRQ( ierr ); 5459d5b6da9SMark F. Adams if( !flag ) emaxs[level] = -1.; 5469d5b6da9SMark F. Adams } 5479d5b6da9SMark F. Adams else emaxs[level] = -1.; 548baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 549*c8b0795cSMark F. Adams if( !Parr[level1] ) { 550*c8b0795cSMark F. Adams if (pc_gamg->verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level); 551*c8b0795cSMark F. Adams break; 552*c8b0795cSMark F. Adams } 553*c8b0795cSMark F. Adams else { 554b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 555b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 556b4fbaa2aSMark F. Adams #endif 557*c8b0795cSMark F. Adams ierr = createLevel( pc, Aarr[level], bs, 5589d5b6da9SMark F. Adams (PetscBool)(level==pc_gamg->Nlevels-2), 559*c8b0795cSMark F. Adams &Parr[level1], &nactivepe, &Aarr[level1] ); 5603530afc2SMark F. Adams CHKERRQ(ierr); 561b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 562b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 563b4fbaa2aSMark F. Adams #endif 5643530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 565737a81a9SMark F. Adams ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); CHKERRQ(ierr); 566*c8b0795cSMark F. Adams if (pc_gamg->verbose){ 567*c8b0795cSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 568*c8b0795cSMark F. Adams mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols, 5692c047e2dSMark F. Adams (int)(info.nz_used/(PetscReal)N),nactivepe); 570*c8b0795cSMark F. Adams } 57181708dccSMark F. Adams /* stop if one node */ 572*c8b0795cSMark F. Adams if( M/pc_gamg->data_cell_cols < 2 ) { 57381708dccSMark F. Adams level++; 57481708dccSMark F. Adams break; 57581708dccSMark F. Adams } 576*c8b0795cSMark F. Adams /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */ 577*c8b0795cSMark F. Adams /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */ 578*c8b0795cSMark F. Adams /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */ 579*c8b0795cSMark F. Adams /* nloceq = Iend-Istart; */ 580*c8b0795cSMark F. Adams /* ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); */ 581*c8b0795cSMark F. Adams /* ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); */ 582*c8b0795cSMark F. Adams /* ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); */ 583*c8b0795cSMark F. Adams /* for(kk=0;kk<nloceq;kk++){ */ 584*c8b0795cSMark F. Adams /* if(data_arr[kk]==0.0) { */ 585*c8b0795cSMark F. Adams /* id = kk + Istart; */ 586*c8b0795cSMark F. Adams /* ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */ 587*c8b0795cSMark F. Adams /* CHKERRQ(ierr); */ 588*c8b0795cSMark F. Adams /* PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */ 589*c8b0795cSMark F. Adams /* } */ 590*c8b0795cSMark F. Adams /* } */ 591*c8b0795cSMark F. Adams /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */ 592*c8b0795cSMark F. Adams /* ierr = VecDestroy( &diag ); CHKERRQ(ierr); */ 593*c8b0795cSMark F. Adams /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */ 594*c8b0795cSMark F. Adams /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */ 595e33ef3b1SMark F. Adams } 596b4fbaa2aSMark F. Adams 597b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 598b4fbaa2aSMark F. Adams ierr = PetscLogStagePop(); CHKERRQ( ierr ); 599b4fbaa2aSMark F. Adams #endif 600*c8b0795cSMark F. Adams } /* levels */ 601*c8b0795cSMark F. Adams 602*c8b0795cSMark F. Adams if( pc_gamg->data ) { 603*c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 604*c8b0795cSMark F. Adams pc_gamg->data = 0; 6055b89ad90SMark F. Adams } 606*c8b0795cSMark F. Adams 6079d5b6da9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 6089d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6095b89ad90SMark F. Adams fine_level = level; 6109d5b6da9SMark F. Adams ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr); 6115b89ad90SMark F. Adams 6129d5b6da9SMark F. Adams if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if */ 613fc4362bfSMark F. Adams /* set default smoothers */ 6149d5b6da9SMark F. Adams for ( lidx = 1, level = pc_gamg->Nlevels-2; 615587fa25dSMark F. Adams lidx <= fine_level; 616587fa25dSMark F. Adams lidx++, level--) { 6175745e0f5SMark F. Adams PetscReal emax, emin; 6185b89ad90SMark F. Adams KSP smoother; PC subpc; 619f6536408SMark F. Adams PetscBool isCheb; 620f6536408SMark F. Adams /* set defaults */ 6219d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 6225b89ad90SMark F. Adams ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 623f6536408SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 624f6536408SMark F. Adams /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */ 6259d5b6da9SMark F. Adams ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr); 626f6536408SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 627f6536408SMark F. Adams /* overide defaults with input parameters */ 628f6536408SMark F. Adams ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr); 629f6536408SMark F. Adams 630f6536408SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 631f6536408SMark F. Adams /* do my own cheby */ 632f6536408SMark F. Adams ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr); 63341139db5SMark F. Adams 634f6536408SMark F. Adams if( isCheb ) { 6359d5b6da9SMark F. Adams ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &isCheb ); CHKERRQ(ierr); 636f6536408SMark F. Adams if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 637587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 638587fa25dSMark F. Adams KSP eksp; Mat Lmat = Aarr[level]; 639fc4362bfSMark F. Adams Vec bb, xx; PC pc; 640f6536408SMark F. Adams const PCType type; 641038e3b61SMark F. Adams 642f6536408SMark F. Adams ierr = PCGetType( subpc, &type ); CHKERRQ(ierr); 6435745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 6445745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 645fc4362bfSMark F. Adams { 646fc4362bfSMark F. Adams PetscRandom rctx; 647fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 648fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 649fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 650fc4362bfSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 6515b89ad90SMark F. Adams } 652fc4362bfSMark F. Adams ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 653038e3b61SMark F. Adams ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 654fc4362bfSMark F. Adams CHKERRQ(ierr); 655fc4362bfSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 656fc4362bfSMark F. Adams 657f6536408SMark F. Adams ierr = KSPAppendOptionsPrefix( eksp, "est_"); CHKERRQ(ierr); 658f6536408SMark F. Adams ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 659f6536408SMark F. Adams 660f6536408SMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 661f6536408SMark F. Adams ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 662f6536408SMark F. Adams ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 663f6536408SMark F. Adams ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */ 664f6536408SMark F. Adams 665fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 6665a9b9e01SMark F. Adams 6675a9b9e01SMark F. Adams /* solve - keep stuff out of logging */ 6685a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 6695a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 670fc4362bfSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 6715a9b9e01SMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 6725a9b9e01SMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 6735a9b9e01SMark F. Adams 674fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 6755a9b9e01SMark F. Adams 676fc4362bfSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 677fc4362bfSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 678fc4362bfSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 679f6536408SMark F. Adams 6809d5b6da9SMark F. Adams if (pc_gamg->verbose) { 68141139db5SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e\n", 68241139db5SMark F. Adams __FUNCT__,emax,emin); 683f6536408SMark F. Adams } 684fc4362bfSMark F. Adams } 685038e3b61SMark F. Adams { 686038e3b61SMark F. Adams PetscInt N1, N0, tt; 687038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 688038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 689f6536408SMark F. Adams /* heuristic - is this crap? */ 690f6536408SMark F. Adams emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); 691038e3b61SMark F. Adams emax *= 1.05; 692038e3b61SMark F. Adams } 693fc4362bfSMark F. Adams ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 694f6536408SMark F. Adams } 6955745e0f5SMark F. Adams } 6965745e0f5SMark F. Adams { 697ecd57171SMark F. Adams /* coarse grid */ 698ecd57171SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 6999d5b6da9SMark F. Adams Mat Lmat = Aarr[pc_gamg->Nlevels-1]; 7009d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 70158471d46SMark F. Adams ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 7025745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 703ecd57171SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 704ecd57171SMark F. Adams ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 705ecd57171SMark F. Adams ierr = PCSetUp( subpc ); CHKERRQ(ierr); 706ecd57171SMark F. Adams ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 707ecd57171SMark F. Adams assert(ii==1); 708ecd57171SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 709ecd57171SMark F. Adams ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 710fc4362bfSMark F. Adams } 711737a81a9SMark F. Adams 712fc4362bfSMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 713ce4cda84SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 7149d5b6da9SMark F. Adams ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); 715ce4cda84SJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 716ce4cda84SJed Brown if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 7175745e0f5SMark F. Adams 71858471d46SMark F. Adams /* set interpolation between the levels, clean up */ 7199d5b6da9SMark F. Adams for (lidx=0,level=pc_gamg->Nlevels-1; 72058471d46SMark F. Adams lidx<fine_level; 72158471d46SMark F. Adams lidx++, level--){ 7229d5b6da9SMark F. Adams ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr); 723587fa25dSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 724587fa25dSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 7255b89ad90SMark F. Adams } 7265745e0f5SMark F. Adams 7275b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 7289d5b6da9SMark F. Adams pc->setupcalled = 0; 7299d5b6da9SMark F. Adams ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 7309d5b6da9SMark F. Adams pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ 73158471d46SMark F. Adams 73258471d46SMark F. Adams { 73358471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 7349d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 73558471d46SMark F. Adams ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 73658471d46SMark F. Adams ierr = KSPSetUp( smoother ); CHKERRQ(ierr); 73758471d46SMark F. Adams } 7385f8cf99dSMark F. Adams } 7395f8cf99dSMark F. Adams else { 7405f8cf99dSMark F. Adams KSP smoother; 7419d5b6da9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__); 7429d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 7435f8cf99dSMark F. Adams ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 7445f8cf99dSMark F. Adams ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 7455f8cf99dSMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 7469d5b6da9SMark F. Adams pc->setupcalled = 0; 7479d5b6da9SMark F. Adams ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 7489d5b6da9SMark F. Adams pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ 7495f8cf99dSMark F. Adams } 7505b89ad90SMark F. Adams PetscFunctionReturn(0); 7515b89ad90SMark F. Adams } 7525b89ad90SMark F. Adams 753eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 7545b89ad90SMark F. Adams /* 7555b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 7565b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 7575b89ad90SMark F. Adams 7585b89ad90SMark F. Adams Input Parameter: 7595b89ad90SMark F. Adams . pc - the preconditioner context 7605b89ad90SMark F. Adams 7615b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 7625b89ad90SMark F. Adams */ 7635b89ad90SMark F. Adams #undef __FUNCT__ 7645b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 7655b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG( PC pc ) 7665b89ad90SMark F. Adams { 7675b89ad90SMark F. Adams PetscErrorCode ierr; 7685b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7695b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 7705b89ad90SMark F. Adams 7715b89ad90SMark F. Adams PetscFunctionBegin; 7725b89ad90SMark F. Adams ierr = PCReset_GAMG( pc );CHKERRQ(ierr); 7735b89ad90SMark F. Adams ierr = PetscFree( pc_gamg );CHKERRQ(ierr); 7745b89ad90SMark F. Adams ierr = PCDestroy_MG( pc );CHKERRQ(ierr); 7755b89ad90SMark F. Adams PetscFunctionReturn(0); 7765b89ad90SMark F. Adams } 7775b89ad90SMark F. Adams 778676e1743SMark F. Adams 779676e1743SMark F. Adams #undef __FUNCT__ 780676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 781676e1743SMark F. Adams /*@ 782676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 783676e1743SMark F. Adams processor reduction. 784676e1743SMark F. Adams 785676e1743SMark F. Adams Not Collective on PC 786676e1743SMark F. Adams 787676e1743SMark F. Adams Input Parameters: 788676e1743SMark F. Adams . pc - the preconditioner context 789676e1743SMark F. Adams 790676e1743SMark F. Adams 791676e1743SMark F. Adams Options Database Key: 792676e1743SMark F. Adams . -pc_gamg_process_eq_limit 793676e1743SMark F. Adams 794676e1743SMark F. Adams Level: intermediate 795676e1743SMark F. Adams 796676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 797676e1743SMark F. Adams 798676e1743SMark F. Adams .seealso: () 799676e1743SMark F. Adams @*/ 800676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 801676e1743SMark F. Adams { 802676e1743SMark F. Adams PetscErrorCode ierr; 803676e1743SMark F. Adams 804676e1743SMark F. Adams PetscFunctionBegin; 805676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 806676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 807676e1743SMark F. Adams PetscFunctionReturn(0); 808676e1743SMark F. Adams } 809676e1743SMark F. Adams 810676e1743SMark F. Adams EXTERN_C_BEGIN 811676e1743SMark F. Adams #undef __FUNCT__ 812676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 813676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 814676e1743SMark F. Adams { 815c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 816c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 817676e1743SMark F. Adams 818676e1743SMark F. Adams PetscFunctionBegin; 8199d5b6da9SMark F. Adams if(n>0) pc_gamg->min_eq_proc = n; 820676e1743SMark F. Adams PetscFunctionReturn(0); 821676e1743SMark F. Adams } 822676e1743SMark F. Adams EXTERN_C_END 823676e1743SMark F. Adams 824676e1743SMark F. Adams #undef __FUNCT__ 825389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 826389730f3SMark F. Adams /*@ 827389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 828389730f3SMark F. Adams 829389730f3SMark F. Adams Collective on PC 830389730f3SMark F. Adams 831389730f3SMark F. Adams Input Parameters: 832389730f3SMark F. Adams . pc - the preconditioner context 833389730f3SMark F. Adams 834389730f3SMark F. Adams 835389730f3SMark F. Adams Options Database Key: 836389730f3SMark F. Adams . -pc_gamg_coarse_eq_limit 837389730f3SMark F. Adams 838389730f3SMark F. Adams Level: intermediate 839389730f3SMark F. Adams 840389730f3SMark F. Adams Concepts: Unstructured multrigrid preconditioner 841389730f3SMark F. Adams 842389730f3SMark F. Adams .seealso: () 843389730f3SMark F. Adams @*/ 844389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 845389730f3SMark F. Adams { 846389730f3SMark F. Adams PetscErrorCode ierr; 847389730f3SMark F. Adams 848389730f3SMark F. Adams PetscFunctionBegin; 849389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 850389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 851389730f3SMark F. Adams PetscFunctionReturn(0); 852389730f3SMark F. Adams } 853389730f3SMark F. Adams 854389730f3SMark F. Adams EXTERN_C_BEGIN 855389730f3SMark F. Adams #undef __FUNCT__ 856389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 857389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 858389730f3SMark F. Adams { 859389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 860389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 861389730f3SMark F. Adams 862389730f3SMark F. Adams PetscFunctionBegin; 8639d5b6da9SMark F. Adams if(n>0) pc_gamg->coarse_eq_limit = n; 864389730f3SMark F. Adams PetscFunctionReturn(0); 865389730f3SMark F. Adams } 866389730f3SMark F. Adams EXTERN_C_END 867389730f3SMark F. Adams 868389730f3SMark F. Adams #undef __FUNCT__ 8698263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 870676e1743SMark F. Adams /*@ 8718263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 872676e1743SMark F. Adams 873676e1743SMark F. Adams Collective on PC 874676e1743SMark F. Adams 875676e1743SMark F. Adams Input Parameters: 876676e1743SMark F. Adams . pc - the preconditioner context 877676e1743SMark F. Adams 878676e1743SMark F. Adams 879676e1743SMark F. Adams Options Database Key: 8808263b398SMark F. Adams . -pc_gamg_repartition 881676e1743SMark F. Adams 882676e1743SMark F. Adams Level: intermediate 883676e1743SMark F. Adams 884676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 885676e1743SMark F. Adams 886676e1743SMark F. Adams .seealso: () 887676e1743SMark F. Adams @*/ 8888263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 889676e1743SMark F. Adams { 890676e1743SMark F. Adams PetscErrorCode ierr; 891676e1743SMark F. Adams 892676e1743SMark F. Adams PetscFunctionBegin; 893676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8948263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 895676e1743SMark F. Adams PetscFunctionReturn(0); 896676e1743SMark F. Adams } 897676e1743SMark F. Adams 898676e1743SMark F. Adams EXTERN_C_BEGIN 899676e1743SMark F. Adams #undef __FUNCT__ 9008263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 9018263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 902676e1743SMark F. Adams { 903c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 904c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 905676e1743SMark F. Adams 906676e1743SMark F. Adams PetscFunctionBegin; 9079d5b6da9SMark F. Adams pc_gamg->repart = n; 908676e1743SMark F. Adams PetscFunctionReturn(0); 909676e1743SMark F. Adams } 910676e1743SMark F. Adams EXTERN_C_END 911676e1743SMark F. Adams 912676e1743SMark F. Adams #undef __FUNCT__ 9134ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 9144ef23d27SMark F. Adams /*@ 9154ef23d27SMark F. Adams PCGAMGSetNlevels - 9164ef23d27SMark F. Adams 9174ef23d27SMark F. Adams Not collective on PC 9184ef23d27SMark F. Adams 9194ef23d27SMark F. Adams Input Parameters: 9204ef23d27SMark F. Adams . pc - the preconditioner context 9214ef23d27SMark F. Adams 9224ef23d27SMark F. Adams 9234ef23d27SMark F. Adams Options Database Key: 9244ef23d27SMark F. Adams . -pc_mg_levels 9254ef23d27SMark F. Adams 9264ef23d27SMark F. Adams Level: intermediate 9274ef23d27SMark F. Adams 9284ef23d27SMark F. Adams Concepts: Unstructured multrigrid preconditioner 9294ef23d27SMark F. Adams 9304ef23d27SMark F. Adams .seealso: () 9314ef23d27SMark F. Adams @*/ 9324ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 9334ef23d27SMark F. Adams { 9344ef23d27SMark F. Adams PetscErrorCode ierr; 9354ef23d27SMark F. Adams 9364ef23d27SMark F. Adams PetscFunctionBegin; 9374ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9384ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 9394ef23d27SMark F. Adams PetscFunctionReturn(0); 9404ef23d27SMark F. Adams } 9414ef23d27SMark F. Adams 9424ef23d27SMark F. Adams EXTERN_C_BEGIN 9434ef23d27SMark F. Adams #undef __FUNCT__ 9444ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 9454ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 9464ef23d27SMark F. Adams { 9474ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9484ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9494ef23d27SMark F. Adams 9504ef23d27SMark F. Adams PetscFunctionBegin; 9519d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 9524ef23d27SMark F. Adams PetscFunctionReturn(0); 9534ef23d27SMark F. Adams } 9544ef23d27SMark F. Adams EXTERN_C_END 9554ef23d27SMark F. Adams 9564ef23d27SMark F. Adams #undef __FUNCT__ 9573542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 9583542efc5SMark F. Adams /*@ 9593542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 9603542efc5SMark F. Adams 9613542efc5SMark F. Adams Not collective on PC 9623542efc5SMark F. Adams 9633542efc5SMark F. Adams Input Parameters: 9643542efc5SMark F. Adams . pc - the preconditioner context 9653542efc5SMark F. Adams 9663542efc5SMark F. Adams 9673542efc5SMark F. Adams Options Database Key: 9683542efc5SMark F. Adams . -pc_gamg_threshold 9693542efc5SMark F. Adams 9703542efc5SMark F. Adams Level: intermediate 9713542efc5SMark F. Adams 9723542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 9733542efc5SMark F. Adams 9743542efc5SMark F. Adams .seealso: () 9753542efc5SMark F. Adams @*/ 9763542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 9773542efc5SMark F. Adams { 9783542efc5SMark F. Adams PetscErrorCode ierr; 9793542efc5SMark F. Adams 9803542efc5SMark F. Adams PetscFunctionBegin; 9813542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9823542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 9833542efc5SMark F. Adams PetscFunctionReturn(0); 9843542efc5SMark F. Adams } 9853542efc5SMark F. Adams 9863542efc5SMark F. Adams EXTERN_C_BEGIN 9873542efc5SMark F. Adams #undef __FUNCT__ 9883542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 9893542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 9903542efc5SMark F. Adams { 991c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 992c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9933542efc5SMark F. Adams 9943542efc5SMark F. Adams PetscFunctionBegin; 9959d5b6da9SMark F. Adams pc_gamg->threshold = n; 9963542efc5SMark F. Adams PetscFunctionReturn(0); 9973542efc5SMark F. Adams } 9983542efc5SMark F. Adams EXTERN_C_END 9993542efc5SMark F. Adams 10003542efc5SMark F. Adams #undef __FUNCT__ 10019d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1002676e1743SMark F. Adams /*@ 10039d5b6da9SMark F. Adams PCGAMGSetType - Set solution method - calls sub create method 1004676e1743SMark F. Adams 1005676e1743SMark F. Adams Collective on PC 1006676e1743SMark F. Adams 1007676e1743SMark F. Adams Input Parameters: 1008676e1743SMark F. Adams . pc - the preconditioner context 1009676e1743SMark F. Adams 1010676e1743SMark F. Adams 1011676e1743SMark F. Adams Options Database Key: 10123542efc5SMark F. Adams . -pc_gamg_type 1013676e1743SMark F. Adams 1014676e1743SMark F. Adams Level: intermediate 1015676e1743SMark F. Adams 1016676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1017676e1743SMark F. Adams 1018676e1743SMark F. Adams .seealso: () 1019676e1743SMark F. Adams @*/ 10200202ef74SSatish Balay PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type ) 1021676e1743SMark F. Adams { 1022676e1743SMark F. Adams PetscErrorCode ierr; 1023676e1743SMark F. Adams 1024676e1743SMark F. Adams PetscFunctionBegin; 1025676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10260202ef74SSatish Balay ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type)); 1027676e1743SMark F. Adams CHKERRQ(ierr); 1028676e1743SMark F. Adams PetscFunctionReturn(0); 1029676e1743SMark F. Adams } 1030676e1743SMark F. Adams 1031676e1743SMark F. Adams EXTERN_C_BEGIN 1032676e1743SMark F. Adams #undef __FUNCT__ 10339d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 10340202ef74SSatish Balay PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type ) 1035676e1743SMark F. Adams { 10369d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 1037676e1743SMark F. Adams 1038676e1743SMark F. Adams PetscFunctionBegin; 10399d5b6da9SMark F. Adams ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r); 10409d5b6da9SMark F. Adams CHKERRQ(ierr); 10419d5b6da9SMark F. Adams 10429d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 10439d5b6da9SMark F. Adams 10449d5b6da9SMark F. Adams /* call sub create method */ 10459d5b6da9SMark F. Adams ierr = (*r)(pc); CHKERRQ(ierr); 10469d5b6da9SMark F. Adams 1047676e1743SMark F. Adams PetscFunctionReturn(0); 1048676e1743SMark F. Adams } 1049676e1743SMark F. Adams EXTERN_C_END 1050676e1743SMark F. Adams 10515b89ad90SMark F. Adams #undef __FUNCT__ 10525b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 10535b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG( PC pc ) 10545b89ad90SMark F. Adams { 1055676e1743SMark F. Adams PetscErrorCode ierr; 1056676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1057676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1058676e1743SMark F. Adams PetscBool flag; 10595b89ad90SMark F. Adams 10605b89ad90SMark F. Adams PetscFunctionBegin; 1061676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1062676e1743SMark F. Adams { 106375b74bdaSMark F. Adams /* -pc_gamg_verbose */ 10649d5b6da9SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 10659d5b6da9SMark F. Adams "none", pc_gamg->verbose, 10669d5b6da9SMark F. Adams &pc_gamg->verbose, PETSC_NULL ); 10679d5b6da9SMark F. Adams CHKERRQ(ierr); 106875b74bdaSMark F. Adams 10698263b398SMark F. Adams /* -pc_gamg_repartition */ 10708263b398SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_repartition", 10718263b398SMark F. Adams "Repartion coarse grids (false)", 10728263b398SMark F. Adams "PCGAMGRepartitioning", 10739d5b6da9SMark F. Adams pc_gamg->repart, 10749d5b6da9SMark F. Adams &pc_gamg->repart, 1075676e1743SMark F. Adams &flag); 1076676e1743SMark F. Adams CHKERRQ(ierr); 1077676e1743SMark F. Adams 1078c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 1079676e1743SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1080676e1743SMark F. Adams "Limit (goal) on number of equations per process on coarse grids", 1081676e1743SMark F. Adams "PCGAMGSetProcEqLim", 10829d5b6da9SMark F. Adams pc_gamg->min_eq_proc, 10839d5b6da9SMark F. Adams &pc_gamg->min_eq_proc, 1084676e1743SMark F. Adams &flag ); 1085676e1743SMark F. Adams CHKERRQ(ierr); 10863542efc5SMark F. Adams 1087389730f3SMark F. Adams /* -pc_gamg_coarse_eq_limit */ 1088389730f3SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1089389730f3SMark F. Adams "Limit on number of equations for the coarse grid", 1090389730f3SMark F. Adams "PCGAMGSetCoarseEqLim", 10919d5b6da9SMark F. Adams pc_gamg->coarse_eq_limit, 10929d5b6da9SMark F. Adams &pc_gamg->coarse_eq_limit, 1093389730f3SMark F. Adams &flag ); 1094389730f3SMark F. Adams CHKERRQ(ierr); 1095389730f3SMark F. Adams 1096c20e4228SMark F. Adams /* -pc_gamg_threshold */ 10973542efc5SMark F. Adams ierr = PetscOptionsReal("-pc_gamg_threshold", 10983542efc5SMark F. Adams "Relative threshold to use for dropping edges in aggregation graph", 10993542efc5SMark F. Adams "PCGAMGSetThreshold", 11009d5b6da9SMark F. Adams pc_gamg->threshold, 11019d5b6da9SMark F. Adams &pc_gamg->threshold, 11023542efc5SMark F. Adams &flag ); 11033542efc5SMark F. Adams CHKERRQ(ierr); 11044ef23d27SMark F. Adams 11054ef23d27SMark F. Adams ierr = PetscOptionsInt("-pc_mg_levels", 11064ef23d27SMark F. Adams "Set number of MG levels", 11074ef23d27SMark F. Adams "PCGAMGSetNlevels", 11089d5b6da9SMark F. Adams pc_gamg->Nlevels, 11099d5b6da9SMark F. Adams &pc_gamg->Nlevels, 11104ef23d27SMark F. Adams &flag ); 1111676e1743SMark F. Adams } 1112676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1113676e1743SMark F. Adams 11145b89ad90SMark F. Adams PetscFunctionReturn(0); 11155b89ad90SMark F. Adams } 11165b89ad90SMark F. Adams 11175b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 11185b89ad90SMark F. Adams /*MC 11199d5b6da9SMark F. Adams PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 11209d5b6da9SMark F. Adams - This is the entry point to GAMG, registered in pcregis.c 11215b89ad90SMark F. Adams 1122280d9858SJed Brown Options Database Keys: 11235b89ad90SMark F. Adams Multigrid options(inherited) 1124280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1125280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1126280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1127280d9858SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 11285b89ad90SMark F. Adams 11295b89ad90SMark F. Adams Level: intermediate 1130280d9858SJed Brown 11315b89ad90SMark F. Adams Concepts: multigrid 11325b89ad90SMark F. Adams 11335b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1134280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 11355b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 11365b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 11375b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 11385b89ad90SMark F. Adams M*/ 11395b89ad90SMark F. Adams EXTERN_C_BEGIN 11405b89ad90SMark F. Adams #undef __FUNCT__ 11415b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 11425b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG( PC pc ) 11435b89ad90SMark F. Adams { 11445b89ad90SMark F. Adams PetscErrorCode ierr; 11455b89ad90SMark F. Adams PC_GAMG *pc_gamg; 11465b89ad90SMark F. Adams PC_MG *mg; 11475ef31b24SMark F. Adams PetscClassId cookie; 11489d5b6da9SMark F. Adams 11495ee9c036SSatish Balay #if defined PETSC_USE_LOG 11509d5b6da9SMark F. Adams static long count = 0; 11515ee9c036SSatish Balay #endif 11525b89ad90SMark F. Adams 11535b89ad90SMark F. Adams PetscFunctionBegin; 11545b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 11555b89ad90SMark F. Adams ierr = PCSetType( pc, PCMG ); CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 11565b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr); 11575b89ad90SMark F. Adams 11585b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 11595b89ad90SMark F. Adams ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr); 11605b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 1161ce4cda84SJed Brown mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 11625b89ad90SMark F. Adams mg->innerctx = pc_gamg; 11635b89ad90SMark F. Adams 11649d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 11659d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 11669d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 11679d5b6da9SMark F. Adams pc_gamg->data = 0; 11685b89ad90SMark F. Adams 11699d5b6da9SMark F. Adams /* register AMG type */ 11709d5b6da9SMark F. Adams if( !GAMGList ){ 11719d5b6da9SMark F. Adams ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 11729d5b6da9SMark F. Adams ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 11739d5b6da9SMark F. Adams } 11749d5b6da9SMark F. Adams 11759d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 11765b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 11775b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 11785b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 11795b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 11805b89ad90SMark F. Adams 11815b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1182676e1743SMark F. Adams "PCGAMGSetProcEqLim_C", 1183676e1743SMark F. Adams "PCGAMGSetProcEqLim_GAMG", 1184676e1743SMark F. Adams PCGAMGSetProcEqLim_GAMG); 1185676e1743SMark F. Adams CHKERRQ(ierr); 1186676e1743SMark F. Adams 1187676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1188389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_C", 1189389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_GAMG", 1190389730f3SMark F. Adams PCGAMGSetCoarseEqLim_GAMG); 1191389730f3SMark F. Adams CHKERRQ(ierr); 1192389730f3SMark F. Adams 1193389730f3SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 11948263b398SMark F. Adams "PCGAMGSetRepartitioning_C", 11958263b398SMark F. Adams "PCGAMGSetRepartitioning_GAMG", 11968263b398SMark F. Adams PCGAMGSetRepartitioning_GAMG); 1197676e1743SMark F. Adams CHKERRQ(ierr); 1198676e1743SMark F. Adams 1199676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 12003542efc5SMark F. Adams "PCGAMGSetThreshold_C", 12013542efc5SMark F. Adams "PCGAMGSetThreshold_GAMG", 12023542efc5SMark F. Adams PCGAMGSetThreshold_GAMG); 12033542efc5SMark F. Adams CHKERRQ(ierr); 12043542efc5SMark F. Adams 12053542efc5SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 12069d5b6da9SMark F. Adams "PCGAMGSetType_C", 12079d5b6da9SMark F. Adams "PCGAMGSetType_GAMG", 12089d5b6da9SMark F. Adams PCGAMGSetType_GAMG); 1209676e1743SMark F. Adams CHKERRQ(ierr); 1210c97e1df0SMark F. Adams 12119d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 12129d5b6da9SMark F. Adams pc_gamg->min_eq_proc = 100; 1213*c8b0795cSMark F. Adams pc_gamg->coarse_eq_limit = 800; 12149d5b6da9SMark F. Adams pc_gamg->threshold = 0.05; 12159d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 12169d5b6da9SMark F. Adams pc_gamg->verbose = 0; 12179d5b6da9SMark F. Adams pc_gamg->emax_id = -1; 12189d5b6da9SMark F. Adams pc_gamg->col_bs_id = -1; 12199d5b6da9SMark F. Adams 12209d5b6da9SMark F. Adams /* instantiate derived type */ 12219d5b6da9SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 12229d5b6da9SMark F. Adams { 12239d5b6da9SMark F. Adams char tname[256] = GAMGAGG; 12249d5b6da9SMark F. Adams ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType", 12259d5b6da9SMark F. Adams GAMGList, tname, tname, sizeof(tname), PETSC_NULL ); 12269d5b6da9SMark F. Adams CHKERRQ(ierr); 12279d5b6da9SMark F. Adams ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr); 12289d5b6da9SMark F. Adams } 12299d5b6da9SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 12309d5b6da9SMark F. Adams 1231b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 1232785cba28SMark F. Adams if( count++ == 0 ) { 12335ef31b24SMark F. Adams PetscClassIdRegister("GAMG Setup",&cookie); 1234b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); 12352a44bfbaSMark F. Adams PetscLogEventRegister(" Graph", cookie, &gamg_setup_events[GRAPH]); 1236*c8b0795cSMark F. Adams /* PetscLogEventRegister(" G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); */ 1237*c8b0795cSMark F. Adams /* PetscLogEventRegister(" G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); */ 1238*c8b0795cSMark F. Adams /* PetscLogEventRegister(" G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); */ 1239b4fbaa2aSMark F. Adams PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); 1240b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); 1241b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); 1242b4fbaa2aSMark F. Adams PetscLogEventRegister(" search&set", cookie, &gamg_setup_events[FIND_V]); 1243*c8b0795cSMark F. Adams PetscLogEventRegister(" SA: colect data", cookie, &gamg_setup_events[SET7]); 1244*c8b0795cSMark F. Adams PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); 1245b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); 1246b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); 1247ed3f9983SMark F. Adams PetscLogEventRegister(" repartition", cookie, &gamg_setup_events[SET12]); 1248ed3f9983SMark F. Adams PetscLogEventRegister(" Invert-Sort", cookie, &gamg_setup_events[SET13]); 1249ed3f9983SMark F. Adams PetscLogEventRegister(" Move A", cookie, &gamg_setup_events[SET14]); 1250ed3f9983SMark F. Adams PetscLogEventRegister(" Move P", cookie, &gamg_setup_events[SET15]); 1251f852f58cSMark F. Adams 1252b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ 1253b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ 1254b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ 1255b4fbaa2aSMark F. Adams /* create timer stages */ 1256b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1257b4fbaa2aSMark F. Adams { 1258b4fbaa2aSMark F. Adams char str[32]; 1259b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d (finest)",0); 1260b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[0]); 1261b4fbaa2aSMark F. Adams PetscInt lidx; 1262b4fbaa2aSMark F. Adams for (lidx=1;lidx<9;lidx++){ 1263b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1264b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[lidx]); 1265b4fbaa2aSMark F. Adams } 1266b4fbaa2aSMark F. Adams } 1267b4fbaa2aSMark F. Adams #endif 1268b4fbaa2aSMark F. Adams } 1269b4fbaa2aSMark F. Adams #endif 12705b89ad90SMark F. Adams PetscFunctionReturn(0); 12715b89ad90SMark F. Adams } 12725b89ad90SMark F. Adams EXTERN_C_END 1273