15b89ad90SMark F. Adams /* 25b89ad90SMark F. Adams GAMG geometric-algebric multiogrid 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*/ 6f96513f1SMatthew G Knepley 7b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 8b4fbaa2aSMark F. Adams PetscLogEvent gamg_setup_events[NUM_SET]; 9b4fbaa2aSMark F. Adams #endif 10b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 11b4fbaa2aSMark F. Adams 12b8fd24d8SMark F. Adams /*#define GAMG_STAGES*/ 13b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 14b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 15b4fbaa2aSMark F. Adams #endif 16f96513f1SMatthew G Knepley 175b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 185b89ad90SMark F. Adams /* 195b89ad90SMark F. Adams PCSetCoordinates_GAMG 205b89ad90SMark F. Adams 215b89ad90SMark F. Adams Input Parameter: 225b89ad90SMark F. Adams . pc - the preconditioner context 235b89ad90SMark F. Adams */ 24a92563c5SMark F. Adams EXTERN_C_BEGIN 255b89ad90SMark F. Adams #undef __FUNCT__ 265b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 27eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords ) 285b89ad90SMark F. Adams { 29eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 305b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 316c237d78SBarry Smith PetscErrorCode ierr; 32d3d6bff4SMark F. Adams PetscInt arrsz,bs,my0,kk,ii,jj,nloc,Iend; 33038e3b61SMark F. Adams Mat Amat = a_pc->pmat; 345b89ad90SMark F. Adams 355b89ad90SMark F. Adams PetscFunctionBegin; 3658471d46SMark F. Adams PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 ); 37038e3b61SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 38d3d6bff4SMark F. Adams ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 39d3d6bff4SMark F. Adams nloc = (Iend-my0)/bs; 40d3d6bff4SMark F. Adams if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 41038e3b61SMark F. Adams 42d3d6bff4SMark F. Adams pc_gamg->m_data_rows = 1; 43f6536408SMark F. Adams if(a_coords==0 && pc_gamg->m_method==0) { 44f6536408SMark F. Adams SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 45f6536408SMark F. Adams } 46470e26f8SMark F. Adams if( pc_gamg->m_method==0 ) pc_gamg->m_data_cols = a_ndm; /* coordinates */ 47038e3b61SMark F. Adams else{ /* SA: null space vectors */ 48d3d6bff4SMark F. Adams if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */ 49d3d6bff4SMark F. Adams else if(a_coords != 0 ) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */ 50d3d6bff4SMark F. Adams else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */ 51d3d6bff4SMark F. Adams pc_gamg->m_data_rows = bs; 52038e3b61SMark F. Adams } 53d3d6bff4SMark F. Adams arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols; 545b89ad90SMark F. Adams 55038e3b61SMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 566e3e101aSMark F. Adams if (pc_gamg->m_data==0 || (pc_gamg->m_data_sz != arrsz)) { 575b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 5833a2b366SJed Brown ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->m_data ); CHKERRQ(ierr); 595b89ad90SMark F. Adams } 60038e3b61SMark F. Adams for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.; 61eb07cef2SMark F. Adams pc_gamg->m_data[arrsz] = -99.; 62038e3b61SMark F. Adams /* copy data in - column oriented */ 63470e26f8SMark F. Adams if( pc_gamg->m_method != 0 ) { 64d3d6bff4SMark F. Adams const PetscInt M = Iend - my0; 65038e3b61SMark F. Adams for(kk=0;kk<nloc;kk++){ 66038e3b61SMark F. Adams PetscReal *data = &pc_gamg->m_data[kk*bs]; 67d3d6bff4SMark F. Adams if( pc_gamg->m_data_cols==1 ) *data = 1.0; 68038e3b61SMark F. Adams else { 69038e3b61SMark F. Adams for(ii=0;ii<bs;ii++) 70038e3b61SMark F. Adams for(jj=0;jj<bs;jj++) 71038e3b61SMark F. Adams if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */ 72038e3b61SMark F. Adams else data[ii*M + jj] = 0.0; 73eb07cef2SMark F. Adams if( a_coords != 0 ) { 74038e3b61SMark F. Adams if( a_ndm == 2 ){ /* rotational modes */ 75038e3b61SMark F. Adams data += 2*M; 76038e3b61SMark F. Adams data[0] = -a_coords[2*kk+1]; 77038e3b61SMark F. Adams data[1] = a_coords[2*kk]; 785b89ad90SMark F. Adams } 79eb07cef2SMark F. Adams else { 80038e3b61SMark F. Adams data += 3*M; 81038e3b61SMark F. Adams data[0] = 0.0; data[M+0] = a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1]; 82038e3b61SMark F. Adams data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = a_coords[3*kk]; 83038e3b61SMark F. Adams data[2] = a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk]; data[2*M+2] = 0.0; 84038e3b61SMark F. Adams } 85eb07cef2SMark F. Adams } 86eb07cef2SMark F. Adams } 87eb07cef2SMark F. Adams } 88eb07cef2SMark F. Adams } 89eb07cef2SMark F. Adams else { 90038e3b61SMark F. Adams for( kk = 0 ; kk < nloc ; kk++ ){ 91038e3b61SMark F. Adams for( ii = 0 ; ii < a_ndm ; ii++ ) { 92038e3b61SMark F. Adams pc_gamg->m_data[ii*nloc + kk] = a_coords[kk*a_ndm + ii]; 93eb07cef2SMark F. Adams } 94eb07cef2SMark F. Adams } 95038e3b61SMark F. Adams } 96eb07cef2SMark F. Adams assert(pc_gamg->m_data[arrsz] == -99.); 97038e3b61SMark F. Adams 985b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 99eb07cef2SMark F. Adams pc_gamg->m_dim = a_ndm; 1005b89ad90SMark F. Adams 1015b89ad90SMark F. Adams PetscFunctionReturn(0); 1025b89ad90SMark F. Adams } 103a92563c5SMark F. Adams EXTERN_C_END 1045b89ad90SMark F. Adams 105d3d6bff4SMark F. Adams 106d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 107d3d6bff4SMark F. Adams #undef __FUNCT__ 108d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 109d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 110d3d6bff4SMark F. Adams { 111d3d6bff4SMark F. Adams PetscErrorCode ierr; 112d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 113d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 114d3d6bff4SMark F. Adams 115d3d6bff4SMark F. Adams PetscFunctionBegin; 11658471d46SMark F. Adams if( pc_gamg->m_data != 0 ) { /* this should not happen, cleaned up in SetUp */ 117d3d6bff4SMark F. Adams ierr = PetscFree(pc_gamg->m_data); CHKERRQ(ierr); 11858471d46SMark F. Adams } 119d3d6bff4SMark F. Adams pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0; 120d3d6bff4SMark F. Adams PetscFunctionReturn(0); 121d3d6bff4SMark F. Adams } 122d3d6bff4SMark F. Adams 1235b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1245b89ad90SMark F. Adams /* 1258263b398SMark F. Adams createLevel 1265b89ad90SMark F. Adams 1275b89ad90SMark F. Adams Input Parameter: 1288263b398SMark F. Adams . a_pc - parameters 1293530afc2SMark F. Adams . a_Amat_fine - matrix on this fine (k) level 130d3d6bff4SMark F. Adams . a_ndata_rows - size of data to move (coarse grid) 131d3d6bff4SMark F. Adams . a_ndata_cols - size of data to move (coarse grid) 1328263b398SMark F. Adams . a_cbs - coarse block size 1338263b398SMark F. Adams . a_isLast - 1343530afc2SMark F. Adams In/Output Parameter: 1353530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 136eb07cef2SMark F. Adams . a_coarse_data - data that need to be moved 137afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 13811e60469SMark F. Adams Output Parameter: 1393530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 1405b89ad90SMark F. Adams */ 1415cb416c2SMark F. Adams 1425b89ad90SMark F. Adams #undef __FUNCT__ 1438263b398SMark F. Adams #define __FUNCT__ "createLevel" 1448263b398SMark F. Adams PetscErrorCode createLevel( const PC a_pc, 1458263b398SMark F. Adams const Mat a_Amat_fine, 1468263b398SMark F. Adams const PetscInt a_ndata_rows, 1478263b398SMark F. Adams const PetscInt a_ndata_cols, 1488263b398SMark F. Adams const PetscInt a_cbs, 1498263b398SMark F. Adams const PetscBool a_isLast, 1503530afc2SMark F. Adams Mat *a_P_inout, 151eb07cef2SMark F. Adams PetscReal **a_coarse_data, 152afc97cdcSMark F. Adams PetscMPIInt *a_nactive_proc, 153389730f3SMark F. Adams Mat *a_Amat_crs 15411e60469SMark F. Adams ) 1555b89ad90SMark F. Adams { 1568263b398SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 157486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1588263b398SMark F. Adams const PetscBool repart = pc_gamg->m_repart; 159389730f3SMark F. Adams const PetscInt min_eq_proc = pc_gamg->m_min_eq_proc, coarse_max = pc_gamg->m_coarse_eq_limit; 1605b89ad90SMark F. Adams PetscErrorCode ierr; 161038e3b61SMark F. Adams Mat Cmat,Pnew,Pold=*a_P_inout; 16211e60469SMark F. Adams IS new_indices,isnum; 1633530afc2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 164fa31c87dSMark F. Adams PetscMPIInt mype,npe,new_npe,nactive; 165fa31c87dSMark F. Adams PetscInt neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0; 1665b89ad90SMark F. Adams 1675b89ad90SMark F. Adams PetscFunctionBegin; 16811e60469SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 16911e60469SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 170f6536408SMark F. Adams 17111e60469SMark F. Adams /* RAP */ 17296434bc1SMark F. Adams #ifdef USE_R 17396434bc1SMark F. Adams /* make R wih brute force for now */ 17496434bc1SMark F. Adams ierr = MatTranspose( Pold, Pnew ); 17596434bc1SMark F. Adams ierr = MatDestroy( &Pold ); CHKERRQ(ierr); 17696434bc1SMark F. Adams ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 17796434bc1SMark F. Adams Pold = Pnew; 17896434bc1SMark F. Adams #else 179038e3b61SMark F. Adams ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 18096434bc1SMark F. Adams #endif 181038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 182038e3b61SMark F. Adams ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); 183038e3b61SMark F. Adams ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0); 184038e3b61SMark F. Adams 185f852f58cSMark F. Adams /* get number of PEs to make active, reduce */ 186f852f58cSMark F. Adams ierr = MatGetSize( Cmat, &neq, &NN ); CHKERRQ(ierr); 187f852f58cSMark F. Adams new_npe = neq/min_eq_proc; /* hardwire min. number of eq/proc */ 188f852f58cSMark F. Adams if( new_npe == 0 || neq < coarse_max ) new_npe = 1; 189f852f58cSMark F. Adams else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */ 1908263b398SMark F. Adams if( a_isLast ) new_npe = 1; 191f852f58cSMark F. Adams 1928263b398SMark F. Adams if( !repart && !(new_npe == 1 && *a_nactive_proc != 1) ) { 19322063be5SMark F. Adams *a_Amat_crs = Cmat; /* output */ 19422063be5SMark F. Adams } 19522063be5SMark F. Adams else { 19622063be5SMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 1975ef31b24SMark F. Adams Mat adj; 19822063be5SMark F. Adams const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols; 199b4fbaa2aSMark F. Adams const PetscInt stride0=ncrs0*a_ndata_rows; 2008263b398SMark F. Adams PetscInt is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx; 201c9a0b8beSMark F. Adams /* create sub communicator */ 20271d60426SBarry Smith MPI_Comm cm; 203afc97cdcSMark F. Adams MPI_Group wg, g2; 204fd3c6afaSMark F. Adams PetscInt *counts,inpe; 205fd3c6afaSMark F. Adams PetscMPIInt *ranks; 20622063be5SMark F. Adams IS isscat; 20722063be5SMark F. Adams PetscScalar *array; 20822063be5SMark F. Adams Vec src_crd, dest_crd; 20922063be5SMark F. Adams PetscReal *data = *a_coarse_data; 21022063be5SMark F. Adams VecScatter vecscat; 211b4fbaa2aSMark F. Adams IS isnewproc; 212e33ef3b1SMark F. Adams 21392a756f0SMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr); 214fd3c6afaSMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 21592a756f0SMark F. Adams 216fd3c6afaSMark F. Adams ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr); 217afc97cdcSMark F. Adams assert(counts[mype]==ncrs0); 218afc97cdcSMark F. Adams /* count real active pes */ 21922063be5SMark F. Adams for( nactive = jj = 0 ; jj < npe ; jj++) { 220afc97cdcSMark F. Adams if( counts[jj] != 0 ) { 22122063be5SMark F. Adams ranks[nactive++] = jj; 222afc97cdcSMark F. Adams } 223afc97cdcSMark F. Adams } 2243303bcf9Sadams 2253303bcf9Sadams if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */ 22622063be5SMark F. Adams 227bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s npe (active): %d --> %d. new npe = %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,nactive,new_npe,neq); 22858471d46SMark F. Adams 22922063be5SMark F. Adams *a_nactive_proc = new_npe; /* output */ 2302f03bc48SMark F. Adams 231afc97cdcSMark F. Adams ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr); 23222063be5SMark F. Adams ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr); 233afc97cdcSMark F. Adams ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr); 234b4fbaa2aSMark F. Adams 235afc97cdcSMark F. Adams ierr = MPI_Group_free( &wg ); CHKERRQ(ierr); 236afc97cdcSMark F. Adams ierr = MPI_Group_free( &g2 ); CHKERRQ(ierr); 237c9a0b8beSMark F. Adams 2385ef31b24SMark F. Adams /* MatPartitioningApply call MatConvert, which is collective */ 239b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 240b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 241b4fbaa2aSMark F. Adams #endif 242038e3b61SMark F. Adams if( a_cbs == 1) { 243038e3b61SMark F. Adams ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 244eb07cef2SMark F. Adams } 245eb07cef2SMark F. Adams else{ 246038e3b61SMark F. Adams /* make a scalar matrix to partition */ 247eb07cef2SMark F. Adams Mat tMat; 24858471d46SMark F. Adams PetscInt ncols,jj,Ii; 249b4fbaa2aSMark F. Adams const PetscScalar *vals; 250b4fbaa2aSMark F. Adams const PetscInt *idx; 25158471d46SMark F. Adams PetscInt *d_nnz; 2525ee9c036SSatish Balay static int llev = 0; 253b4fbaa2aSMark F. Adams 25458471d46SMark F. Adams ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); 25558471d46SMark F. Adams for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) { 25658471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 25758471d46SMark F. Adams d_nnz[jj] = ncols/a_cbs; 25858471d46SMark F. Adams if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0; 25958471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 26058471d46SMark F. Adams } 2616876a03eSMark F. Adams 262eb07cef2SMark F. Adams ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 263eb07cef2SMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 26458471d46SMark F. Adams 0, d_nnz, 0, d_nnz, 265eb07cef2SMark F. Adams &tMat ); 2666876a03eSMark F. Adams CHKERRQ(ierr); 26758471d46SMark F. Adams ierr = PetscFree( d_nnz ); CHKERRQ(ierr); 268eb07cef2SMark F. Adams 26922063be5SMark F. Adams for ( ii = Istart0; ii < Iend0; ii++ ) { 27022063be5SMark F. Adams PetscInt dest_row = ii/a_cbs; 27122063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 272eb07cef2SMark F. Adams for( jj = 0 ; jj < ncols ; jj++ ){ 273038e3b61SMark F. Adams PetscInt dest_col = idx[jj]/a_cbs; 274eb07cef2SMark F. Adams PetscScalar v = 1.0; 275eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 276eb07cef2SMark F. Adams } 27722063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 278eb07cef2SMark F. Adams } 279eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 280eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 281eb07cef2SMark F. Adams 282b4fbaa2aSMark F. Adams if( llev++ == -1 ) { 283b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 284b4fbaa2aSMark F. Adams sprintf(fname,"part_mat_%d.mat",llev); 285b4fbaa2aSMark F. Adams PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 286b4fbaa2aSMark F. Adams ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 287b4fbaa2aSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 288b4fbaa2aSMark F. Adams } 289b4fbaa2aSMark F. Adams 290eb07cef2SMark F. Adams ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 291eb07cef2SMark F. Adams 292eb07cef2SMark F. Adams ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 293eb07cef2SMark F. Adams } 294f150b916SMark F. Adams 295afc97cdcSMark F. Adams if( ncrs0 != 0 ){ 296b4fbaa2aSMark F. Adams const PetscInt *is_idx; 297b4fbaa2aSMark F. Adams MatPartitioning mpart; 298a4b7d37bSMark F. Adams IS proc_is; 2995ef31b24SMark F. Adams /* hack to fix global data that pmetis.c uses in 'adj' */ 30022063be5SMark F. Adams for( nactive = jj = 0 ; jj < npe ; jj++ ) { 301afc97cdcSMark F. Adams if( counts[jj] != 0 ) { 30222063be5SMark F. Adams adj->rmap->range[nactive++] = adj->rmap->range[jj]; 303afc97cdcSMark F. Adams } 3045ef31b24SMark F. Adams } 30522063be5SMark F. Adams adj->rmap->range[nactive] = adj->rmap->range[npe]; 3062f03bc48SMark F. Adams 3078263b398SMark F. Adams if( new_npe == 1 ) { 3088263b398SMark F. Adams ierr = MatGetLocalSize( adj, &is_sz, &ii ); CHKERRQ(ierr); 309a4b7d37bSMark F. Adams ierr = ISCreateStride( cm, is_sz, 0, 0, &proc_is ); CHKERRQ(ierr); 31059a0be82SJed Brown } else { 31159a0be82SJed Brown char prefix[256]; 31259a0be82SJed Brown const char *pcpre; 31371d60426SBarry Smith ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr); 3145ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 31559a0be82SJed Brown ierr = PCGetOptionsPrefix(a_pc,&pcpre);CHKERRQ(ierr); 31659a0be82SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr); 31759a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 31811e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 3195ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 320a4b7d37bSMark F. Adams ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr); 32111e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 3228263b398SMark F. Adams } 3235ef31b24SMark F. Adams /* collect IS info */ 324a4b7d37bSMark F. Adams ierr = ISGetLocalSize( proc_is, &is_sz ); CHKERRQ(ierr); 3258263b398SMark F. Adams ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr); 326a4b7d37bSMark F. Adams ierr = ISGetIndices( proc_is, &is_idx ); CHKERRQ(ierr); 327b4fbaa2aSMark F. Adams /* spread partitioning across machine - best way ??? */ 328ecd57171SMark F. Adams NN = 1; /*npe/new_npe;*/ 329b4fbaa2aSMark F. Adams for( kk = jj = 0 ; kk < is_sz ; kk++ ){ 330afc97cdcSMark F. Adams for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) { 3318263b398SMark F. Adams newproc_idx[jj] = is_idx[kk] * NN; /* distribution */ 332eb07cef2SMark F. Adams } 3335ef31b24SMark F. Adams } 334a4b7d37bSMark F. Adams ierr = ISRestoreIndices( proc_is, &is_idx ); CHKERRQ(ierr); 335a4b7d37bSMark F. Adams ierr = ISDestroy( &proc_is ); CHKERRQ(ierr); 33671d60426SBarry Smith ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 337b4fbaa2aSMark F. Adams 338b4fbaa2aSMark F. Adams is_sz *= a_cbs; 3395ef31b24SMark F. Adams } 3405ef31b24SMark F. Adams else{ 3418263b398SMark F. Adams newproc_idx = 0; 3425ef31b24SMark F. Adams is_sz = 0; 3435ef31b24SMark F. Adams } 344f150b916SMark F. Adams 3455ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 3468263b398SMark F. Adams ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc ); 3478263b398SMark F. Adams if( newproc_idx != 0 ) { 3488263b398SMark F. Adams ierr = PetscFree( newproc_idx ); CHKERRQ(ierr); 3495ef31b24SMark F. Adams } 350e33ef3b1SMark F. Adams 35111e60469SMark F. Adams /* 35211e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 35311e60469SMark F. Adams */ 35411e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 35511e60469SMark F. Adams /* 35611e60469SMark F. Adams Determine how many elements are assigned to each processor 35711e60469SMark F. Adams */ 358fd3c6afaSMark F. Adams inpe = npe; 359fd3c6afaSMark F. Adams ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr); 36011e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 361038e3b61SMark F. Adams ncrs_new = counts[mype]/a_cbs; 36222063be5SMark F. Adams strideNew = ncrs_new*a_ndata_rows; 363b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 364b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 365b4fbaa2aSMark F. Adams #endif 36622063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 36711e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 368d3d6bff4SMark F. Adams ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 369470e26f8SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ 37011e60469SMark F. Adams /* 371d3d6bff4SMark F. Adams There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having 372d3d6bff4SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'a_cbs'. 37311e60469SMark F. Adams */ 37492a756f0SMark F. Adams ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 37511e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 376038e3b61SMark F. Adams for(ii=0,jj=0; ii<ncrs0 ; ii++) { 377d3d6bff4SMark F. Adams PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */ 378d3d6bff4SMark F. Adams for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk; 37911e60469SMark F. Adams } 380038e3b61SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 381d3d6bff4SMark F. Adams ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 38211e60469SMark F. Adams CHKERRQ(ierr); 38392a756f0SMark F. Adams ierr = PetscFree( tidx ); CHKERRQ(ierr); 38411e60469SMark F. Adams /* 38511e60469SMark F. Adams Create a vector to contain the original vertex information for each element 38611e60469SMark F. Adams */ 387d3d6bff4SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 388d3d6bff4SMark F. Adams for( jj=0; jj<a_ndata_cols ; jj++ ) { 389d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs0 ; ii++) { 390d3d6bff4SMark F. Adams for( kk=0; kk<a_ndata_rows ; kk++ ) { 391d3d6bff4SMark F. Adams PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj; 392676e1743SMark F. Adams PetscScalar tt = (PetscScalar)data[ix]; 393676e1743SMark F. Adams ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES ); CHKERRQ(ierr); 394d3d6bff4SMark F. Adams } 395038e3b61SMark F. Adams } 396eb07cef2SMark F. Adams } 397eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 398eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 39911e60469SMark F. Adams /* 40011e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 40111e60469SMark F. Adams to the correct processor 40211e60469SMark F. Adams */ 40311e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 40411e60469SMark F. Adams CHKERRQ(ierr); 40511e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 40611e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40711e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40811e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 40911e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 41011e60469SMark F. Adams /* 41111e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 41211e60469SMark F. Adams */ 413eb07cef2SMark F. Adams ierr = PetscFree( *a_coarse_data ); CHKERRQ(ierr); 414d3d6bff4SMark F. Adams ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data ); CHKERRQ(ierr); 415eb07cef2SMark F. Adams 41611e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 417eb07cef2SMark F. Adams data = *a_coarse_data; 418d3d6bff4SMark F. Adams for( jj=0; jj<a_ndata_cols ; jj++ ) { 419d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs_new ; ii++) { 420d3d6bff4SMark F. Adams for( kk=0; kk<a_ndata_rows ; kk++ ) { 421d3d6bff4SMark F. Adams PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj; 422d3d6bff4SMark F. Adams data[ix] = PetscRealPart(array[jx]); 423d3d6bff4SMark F. Adams array[jx] = 1.e300; 424d3d6bff4SMark F. Adams } 425038e3b61SMark F. Adams } 426038e3b61SMark F. Adams } 42711e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 42811e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 42911e60469SMark F. Adams /* 43011e60469SMark F. Adams Invert for MatGetSubMatrix 43111e60469SMark F. Adams */ 432038e3b61SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr); 43311e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 43411e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 43511e60469SMark F. Adams /* A_crs output */ 436038e3b61SMark F. Adams ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 43711e60469SMark F. Adams CHKERRQ(ierr); 438eb07cef2SMark F. Adams 439038e3b61SMark F. Adams ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 440e33ef3b1SMark F. Adams Cmat = *a_Amat_crs; /* output */ 441038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 442eb07cef2SMark F. Adams 44311e60469SMark F. Adams /* prolongator */ 44411e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 44511e60469SMark F. Adams { 44611e60469SMark F. Adams IS findices; 44711e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 44896434bc1SMark F. Adams #ifdef USE_R 44996434bc1SMark F. Adams ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew ); 45096434bc1SMark F. Adams #else 45111e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 45296434bc1SMark F. Adams #endif 45311e60469SMark F. Adams CHKERRQ(ierr); 454d61a3a7dSMark F. Adams 45511e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 45611e60469SMark F. Adams } 457d61a3a7dSMark F. Adams 4583530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 4593530afc2SMark F. Adams *a_P_inout = Pnew; /* output */ 46092a756f0SMark F. Adams 46111e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 46292a756f0SMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 46392a756f0SMark F. Adams ierr = PetscFree( ranks ); CHKERRQ(ierr); 464e33ef3b1SMark F. Adams } 4655b89ad90SMark F. Adams 4665b89ad90SMark F. Adams PetscFunctionReturn(0); 4675b89ad90SMark F. Adams } 4685b89ad90SMark F. Adams 4695b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4705b89ad90SMark F. Adams /* 4715b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4725b89ad90SMark F. Adams by setting data structures and options. 4735b89ad90SMark F. Adams 4745b89ad90SMark F. Adams Input Parameter: 4755b89ad90SMark F. Adams . pc - the preconditioner context 4765b89ad90SMark F. Adams 4775b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 4785b89ad90SMark F. Adams 4795b89ad90SMark F. Adams Notes: 4805b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 4815b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 4825b89ad90SMark F. Adams */ 4835b89ad90SMark F. Adams #undef __FUNCT__ 4845b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 485eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc ) 4865b89ad90SMark F. Adams { 4875b89ad90SMark F. Adams PetscErrorCode ierr; 488eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 4895b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 49058471d46SMark F. Adams PC_MG_Levels **mglevels = mg->levels; 491eb07cef2SMark F. Adams Mat Amat = a_pc->mat, Pmat = a_pc->pmat; 492d3d6bff4SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend; 493eb07cef2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_pc)->comm; 4943530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 495038e3b61SMark F. Adams PetscBool isOK; 496587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 497587fa25dSMark F. Adams PetscReal *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS]; 498737a81a9SMark F. Adams MatInfo info; 4995ef31b24SMark F. Adams 5005b89ad90SMark F. Adams PetscFunctionBegin; 50103a628feSMark F. Adams pc_gamg->m_count++; 502fca9ed99SMark F. Adams 50358471d46SMark F. Adams if( a_pc->setupcalled > 0 ) { 50403a628feSMark F. Adams /* just do Galerkin grids */ 50558471d46SMark F. Adams Mat B,dA,dB; 50658471d46SMark F. Adams 50758471d46SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 50858471d46SMark F. Adams ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 50958471d46SMark F. Adams 51058471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 51158471d46SMark F. Adams ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 51258471d46SMark F. Adams /* (re)set to get dirty flag */ 51358471d46SMark F. Adams ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 51458471d46SMark F. Adams ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr); 51558471d46SMark F. Adams 51658471d46SMark F. Adams for (level=pc_gamg->m_Nlevels-2; level>-1; level--) { 51758471d46SMark F. Adams ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 51803a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 51903a628feSMark F. Adams if( pc_gamg->m_count == 2 ) { 52003a628feSMark F. Adams ierr = MatDestroy( &B ); CHKERRQ(ierr); 52103a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 52203a628feSMark F. Adams mglevels[level]->A = B; 52303a628feSMark F. Adams } 52403a628feSMark F. Adams else { 52558471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 52603a628feSMark F. Adams } 52758471d46SMark F. Adams ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 52858471d46SMark F. Adams dB = B; 52958471d46SMark F. Adams /* setup KSP/PC */ 53058471d46SMark F. Adams ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr); 53158471d46SMark F. Adams } 53258471d46SMark F. Adams 533f6536408SMark F. Adams #define PRINT_MATS PETSC_FALSE 53458471d46SMark F. Adams /* plot levels - A */ 53558471d46SMark F. Adams if( PRINT_MATS ) { 53658471d46SMark F. Adams for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){ 53758471d46SMark F. Adams PetscViewer viewer; 53858471d46SMark F. Adams char fname[32]; KSP smoother; Mat Tmat, TTm; 53958471d46SMark F. Adams ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 54058471d46SMark F. Adams ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr); 54103a628feSMark F. Adams sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level); 54258471d46SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 54358471d46SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 54458471d46SMark F. Adams ierr = MatView( Tmat, viewer ); CHKERRQ(ierr); 54558471d46SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 54658471d46SMark F. Adams } 54758471d46SMark F. Adams } 54858471d46SMark F. Adams 54903a628feSMark F. Adams a_pc->setupcalled = 2; 55003a628feSMark F. Adams 55158471d46SMark F. Adams PetscFunctionReturn(0); 552eb07cef2SMark F. Adams } 553f6536408SMark F. Adams 554baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 555baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 5565b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 5575b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 5583530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 559eb07cef2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 560eb07cef2SMark F. Adams nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 561eb07cef2SMark F. Adams 562038e3b61SMark F. Adams /* get data of not around */ 563038e3b61SMark F. Adams if( pc_gamg->m_data == 0 && nloc > 0 ) { 564038e3b61SMark F. Adams ierr = PCSetCoordinates_GAMG( a_pc, -1, 0 ); CHKERRQ( ierr ); 565eb07cef2SMark F. Adams } 566eb07cef2SMark F. Adams data = pc_gamg->m_data; 567038e3b61SMark F. Adams 568eb07cef2SMark F. Adams /* Get A_i and R_i */ 569737a81a9SMark F. Adams ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 570bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 5712c047e2dSMark F. Adams mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols, 5722c047e2dSMark F. Adams (int)(info.nz_used/(PetscReal)N),npe); 5738f4b7eb5SMark F. Adams for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 5744ef23d27SMark F. Adams level < (pc_gamg->m_Nlevels-1) && (level==0 || M>pc_gamg->m_coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 5750205a208SMark F. Adams level++ ){ 5765b89ad90SMark F. Adams level1 = level + 1; 577b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 578b4fbaa2aSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 579b4fbaa2aSMark F. Adams #endif 580b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 581b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 582b4fbaa2aSMark F. Adams #endif 583389730f3SMark F. Adams ierr = createProlongation(Aarr[level], data, level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level], pc_gamg ); 5845b89ad90SMark F. Adams CHKERRQ(ierr); 585d3d6bff4SMark F. Adams ierr = PetscFree( data ); CHKERRQ( ierr ); 586b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 587b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 5885c7a89a6SBarry Smith 589b4fbaa2aSMark F. Adams #endif 590baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 591baa4e9faSMark F. Adams if( isOK ) { 592b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 593b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 594b4fbaa2aSMark F. Adams #endif 5958263b398SMark F. Adams ierr = createLevel( a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs, 5968263b398SMark F. Adams (PetscBool)(level == pc_gamg->m_Nlevels-2), 597389730f3SMark F. Adams &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] ); 5983530afc2SMark F. Adams CHKERRQ(ierr); 599b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 600b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 601b4fbaa2aSMark F. Adams #endif 6023530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 603737a81a9SMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 604bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 6052c047e2dSMark F. Adams mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols, 6062c047e2dSMark F. Adams (int)(info.nz_used/(PetscReal)N),nactivepe); 607e33ef3b1SMark F. Adams /* coarse grids with SA can have zero row/cols from singleton aggregates */ 60858471d46SMark F. Adams /* aggregation method should gaurrentee this does not happen! */ 609737a81a9SMark F. Adams 610bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) { 611785cba28SMark F. Adams Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; 612785cba28SMark F. Adams v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */ 613e33ef3b1SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); 614785cba28SMark F. Adams nloceq = Iend-Istart; 615e33ef3b1SMark F. Adams ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); 616e33ef3b1SMark F. Adams ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); 617e33ef3b1SMark F. Adams ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); 618785cba28SMark F. Adams for(kk=0;kk<nloceq;kk++){ 619e33ef3b1SMark F. Adams if(data_arr[kk]==0.0) { 620785cba28SMark F. Adams id = kk + Istart; 621e33ef3b1SMark F. Adams ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); 622e33ef3b1SMark F. Adams CHKERRQ(ierr); 623ecd57171SMark F. Adams PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); 624e33ef3b1SMark F. Adams } 625e33ef3b1SMark F. Adams } 626e33ef3b1SMark F. Adams ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); 627e33ef3b1SMark F. Adams ierr = VecDestroy( &diag ); CHKERRQ(ierr); 628e33ef3b1SMark F. Adams ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 629e33ef3b1SMark F. Adams ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 630e33ef3b1SMark F. Adams } 631486a8d0bSJed Brown } else { 632be544d3cSMark F. Adams coarse_data = 0; 633baa4e9faSMark F. Adams break; 634baa4e9faSMark F. Adams } 635eb07cef2SMark F. Adams data = coarse_data; 636b4fbaa2aSMark F. Adams 637b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 638b4fbaa2aSMark F. Adams ierr = PetscLogStagePop(); CHKERRQ( ierr ); 639b4fbaa2aSMark F. Adams #endif 6405b89ad90SMark F. Adams } 641be544d3cSMark F. Adams if( coarse_data ) { 642eb07cef2SMark F. Adams ierr = PetscFree( coarse_data ); CHKERRQ( ierr ); 643be544d3cSMark F. Adams } 644bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 6455b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 6465b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 6475b89ad90SMark F. Adams fine_level = level; 648eb07cef2SMark F. Adams ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 6495b89ad90SMark F. Adams 650fc4362bfSMark F. Adams /* set default smoothers */ 651587fa25dSMark F. Adams for ( lidx=1, level = pc_gamg->m_Nlevels-2; 652587fa25dSMark F. Adams lidx <= fine_level; 653587fa25dSMark F. Adams lidx++, level--) { 6545745e0f5SMark F. Adams PetscReal emax, emin; 6555b89ad90SMark F. Adams KSP smoother; PC subpc; 656f6536408SMark F. Adams PetscBool isCheb; 657f6536408SMark F. Adams /* set defaults */ 658587fa25dSMark F. Adams ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 6595b89ad90SMark F. Adams ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 660f6536408SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 661f6536408SMark F. Adams /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */ 662f6536408SMark F. Adams ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); 663f6536408SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 664f6536408SMark F. Adams /* overide defaults with input parameters */ 665f6536408SMark F. Adams ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr); 666f6536408SMark F. Adams 667f6536408SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 668f6536408SMark F. Adams /* do my own cheby */ 669f6536408SMark F. Adams ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr); 670f6536408SMark F. Adams if( isCheb ) { 671f6536408SMark F. Adams ierr = PetscTypeCompare( (PetscObject)subpc, PETSC_GAMG_SMOOTHER, &isCheb ); CHKERRQ(ierr); 672f6536408SMark F. Adams if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 673587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 674587fa25dSMark F. Adams KSP eksp; Mat Lmat = Aarr[level]; 675fc4362bfSMark F. Adams Vec bb, xx; PC pc; 676f6536408SMark F. Adams const PCType type; 677038e3b61SMark F. Adams 678f6536408SMark F. Adams ierr = PCGetType( subpc, &type ); CHKERRQ(ierr); 6795745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 6805745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 681fc4362bfSMark F. Adams { 682fc4362bfSMark F. Adams PetscRandom rctx; 683fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 684fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 685fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 686fc4362bfSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 6875b89ad90SMark F. Adams } 688fc4362bfSMark F. Adams ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 689e33ef3b1SMark F. Adams ierr = KSPSetType( eksp, KSPCG ); CHKERRQ(ierr); 690038e3b61SMark F. Adams ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 691fc4362bfSMark F. Adams CHKERRQ(ierr); 692fc4362bfSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 693fc4362bfSMark F. Adams 694f6536408SMark F. Adams ierr = KSPAppendOptionsPrefix( eksp, "est_"); CHKERRQ(ierr); 695f6536408SMark F. Adams ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 696f6536408SMark F. Adams 697f6536408SMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 698f6536408SMark F. Adams ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 699f6536408SMark F. Adams ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 700f6536408SMark F. Adams ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */ 701f6536408SMark F. Adams 702fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 703fc4362bfSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 704fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 705fc4362bfSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 706fc4362bfSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 707fc4362bfSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 708f6536408SMark F. Adams 709f6536408SMark F. Adams if (pc_gamg->m_verbose) { 710f6536408SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e PC=%s\n", 711f6536408SMark F. Adams __FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER); 712f6536408SMark F. Adams } 713fc4362bfSMark F. Adams } 714038e3b61SMark F. Adams { 715038e3b61SMark F. Adams PetscInt N1, N0, tt; 716038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 717038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 718f6536408SMark F. Adams /* heuristic - is this crap? */ 719f6536408SMark F. Adams emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); 720038e3b61SMark F. Adams emax *= 1.05; 721038e3b61SMark F. Adams } 722fc4362bfSMark F. Adams ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 723f6536408SMark F. Adams } 7245745e0f5SMark F. Adams } 7255745e0f5SMark F. Adams { 726ecd57171SMark F. Adams /* coarse grid */ 727ecd57171SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 7285745e0f5SMark F. Adams Mat Lmat = Aarr[pc_gamg->m_Nlevels-1]; 729eb07cef2SMark F. Adams ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 73058471d46SMark F. Adams ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 7315745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 732ecd57171SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 733ecd57171SMark F. Adams ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 734ecd57171SMark F. Adams ierr = PCSetUp( subpc ); CHKERRQ(ierr); 735ecd57171SMark F. Adams ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 736ecd57171SMark F. Adams assert(ii==1); 737ecd57171SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 738ecd57171SMark F. Adams ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 739fc4362bfSMark F. Adams } 740737a81a9SMark F. Adams 741fc4362bfSMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 742eb07cef2SMark F. Adams ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr); 7435b89ad90SMark F. Adams { 7445b89ad90SMark F. Adams PetscBool galerkin; 745eb07cef2SMark F. Adams ierr = PCMGGetGalerkin( a_pc, &galerkin); CHKERRQ(ierr); 7465b89ad90SMark F. Adams if(galerkin){ 7475b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 7485b89ad90SMark F. Adams } 7495b89ad90SMark F. Adams } 7505745e0f5SMark F. Adams 75158471d46SMark F. Adams /* plot levels - R/P */ 75258471d46SMark F. Adams if( PRINT_MATS ) { 75358471d46SMark F. Adams for (level=pc_gamg->m_Nlevels-1;level>0;level--){ 75458471d46SMark F. Adams PetscViewer viewer; 75558471d46SMark F. Adams char fname[32]; 75603a628feSMark F. Adams sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level); 75711e60469SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 7585b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 759038e3b61SMark F. Adams ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr); 7605b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 76103a628feSMark F. Adams sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level); 762e33ef3b1SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 763e33ef3b1SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 764e33ef3b1SMark F. Adams ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr); 765e33ef3b1SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 7665b89ad90SMark F. Adams } 76758471d46SMark F. Adams } 76858471d46SMark F. Adams 76958471d46SMark F. Adams /* set interpolation between the levels, clean up */ 77058471d46SMark F. Adams for (lidx=0,level=pc_gamg->m_Nlevels-1; 77158471d46SMark F. Adams lidx<fine_level; 77258471d46SMark F. Adams lidx++, level--){ 77358471d46SMark F. Adams ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr); 774587fa25dSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 775587fa25dSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 7765b89ad90SMark F. Adams } 7775745e0f5SMark F. Adams 7785b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 779eb07cef2SMark F. Adams a_pc->setupcalled = 0; 780eb07cef2SMark F. Adams ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr ); 78103a628feSMark F. Adams a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ 78258471d46SMark F. Adams 78358471d46SMark F. Adams { 78458471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 78558471d46SMark F. Adams ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 78658471d46SMark F. Adams ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 78758471d46SMark F. Adams ierr = KSPSetUp( smoother ); CHKERRQ(ierr); 78858471d46SMark F. Adams } 789e33ef3b1SMark F. Adams 7905b89ad90SMark F. Adams PetscFunctionReturn(0); 7915b89ad90SMark F. Adams } 7925b89ad90SMark F. Adams 793eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 7945b89ad90SMark F. Adams /* 7955b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 7965b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 7975b89ad90SMark F. Adams 7985b89ad90SMark F. Adams Input Parameter: 7995b89ad90SMark F. Adams . pc - the preconditioner context 8005b89ad90SMark F. Adams 8015b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 8025b89ad90SMark F. Adams */ 8035b89ad90SMark F. Adams #undef __FUNCT__ 8045b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 8055b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 8065b89ad90SMark F. Adams { 8075b89ad90SMark F. Adams PetscErrorCode ierr; 8085b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 8095b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 8105b89ad90SMark F. Adams 8115b89ad90SMark F. Adams PetscFunctionBegin; 8125b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 8135b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 8145b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 8155b89ad90SMark F. Adams PetscFunctionReturn(0); 8165b89ad90SMark F. Adams } 8175b89ad90SMark F. Adams 818676e1743SMark F. Adams 819676e1743SMark F. Adams #undef __FUNCT__ 820676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 821676e1743SMark F. Adams /*@ 822676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 823676e1743SMark F. Adams processor reduction. 824676e1743SMark F. Adams 825676e1743SMark F. Adams Not Collective on PC 826676e1743SMark F. Adams 827676e1743SMark F. Adams Input Parameters: 828676e1743SMark F. Adams . pc - the preconditioner context 829676e1743SMark F. Adams 830676e1743SMark F. Adams 831676e1743SMark F. Adams Options Database Key: 832676e1743SMark F. Adams . -pc_gamg_process_eq_limit 833676e1743SMark F. Adams 834676e1743SMark F. Adams Level: intermediate 835676e1743SMark F. Adams 836676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 837676e1743SMark F. Adams 838676e1743SMark F. Adams .seealso: () 839676e1743SMark F. Adams @*/ 840676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 841676e1743SMark F. Adams { 842676e1743SMark F. Adams PetscErrorCode ierr; 843676e1743SMark F. Adams 844676e1743SMark F. Adams PetscFunctionBegin; 845676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 846676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 847676e1743SMark F. Adams PetscFunctionReturn(0); 848676e1743SMark F. Adams } 849676e1743SMark F. Adams 850676e1743SMark F. Adams EXTERN_C_BEGIN 851676e1743SMark F. Adams #undef __FUNCT__ 852676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 853676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 854676e1743SMark F. Adams { 855c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 856c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 857676e1743SMark F. Adams 858676e1743SMark F. Adams PetscFunctionBegin; 859c20e4228SMark F. Adams if(n>0) pc_gamg->m_min_eq_proc = n; 860676e1743SMark F. Adams PetscFunctionReturn(0); 861676e1743SMark F. Adams } 862676e1743SMark F. Adams EXTERN_C_END 863676e1743SMark F. Adams 864676e1743SMark F. Adams #undef __FUNCT__ 865389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 866389730f3SMark F. Adams /*@ 867389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 868389730f3SMark F. Adams 869389730f3SMark F. Adams Collective on PC 870389730f3SMark F. Adams 871389730f3SMark F. Adams Input Parameters: 872389730f3SMark F. Adams . pc - the preconditioner context 873389730f3SMark F. Adams 874389730f3SMark F. Adams 875389730f3SMark F. Adams Options Database Key: 876389730f3SMark F. Adams . -pc_gamg_coarse_eq_limit 877389730f3SMark F. Adams 878389730f3SMark F. Adams Level: intermediate 879389730f3SMark F. Adams 880389730f3SMark F. Adams Concepts: Unstructured multrigrid preconditioner 881389730f3SMark F. Adams 882389730f3SMark F. Adams .seealso: () 883389730f3SMark F. Adams @*/ 884389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 885389730f3SMark F. Adams { 886389730f3SMark F. Adams PetscErrorCode ierr; 887389730f3SMark F. Adams 888389730f3SMark F. Adams PetscFunctionBegin; 889389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 890389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 891389730f3SMark F. Adams PetscFunctionReturn(0); 892389730f3SMark F. Adams } 893389730f3SMark F. Adams 894389730f3SMark F. Adams EXTERN_C_BEGIN 895389730f3SMark F. Adams #undef __FUNCT__ 896389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 897389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 898389730f3SMark F. Adams { 899389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 900389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 901389730f3SMark F. Adams 902389730f3SMark F. Adams PetscFunctionBegin; 903389730f3SMark F. Adams if(n>0) pc_gamg->m_coarse_eq_limit = n; 904389730f3SMark F. Adams PetscFunctionReturn(0); 905389730f3SMark F. Adams } 906389730f3SMark F. Adams EXTERN_C_END 907389730f3SMark F. Adams 908389730f3SMark F. Adams #undef __FUNCT__ 9098263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 910676e1743SMark F. Adams /*@ 9118263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 912676e1743SMark F. Adams 913676e1743SMark F. Adams Collective on PC 914676e1743SMark F. Adams 915676e1743SMark F. Adams Input Parameters: 916676e1743SMark F. Adams . pc - the preconditioner context 917676e1743SMark F. Adams 918676e1743SMark F. Adams 919676e1743SMark F. Adams Options Database Key: 9208263b398SMark F. Adams . -pc_gamg_repartition 921676e1743SMark F. Adams 922676e1743SMark F. Adams Level: intermediate 923676e1743SMark F. Adams 924676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 925676e1743SMark F. Adams 926676e1743SMark F. Adams .seealso: () 927676e1743SMark F. Adams @*/ 9288263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 929676e1743SMark F. Adams { 930676e1743SMark F. Adams PetscErrorCode ierr; 931676e1743SMark F. Adams 932676e1743SMark F. Adams PetscFunctionBegin; 933676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9348263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 935676e1743SMark F. Adams PetscFunctionReturn(0); 936676e1743SMark F. Adams } 937676e1743SMark F. Adams 938676e1743SMark F. Adams EXTERN_C_BEGIN 939676e1743SMark F. Adams #undef __FUNCT__ 9408263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 9418263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 942676e1743SMark F. Adams { 943c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 944c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 945676e1743SMark F. Adams 946676e1743SMark F. Adams PetscFunctionBegin; 9478263b398SMark F. Adams pc_gamg->m_repart = n; 948676e1743SMark F. Adams PetscFunctionReturn(0); 949676e1743SMark F. Adams } 950676e1743SMark F. Adams EXTERN_C_END 951676e1743SMark F. Adams 952676e1743SMark F. Adams #undef __FUNCT__ 9534ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 9544ef23d27SMark F. Adams /*@ 9554ef23d27SMark F. Adams PCGAMGSetNlevels - 9564ef23d27SMark F. Adams 9574ef23d27SMark F. Adams Not collective on PC 9584ef23d27SMark F. Adams 9594ef23d27SMark F. Adams Input Parameters: 9604ef23d27SMark F. Adams . pc - the preconditioner context 9614ef23d27SMark F. Adams 9624ef23d27SMark F. Adams 9634ef23d27SMark F. Adams Options Database Key: 9644ef23d27SMark F. Adams . -pc_mg_levels 9654ef23d27SMark F. Adams 9664ef23d27SMark F. Adams Level: intermediate 9674ef23d27SMark F. Adams 9684ef23d27SMark F. Adams Concepts: Unstructured multrigrid preconditioner 9694ef23d27SMark F. Adams 9704ef23d27SMark F. Adams .seealso: () 9714ef23d27SMark F. Adams @*/ 9724ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 9734ef23d27SMark F. Adams { 9744ef23d27SMark F. Adams PetscErrorCode ierr; 9754ef23d27SMark F. Adams 9764ef23d27SMark F. Adams PetscFunctionBegin; 9774ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9784ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 9794ef23d27SMark F. Adams PetscFunctionReturn(0); 9804ef23d27SMark F. Adams } 9814ef23d27SMark F. Adams 9824ef23d27SMark F. Adams EXTERN_C_BEGIN 9834ef23d27SMark F. Adams #undef __FUNCT__ 9844ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 9854ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 9864ef23d27SMark F. Adams { 9874ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9884ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9894ef23d27SMark F. Adams 9904ef23d27SMark F. Adams PetscFunctionBegin; 9914ef23d27SMark F. Adams pc_gamg->m_Nlevels = n; 9924ef23d27SMark F. Adams PetscFunctionReturn(0); 9934ef23d27SMark F. Adams } 9944ef23d27SMark F. Adams EXTERN_C_END 9954ef23d27SMark F. Adams 9964ef23d27SMark F. Adams #undef __FUNCT__ 9973542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 9983542efc5SMark F. Adams /*@ 9993542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 10003542efc5SMark F. Adams 10013542efc5SMark F. Adams Not collective on PC 10023542efc5SMark F. Adams 10033542efc5SMark F. Adams Input Parameters: 10043542efc5SMark F. Adams . pc - the preconditioner context 10053542efc5SMark F. Adams 10063542efc5SMark F. Adams 10073542efc5SMark F. Adams Options Database Key: 10083542efc5SMark F. Adams . -pc_gamg_threshold 10093542efc5SMark F. Adams 10103542efc5SMark F. Adams Level: intermediate 10113542efc5SMark F. Adams 10123542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 10133542efc5SMark F. Adams 10143542efc5SMark F. Adams .seealso: () 10153542efc5SMark F. Adams @*/ 10163542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 10173542efc5SMark F. Adams { 10183542efc5SMark F. Adams PetscErrorCode ierr; 10193542efc5SMark F. Adams 10203542efc5SMark F. Adams PetscFunctionBegin; 10213542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10223542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 10233542efc5SMark F. Adams PetscFunctionReturn(0); 10243542efc5SMark F. Adams } 10253542efc5SMark F. Adams 10263542efc5SMark F. Adams EXTERN_C_BEGIN 10273542efc5SMark F. Adams #undef __FUNCT__ 10283542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 10293542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 10303542efc5SMark F. Adams { 1031c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1032c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10333542efc5SMark F. Adams 10343542efc5SMark F. Adams PetscFunctionBegin; 1035c20e4228SMark F. Adams pc_gamg->m_threshold = n; 10363542efc5SMark F. Adams PetscFunctionReturn(0); 10373542efc5SMark F. Adams } 10383542efc5SMark F. Adams EXTERN_C_END 10393542efc5SMark F. Adams 10403542efc5SMark F. Adams #undef __FUNCT__ 1041676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType" 1042676e1743SMark F. Adams /*@ 1043676e1743SMark F. Adams PCGAMGSetSolverType - Set solution method. 1044676e1743SMark F. Adams 1045676e1743SMark F. Adams Collective on PC 1046676e1743SMark F. Adams 1047676e1743SMark F. Adams Input Parameters: 1048676e1743SMark F. Adams . pc - the preconditioner context 1049676e1743SMark F. Adams 1050676e1743SMark F. Adams 1051676e1743SMark F. Adams Options Database Key: 10523542efc5SMark F. Adams . -pc_gamg_type 1053676e1743SMark F. Adams 1054676e1743SMark F. Adams Level: intermediate 1055676e1743SMark F. Adams 1056676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1057676e1743SMark F. Adams 1058676e1743SMark F. Adams .seealso: () 1059676e1743SMark F. Adams @*/ 1060676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz ) 1061676e1743SMark F. Adams { 1062676e1743SMark F. Adams PetscErrorCode ierr; 1063676e1743SMark F. Adams 1064676e1743SMark F. Adams PetscFunctionBegin; 1065676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1066676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz)); 1067676e1743SMark F. Adams CHKERRQ(ierr); 1068676e1743SMark F. Adams PetscFunctionReturn(0); 1069676e1743SMark F. Adams } 1070676e1743SMark F. Adams 1071676e1743SMark F. Adams EXTERN_C_BEGIN 1072676e1743SMark F. Adams #undef __FUNCT__ 1073676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG" 1074676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz ) 1075676e1743SMark F. Adams { 1076676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1077676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1078676e1743SMark F. Adams 1079676e1743SMark F. Adams PetscFunctionBegin; 1080676e1743SMark F. Adams if(sz < 64) strcpy(pc_gamg->m_type,str); 1081676e1743SMark F. Adams PetscFunctionReturn(0); 1082676e1743SMark F. Adams } 1083676e1743SMark F. Adams EXTERN_C_END 1084676e1743SMark F. Adams 10855b89ad90SMark F. Adams #undef __FUNCT__ 10865b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 10875b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 10885b89ad90SMark F. Adams { 1089676e1743SMark F. Adams PetscErrorCode ierr; 1090676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1091676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1092676e1743SMark F. Adams PetscBool flag; 10935b89ad90SMark F. Adams 10945b89ad90SMark F. Adams PetscFunctionBegin; 1095676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1096676e1743SMark F. Adams { 1097676e1743SMark F. Adams ierr = PetscOptionsString("-pc_gamg_type", 1098470e26f8SMark F. Adams "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)", 1099676e1743SMark F. Adams "PCGAMGSetSolverType", 1100676e1743SMark F. Adams pc_gamg->m_type, 1101676e1743SMark F. Adams pc_gamg->m_type, 1102676e1743SMark F. Adams 64, 1103676e1743SMark F. Adams &flag ); 1104676e1743SMark F. Adams CHKERRQ(ierr); 1105*d8c859f0SMark F. Adams if( flag && pc_gamg->m_data != 0 ) { 1106*d8c859f0SMark F. Adams if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) || 1107*d8c859f0SMark F. Adams (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) || 1108*d8c859f0SMark F. Adams (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) { 1109*d8c859f0SMark F. Adams SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)"); 1110*d8c859f0SMark F. Adams } 1111*d8c859f0SMark F. Adams } 1112bed7c2b9SMark F. Adams 1113bed7c2b9SMark F. Adams /* -pc_gamg_verbose */ 1114bed7c2b9SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none",pc_gamg->m_verbose,&pc_gamg->m_verbose,PETSC_NULL);CHKERRQ(ierr); 11152a44bfbaSMark F. Adams 1116f6536408SMark F. Adams pc_gamg->m_method = 1; /* default to plane aggregation */ 1117f6536408SMark F. Adams if (flag ) { 1118f6536408SMark F. Adams if( strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2; 1119f6536408SMark F. Adams else if( strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1; 1120f6536408SMark F. Adams else if( strcmp(pc_gamg->m_type,"geo") == 0) pc_gamg->m_method = 0; 1121f6536408SMark F. Adams else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid gamg type: %s",pc_gamg->m_type); 1122f6536408SMark F. Adams } 11238263b398SMark F. Adams /* -pc_gamg_repartition */ 11248263b398SMark F. Adams pc_gamg->m_repart = PETSC_FALSE; 11258263b398SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_repartition", 11268263b398SMark F. Adams "Repartion coarse grids (false)", 11278263b398SMark F. Adams "PCGAMGRepartitioning", 11288263b398SMark F. Adams pc_gamg->m_repart, 11298263b398SMark F. Adams &pc_gamg->m_repart, 1130676e1743SMark F. Adams &flag); 1131676e1743SMark F. Adams CHKERRQ(ierr); 1132676e1743SMark F. Adams 1133c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 1134c20e4228SMark F. Adams pc_gamg->m_min_eq_proc = 600; 1135676e1743SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1136676e1743SMark F. Adams "Limit (goal) on number of equations per process on coarse grids", 1137676e1743SMark F. Adams "PCGAMGSetProcEqLim", 1138c20e4228SMark F. Adams pc_gamg->m_min_eq_proc, 1139c20e4228SMark F. Adams &pc_gamg->m_min_eq_proc, 1140676e1743SMark F. Adams &flag ); 1141676e1743SMark F. Adams CHKERRQ(ierr); 11423542efc5SMark F. Adams 1143389730f3SMark F. Adams /* -pc_gamg_coarse_eq_limit */ 1144389730f3SMark F. Adams pc_gamg->m_coarse_eq_limit = 1500; 1145389730f3SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1146389730f3SMark F. Adams "Limit on number of equations for the coarse grid", 1147389730f3SMark F. Adams "PCGAMGSetCoarseEqLim", 1148389730f3SMark F. Adams pc_gamg->m_coarse_eq_limit, 1149389730f3SMark F. Adams &pc_gamg->m_coarse_eq_limit, 1150389730f3SMark F. Adams &flag ); 1151389730f3SMark F. Adams CHKERRQ(ierr); 1152389730f3SMark F. Adams 1153c20e4228SMark F. Adams /* -pc_gamg_threshold */ 1154c20e4228SMark F. Adams pc_gamg->m_threshold = 0.05; 11553542efc5SMark F. Adams ierr = PetscOptionsReal("-pc_gamg_threshold", 11563542efc5SMark F. Adams "Relative threshold to use for dropping edges in aggregation graph", 11573542efc5SMark F. Adams "PCGAMGSetThreshold", 1158c20e4228SMark F. Adams pc_gamg->m_threshold, 1159c20e4228SMark F. Adams &pc_gamg->m_threshold, 11603542efc5SMark F. Adams &flag ); 11613542efc5SMark F. Adams CHKERRQ(ierr); 11624ef23d27SMark F. Adams 11634ef23d27SMark F. Adams pc_gamg->m_Nlevels = GAMG_MAXLEVELS; 11644ef23d27SMark F. Adams ierr = PetscOptionsInt("-pc_mg_levels", 11654ef23d27SMark F. Adams "Set number of MG levels", 11664ef23d27SMark F. Adams "PCGAMGSetNlevels", 11674ef23d27SMark F. Adams pc_gamg->m_Nlevels, 11684ef23d27SMark F. Adams &pc_gamg->m_Nlevels, 11694ef23d27SMark F. Adams &flag ); 1170676e1743SMark F. Adams } 1171676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1172676e1743SMark F. Adams 11735b89ad90SMark F. Adams PetscFunctionReturn(0); 11745b89ad90SMark F. Adams } 11755b89ad90SMark F. Adams 11765b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 11775b89ad90SMark F. Adams /*MC 1178dc76db92SMark F. Adams PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two 1179dc76db92SMark F. Adams AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each 1180dc76db92SMark F. Adams vertex, and 2) smoothed aggregation. Smoothed aggregation (SA) can work without coordinates but it 1181dc76db92SMark F. Adams will generate some common non-trivial null spaces if coordinates are provided. The input fine grid matrix 1182dc76db92SMark F. Adams must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly. 1183dc76db92SMark F. Adams SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational 1184dc76db92SMark F. Adams modes, when coordinates are provide in 2D and 3D. 11855b89ad90SMark F. Adams 1186280d9858SJed Brown Options Database Keys: 11875b89ad90SMark F. Adams Multigrid options(inherited) 1188280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1189280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1190280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1191280d9858SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 11925b89ad90SMark F. Adams 11935b89ad90SMark F. Adams Level: intermediate 1194280d9858SJed Brown 11955b89ad90SMark F. Adams Concepts: multigrid 11965b89ad90SMark F. Adams 11975b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1198280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 11995b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 12005b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 12015b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 12025b89ad90SMark F. Adams M*/ 12035b89ad90SMark F. Adams 12045b89ad90SMark F. Adams EXTERN_C_BEGIN 12055b89ad90SMark F. Adams #undef __FUNCT__ 12065b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 12075b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 12085b89ad90SMark F. Adams { 12095b89ad90SMark F. Adams PetscErrorCode ierr; 12105b89ad90SMark F. Adams PC_GAMG *pc_gamg; 12115b89ad90SMark F. Adams PC_MG *mg; 12125ef31b24SMark F. Adams PetscClassId cookie; 12135ee9c036SSatish Balay #if defined PETSC_USE_LOG 12145ee9c036SSatish Balay static int count = 0; 12155ee9c036SSatish Balay #endif 12165b89ad90SMark F. Adams 12175b89ad90SMark F. Adams PetscFunctionBegin; 12185b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 12195b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 12205b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 12215b89ad90SMark F. Adams 12225b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 12235b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 122403a628feSMark F. Adams pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0; 12255b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 12265b89ad90SMark F. Adams mg->innerctx = pc_gamg; 12275b89ad90SMark F. Adams 12285b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 12295b89ad90SMark F. Adams 12305b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 12315b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 12325b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 12335b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 12345b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 12355b89ad90SMark F. Adams 12365b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 12375b89ad90SMark F. Adams "PCSetCoordinates_C", 12385b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 1239c97e1df0SMark F. Adams PCSetCoordinates_GAMG); 1240c97e1df0SMark F. Adams CHKERRQ(ierr); 1241c97e1df0SMark F. Adams 1242676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1243676e1743SMark F. Adams "PCGAMGSetProcEqLim_C", 1244676e1743SMark F. Adams "PCGAMGSetProcEqLim_GAMG", 1245676e1743SMark F. Adams PCGAMGSetProcEqLim_GAMG); 1246676e1743SMark F. Adams CHKERRQ(ierr); 1247676e1743SMark F. Adams 1248676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1249389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_C", 1250389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_GAMG", 1251389730f3SMark F. Adams PCGAMGSetCoarseEqLim_GAMG); 1252389730f3SMark F. Adams CHKERRQ(ierr); 1253389730f3SMark F. Adams 1254389730f3SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 12558263b398SMark F. Adams "PCGAMGSetRepartitioning_C", 12568263b398SMark F. Adams "PCGAMGSetRepartitioning_GAMG", 12578263b398SMark F. Adams PCGAMGSetRepartitioning_GAMG); 1258676e1743SMark F. Adams CHKERRQ(ierr); 1259676e1743SMark F. Adams 1260676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 12613542efc5SMark F. Adams "PCGAMGSetThreshold_C", 12623542efc5SMark F. Adams "PCGAMGSetThreshold_GAMG", 12633542efc5SMark F. Adams PCGAMGSetThreshold_GAMG); 12643542efc5SMark F. Adams CHKERRQ(ierr); 12653542efc5SMark F. Adams 12663542efc5SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1267676e1743SMark F. Adams "PCGAMGSetSolverType_C", 1268676e1743SMark F. Adams "PCGAMGSetSolverType_GAMG", 1269676e1743SMark F. Adams PCGAMGSetSolverType_GAMG); 1270676e1743SMark F. Adams CHKERRQ(ierr); 1271c97e1df0SMark F. Adams 1272b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 1273785cba28SMark F. Adams if( count++ == 0 ) { 12745ef31b24SMark F. Adams PetscClassIdRegister("GAMG Setup",&cookie); 1275b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); 12762a44bfbaSMark F. Adams PetscLogEventRegister(" Graph", cookie, &gamg_setup_events[GRAPH]); 12772a44bfbaSMark F. Adams PetscLogEventRegister(" G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); 12782a44bfbaSMark F. Adams PetscLogEventRegister(" G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); 12792a44bfbaSMark F. Adams PetscLogEventRegister(" G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); 1280b4fbaa2aSMark F. Adams PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); 1281b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); 1282b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); 1283b4fbaa2aSMark F. Adams PetscLogEventRegister(" search&set", cookie, &gamg_setup_events[FIND_V]); 1284b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: init", cookie, &gamg_setup_events[SET7]); 1285b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); */ 1286b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); 1287b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); 1288b4fbaa2aSMark F. Adams PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_events[SET12]); 1289f852f58cSMark F. Adams 1290fca9ed99SMark F. Adams /* PetscLogEventRegister(" PL 1", cookie, &gamg_setup_events[SET13]); */ 1291fca9ed99SMark F. Adams /* PetscLogEventRegister(" PL 2", cookie, &gamg_setup_events[SET14]); */ 1292fca9ed99SMark F. Adams /* PetscLogEventRegister(" PL 3", cookie, &gamg_setup_events[SET15]); */ 1293f852f58cSMark F. Adams 1294b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ 1295b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ 1296b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ 12975ef31b24SMark F. Adams 1298b4fbaa2aSMark F. Adams /* create timer stages */ 1299b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1300b4fbaa2aSMark F. Adams { 1301b4fbaa2aSMark F. Adams char str[32]; 1302b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d (finest)",0); 1303b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[0]); 1304b4fbaa2aSMark F. Adams PetscInt lidx; 1305b4fbaa2aSMark F. Adams for (lidx=1;lidx<9;lidx++){ 1306b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1307b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[lidx]); 1308b4fbaa2aSMark F. Adams } 1309b4fbaa2aSMark F. Adams } 1310b4fbaa2aSMark F. Adams #endif 1311b4fbaa2aSMark F. Adams } 1312b4fbaa2aSMark F. Adams #endif 13135b89ad90SMark F. Adams PetscFunctionReturn(0); 13145b89ad90SMark F. Adams } 13155b89ad90SMark F. Adams EXTERN_C_END 1316