15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4b45d2f2cSJed Brown #include "petsc-private/matimpl.h" 5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h> 7f96513f1SMatthew G Knepley 80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 10b4fbaa2aSMark F. Adams #endif 110cbbd2e1SMark F. Adams 120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 130cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG; 140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO; 150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG; 20a2f3521dSMark F. Adams PetscLogEvent PC_GAMGKKTProl_AGG; 210cbbd2e1SMark F. Adams #endif 220cbbd2e1SMark F. Adams 23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 24b4fbaa2aSMark F. Adams 25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 28b4fbaa2aSMark F. Adams #endif 29f96513f1SMatthew G Knepley 309d5b6da9SMark F. Adams static PetscFList GAMGList = 0; 319d5b6da9SMark F. Adams 32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 33d3d6bff4SMark F. Adams #undef __FUNCT__ 34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 36d3d6bff4SMark F. Adams { 37d3d6bff4SMark F. Adams PetscErrorCode ierr; 38d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 39d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 40d3d6bff4SMark F. Adams 41d3d6bff4SMark F. Adams PetscFunctionBegin; 42a2f3521dSMark F. Adams if( pc_gamg->data ) { /* this should not happen, cleaned up in SetUp */ 439d5b6da9SMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 4458471d46SMark F. Adams } 45a2f3521dSMark F. Adams pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0; 46a2f3521dSMark F. Adams PetscFunctionReturn(0); 47a2f3521dSMark F. Adams } 48a2f3521dSMark F. Adams 49a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */ 50a2f3521dSMark F. Adams typedef struct{ 51a2f3521dSMark F. Adams Mat A11,A21,A12,Amat; 52a2f3521dSMark F. Adams IS prim_is,constr_is; 53a2f3521dSMark F. Adams }GAMGKKTMat; 54a2f3521dSMark F. Adams 55a2f3521dSMark F. Adams #undef __FUNCT__ 56a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate" 57a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate( Mat A, PetscBool iskkt, GAMGKKTMat *out ) 58a2f3521dSMark F. Adams { 59a2f3521dSMark F. Adams PetscFunctionBegin; 60a2f3521dSMark F. Adams out->Amat = A; 61a2f3521dSMark F. Adams if( iskkt ){ 62a2f3521dSMark F. Adams PetscErrorCode ierr; 63a2f3521dSMark F. Adams IS is_constraint, is_prime; 64a2f3521dSMark F. Adams PetscInt nmin,nmax; 65a2f3521dSMark F. Adams 66a2f3521dSMark F. Adams ierr = MatGetOwnershipRange( A, &nmin, &nmax ); CHKERRQ(ierr); 67a2f3521dSMark F. Adams ierr = MatFindZeroDiagonals( A, &is_constraint ); CHKERRQ(ierr); 68a2f3521dSMark F. Adams ierr = ISComplement( is_constraint, nmin, nmax, &is_prime ); CHKERRQ(ierr); 69a2f3521dSMark F. Adams out->prim_is = is_prime; 70a2f3521dSMark F. Adams out->constr_is = is_constraint; 71a2f3521dSMark F. Adams 72a2f3521dSMark F. Adams ierr = MatGetSubMatrix( A, is_prime, is_prime, MAT_INITIAL_MATRIX, &out->A11); CHKERRQ(ierr); 73a2f3521dSMark F. Adams ierr = MatGetSubMatrix( A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12); CHKERRQ(ierr); 74a2f3521dSMark F. Adams ierr = MatGetSubMatrix( A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21); CHKERRQ(ierr); 759d1b15c3SMark F. Adams PetscPrintf(((PetscObject)A)->comm,"[%d]%s N=%d N_11=%d\n",0,__FUNCT__,A->rmap->N,out->A11->rmap->N); 76a2f3521dSMark F. Adams } 77a2f3521dSMark F. Adams else { 78a2f3521dSMark F. Adams out->A11 = A; 79a2f3521dSMark F. Adams out->A21 = PETSC_NULL; 80a2f3521dSMark F. Adams out->A12 = PETSC_NULL; 81a2f3521dSMark F. Adams out->prim_is = PETSC_NULL; 82a2f3521dSMark F. Adams out->constr_is = PETSC_NULL; 83a2f3521dSMark F. Adams } 84a2f3521dSMark F. Adams PetscFunctionReturn(0); 85a2f3521dSMark F. Adams } 86a2f3521dSMark F. Adams 87a2f3521dSMark F. Adams #undef __FUNCT__ 88a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy" 89a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy( GAMGKKTMat *mat ) 90a2f3521dSMark F. Adams { 91a2f3521dSMark F. Adams PetscErrorCode ierr; 92a2f3521dSMark F. Adams 93a2f3521dSMark F. Adams PetscFunctionBegin; 94a2f3521dSMark F. Adams if( mat->A11 && mat->A11 != mat->Amat ) { 95a2f3521dSMark F. Adams ierr = MatDestroy( &mat->A11 ); CHKERRQ(ierr); 96a2f3521dSMark F. Adams } 97a2f3521dSMark F. Adams ierr = MatDestroy( &mat->A21 ); CHKERRQ(ierr); 98a2f3521dSMark F. Adams ierr = MatDestroy( &mat->A12 ); CHKERRQ(ierr); 99a2f3521dSMark F. Adams 100a2f3521dSMark F. Adams ierr = ISDestroy( &mat->prim_is ); CHKERRQ(ierr); 101a2f3521dSMark F. Adams ierr = ISDestroy( &mat->constr_is ); CHKERRQ(ierr); 102a2f3521dSMark F. Adams 103d3d6bff4SMark F. Adams PetscFunctionReturn(0); 104d3d6bff4SMark F. Adams } 105d3d6bff4SMark F. Adams 1065b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1075b89ad90SMark F. Adams /* 108a147abb0SMark F. Adams createLevel: create coarse op with RAP. repartition and/or reduce number 109a147abb0SMark F. Adams of active processors. 1105b89ad90SMark F. Adams 1115b89ad90SMark F. Adams Input Parameter: 112a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 113a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 1149d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 1159d5b6da9SMark F. Adams . cbs - coarse block size 1169d5b6da9SMark F. Adams . isLast - 117a2f3521dSMark F. Adams . stokes - 1183530afc2SMark F. Adams In/Output Parameter: 119a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 120afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 12111e60469SMark F. Adams Output Parameter: 1223530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 1235b89ad90SMark F. Adams */ 1245cb416c2SMark F. Adams 1255b89ad90SMark F. Adams #undef __FUNCT__ 1268263b398SMark F. Adams #define __FUNCT__ "createLevel" 1270cbbd2e1SMark F. Adams static PetscErrorCode createLevel( const PC pc, 1289d5b6da9SMark F. Adams const Mat Amat_fine, 1299d5b6da9SMark F. Adams const PetscInt cbs, 1309d5b6da9SMark F. Adams const PetscBool isLast, 131a2f3521dSMark F. Adams const PetscBool stokes, 1323530afc2SMark F. Adams Mat *a_P_inout, 133a2f3521dSMark F. Adams Mat *a_Amat_crs, 134a2f3521dSMark F. Adams PetscMPIInt *a_nactive_proc 13511e60469SMark F. Adams ) 1365b89ad90SMark F. Adams { 137a2f3521dSMark F. Adams PetscErrorCode ierr; 1389d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 139486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1409d5b6da9SMark F. Adams const PetscBool repart = pc_gamg->repart; 1419d5b6da9SMark F. Adams const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 142a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 1439d5b6da9SMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat_fine)->comm; 1445a9b9e01SMark F. Adams PetscMPIInt mype,npe,new_npe,nactive=*a_nactive_proc; 145a2f3521dSMark F. Adams PetscInt ncrs_eq,ncrs_prim; 1465b89ad90SMark F. Adams 1475b89ad90SMark F. Adams PetscFunctionBegin; 14811e60469SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 14911e60469SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 15011e60469SMark F. Adams /* RAP */ 1519d5b6da9SMark F. Adams ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 152038e3b61SMark F. Adams 153a2f3521dSMark F. Adams /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/ 154a2f3521dSMark F. Adams ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 155a2f3521dSMark F. Adams assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0); 156a2f3521dSMark F. Adams ierr = MatGetLocalSize( Cmat, &ncrs_eq, PETSC_NULL ); CHKERRQ(ierr); 157a2f3521dSMark F. Adams 158a2f3521dSMark F. Adams /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */ 159a2f3521dSMark F. Adams { 160a2f3521dSMark F. Adams PetscInt ncrs_eq_glob,ncrs_eq_ave; 161a2f3521dSMark F. Adams ierr = MatGetSize( Cmat, &ncrs_eq_glob, PETSC_NULL ); CHKERRQ(ierr); 162a2f3521dSMark F. Adams ncrs_eq_ave = ncrs_eq_glob/npe; 163a2f3521dSMark F. Adams new_npe = (PetscMPIInt)((float)ncrs_eq_ave/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 164a2f3521dSMark F. Adams if( new_npe == 0 || ncrs_eq_ave < coarse_max ) new_npe = 1; 1655a9b9e01SMark F. Adams else if ( new_npe >= nactive ) new_npe = nactive; /* no change, rare */ 1669d5b6da9SMark F. Adams if( isLast ) new_npe = 1; 167a2f3521dSMark F. Adams } 168f852f58cSMark F. Adams 1695a9b9e01SMark F. Adams if( !repart && new_npe==nactive ) { 170a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 17122063be5SMark F. Adams } 17222063be5SMark F. Adams else { 173a2f3521dSMark F. Adams const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols; 174a2f3521dSMark F. Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old; 175a2f3521dSMark F. Adams IS is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices; 1765a9b9e01SMark F. Adams VecScatter vecscat; 17722063be5SMark F. Adams PetscScalar *array; 17822063be5SMark F. Adams Vec src_crd, dest_crd; 179e33ef3b1SMark F. Adams 180a2f3521dSMark F. Adams nloc_old = ncrs_eq/cbs; assert(ncrs_eq%cbs==0); 1810cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1820cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 183b4fbaa2aSMark F. Adams #endif 184a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 185a2f3521dSMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 186a2f3521dSMark F. Adams if( repart && !stokes ) { 187a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1885a9b9e01SMark F. Adams Mat adj; 1895a9b9e01SMark F. Adams 190a2f3521dSMark F. Adams if (pc_gamg->verbose>0) { 191a2f3521dSMark F. Adams if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq); 192a2f3521dSMark F. Adams else { 193a2f3521dSMark F. Adams PetscInt n; 194a2f3521dSMark F. Adams ierr = MPI_Allreduce( &ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm );CHKERRQ(ierr); 195a2f3521dSMark F. Adams PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n); 196a2f3521dSMark F. Adams } 197a2f3521dSMark F. Adams } 1985a9b9e01SMark F. Adams 199a2f3521dSMark F. Adams /* get 'adj' */ 2009d5b6da9SMark F. Adams if( cbs == 1 ) { 201038e3b61SMark F. Adams ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 202eb07cef2SMark F. Adams } 203eb07cef2SMark F. Adams else{ 204a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 205eb07cef2SMark F. Adams Mat tMat; 206a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 207b4fbaa2aSMark F. Adams const PetscScalar *vals; 208b4fbaa2aSMark F. Adams const PetscInt *idx; 209a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 2109057884aSMark F. Adams static PetscInt llev = 0; 211b4fbaa2aSMark F. Adams 212a2f3521dSMark F. Adams ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); 213a2f3521dSMark F. Adams ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr); 214a2f3521dSMark F. Adams ierr = MatGetOwnershipRange( Cmat, &Istart_crs, &Iend_crs ); CHKERRQ(ierr); 215a2f3521dSMark F. Adams ierr = MatGetSize( Cmat, &M, &N );CHKERRQ(ierr); 216a2f3521dSMark F. Adams for ( Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cbs, jj++ ) { 21758471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 2189d5b6da9SMark F. Adams d_nnz[jj] = ncols/cbs; 2199d5b6da9SMark F. Adams o_nnz[jj] = ncols/cbs; 22058471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 221a2f3521dSMark F. Adams if( d_nnz[jj] > ncrs_prim ) d_nnz[jj] = ncrs_prim; 222a2f3521dSMark F. Adams if( o_nnz[jj] > (M/cbs-ncrs_prim) ) o_nnz[jj] = M/cbs-ncrs_prim; 22358471d46SMark F. Adams } 2246876a03eSMark F. Adams 225a2f3521dSMark F. Adams ierr = MatCreate( wcomm, &tMat ); CHKERRQ(ierr); 226a2f3521dSMark F. Adams ierr = MatSetSizes( tMat, ncrs_prim, ncrs_prim, 227a2f3521dSMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE ); 2286876a03eSMark F. Adams CHKERRQ(ierr); 229a2f3521dSMark F. Adams ierr = MatSetType(tMat,MATAIJ); CHKERRQ(ierr); 230a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 231a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 23258471d46SMark F. Adams ierr = PetscFree( d_nnz ); CHKERRQ(ierr); 2335f8cf99dSMark F. Adams ierr = PetscFree( o_nnz ); CHKERRQ(ierr); 234eb07cef2SMark F. Adams 235a2f3521dSMark F. Adams for ( ii = Istart_crs; ii < Iend_crs; ii++ ) { 2369d5b6da9SMark F. Adams PetscInt dest_row = ii/cbs; 23722063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 238eb07cef2SMark F. Adams for( jj = 0 ; jj < ncols ; jj++ ){ 2399d5b6da9SMark F. Adams PetscInt dest_col = idx[jj]/cbs; 240eb07cef2SMark F. Adams PetscScalar v = 1.0; 241eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 242eb07cef2SMark F. Adams } 24322063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 244eb07cef2SMark F. Adams } 245eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 246eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 247eb07cef2SMark F. Adams 248b4fbaa2aSMark F. Adams if( llev++ == -1 ) { 249b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2504704b1adSJed Brown ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr); 251b4fbaa2aSMark F. Adams PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 252b4fbaa2aSMark F. Adams ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 253b4fbaa2aSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 254b4fbaa2aSMark F. Adams } 255b4fbaa2aSMark F. Adams 256eb07cef2SMark F. Adams ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 257eb07cef2SMark F. Adams 258eb07cef2SMark F. Adams ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 259a2f3521dSMark F. Adams } /* create 'adj' */ 260f150b916SMark F. Adams 261a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2625a9b9e01SMark F. Adams char prefix[256]; 2635a9b9e01SMark F. Adams const char *pcpre; 264b4fbaa2aSMark F. Adams const PetscInt *is_idx; 265b4fbaa2aSMark F. Adams MatPartitioning mpart; 266a4b7d37bSMark F. Adams IS proc_is; 267a2f3521dSMark F. Adams PetscInt targetPE; 2682f03bc48SMark F. Adams 2695a9b9e01SMark F. Adams ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr); 2705ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 2719d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr); 27259a0be82SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr); 27359a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 27411e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 2755ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 276a4b7d37bSMark F. Adams ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr); 27711e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 2785a9b9e01SMark F. Adams 2795ef31b24SMark F. Adams /* collect IS info */ 280a2f3521dSMark F. Adams ierr = PetscMalloc( ncrs_eq*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr); 281a4b7d37bSMark F. Adams ierr = ISGetIndices( proc_is, &is_idx ); CHKERRQ(ierr); 282a2f3521dSMark F. Adams targetPE = 1; /* bring to "front" of machine */ 283a2f3521dSMark F. Adams /*targetPE = npe/new_npe;*/ /* spread partitioning across machine */ 284a2f3521dSMark F. Adams for( kk = jj = 0 ; kk < nloc_old ; kk++ ){ 2859d5b6da9SMark F. Adams for( ii = 0 ; ii < cbs ; ii++, jj++ ){ 286a2f3521dSMark F. Adams newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 287eb07cef2SMark F. Adams } 2885ef31b24SMark F. Adams } 289a4b7d37bSMark F. Adams ierr = ISRestoreIndices( proc_is, &is_idx ); CHKERRQ(ierr); 290a4b7d37bSMark F. Adams ierr = ISDestroy( &proc_is ); CHKERRQ(ierr); 2915ef31b24SMark F. Adams } 2925ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 2935a9b9e01SMark F. Adams 294a2f3521dSMark F. Adams ierr = ISCreateGeneral( wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc ); 2955a9b9e01SMark F. Adams CHKERRQ(ierr); 2968263b398SMark F. Adams if( newproc_idx != 0 ) { 2978263b398SMark F. Adams ierr = PetscFree( newproc_idx ); CHKERRQ(ierr); 2985ef31b24SMark F. Adams } 299a2f3521dSMark F. Adams } /* repartitioning */ 300a2f3521dSMark F. Adams else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 301a2f3521dSMark F. Adams 302a2f3521dSMark F. Adams PetscInt rfactor,targetPE; 3035a9b9e01SMark F. Adams /* find factor */ 3045a9b9e01SMark F. Adams if( new_npe == 1 ) rfactor = npe; /* easy */ 3055a9b9e01SMark F. Adams else { 3065a9b9e01SMark F. Adams PetscReal best_fact = 0.; 3075a9b9e01SMark F. Adams jj = -1; 3085a9b9e01SMark F. Adams for( kk = 1 ; kk <= npe ; kk++ ){ 3095a9b9e01SMark F. Adams if( npe%kk==0 ) { /* a candidate */ 3105a9b9e01SMark F. Adams PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe; 3115a9b9e01SMark F. Adams if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 3125a9b9e01SMark F. Adams if( fact > best_fact ) { 3135a9b9e01SMark F. Adams best_fact = fact; jj = kk; 3145a9b9e01SMark F. Adams } 3155a9b9e01SMark F. Adams } 3165a9b9e01SMark F. Adams } 3175a9b9e01SMark F. Adams if( jj != -1 ) rfactor = jj; 318a2f3521dSMark F. Adams else rfactor = 1; /* does this happen .. a prime */ 3195a9b9e01SMark F. Adams } 3205a9b9e01SMark F. Adams new_npe = npe/rfactor; 3215a9b9e01SMark F. Adams 3225a9b9e01SMark F. Adams if( new_npe==nactive ) { 323a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 3245a9b9e01SMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 325a2f3521dSMark F. Adams if (pc_gamg->verbose>0){ 326a2f3521dSMark F. Adams PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq); 327a2f3521dSMark F. Adams } 3285a9b9e01SMark F. Adams PetscFunctionReturn(0); 3295a9b9e01SMark F. Adams } 3305a9b9e01SMark F. Adams 331a2f3521dSMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq); 3325a9b9e01SMark F. Adams targetPE = mype/rfactor; 333a2f3521dSMark F. Adams ierr = ISCreateStride( wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc ); CHKERRQ(ierr); 3345a9b9e01SMark F. Adams 335a2f3521dSMark F. Adams if( stokes ) { 336a2f3521dSMark F. Adams ierr = ISCreateStride( wcomm, ncrs_prim*cbs, targetPE, 0, &is_eq_newproc_prim ); CHKERRQ(ierr); 3375a9b9e01SMark F. Adams } 338a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 339e33ef3b1SMark F. Adams 34011e60469SMark F. Adams /* 341a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 34211e60469SMark F. Adams */ 343a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering( is_eq_newproc, &is_eq_num ); CHKERRQ(ierr); 344a2f3521dSMark F. Adams if( stokes ) { 345a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering( is_eq_newproc_prim, &is_eq_num_prim ); CHKERRQ(ierr); 346a2f3521dSMark F. Adams } 347a2f3521dSMark F. Adams else is_eq_num_prim = is_eq_num; 34811e60469SMark F. Adams /* 349a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 35011e60469SMark F. Adams */ 351a2f3521dSMark F. Adams ierr = ISPartitioningCount( is_eq_newproc, npe, counts ); CHKERRQ(ierr); 352a2f3521dSMark F. Adams ncrs_eq_new = counts[mype]; 353a2f3521dSMark F. Adams ierr = ISDestroy( &is_eq_newproc ); CHKERRQ(ierr); 354a2f3521dSMark F. Adams if( stokes ) { 355a2f3521dSMark F. Adams ierr = ISPartitioningCount( is_eq_newproc_prim, npe, counts ); CHKERRQ(ierr); 356a2f3521dSMark F. Adams ierr = ISDestroy( &is_eq_newproc_prim ); CHKERRQ(ierr); 357a2f3521dSMark F. Adams ncrs_prim_new = counts[mype]/cbs; /* nodes */ 358a2f3521dSMark F. Adams } 359a2f3521dSMark F. Adams else ncrs_prim_new = ncrs_eq_new/cbs; /* eqs */ 360a2f3521dSMark F. Adams 361a2f3521dSMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 3620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3630cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 364b4fbaa2aSMark F. Adams #endif 365a2f3521dSMark F. Adams 366a2f3521dSMark F. Adams /* move data (for primal equations only) */ 36722063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 36811e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 369a2f3521dSMark F. Adams ierr = VecSetSizes( dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE ); CHKERRQ(ierr); 370470e26f8SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ 37111e60469SMark F. Adams /* 3729d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 3739d5b6da9SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cbs'. 37411e60469SMark F. Adams */ 375a2f3521dSMark F. Adams ierr = PetscMalloc( (ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 376a2f3521dSMark F. Adams ierr = ISGetIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr); 377a2f3521dSMark F. Adams for(ii=0,jj=0; ii<ncrs_prim ; ii++) { 3789d5b6da9SMark F. Adams PetscInt id = idx[ii*cbs]/cbs; /* get node back */ 379a2f3521dSMark F. Adams for( kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 38011e60469SMark F. Adams } 381a2f3521dSMark F. Adams ierr = ISRestoreIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr); 382a2f3521dSMark F. Adams ierr = ISCreateGeneral( wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat ); 38311e60469SMark F. Adams CHKERRQ(ierr); 38492a756f0SMark F. Adams ierr = PetscFree( tidx ); CHKERRQ(ierr); 38511e60469SMark F. Adams /* 38611e60469SMark F. Adams Create a vector to contain the original vertex information for each element 38711e60469SMark F. Adams */ 388a2f3521dSMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd ); CHKERRQ(ierr); 3899d5b6da9SMark F. Adams for( jj=0; jj<ndata_cols ; jj++ ) { 390a2f3521dSMark F. Adams const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows; 391a2f3521dSMark F. Adams for( ii=0 ; ii<ncrs_prim ; ii++) { 3929d5b6da9SMark F. Adams for( kk=0; kk<ndata_rows ; kk++ ) { 393a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 394c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 395676e1743SMark F. Adams ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES ); CHKERRQ(ierr); 396d3d6bff4SMark F. Adams } 397038e3b61SMark F. Adams } 398eb07cef2SMark F. Adams } 399eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 400eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 40111e60469SMark F. Adams /* 40211e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 40311e60469SMark F. Adams to the correct processor 40411e60469SMark F. Adams */ 40511e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 40611e60469SMark F. Adams CHKERRQ(ierr); 40711e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 40811e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40911e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 41011e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 41111e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 41211e60469SMark F. Adams /* 41311e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 41411e60469SMark F. Adams */ 415c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 416a2f3521dSMark F. Adams ierr = PetscMalloc( node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 417a2f3521dSMark F. Adams pc_gamg->data_sz = node_data_sz*ncrs_prim_new; 418a2f3521dSMark F. Adams strideNew = ncrs_prim_new*ndata_rows; 41911e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 4209d5b6da9SMark F. Adams for( jj=0; jj<ndata_cols ; jj++ ) { 421a2f3521dSMark F. Adams for( ii=0 ; ii<ncrs_prim_new ; ii++) { 4229d5b6da9SMark F. Adams for( kk=0; kk<ndata_rows ; kk++ ) { 423a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 424c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 425d3d6bff4SMark F. Adams } 426038e3b61SMark F. Adams } 427038e3b61SMark F. Adams } 42811e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 42911e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 430a2f3521dSMark F. Adams 431a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 4320cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4330cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 434ed3f9983SMark F. Adams #endif 435a2f3521dSMark F. Adams 43611e60469SMark F. Adams /* 43711e60469SMark F. Adams Invert for MatGetSubMatrix 43811e60469SMark F. Adams */ 439a2f3521dSMark F. Adams ierr = ISInvertPermutation( is_eq_num, ncrs_eq_new, &new_eq_indices ); CHKERRQ(ierr); 440a2f3521dSMark F. Adams ierr = ISSort( new_eq_indices ); CHKERRQ(ierr); /* is this needed? */ 441a2f3521dSMark F. Adams if(is_eq_num != is_eq_num_prim) { 442a2f3521dSMark F. Adams ierr = ISDestroy( &is_eq_num_prim ); CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 443a2f3521dSMark F. Adams } 444a2f3521dSMark F. Adams ierr = ISDestroy( &is_eq_num ); CHKERRQ(ierr); 4450cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4460cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 4470cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 448ed3f9983SMark F. Adams #endif 449a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 450a2f3521dSMark F. Adams { 451a2f3521dSMark F. Adams Mat mat; 452a2f3521dSMark F. Adams ierr = MatGetSubMatrix( Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat ); 45311e60469SMark F. Adams CHKERRQ(ierr); 454a2f3521dSMark F. Adams *a_Amat_crs = mat; 455a2f3521dSMark F. Adams } 456038e3b61SMark F. Adams ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 457a2f3521dSMark F. Adams 4580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4590cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 460ed3f9983SMark F. Adams #endif 46111e60469SMark F. Adams /* prolongator */ 46211e60469SMark F. Adams { 46311e60469SMark F. Adams IS findices; 464a2f3521dSMark F. Adams PetscInt Istart,Iend; 465a2f3521dSMark F. Adams Mat Pnew; 466a2f3521dSMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 4670cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4680cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 469ed3f9983SMark F. Adams #endif 47011e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 471a2f3521dSMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew ); 47211e60469SMark F. Adams CHKERRQ(ierr); 47311e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 4740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4750cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 476ed3f9983SMark F. Adams #endif 4773530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 478a2f3521dSMark F. Adams 479a2f3521dSMark F. Adams /* output - repartitioned */ 480a2f3521dSMark F. Adams *a_P_inout = Pnew; 481e33ef3b1SMark F. Adams } 482a2f3521dSMark F. Adams ierr = ISDestroy( &new_eq_indices ); CHKERRQ(ierr); 4835b89ad90SMark F. Adams 4845a9b9e01SMark F. Adams *a_nactive_proc = new_npe; /* output */ 485a2f3521dSMark F. Adams } 4865a9b9e01SMark F. Adams 487a2f3521dSMark F. Adams /* outout matrix data */ 488c8b0795cSMark F. Adams if( !PETSC_TRUE ) { 489c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 490c8b0795cSMark F. Adams if(llev==0) { 491c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 492c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 493c8b0795cSMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 494c8b0795cSMark F. Adams ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr); 495c8b0795cSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 496c8b0795cSMark F. Adams } 497c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 498c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 499c8b0795cSMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 500c8b0795cSMark F. Adams ierr = MatView(Cmat, viewer ); CHKERRQ(ierr); 501c8b0795cSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 502c8b0795cSMark F. Adams } 503c8b0795cSMark F. Adams 5045b89ad90SMark F. Adams PetscFunctionReturn(0); 5055b89ad90SMark F. Adams } 5065b89ad90SMark F. Adams 5075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 5085b89ad90SMark F. Adams /* 5095b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 5105b89ad90SMark F. Adams by setting data structures and options. 5115b89ad90SMark F. Adams 5125b89ad90SMark F. Adams Input Parameter: 5135b89ad90SMark F. Adams . pc - the preconditioner context 5145b89ad90SMark F. Adams 5155b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 5165b89ad90SMark F. Adams 5175b89ad90SMark F. Adams Notes: 5185b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 5195b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 5205b89ad90SMark F. Adams */ 5215b89ad90SMark F. Adams #undef __FUNCT__ 5225b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 5239d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc ) 5245b89ad90SMark F. Adams { 5255b89ad90SMark F. Adams PetscErrorCode ierr; 5269d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 5275b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 5282adcac29SMark F. Adams Mat Pmat = pc->pmat; 529a2f3521dSMark F. Adams PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 5309d5b6da9SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 5313530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 532587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 533c8b0795cSMark F. Adams PetscReal emaxs[GAMG_MAXLEVELS]; 534a94c3b12SMark F. Adams IS *ASMLocalIDsArr[GAMG_MAXLEVELS],removedEqs[GAMG_MAXLEVELS]; 535a94c3b12SMark F. Adams PetscInt level_bs[GAMG_MAXLEVELS]; 536a2f3521dSMark F. Adams GAMGKKTMat kktMatsArr[GAMG_MAXLEVELS]; 537a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 538737a81a9SMark F. Adams MatInfo info; 539302f38e8SMark F. Adams PetscBool stokes = PETSC_FALSE; 5405ef31b24SMark F. Adams 5415b89ad90SMark F. Adams PetscFunctionBegin; 54284d3f75bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 54384d3f75bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 544a2f3521dSMark F. Adams if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled); 54584d3f75bSMark F. Adams if( pc_gamg->setup_count++ > 0 ) { 54684d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 54703a628feSMark F. Adams /* just do Galerkin grids */ 54858471d46SMark F. Adams Mat B,dA,dB; 549d5280255SMark F. Adams assert(pc->setupcalled); 55058471d46SMark F. Adams 5519d5b6da9SMark F. Adams if( pc_gamg->Nlevels > 1 ) { 55258471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5539d5b6da9SMark F. Adams ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 55458471d46SMark F. Adams /* (re)set to get dirty flag */ 5559d5b6da9SMark F. Adams ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 55658471d46SMark F. Adams 5579d5b6da9SMark F. Adams for (level=pc_gamg->Nlevels-2; level>-1; level--) { 55803a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 55984d3f75bSMark F. Adams if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) { 56003a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 561084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 56203a628feSMark F. Adams mglevels[level]->A = B; 56303a628feSMark F. Adams } 56403a628feSMark F. Adams else { 565084a8fe3SJed Brown ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 56658471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 56703a628feSMark F. Adams } 56858471d46SMark F. Adams ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 56958471d46SMark F. Adams dB = B; 57058471d46SMark F. Adams } 5715f8cf99dSMark F. Adams } 572d5280255SMark F. Adams 573d5280255SMark F. Adams ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 574d5280255SMark F. Adams 575d5280255SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 576d5280255SMark F. Adams ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 577d5280255SMark F. Adams 57858471d46SMark F. Adams PetscFunctionReturn(0); 579eb07cef2SMark F. Adams } 58084d3f75bSMark F. Adams assert(pc->setupcalled == 0); 581f6536408SMark F. Adams 582302f38e8SMark F. Adams ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 583eb07cef2SMark F. Adams 5849d1b15c3SMark F. Adams ierr = GAMGKKTMatCreate( Pmat, stokes, &kktMatsArr[0] ); CHKERRQ(ierr); 5859d1b15c3SMark F. Adams 586302f38e8SMark F. Adams if( pc_gamg->data == 0 ) { 5879d5b6da9SMark F. Adams if( !pc_gamg->createdefaultdata ){ 588a2f3521dSMark F. Adams SETERRQ(wcomm,PETSC_ERR_LIB,"'createdefaultdata' not set(?) need to support NULL data"); 589302f38e8SMark F. Adams } 590302f38e8SMark F. Adams if( stokes ) { 591302f38e8SMark F. Adams SETERRQ(wcomm,PETSC_ERR_LIB,"Need data (eg, PCSetCoordinates) for Stokes problems"); 592eb07cef2SMark F. Adams } 5939d1b15c3SMark F. Adams ierr = pc_gamg->createdefaultdata( pc, kktMatsArr[0].A11 ); CHKERRQ(ierr); 5949d5b6da9SMark F. Adams } 595038e3b61SMark F. Adams 596302f38e8SMark F. Adams /* get basic dims */ 597302f38e8SMark F. Adams if( stokes ) { 598a2f3521dSMark F. Adams bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */ 599302f38e8SMark F. Adams } 600302f38e8SMark F. Adams else { 601302f38e8SMark F. Adams ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr); 602302f38e8SMark F. Adams } 603a2f3521dSMark F. Adams 604a2f3521dSMark F. Adams ierr = MatGetSize( Pmat, &M, &qq );CHKERRQ(ierr); 605c8b0795cSMark F. Adams if (pc_gamg->verbose) { 6062adcac29SMark F. Adams if(pc_gamg->verbose==1) ierr = MatGetInfo(Pmat,MAT_LOCAL,&info); 6072adcac29SMark F. Adams else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); 608b2a4f308SMark F. Adams CHKERRQ(ierr); 609b2a4f308SMark F. Adams nnz0 = info.nz_used; 610b2a4f308SMark F. Adams nnztot = info.nz_used; 611b43b03e9SMark F. Adams PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 612a2f3521dSMark F. Adams mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 613a2f3521dSMark F. Adams (int)(nnz0/(PetscReal)M),npe); 614c8b0795cSMark F. Adams } 61584d3f75bSMark F. Adams 616a2f3521dSMark F. Adams /* Get A_i and R_i */ 6178f4b7eb5SMark F. Adams for ( level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */ 6189d5b6da9SMark F. Adams level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 6190205a208SMark F. Adams level++ ){ 6205b89ad90SMark F. Adams level1 = level + 1; 6210cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6220cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr); 623a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 624a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 625b4fbaa2aSMark F. Adams #endif 626a2f3521dSMark F. Adams #endif 627a2f3521dSMark F. Adams /* deal with Stokes, get sub matrices */ 6289d1b15c3SMark F. Adams if( level > 0 ) { 629a2f3521dSMark F. Adams ierr = GAMGKKTMatCreate( Aarr[level], stokes, &kktMatsArr[level] ); CHKERRQ(ierr); 6309d1b15c3SMark F. Adams } 631c8b0795cSMark F. Adams { /* construct prolongator */ 632725b86d8SJed Brown Mat Gmat; 6330cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 634a2f3521dSMark F. Adams Mat Prol11,Prol22; 635c8b0795cSMark F. Adams 636a94c3b12SMark F. Adams level_bs[level] = bs; 637a2f3521dSMark F. Adams ierr = pc_gamg->graph( pc,kktMatsArr[level].A11, &Gmat ); CHKERRQ(ierr); 6380cbbd2e1SMark F. Adams ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr); 639a2f3521dSMark F. Adams ierr = pc_gamg->prolongator( pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11 ); CHKERRQ(ierr); 640c8b0795cSMark F. Adams 641a2f3521dSMark F. Adams /* could have failed to create new level */ 642a2f3521dSMark F. Adams if( Prol11 ){ 6439d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6449d1b15c3SMark F. Adams ierr = MatGetBlockSizes( Prol11, PETSC_NULL, &bs ); CHKERRQ(ierr); 645a2f3521dSMark F. Adams 646a2f3521dSMark F. Adams if( stokes ) { 647a2f3521dSMark F. Adams if(!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method."); 648a2f3521dSMark F. Adams /* R A12 == (T = A21 P)'; G = T' T; coarsen G; form plain agg with G */ 649a2f3521dSMark F. Adams ierr = pc_gamg->formkktprol( pc, Prol11, kktMatsArr[level].A21, &Prol22 ); CHKERRQ(ierr); 650c8b0795cSMark F. Adams } 651c8b0795cSMark F. Adams 652a2f3521dSMark F. Adams if( pc_gamg->optprol ){ 653c8b0795cSMark F. Adams /* smooth */ 654a2f3521dSMark F. Adams ierr = pc_gamg->optprol( pc, kktMatsArr[level].A11, &Prol11 ); CHKERRQ(ierr); 655c8b0795cSMark F. Adams } 656c8b0795cSMark F. Adams 657a2f3521dSMark F. Adams if( stokes ) { 658a2f3521dSMark F. Adams IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is}; 659a2f3521dSMark F. Adams Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 }; 660a2f3521dSMark F. Adams ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1] ); CHKERRQ(ierr); 661a2f3521dSMark F. Adams } 662a2f3521dSMark F. Adams else { 663a2f3521dSMark F. Adams Parr[level1] = Prol11; 664a2f3521dSMark F. Adams } 665a2f3521dSMark F. Adams } 666a2f3521dSMark F. Adams else Parr[level1] = PETSC_NULL; 667ffc955d6SMark F. Adams 668ffc955d6SMark F. Adams if ( pc_gamg->use_aggs_in_gasm ) { 669a94c3b12SMark F. Adams ierr = PetscCDGetASMBlocks(agg_lists, level_bs[level], &nASMBlocksArr[level], &ASMLocalIDsArr[level] ); CHKERRQ(ierr); 670ffc955d6SMark F. Adams } 671ffc955d6SMark F. Adams 672a94c3b12SMark F. Adams ierr = PetscCDGetRemovedIS( agg_lists, &removedEqs[level] ); CHKERRQ(ierr); 673a94c3b12SMark F. Adams 674a2f3521dSMark F. Adams ierr = MatDestroy( &Gmat ); CHKERRQ(ierr); 67541b27cdeSMark F. Adams ierr = PetscCDDestroy( agg_lists ); CHKERRQ(ierr); 676a2f3521dSMark F. Adams } /* construct prolongator scope */ 6770cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6780cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 679c8b0795cSMark F. Adams #endif 6809d5b6da9SMark F. Adams /* cache eigen estimate */ 6819d5b6da9SMark F. Adams if( pc_gamg->emax_id != -1 ){ 6829d5b6da9SMark F. Adams PetscBool flag; 683a2f3521dSMark F. Adams ierr = PetscObjectComposedDataGetReal( (PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag ); 6849d5b6da9SMark F. Adams CHKERRQ( ierr ); 6859d5b6da9SMark F. Adams if( !flag ) emaxs[level] = -1.; 6869d5b6da9SMark F. Adams } 6879d5b6da9SMark F. Adams else emaxs[level] = -1.; 6882adcac29SMark F. Adams if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 689c8b0795cSMark F. Adams if( !Parr[level1] ) { 690b43b03e9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level); 691c8b0795cSMark F. Adams break; 692c8b0795cSMark F. Adams } 6930cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6940cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 695b4fbaa2aSMark F. Adams #endif 696a2f3521dSMark F. Adams 697a2f3521dSMark F. Adams ierr = createLevel( pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2), 698a2f3521dSMark F. Adams stokes, &Parr[level1], &Aarr[level1], &nactivepe ); 6993530afc2SMark F. Adams CHKERRQ(ierr); 700a2f3521dSMark F. Adams 7010cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7020cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 703b4fbaa2aSMark F. Adams #endif 704a2f3521dSMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &qq );CHKERRQ(ierr); 705a2f3521dSMark F. Adams 706a2f3521dSMark F. Adams if (pc_gamg->verbose > 0){ 7070cbbd2e1SMark F. Adams PetscInt NN = M; 7080cbbd2e1SMark F. Adams if(pc_gamg->verbose==1) { 709a2f3521dSMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); CHKERRQ(ierr); 710a2f3521dSMark F. Adams ierr = MatGetLocalSize( Aarr[level1], &NN, &qq ); 7110cbbd2e1SMark F. Adams } 7120cbbd2e1SMark F. Adams else ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); 713a2f3521dSMark F. Adams 7140cbbd2e1SMark F. Adams CHKERRQ(ierr); 7150cbbd2e1SMark F. Adams nnztot += info.nz_used; 716b43b03e9SMark F. Adams PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 717a2f3521dSMark F. Adams mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols, 7180cbbd2e1SMark F. Adams (int)(info.nz_used/(PetscReal)NN), nactivepe ); 7190cbbd2e1SMark F. Adams CHKERRQ(ierr); 720c8b0795cSMark F. Adams } 721a2f3521dSMark F. Adams 722a2f3521dSMark F. Adams /* stop if one node -- could pull back for singular problems */ 723c8b0795cSMark F. Adams if( M/pc_gamg->data_cell_cols < 2 ) { 72481708dccSMark F. Adams level++; 72581708dccSMark F. Adams break; 72681708dccSMark F. Adams } 7270cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 728b4fbaa2aSMark F. Adams ierr = PetscLogStagePop(); CHKERRQ( ierr ); 729b4fbaa2aSMark F. Adams #endif 730c8b0795cSMark F. Adams } /* levels */ 731c8b0795cSMark F. Adams 732c8b0795cSMark F. Adams if( pc_gamg->data ) { 733c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 734a2f3521dSMark F. Adams pc_gamg->data = PETSC_NULL; 7355b89ad90SMark F. Adams } 736c8b0795cSMark F. Adams 7370cbbd2e1SMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 7389d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 7395b89ad90SMark F. Adams fine_level = level; 7409d5b6da9SMark F. Adams ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr); 7415b89ad90SMark F. Adams 74284d3f75bSMark F. Adams /* simple setup */ 74384d3f75bSMark F. Adams if( !PETSC_TRUE ){ 74484d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 74584d3f75bSMark F. Adams for (lidx=0,level=pc_gamg->Nlevels-1; 74684d3f75bSMark F. Adams lidx<fine_level; 74784d3f75bSMark F. Adams lidx++, level--){ 74884d3f75bSMark F. Adams ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr); 74984d3f75bSMark F. Adams ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr); 75084d3f75bSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 75184d3f75bSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 75284d3f75bSMark F. Adams } 75384d3f75bSMark F. Adams ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 75484d3f75bSMark F. Adams 75584d3f75bSMark F. Adams ierr = PCSetUp_MG( pc ); CHKERRQ( ierr ); 75684d3f75bSMark F. Adams } 75784d3f75bSMark F. Adams else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */ 758d5280255SMark F. Adams /* set default smoothers & set operators */ 7599d5b6da9SMark F. Adams for ( lidx = 1, level = pc_gamg->Nlevels-2; 760587fa25dSMark F. Adams lidx <= fine_level; 761587fa25dSMark F. Adams lidx++, level--) { 762ffc955d6SMark F. Adams KSP smoother; 763ffc955d6SMark F. Adams PC subpc; 764a2f3521dSMark F. Adams 7659d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 766f6536408SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 767ffc955d6SMark F. Adams 768a2f3521dSMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 769a2f3521dSMark F. Adams /* set ops */ 770a2f3521dSMark F. Adams ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 771a2f3521dSMark F. Adams ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr); 772a2f3521dSMark F. Adams 773a2f3521dSMark F. Adams /* create field split PC, get subsmoother */ 774a2f3521dSMark F. Adams if( stokes ) { 775a2f3521dSMark F. Adams KSP *ksps; 776a2f3521dSMark F. Adams PetscInt nn; 777a2f3521dSMark F. Adams ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is); CHKERRQ(ierr); 778a2f3521dSMark F. Adams ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is); CHKERRQ(ierr); 779a2f3521dSMark F. Adams ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr); 780a2f3521dSMark F. Adams smoother = ksps[0]; 781a2f3521dSMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 782a2f3521dSMark F. Adams ierr = PetscFree( ksps ); CHKERRQ(ierr); 783a2f3521dSMark F. Adams } 784a2f3521dSMark F. Adams ierr = GAMGKKTMatDestroy( &kktMatsArr[level] ); CHKERRQ(ierr); 785a2f3521dSMark F. Adams 786a2f3521dSMark F. Adams /* set defaults */ 787a2f3521dSMark F. Adams ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 788a2f3521dSMark F. Adams 789ffc955d6SMark F. Adams /* override defaults and command line args (!) */ 790ffc955d6SMark F. Adams if ( pc_gamg->use_aggs_in_gasm ) { 7912d3561bbSSatish Balay PetscInt sz; 7922d3561bbSSatish Balay IS *is; 793a2f3521dSMark F. Adams 7942d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7952d3561bbSSatish Balay is = ASMLocalIDsArr[level]; 796ffc955d6SMark F. Adams ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr); 797ffc955d6SMark F. Adams if(sz==0){ 798ffc955d6SMark F. Adams IS is; 799ffc955d6SMark F. Adams PetscInt my0,kk; 800ffc955d6SMark F. Adams ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr); 801ffc955d6SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr); 802*06b43e7eSDmitry Karpeev ierr = PCGASMSetSubdomains( subpc, 1, &is, PETSC_NULL ); CHKERRQ(ierr); 803a94c3b12SMark F. Adams ierr = ISDestroy( &is ); CHKERRQ(ierr); 804ffc955d6SMark F. Adams } 805ffc955d6SMark F. Adams else { 806a94c3b12SMark F. Adams PetscInt kk; 807*06b43e7eSDmitry Karpeev ierr = PCGASMSetSubdomains( subpc, sz, is, PETSC_NULL ); CHKERRQ(ierr); 808a94c3b12SMark F. Adams for(kk=0;kk<sz;kk++){ 809a94c3b12SMark F. Adams ierr = ISDestroy( &is[kk] ); CHKERRQ(ierr); 810a94c3b12SMark F. Adams } 811ffc955d6SMark F. Adams ierr = PetscFree( is ); CHKERRQ(ierr); 812ffc955d6SMark F. Adams } 813ffc955d6SMark F. Adams ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr); 814ffc955d6SMark F. Adams 815ffc955d6SMark F. Adams ASMLocalIDsArr[level] = PETSC_NULL; 816ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 817ffc955d6SMark F. Adams ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr); 818ffc955d6SMark F. Adams } 819ffc955d6SMark F. Adams else { 8209d5b6da9SMark F. Adams ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr); 821ffc955d6SMark F. Adams } 822d5280255SMark F. Adams } 823d5280255SMark F. Adams { 824d5280255SMark F. Adams /* coarse grid */ 825d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 826d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 827d5280255SMark F. Adams ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 828d5280255SMark F. Adams ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 829d5280255SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 830d5280255SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 831d5280255SMark F. Adams ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 832d5280255SMark F. Adams ierr = PCSetUp( subpc ); CHKERRQ(ierr); 833d5280255SMark F. Adams ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); assert(ii==1); 834d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 835d5280255SMark F. Adams ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 836d5280255SMark F. Adams } 837d5280255SMark F. Adams 838d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 839d5280255SMark F. Adams ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr); 840d5280255SMark F. Adams ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr); 841d5280255SMark F. Adams ierr = PetscOptionsEnd(); CHKERRQ(ierr); 842d5280255SMark F. Adams if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 843d5280255SMark F. Adams 844d5280255SMark F. Adams /* create cheby smoothers */ 845d5280255SMark F. Adams for ( lidx = 1, level = pc_gamg->Nlevels-2; 846d5280255SMark F. Adams lidx <= fine_level; 847d5280255SMark F. Adams lidx++, level--) { 848d5280255SMark F. Adams KSP smoother; 849ffc955d6SMark F. Adams PetscBool flag; 850d5280255SMark F. Adams PC subpc; 851d5280255SMark F. Adams 852ffc955d6SMark F. Adams ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 853a2f3521dSMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 854a2f3521dSMark F. Adams 855a2f3521dSMark F. Adams /* create field split PC, get subsmoother */ 856a2f3521dSMark F. Adams if( stokes ) { 857a2f3521dSMark F. Adams KSP *ksps; 858a2f3521dSMark F. Adams PetscInt nn; 859a2f3521dSMark F. Adams ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr); 860a2f3521dSMark F. Adams smoother = ksps[0]; 861a2f3521dSMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 862a2f3521dSMark F. Adams ierr = PetscFree( ksps ); CHKERRQ(ierr); 863a2f3521dSMark F. Adams } 864ffc955d6SMark F. Adams 865ffc955d6SMark F. Adams /* do my own cheby */ 866251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &flag ); CHKERRQ(ierr); 867ffc955d6SMark F. Adams if( flag ) { 868ffc955d6SMark F. Adams PetscReal emax, emin; 869251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr); 870ffc955d6SMark F. Adams if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 871587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 872587fa25dSMark F. Adams KSP eksp; Mat Lmat = Aarr[level]; 873b2a4f308SMark F. Adams Vec bb, xx; 874038e3b61SMark F. Adams 8755745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 8765745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 877fc4362bfSMark F. Adams { 878fc4362bfSMark F. Adams PetscRandom rctx; 879fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 880fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 881fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 882fc4362bfSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 8835b89ad90SMark F. Adams } 884a2f3521dSMark F. Adams 885a94c3b12SMark F. Adams if( removedEqs[level] ) { 886302f38e8SMark F. Adams /* being very careful - zeroing out BC rows (this is not done in agg.c estimates) */ 887a94c3b12SMark F. Adams PetscScalar *zeros; 888bea1fdf0SMark F. Adams PetscInt ii,jj, *idx_bs, sz, bs=level_bs[level]; 889bea1fdf0SMark F. Adams const PetscInt *idx; 890a94c3b12SMark F. Adams ierr = ISGetLocalSize( removedEqs[level], &sz ); CHKERRQ(ierr); 891a94c3b12SMark F. Adams ierr = PetscMalloc( bs*sz*sizeof(PetscScalar), &zeros ); CHKERRQ(ierr); 892a94c3b12SMark F. Adams for(ii=0;ii<bs*sz;ii++) zeros[ii] = 0.; 893a94c3b12SMark F. Adams ierr = PetscMalloc( bs*sz*sizeof(PetscInt), &idx_bs ); CHKERRQ(ierr); 894a94c3b12SMark F. Adams ierr = ISGetIndices( removedEqs[level], &idx); CHKERRQ(ierr); 895a94c3b12SMark F. Adams for(ii=0;ii<sz;ii++) { 896a94c3b12SMark F. Adams for(jj=0;jj<bs;jj++) { 897a94c3b12SMark F. Adams idx_bs[ii] = bs*idx[ii]+jj; 898a94c3b12SMark F. Adams } 899a94c3b12SMark F. Adams } 900a94c3b12SMark F. Adams ierr = ISRestoreIndices( removedEqs[level], &idx ); CHKERRQ(ierr); 901a94c3b12SMark F. Adams if( sz > 0 ) { 902a94c3b12SMark F. Adams ierr = VecSetValues( bb, sz, idx_bs, zeros, INSERT_VALUES ); CHKERRQ(ierr); 903a94c3b12SMark F. Adams } 904a94c3b12SMark F. Adams ierr = PetscFree( idx_bs ); CHKERRQ(ierr); 905a94c3b12SMark F. Adams ierr = PetscFree( zeros ); CHKERRQ(ierr); 906a94c3b12SMark F. Adams ierr = VecAssemblyBegin(bb); CHKERRQ(ierr); 907a94c3b12SMark F. Adams ierr = VecAssemblyEnd(bb); CHKERRQ(ierr); 908a94c3b12SMark F. Adams } 909fc4362bfSMark F. Adams ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr); 910038e3b61SMark F. Adams ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 911fc4362bfSMark F. Adams CHKERRQ(ierr); 912fc4362bfSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 9131a166f3bSJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 9141a166f3bSJed Brown ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_"); CHKERRQ(ierr); 915f6536408SMark F. Adams ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 916f6536408SMark F. Adams 917f6536408SMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 918f6536408SMark F. Adams ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 919fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 9205a9b9e01SMark F. Adams 921d3d0db20SJed Brown /* set PC type to be same as smoother */ 922ffc955d6SMark F. Adams ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr ); 923b2a4f308SMark F. Adams 9245a9b9e01SMark F. Adams /* solve - keep stuff out of logging */ 9255a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 9265a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 927fc4362bfSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 9285a9b9e01SMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 9295a9b9e01SMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 9305a9b9e01SMark F. Adams 931fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 9325a9b9e01SMark F. Adams 933fc4362bfSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 934fc4362bfSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 935fc4362bfSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 936f6536408SMark F. Adams 937ffc955d6SMark F. Adams if( pc_gamg->verbose > 0 ) { 938a94c3b12SMark F. Adams PetscInt N1, tt; 939a94c3b12SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 940a94c3b12SMark F. Adams PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 941f6536408SMark F. Adams } 942fc4362bfSMark F. Adams } 943038e3b61SMark F. Adams { 944038e3b61SMark F. Adams PetscInt N1, N0, tt; 945038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 946038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 947f6536408SMark F. Adams /* heuristic - is this crap? */ 948f6536408SMark F. Adams emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); 949038e3b61SMark F. Adams emax *= 1.05; 950038e3b61SMark F. Adams } 951fc4362bfSMark F. Adams ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 952ffc955d6SMark F. Adams } /* setup checby flag */ 953a2f3521dSMark F. Adams 954a2f3521dSMark F. Adams if( removedEqs[level] ) { 955a2f3521dSMark F. Adams ierr = ISDestroy( &removedEqs[level] ); CHKERRQ(ierr); 956a2f3521dSMark F. Adams } 957ffc955d6SMark F. Adams } /* non-coarse levels */ 958737a81a9SMark F. Adams 959d5280255SMark F. Adams /* clean up */ 960d5280255SMark F. Adams for(level=1;level<pc_gamg->Nlevels;level++){ 961587fa25dSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 962587fa25dSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 9635b89ad90SMark F. Adams } 9640cbbd2e1SMark F. Adams 9650cbbd2e1SMark F. Adams ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 9660cbbd2e1SMark F. Adams 967a94c3b12SMark F. Adams if( PETSC_FALSE ){ 96858471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 9699d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 97058471d46SMark F. Adams ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 97158471d46SMark F. Adams } 9725f8cf99dSMark F. Adams } 9735f8cf99dSMark F. Adams else { 9745f8cf99dSMark F. Adams KSP smoother; 975b43b03e9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__); 9769d5b6da9SMark F. Adams ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 9775f8cf99dSMark F. Adams ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 9785f8cf99dSMark F. Adams ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 9799d5b6da9SMark F. Adams ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 9805f8cf99dSMark F. Adams } 98184d3f75bSMark F. Adams 9825b89ad90SMark F. Adams PetscFunctionReturn(0); 9835b89ad90SMark F. Adams } 9845b89ad90SMark F. Adams 985eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 9865b89ad90SMark F. Adams /* 9875b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 9885b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 9895b89ad90SMark F. Adams 9905b89ad90SMark F. Adams Input Parameter: 9915b89ad90SMark F. Adams . pc - the preconditioner context 9925b89ad90SMark F. Adams 9935b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 9945b89ad90SMark F. Adams */ 9955b89ad90SMark F. Adams #undef __FUNCT__ 9965b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 9975b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG( PC pc ) 9985b89ad90SMark F. Adams { 9995b89ad90SMark F. Adams PetscErrorCode ierr; 10005b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 10015b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 10025b89ad90SMark F. Adams 10035b89ad90SMark F. Adams PetscFunctionBegin; 10045b89ad90SMark F. Adams ierr = PCReset_GAMG( pc );CHKERRQ(ierr); 10055b89ad90SMark F. Adams ierr = PetscFree( pc_gamg );CHKERRQ(ierr); 10065b89ad90SMark F. Adams ierr = PCDestroy_MG( pc );CHKERRQ(ierr); 10075b89ad90SMark F. Adams PetscFunctionReturn(0); 10085b89ad90SMark F. Adams } 10095b89ad90SMark F. Adams 1010676e1743SMark F. Adams 1011676e1743SMark F. Adams #undef __FUNCT__ 1012676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 1013676e1743SMark F. Adams /*@ 1014676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 1015676e1743SMark F. Adams processor reduction. 1016676e1743SMark F. Adams 1017676e1743SMark F. Adams Not Collective on PC 1018676e1743SMark F. Adams 1019676e1743SMark F. Adams Input Parameters: 1020676e1743SMark F. Adams . pc - the preconditioner context 1021676e1743SMark F. Adams 1022676e1743SMark F. Adams 1023676e1743SMark F. Adams Options Database Key: 1024676e1743SMark F. Adams . -pc_gamg_process_eq_limit 1025676e1743SMark F. Adams 1026676e1743SMark F. Adams Level: intermediate 1027676e1743SMark F. Adams 1028676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1029676e1743SMark F. Adams 1030676e1743SMark F. Adams .seealso: () 1031676e1743SMark F. Adams @*/ 1032676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 1033676e1743SMark F. Adams { 1034676e1743SMark F. Adams PetscErrorCode ierr; 1035676e1743SMark F. Adams 1036676e1743SMark F. Adams PetscFunctionBegin; 1037676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1038676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1039676e1743SMark F. Adams PetscFunctionReturn(0); 1040676e1743SMark F. Adams } 1041676e1743SMark F. Adams 1042676e1743SMark F. Adams EXTERN_C_BEGIN 1043676e1743SMark F. Adams #undef __FUNCT__ 1044676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 1045676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 1046676e1743SMark F. Adams { 1047c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1048c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1049676e1743SMark F. Adams 1050676e1743SMark F. Adams PetscFunctionBegin; 10519d5b6da9SMark F. Adams if(n>0) pc_gamg->min_eq_proc = n; 1052676e1743SMark F. Adams PetscFunctionReturn(0); 1053676e1743SMark F. Adams } 1054676e1743SMark F. Adams EXTERN_C_END 1055676e1743SMark F. Adams 1056676e1743SMark F. Adams #undef __FUNCT__ 1057389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 1058389730f3SMark F. Adams /*@ 1059389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 1060389730f3SMark F. Adams 1061389730f3SMark F. Adams Collective on PC 1062389730f3SMark F. Adams 1063389730f3SMark F. Adams Input Parameters: 1064389730f3SMark F. Adams . pc - the preconditioner context 1065389730f3SMark F. Adams 1066389730f3SMark F. Adams 1067389730f3SMark F. Adams Options Database Key: 1068389730f3SMark F. Adams . -pc_gamg_coarse_eq_limit 1069389730f3SMark F. Adams 1070389730f3SMark F. Adams Level: intermediate 1071389730f3SMark F. Adams 1072389730f3SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1073389730f3SMark F. Adams 1074389730f3SMark F. Adams .seealso: () 1075389730f3SMark F. Adams @*/ 1076389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 1077389730f3SMark F. Adams { 1078389730f3SMark F. Adams PetscErrorCode ierr; 1079389730f3SMark F. Adams 1080389730f3SMark F. Adams PetscFunctionBegin; 1081389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1082389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1083389730f3SMark F. Adams PetscFunctionReturn(0); 1084389730f3SMark F. Adams } 1085389730f3SMark F. Adams 1086389730f3SMark F. Adams EXTERN_C_BEGIN 1087389730f3SMark F. Adams #undef __FUNCT__ 1088389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 1089389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1090389730f3SMark F. Adams { 1091389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1092389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1093389730f3SMark F. Adams 1094389730f3SMark F. Adams PetscFunctionBegin; 10959d5b6da9SMark F. Adams if(n>0) pc_gamg->coarse_eq_limit = n; 1096389730f3SMark F. Adams PetscFunctionReturn(0); 1097389730f3SMark F. Adams } 1098389730f3SMark F. Adams EXTERN_C_END 1099389730f3SMark F. Adams 1100389730f3SMark F. Adams #undef __FUNCT__ 11018263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 1102676e1743SMark F. Adams /*@ 11038263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 1104676e1743SMark F. Adams 1105676e1743SMark F. Adams Collective on PC 1106676e1743SMark F. Adams 1107676e1743SMark F. Adams Input Parameters: 1108676e1743SMark F. Adams . pc - the preconditioner context 1109676e1743SMark F. Adams 1110676e1743SMark F. Adams 1111676e1743SMark F. Adams Options Database Key: 11128263b398SMark F. Adams . -pc_gamg_repartition 1113676e1743SMark F. Adams 1114676e1743SMark F. Adams Level: intermediate 1115676e1743SMark F. Adams 1116676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1117676e1743SMark F. Adams 1118676e1743SMark F. Adams .seealso: () 1119676e1743SMark F. Adams @*/ 11208263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1121676e1743SMark F. Adams { 1122676e1743SMark F. Adams PetscErrorCode ierr; 1123676e1743SMark F. Adams 1124676e1743SMark F. Adams PetscFunctionBegin; 1125676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11268263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1127676e1743SMark F. Adams PetscFunctionReturn(0); 1128676e1743SMark F. Adams } 1129676e1743SMark F. Adams 1130676e1743SMark F. Adams EXTERN_C_BEGIN 1131676e1743SMark F. Adams #undef __FUNCT__ 11328263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 11338263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 1134676e1743SMark F. Adams { 1135c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1136c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1137676e1743SMark F. Adams 1138676e1743SMark F. Adams PetscFunctionBegin; 11399d5b6da9SMark F. Adams pc_gamg->repart = n; 1140676e1743SMark F. Adams PetscFunctionReturn(0); 1141676e1743SMark F. Adams } 1142676e1743SMark F. Adams EXTERN_C_END 1143676e1743SMark F. Adams 1144676e1743SMark F. Adams #undef __FUNCT__ 1145ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs" 1146ffc955d6SMark F. Adams /*@ 1147ffc955d6SMark F. Adams PCGAMGSetUseASMAggs - 1148ffc955d6SMark F. Adams 1149ffc955d6SMark F. Adams Collective on PC 1150ffc955d6SMark F. Adams 1151ffc955d6SMark F. Adams Input Parameters: 1152ffc955d6SMark F. Adams . pc - the preconditioner context 1153ffc955d6SMark F. Adams 1154ffc955d6SMark F. Adams 1155ffc955d6SMark F. Adams Options Database Key: 1156ffc955d6SMark F. Adams . -pc_gamg_use_agg_gasm 1157ffc955d6SMark F. Adams 1158ffc955d6SMark F. Adams Level: intermediate 1159ffc955d6SMark F. Adams 1160ffc955d6SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1161ffc955d6SMark F. Adams 1162ffc955d6SMark F. Adams .seealso: () 1163ffc955d6SMark F. Adams @*/ 1164ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1165ffc955d6SMark F. Adams { 1166ffc955d6SMark F. Adams PetscErrorCode ierr; 1167ffc955d6SMark F. Adams 1168ffc955d6SMark F. Adams PetscFunctionBegin; 1169ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1170ffc955d6SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1171ffc955d6SMark F. Adams PetscFunctionReturn(0); 1172ffc955d6SMark F. Adams } 1173ffc955d6SMark F. Adams 1174ffc955d6SMark F. Adams EXTERN_C_BEGIN 1175ffc955d6SMark F. Adams #undef __FUNCT__ 1176ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 1177ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1178ffc955d6SMark F. Adams { 1179ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1180ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1181ffc955d6SMark F. Adams 1182ffc955d6SMark F. Adams PetscFunctionBegin; 1183ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = n; 1184ffc955d6SMark F. Adams PetscFunctionReturn(0); 1185ffc955d6SMark F. Adams } 1186ffc955d6SMark F. Adams EXTERN_C_END 1187ffc955d6SMark F. Adams 1188ffc955d6SMark F. Adams #undef __FUNCT__ 11894ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 11904ef23d27SMark F. Adams /*@ 11914ef23d27SMark F. Adams PCGAMGSetNlevels - 11924ef23d27SMark F. Adams 11934ef23d27SMark F. Adams Not collective on PC 11944ef23d27SMark F. Adams 11954ef23d27SMark F. Adams Input Parameters: 11964ef23d27SMark F. Adams . pc - the preconditioner context 11974ef23d27SMark F. Adams 11984ef23d27SMark F. Adams 11994ef23d27SMark F. Adams Options Database Key: 12004ef23d27SMark F. Adams . -pc_mg_levels 12014ef23d27SMark F. Adams 12024ef23d27SMark F. Adams Level: intermediate 12034ef23d27SMark F. Adams 12044ef23d27SMark F. Adams Concepts: Unstructured multrigrid preconditioner 12054ef23d27SMark F. Adams 12064ef23d27SMark F. Adams .seealso: () 12074ef23d27SMark F. Adams @*/ 12084ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 12094ef23d27SMark F. Adams { 12104ef23d27SMark F. Adams PetscErrorCode ierr; 12114ef23d27SMark F. Adams 12124ef23d27SMark F. Adams PetscFunctionBegin; 12134ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12144ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 12154ef23d27SMark F. Adams PetscFunctionReturn(0); 12164ef23d27SMark F. Adams } 12174ef23d27SMark F. Adams 12184ef23d27SMark F. Adams EXTERN_C_BEGIN 12194ef23d27SMark F. Adams #undef __FUNCT__ 12204ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 12214ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 12224ef23d27SMark F. Adams { 12234ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 12244ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12254ef23d27SMark F. Adams 12264ef23d27SMark F. Adams PetscFunctionBegin; 12279d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12284ef23d27SMark F. Adams PetscFunctionReturn(0); 12294ef23d27SMark F. Adams } 12304ef23d27SMark F. Adams EXTERN_C_END 12314ef23d27SMark F. Adams 12324ef23d27SMark F. Adams #undef __FUNCT__ 12333542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 12343542efc5SMark F. Adams /*@ 12353542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12363542efc5SMark F. Adams 12373542efc5SMark F. Adams Not collective on PC 12383542efc5SMark F. Adams 12393542efc5SMark F. Adams Input Parameters: 12403542efc5SMark F. Adams . pc - the preconditioner context 12413542efc5SMark F. Adams 12423542efc5SMark F. Adams 12433542efc5SMark F. Adams Options Database Key: 12443542efc5SMark F. Adams . -pc_gamg_threshold 12453542efc5SMark F. Adams 12463542efc5SMark F. Adams Level: intermediate 12473542efc5SMark F. Adams 12483542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 12493542efc5SMark F. Adams 12503542efc5SMark F. Adams .seealso: () 12513542efc5SMark F. Adams @*/ 12523542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 12533542efc5SMark F. Adams { 12543542efc5SMark F. Adams PetscErrorCode ierr; 12553542efc5SMark F. Adams 12563542efc5SMark F. Adams PetscFunctionBegin; 12573542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12583542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 12593542efc5SMark F. Adams PetscFunctionReturn(0); 12603542efc5SMark F. Adams } 12613542efc5SMark F. Adams 12623542efc5SMark F. Adams EXTERN_C_BEGIN 12633542efc5SMark F. Adams #undef __FUNCT__ 12643542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 12653542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 12663542efc5SMark F. Adams { 1267c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1268c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12693542efc5SMark F. Adams 12703542efc5SMark F. Adams PetscFunctionBegin; 12719d5b6da9SMark F. Adams pc_gamg->threshold = n; 12723542efc5SMark F. Adams PetscFunctionReturn(0); 12733542efc5SMark F. Adams } 12743542efc5SMark F. Adams EXTERN_C_END 12753542efc5SMark F. Adams 12763542efc5SMark F. Adams #undef __FUNCT__ 12779d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1278676e1743SMark F. Adams /*@ 12799d5b6da9SMark F. Adams PCGAMGSetType - Set solution method - calls sub create method 1280676e1743SMark F. Adams 1281676e1743SMark F. Adams Collective on PC 1282676e1743SMark F. Adams 1283676e1743SMark F. Adams Input Parameters: 1284676e1743SMark F. Adams . pc - the preconditioner context 1285676e1743SMark F. Adams 1286676e1743SMark F. Adams 1287676e1743SMark F. Adams Options Database Key: 12883542efc5SMark F. Adams . -pc_gamg_type 1289676e1743SMark F. Adams 1290676e1743SMark F. Adams Level: intermediate 1291676e1743SMark F. Adams 1292676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1293676e1743SMark F. Adams 1294676e1743SMark F. Adams .seealso: () 1295676e1743SMark F. Adams @*/ 12960202ef74SSatish Balay PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type ) 1297676e1743SMark F. Adams { 1298676e1743SMark F. Adams PetscErrorCode ierr; 1299676e1743SMark F. Adams 1300676e1743SMark F. Adams PetscFunctionBegin; 1301676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13020202ef74SSatish Balay ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type)); 1303676e1743SMark F. Adams CHKERRQ(ierr); 1304676e1743SMark F. Adams PetscFunctionReturn(0); 1305676e1743SMark F. Adams } 1306676e1743SMark F. Adams 1307676e1743SMark F. Adams EXTERN_C_BEGIN 1308676e1743SMark F. Adams #undef __FUNCT__ 13099d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 13100202ef74SSatish Balay PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type ) 1311676e1743SMark F. Adams { 13129d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 1313676e1743SMark F. Adams 1314676e1743SMark F. Adams PetscFunctionBegin; 13159d5b6da9SMark F. Adams ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r); 13169d5b6da9SMark F. Adams CHKERRQ(ierr); 13179d5b6da9SMark F. Adams 13189d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 13199d5b6da9SMark F. Adams 13209d5b6da9SMark F. Adams /* call sub create method */ 13219d5b6da9SMark F. Adams ierr = (*r)(pc); CHKERRQ(ierr); 13229d5b6da9SMark F. Adams 1323676e1743SMark F. Adams PetscFunctionReturn(0); 1324676e1743SMark F. Adams } 1325676e1743SMark F. Adams EXTERN_C_END 1326676e1743SMark F. Adams 13275b89ad90SMark F. Adams #undef __FUNCT__ 13285b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 13295b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG( PC pc ) 13305b89ad90SMark F. Adams { 1331676e1743SMark F. Adams PetscErrorCode ierr; 1332676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1333676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1334676e1743SMark F. Adams PetscBool flag; 1335b43b03e9SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 13365b89ad90SMark F. Adams 13375b89ad90SMark F. Adams PetscFunctionBegin; 1338676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1339676e1743SMark F. Adams { 134075b74bdaSMark F. Adams /* -pc_gamg_verbose */ 13419d5b6da9SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 13429d5b6da9SMark F. Adams "none", pc_gamg->verbose, 13439d5b6da9SMark F. Adams &pc_gamg->verbose, PETSC_NULL ); 13449d5b6da9SMark F. Adams CHKERRQ(ierr); 134575b74bdaSMark F. Adams 13468263b398SMark F. Adams /* -pc_gamg_repartition */ 13478263b398SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_repartition", 13488263b398SMark F. Adams "Repartion coarse grids (false)", 13498263b398SMark F. Adams "PCGAMGRepartitioning", 13509d5b6da9SMark F. Adams pc_gamg->repart, 13519d5b6da9SMark F. Adams &pc_gamg->repart, 1352676e1743SMark F. Adams &flag); 1353676e1743SMark F. Adams CHKERRQ(ierr); 1354676e1743SMark F. Adams 1355ffc955d6SMark F. Adams /* -pc_gamg_use_agg_gasm */ 1356ffc955d6SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm", 1357ffc955d6SMark F. Adams "Use aggregation agragates for GASM smoother (false)", 1358ffc955d6SMark F. Adams "PCGAMGUseASMAggs", 1359ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm, 1360ffc955d6SMark F. Adams &pc_gamg->use_aggs_in_gasm, 1361ffc955d6SMark F. Adams &flag); 1362ffc955d6SMark F. Adams CHKERRQ(ierr); 1363ffc955d6SMark F. Adams 1364c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 1365676e1743SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1366676e1743SMark F. Adams "Limit (goal) on number of equations per process on coarse grids", 1367676e1743SMark F. Adams "PCGAMGSetProcEqLim", 13689d5b6da9SMark F. Adams pc_gamg->min_eq_proc, 13699d5b6da9SMark F. Adams &pc_gamg->min_eq_proc, 1370676e1743SMark F. Adams &flag ); 1371676e1743SMark F. Adams CHKERRQ(ierr); 13723542efc5SMark F. Adams 1373389730f3SMark F. Adams /* -pc_gamg_coarse_eq_limit */ 1374389730f3SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1375389730f3SMark F. Adams "Limit on number of equations for the coarse grid", 1376389730f3SMark F. Adams "PCGAMGSetCoarseEqLim", 13779d5b6da9SMark F. Adams pc_gamg->coarse_eq_limit, 13789d5b6da9SMark F. Adams &pc_gamg->coarse_eq_limit, 1379389730f3SMark F. Adams &flag ); 1380389730f3SMark F. Adams CHKERRQ(ierr); 1381389730f3SMark F. Adams 1382c20e4228SMark F. Adams /* -pc_gamg_threshold */ 13833542efc5SMark F. Adams ierr = PetscOptionsReal("-pc_gamg_threshold", 13843542efc5SMark F. Adams "Relative threshold to use for dropping edges in aggregation graph", 13853542efc5SMark F. Adams "PCGAMGSetThreshold", 13869d5b6da9SMark F. Adams pc_gamg->threshold, 13879d5b6da9SMark F. Adams &pc_gamg->threshold, 13883542efc5SMark F. Adams &flag ); 13893542efc5SMark F. Adams CHKERRQ(ierr); 1390b43b03e9SMark F. Adams if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold); 13914ef23d27SMark F. Adams 13924ef23d27SMark F. Adams ierr = PetscOptionsInt("-pc_mg_levels", 13934ef23d27SMark F. Adams "Set number of MG levels", 13944ef23d27SMark F. Adams "PCGAMGSetNlevels", 13959d5b6da9SMark F. Adams pc_gamg->Nlevels, 13969d5b6da9SMark F. Adams &pc_gamg->Nlevels, 13974ef23d27SMark F. Adams &flag ); 1398676e1743SMark F. Adams } 1399676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1400676e1743SMark F. Adams 14015b89ad90SMark F. Adams PetscFunctionReturn(0); 14025b89ad90SMark F. Adams } 14035b89ad90SMark F. Adams 14045b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 14055b89ad90SMark F. Adams /*MC 1406856830bfSJed Brown PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 14079d5b6da9SMark F. Adams - This is the entry point to GAMG, registered in pcregis.c 14085b89ad90SMark F. Adams 1409280d9858SJed Brown Options Database Keys: 14105b89ad90SMark F. Adams Multigrid options(inherited) 1411280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1412280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1413280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1414280d9858SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 14155b89ad90SMark F. Adams 14165b89ad90SMark F. Adams Level: intermediate 1417280d9858SJed Brown 14185b89ad90SMark F. Adams Concepts: multigrid 14195b89ad90SMark F. Adams 14205b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1421280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 14225b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 14235b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 14245b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 14255b89ad90SMark F. Adams M*/ 14265b89ad90SMark F. Adams EXTERN_C_BEGIN 14275b89ad90SMark F. Adams #undef __FUNCT__ 14285b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 14295b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG( PC pc ) 14305b89ad90SMark F. Adams { 14315b89ad90SMark F. Adams PetscErrorCode ierr; 14325b89ad90SMark F. Adams PC_GAMG *pc_gamg; 14335b89ad90SMark F. Adams PC_MG *mg; 14340cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 14359d5b6da9SMark F. Adams static long count = 0; 14365ee9c036SSatish Balay #endif 14375b89ad90SMark F. Adams 14385b89ad90SMark F. Adams PetscFunctionBegin; 143960a6d8f8SMark F. Adams 14405b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 14415b89ad90SMark F. Adams ierr = PCSetType( pc, PCMG ); CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 14425b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr); 14435b89ad90SMark F. Adams 14445b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 14455b89ad90SMark F. Adams ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr); 14465b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 1447ce4cda84SJed Brown mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 14485b89ad90SMark F. Adams mg->innerctx = pc_gamg; 14495b89ad90SMark F. Adams 14509d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 14519d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 14529d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 14539d5b6da9SMark F. Adams pc_gamg->data = 0; 14545b89ad90SMark F. Adams 14559d5b6da9SMark F. Adams /* register AMG type */ 14569d5b6da9SMark F. Adams if( !GAMGList ){ 14579d5b6da9SMark F. Adams ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 14589d5b6da9SMark F. Adams ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 14599d5b6da9SMark F. Adams } 14609d5b6da9SMark F. Adams 14619d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 14625b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 14635b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 14645b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 14655b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 14665b89ad90SMark F. Adams 14675b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1468676e1743SMark F. Adams "PCGAMGSetProcEqLim_C", 1469676e1743SMark F. Adams "PCGAMGSetProcEqLim_GAMG", 1470676e1743SMark F. Adams PCGAMGSetProcEqLim_GAMG); 1471676e1743SMark F. Adams CHKERRQ(ierr); 1472676e1743SMark F. Adams 1473676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1474389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_C", 1475389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_GAMG", 1476389730f3SMark F. Adams PCGAMGSetCoarseEqLim_GAMG); 1477389730f3SMark F. Adams CHKERRQ(ierr); 1478389730f3SMark F. Adams 1479389730f3SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 14808263b398SMark F. Adams "PCGAMGSetRepartitioning_C", 14818263b398SMark F. Adams "PCGAMGSetRepartitioning_GAMG", 14828263b398SMark F. Adams PCGAMGSetRepartitioning_GAMG); 1483676e1743SMark F. Adams CHKERRQ(ierr); 1484676e1743SMark F. Adams 1485676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1486ffc955d6SMark F. Adams "PCGAMGSetUseASMAggs_C", 1487ffc955d6SMark F. Adams "PCGAMGSetUseASMAggs_GAMG", 1488ffc955d6SMark F. Adams PCGAMGSetUseASMAggs_GAMG); 1489ffc955d6SMark F. Adams CHKERRQ(ierr); 1490ffc955d6SMark F. Adams 1491ffc955d6SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 14923542efc5SMark F. Adams "PCGAMGSetThreshold_C", 14933542efc5SMark F. Adams "PCGAMGSetThreshold_GAMG", 14943542efc5SMark F. Adams PCGAMGSetThreshold_GAMG); 14953542efc5SMark F. Adams CHKERRQ(ierr); 14963542efc5SMark F. Adams 14973542efc5SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 14989d5b6da9SMark F. Adams "PCGAMGSetType_C", 14999d5b6da9SMark F. Adams "PCGAMGSetType_GAMG", 15009d5b6da9SMark F. Adams PCGAMGSetType_GAMG); 1501676e1743SMark F. Adams CHKERRQ(ierr); 1502c97e1df0SMark F. Adams 15039d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1504ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 15059d5b6da9SMark F. Adams pc_gamg->min_eq_proc = 100; 1506c8b0795cSMark F. Adams pc_gamg->coarse_eq_limit = 800; 1507a2f3521dSMark F. Adams pc_gamg->threshold = 0.001; 15089d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 15099d5b6da9SMark F. Adams pc_gamg->verbose = 0; 15109d5b6da9SMark F. Adams pc_gamg->emax_id = -1; 15119d5b6da9SMark F. Adams 15120cbbd2e1SMark F. Adams /* private events */ 15130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1514785cba28SMark F. Adams if( count++ == 0 ) { 15150cbbd2e1SMark F. Adams PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]); 15160cbbd2e1SMark F. Adams PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]); 15170cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 15180cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 15190cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 15200cbbd2e1SMark F. Adams PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]); 15210cbbd2e1SMark F. Adams PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]); 15220cbbd2e1SMark F. Adams PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]); 15230cbbd2e1SMark F. Adams PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]); 152460a6d8f8SMark F. Adams PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]); 15250cbbd2e1SMark F. Adams PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]); 15260cbbd2e1SMark F. Adams PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]); 15270cbbd2e1SMark F. Adams PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]); 15280cbbd2e1SMark F. Adams PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]); 15290cbbd2e1SMark F. Adams PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]); 15300cbbd2e1SMark F. Adams PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]); 15310cbbd2e1SMark F. Adams PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]); 1532f852f58cSMark F. Adams 15330cbbd2e1SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 15340cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 15350cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1536b4fbaa2aSMark F. Adams /* create timer stages */ 1537b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1538b4fbaa2aSMark F. Adams { 1539b4fbaa2aSMark F. Adams char str[32]; 1540b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d (finest)",0); 1541b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[0]); 1542b4fbaa2aSMark F. Adams PetscInt lidx; 1543b4fbaa2aSMark F. Adams for (lidx=1;lidx<9;lidx++){ 1544b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1545b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[lidx]); 1546b4fbaa2aSMark F. Adams } 1547b4fbaa2aSMark F. Adams } 1548b4fbaa2aSMark F. Adams #endif 1549b4fbaa2aSMark F. Adams } 1550b4fbaa2aSMark F. Adams #endif 15510cbbd2e1SMark F. Adams /* general events */ 15520cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 155360a6d8f8SMark F. Adams PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG); 15540cbbd2e1SMark F. Adams PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO); 15550cbbd2e1SMark F. Adams PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG); 15560cbbd2e1SMark F. Adams PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO); 15570cbbd2e1SMark F. Adams PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG); 15580cbbd2e1SMark F. Adams PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO); 155960a6d8f8SMark F. Adams PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG); 1560a2f3521dSMark F. Adams PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG); 15610cbbd2e1SMark F. Adams #endif 15620cbbd2e1SMark F. Adams 156360a6d8f8SMark F. Adams /* instantiate derived type */ 156460a6d8f8SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 156560a6d8f8SMark F. Adams { 156660a6d8f8SMark F. Adams char tname[256] = GAMGAGG; 156760a6d8f8SMark F. Adams ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType", 156860a6d8f8SMark F. Adams GAMGList, tname, tname, sizeof(tname), PETSC_NULL ); 156960a6d8f8SMark F. Adams CHKERRQ(ierr); 157060a6d8f8SMark F. Adams ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr); 157160a6d8f8SMark F. Adams } 157260a6d8f8SMark F. Adams ierr = PetscOptionsTail(); CHKERRQ(ierr); 157360a6d8f8SMark F. Adams 15745b89ad90SMark F. Adams PetscFunctionReturn(0); 15755b89ad90SMark F. Adams } 15765b89ad90SMark F. Adams EXTERN_C_END 1577