15b89ad90SMark F. Adams /* 25b89ad90SMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4313a3e24SSatish Balay #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 5313a3e24SSatish Balay #include "private/matimpl.h" 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 /* Private context for the GAMG preconditioner */ 185b89ad90SMark F. Adams typedef struct gamg_TAG { 195b89ad90SMark F. Adams PetscInt m_dim; 205b89ad90SMark F. Adams PetscInt m_Nlevels; 215b89ad90SMark F. Adams PetscInt m_data_sz; 22d3d6bff4SMark F. Adams PetscInt m_data_rows; 23d3d6bff4SMark F. Adams PetscInt m_data_cols; 2403a628feSMark F. Adams PetscInt m_count; 252c047e2dSMark F. Adams PetscInt m_method; /* 0: geomg; 1: plain agg 'pa'; 2: smoothed agg 'sa' */ 265b89ad90SMark F. Adams PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 27676e1743SMark F. Adams char m_type[64]; 28c20e4228SMark F. Adams PetscBool m_avoid_repart; 29c20e4228SMark F. Adams PetscInt m_min_eq_proc; 30c20e4228SMark F. Adams PetscReal m_threshold; 31bed7c2b9SMark F. Adams PetscBool m_verbose; 325b89ad90SMark F. Adams } PC_GAMG; 335b89ad90SMark F. Adams 345b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 355b89ad90SMark F. Adams /* 365b89ad90SMark F. Adams PCSetCoordinates_GAMG 375b89ad90SMark F. Adams 385b89ad90SMark F. Adams Input Parameter: 395b89ad90SMark F. Adams . pc - the preconditioner context 405b89ad90SMark F. Adams */ 41a92563c5SMark F. Adams EXTERN_C_BEGIN 425b89ad90SMark F. Adams #undef __FUNCT__ 435b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 44eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords ) 455b89ad90SMark F. Adams { 46eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 475b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 486c237d78SBarry Smith PetscErrorCode ierr; 49d3d6bff4SMark F. Adams PetscInt arrsz,bs,my0,kk,ii,jj,nloc,Iend; 50038e3b61SMark F. Adams Mat Amat = a_pc->pmat; 515b89ad90SMark F. Adams 525b89ad90SMark F. Adams PetscFunctionBegin; 5358471d46SMark F. Adams PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 ); 54038e3b61SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 55d3d6bff4SMark F. Adams ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 56d3d6bff4SMark F. Adams nloc = (Iend-my0)/bs; 57d3d6bff4SMark F. Adams if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 58038e3b61SMark F. Adams 59d3d6bff4SMark F. Adams pc_gamg->m_data_rows = 1; 602a44bfbaSMark F. Adams if(a_coords==0 && pc_gamg->m_method==0) pc_gamg->m_method = 2; /* use SA if no coords */ 61470e26f8SMark F. Adams if( pc_gamg->m_method==0 ) pc_gamg->m_data_cols = a_ndm; /* coordinates */ 62038e3b61SMark F. Adams else{ /* SA: null space vectors */ 63d3d6bff4SMark F. Adams if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */ 64d3d6bff4SMark F. Adams else if(a_coords != 0 ) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */ 65d3d6bff4SMark F. Adams else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */ 66d3d6bff4SMark F. Adams pc_gamg->m_data_rows = bs; 67038e3b61SMark F. Adams } 68d3d6bff4SMark F. Adams arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols; 695b89ad90SMark F. Adams 70038e3b61SMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 716e3e101aSMark F. Adams if (pc_gamg->m_data==0 || (pc_gamg->m_data_sz != arrsz)) { 725b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 7333a2b366SJed Brown ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->m_data ); CHKERRQ(ierr); 745b89ad90SMark F. Adams } 75038e3b61SMark F. Adams for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.; 76eb07cef2SMark F. Adams pc_gamg->m_data[arrsz] = -99.; 77038e3b61SMark F. Adams /* copy data in - column oriented */ 78470e26f8SMark F. Adams if( pc_gamg->m_method != 0 ) { 79d3d6bff4SMark F. Adams const PetscInt M = Iend - my0; 80038e3b61SMark F. Adams for(kk=0;kk<nloc;kk++){ 81038e3b61SMark F. Adams PetscReal *data = &pc_gamg->m_data[kk*bs]; 82d3d6bff4SMark F. Adams if( pc_gamg->m_data_cols==1 ) *data = 1.0; 83038e3b61SMark F. Adams else { 84038e3b61SMark F. Adams for(ii=0;ii<bs;ii++) 85038e3b61SMark F. Adams for(jj=0;jj<bs;jj++) 86038e3b61SMark F. Adams if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */ 87038e3b61SMark F. Adams else data[ii*M + jj] = 0.0; 88eb07cef2SMark F. Adams if( a_coords != 0 ) { 89038e3b61SMark F. Adams if( a_ndm == 2 ){ /* rotational modes */ 90038e3b61SMark F. Adams data += 2*M; 91038e3b61SMark F. Adams data[0] = -a_coords[2*kk+1]; 92038e3b61SMark F. Adams data[1] = a_coords[2*kk]; 935b89ad90SMark F. Adams } 94eb07cef2SMark F. Adams else { 95038e3b61SMark F. Adams data += 3*M; 96038e3b61SMark F. Adams data[0] = 0.0; data[M+0] = a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1]; 97038e3b61SMark F. Adams data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = a_coords[3*kk]; 98038e3b61SMark F. Adams data[2] = a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk]; data[2*M+2] = 0.0; 99038e3b61SMark F. Adams } 100eb07cef2SMark F. Adams } 101eb07cef2SMark F. Adams } 102eb07cef2SMark F. Adams } 103eb07cef2SMark F. Adams } 104eb07cef2SMark F. Adams else { 105038e3b61SMark F. Adams for( kk = 0 ; kk < nloc ; kk++ ){ 106038e3b61SMark F. Adams for( ii = 0 ; ii < a_ndm ; ii++ ) { 107038e3b61SMark F. Adams pc_gamg->m_data[ii*nloc + kk] = a_coords[kk*a_ndm + ii]; 108eb07cef2SMark F. Adams } 109eb07cef2SMark F. Adams } 110038e3b61SMark F. Adams } 111eb07cef2SMark F. Adams assert(pc_gamg->m_data[arrsz] == -99.); 112038e3b61SMark F. Adams 1135b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 114eb07cef2SMark F. Adams pc_gamg->m_dim = a_ndm; 1155b89ad90SMark F. Adams 1165b89ad90SMark F. Adams PetscFunctionReturn(0); 1175b89ad90SMark F. Adams } 118a92563c5SMark F. Adams EXTERN_C_END 1195b89ad90SMark F. Adams 120d3d6bff4SMark F. Adams 121d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 122d3d6bff4SMark F. Adams #undef __FUNCT__ 123d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 124d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 125d3d6bff4SMark F. Adams { 126d3d6bff4SMark F. Adams PetscErrorCode ierr; 127d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 128d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 129d3d6bff4SMark F. Adams 130d3d6bff4SMark F. Adams PetscFunctionBegin; 13158471d46SMark F. Adams if( pc_gamg->m_data != 0 ) { /* this should not happen, cleaned up in SetUp */ 132d3d6bff4SMark F. Adams ierr = PetscFree(pc_gamg->m_data); CHKERRQ(ierr); 13358471d46SMark F. Adams } 134d3d6bff4SMark F. Adams pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0; 135d3d6bff4SMark F. Adams PetscFunctionReturn(0); 136d3d6bff4SMark F. Adams } 137d3d6bff4SMark F. Adams 1385b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1395b89ad90SMark F. Adams /* 140486a8d0bSJed Brown PCGAMGPartitionLevel 1415b89ad90SMark F. Adams 1425b89ad90SMark F. Adams Input Parameter: 1433530afc2SMark F. Adams . a_Amat_fine - matrix on this fine (k) level 144d3d6bff4SMark F. Adams . a_ndata_rows - size of data to move (coarse grid) 145d3d6bff4SMark F. Adams . a_ndata_cols - size of data to move (coarse grid) 1463530afc2SMark F. Adams In/Output Parameter: 1473530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 148eb07cef2SMark F. Adams . a_coarse_data - data that need to be moved 149afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 15011e60469SMark F. Adams Output Parameter: 1513530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 152c20e4228SMark F. Adams . a_avoid_repart - 153c20e4228SMark F. Adams . a_min_eq_proc - goal for number of eqs. to have per proc. after proc reduction 1545b89ad90SMark F. Adams */ 1555cb416c2SMark F. Adams 1565b89ad90SMark F. Adams #undef __FUNCT__ 157486a8d0bSJed Brown #define __FUNCT__ "PCGAMGPartitionLevel" 158486a8d0bSJed Brown PetscErrorCode PCGAMGPartitionLevel(PC pc, Mat a_Amat_fine, 159d3d6bff4SMark F. Adams PetscInt a_ndata_rows, 160d3d6bff4SMark F. Adams PetscInt a_ndata_cols, 161038e3b61SMark F. Adams PetscInt a_cbs, 1623530afc2SMark F. Adams Mat *a_P_inout, 163eb07cef2SMark F. Adams PetscReal **a_coarse_data, 164afc97cdcSMark F. Adams PetscMPIInt *a_nactive_proc, 165c20e4228SMark F. Adams Mat *a_Amat_crs, 166c20e4228SMark F. Adams PetscBool a_avoid_repart, 167c20e4228SMark F. Adams PetscInt a_min_eq_proc 16811e60469SMark F. Adams ) 1695b89ad90SMark F. Adams { 170486a8d0bSJed Brown PC_MG *mg = (PC_MG*)pc->data; 171486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1725b89ad90SMark F. Adams PetscErrorCode ierr; 173038e3b61SMark F. Adams Mat Cmat,Pnew,Pold=*a_P_inout; 17411e60469SMark F. Adams IS new_indices,isnum; 1753530afc2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 176fa31c87dSMark F. Adams PetscMPIInt mype,npe,new_npe,nactive; 177fa31c87dSMark F. Adams PetscInt neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0; 1785b89ad90SMark F. Adams 1795b89ad90SMark F. Adams PetscFunctionBegin; 18011e60469SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 18111e60469SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 18211e60469SMark F. Adams /* RAP */ 18396434bc1SMark F. Adams #ifdef USE_R 18496434bc1SMark F. Adams /* make R wih brute force for now */ 18596434bc1SMark F. Adams ierr = MatTranspose( Pold, Pnew ); 18696434bc1SMark F. Adams ierr = MatDestroy( &Pold ); CHKERRQ(ierr); 18796434bc1SMark F. Adams ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 18896434bc1SMark F. Adams Pold = Pnew; 18996434bc1SMark F. Adams #else 190038e3b61SMark F. Adams ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 19196434bc1SMark F. Adams #endif 192038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 193038e3b61SMark F. Adams ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); 194038e3b61SMark F. Adams ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0); 195038e3b61SMark F. Adams 196c20e4228SMark F. Adams if( a_avoid_repart) { 19722063be5SMark F. Adams *a_Amat_crs = Cmat; /* output */ 19822063be5SMark F. Adams } 19922063be5SMark F. Adams else { 20022063be5SMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 2015ef31b24SMark F. Adams Mat adj; 20222063be5SMark F. Adams const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols; 203b4fbaa2aSMark F. Adams const PetscInt stride0=ncrs0*a_ndata_rows; 20492a756f0SMark F. Adams PetscInt is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx; 205c9a0b8beSMark F. Adams /* create sub communicator */ 20671d60426SBarry Smith MPI_Comm cm; 207afc97cdcSMark F. Adams MPI_Group wg, g2; 208fd3c6afaSMark F. Adams PetscInt *counts,inpe; 209fd3c6afaSMark F. Adams PetscMPIInt *ranks; 21022063be5SMark F. Adams IS isscat; 21122063be5SMark F. Adams PetscScalar *array; 21222063be5SMark F. Adams Vec src_crd, dest_crd; 21322063be5SMark F. Adams PetscReal *data = *a_coarse_data; 21422063be5SMark F. Adams VecScatter vecscat; 215b4fbaa2aSMark F. Adams IS isnewproc; 216e33ef3b1SMark F. Adams 21722063be5SMark F. Adams /* get number of PEs to make active, reduce */ 21822063be5SMark F. Adams ierr = MatGetSize( Cmat, &neq, &NN ); CHKERRQ(ierr); 219c20e4228SMark F. Adams new_npe = neq/a_min_eq_proc; /* hardwire min. number of eq/proc */ 220c20e4228SMark F. Adams if( new_npe == 0 || neq < 2*a_min_eq_proc ) new_npe = 1; 22122063be5SMark F. Adams else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */ 22222063be5SMark F. Adams 22392a756f0SMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr); 224fd3c6afaSMark F. Adams ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 22592a756f0SMark F. Adams 226fd3c6afaSMark F. Adams ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr); 227afc97cdcSMark F. Adams assert(counts[mype]==ncrs0); 228afc97cdcSMark F. Adams /* count real active pes */ 22922063be5SMark F. Adams for( nactive = jj = 0 ; jj < npe ; jj++) { 230afc97cdcSMark F. Adams if( counts[jj] != 0 ) { 23122063be5SMark F. Adams ranks[nactive++] = jj; 232afc97cdcSMark F. Adams } 233afc97cdcSMark F. Adams } 2343303bcf9Sadams 2353303bcf9Sadams if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */ 23622063be5SMark F. Adams 237bed7c2b9SMark 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); 23858471d46SMark F. Adams 23922063be5SMark F. Adams *a_nactive_proc = new_npe; /* output */ 2402f03bc48SMark F. Adams 241afc97cdcSMark F. Adams ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr); 24222063be5SMark F. Adams ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr); 243afc97cdcSMark F. Adams ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr); 244b4fbaa2aSMark F. Adams 245afc97cdcSMark F. Adams ierr = MPI_Group_free( &wg ); CHKERRQ(ierr); 246afc97cdcSMark F. Adams ierr = MPI_Group_free( &g2 ); CHKERRQ(ierr); 247c9a0b8beSMark F. Adams 2485ef31b24SMark F. Adams /* MatPartitioningApply call MatConvert, which is collective */ 249b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 250b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 251b4fbaa2aSMark F. Adams #endif 252038e3b61SMark F. Adams if( a_cbs == 1) { 253038e3b61SMark F. Adams ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 254eb07cef2SMark F. Adams } 255eb07cef2SMark F. Adams else{ 256038e3b61SMark F. Adams /* make a scalar matrix to partition */ 257eb07cef2SMark F. Adams Mat tMat; 25858471d46SMark F. Adams PetscInt ncols,jj,Ii; 259b4fbaa2aSMark F. Adams const PetscScalar *vals; 260b4fbaa2aSMark F. Adams const PetscInt *idx; 26158471d46SMark F. Adams PetscInt *d_nnz; 2625ee9c036SSatish Balay static int llev = 0; 263b4fbaa2aSMark F. Adams 26458471d46SMark F. Adams ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); 26558471d46SMark F. Adams for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) { 26658471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 26758471d46SMark F. Adams d_nnz[jj] = ncols/a_cbs; 26858471d46SMark F. Adams if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0; 26958471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 27058471d46SMark F. Adams } 2716876a03eSMark F. Adams 272eb07cef2SMark F. Adams ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 273eb07cef2SMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 27458471d46SMark F. Adams 0, d_nnz, 0, d_nnz, 275eb07cef2SMark F. Adams &tMat ); 2766876a03eSMark F. Adams CHKERRQ(ierr); 27758471d46SMark F. Adams ierr = PetscFree( d_nnz ); CHKERRQ(ierr); 278eb07cef2SMark F. Adams 27922063be5SMark F. Adams for ( ii = Istart0; ii < Iend0; ii++ ) { 28022063be5SMark F. Adams PetscInt dest_row = ii/a_cbs; 28122063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 282eb07cef2SMark F. Adams for( jj = 0 ; jj < ncols ; jj++ ){ 283038e3b61SMark F. Adams PetscInt dest_col = idx[jj]/a_cbs; 284eb07cef2SMark F. Adams PetscScalar v = 1.0; 285eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 286eb07cef2SMark F. Adams } 28722063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 288eb07cef2SMark F. Adams } 289eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 290eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 291eb07cef2SMark F. Adams 292b4fbaa2aSMark F. Adams if( llev++ == -1 ) { 293b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 294b4fbaa2aSMark F. Adams sprintf(fname,"part_mat_%d.mat",llev); 295b4fbaa2aSMark F. Adams PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 296b4fbaa2aSMark F. Adams ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 297b4fbaa2aSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 298b4fbaa2aSMark F. Adams } 299b4fbaa2aSMark F. Adams 300eb07cef2SMark F. Adams ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 301eb07cef2SMark F. Adams 302eb07cef2SMark F. Adams ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 303eb07cef2SMark F. Adams } 304afc97cdcSMark F. Adams if( ncrs0 != 0 ){ 305b4fbaa2aSMark F. Adams const PetscInt *is_idx; 306b4fbaa2aSMark F. Adams MatPartitioning mpart; 3075ef31b24SMark F. Adams /* hack to fix global data that pmetis.c uses in 'adj' */ 30822063be5SMark F. Adams for( nactive = jj = 0 ; jj < npe ; jj++ ) { 309afc97cdcSMark F. Adams if( counts[jj] != 0 ) { 31022063be5SMark F. Adams adj->rmap->range[nactive++] = adj->rmap->range[jj]; 311afc97cdcSMark F. Adams } 3125ef31b24SMark F. Adams } 31322063be5SMark F. Adams adj->rmap->range[nactive] = adj->rmap->range[npe]; 3142f03bc48SMark F. Adams 31571d60426SBarry Smith ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr); 3165ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 31711e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 3185ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 31911e60469SMark F. Adams ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 32011e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 32122063be5SMark F. Adams 3225ef31b24SMark F. Adams /* collect IS info */ 3235ef31b24SMark F. Adams ierr = ISGetLocalSize( isnewproc, &is_sz ); CHKERRQ(ierr); 324038e3b61SMark F. Adams ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr); 3255ef31b24SMark F. Adams ierr = ISGetIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 326b4fbaa2aSMark F. Adams /* spread partitioning across machine - best way ??? */ 327ecd57171SMark F. Adams NN = 1; /*npe/new_npe;*/ 328b4fbaa2aSMark F. Adams for( kk = jj = 0 ; kk < is_sz ; kk++ ){ 329afc97cdcSMark F. Adams for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) { 33022063be5SMark F. Adams isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */ 331eb07cef2SMark F. Adams } 3325ef31b24SMark F. Adams } 3335ef31b24SMark F. Adams ierr = ISRestoreIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 3345ef31b24SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 33571d60426SBarry Smith ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 336b4fbaa2aSMark F. Adams 337b4fbaa2aSMark F. Adams is_sz *= a_cbs; 3385ef31b24SMark F. Adams } 3395ef31b24SMark F. Adams else{ 3405ef31b24SMark F. Adams isnewproc_idx = 0; 3415ef31b24SMark F. Adams is_sz = 0; 3425ef31b24SMark F. Adams } 3435ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 3445ef31b24SMark F. Adams ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc ); 345afc97cdcSMark F. Adams if( isnewproc_idx != 0 ) { 3465ef31b24SMark F. Adams ierr = PetscFree( isnewproc_idx ); CHKERRQ(ierr); 3475ef31b24SMark F. Adams } 348e33ef3b1SMark F. Adams 34911e60469SMark F. Adams /* 35011e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 35111e60469SMark F. Adams */ 35211e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 35311e60469SMark F. Adams /* 35411e60469SMark F. Adams Determine how many elements are assigned to each processor 35511e60469SMark F. Adams */ 356fd3c6afaSMark F. Adams inpe = npe; 357fd3c6afaSMark F. Adams ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr); 35811e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 359038e3b61SMark F. Adams ncrs_new = counts[mype]/a_cbs; 36022063be5SMark F. Adams strideNew = ncrs_new*a_ndata_rows; 361b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 362b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 363b4fbaa2aSMark F. Adams #endif 36422063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 36511e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 366d3d6bff4SMark F. Adams ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 367470e26f8SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ 36811e60469SMark F. Adams /* 369d3d6bff4SMark F. Adams There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having 370d3d6bff4SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'a_cbs'. 37111e60469SMark F. Adams */ 37292a756f0SMark F. Adams ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 37311e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 374038e3b61SMark F. Adams for(ii=0,jj=0; ii<ncrs0 ; ii++) { 375d3d6bff4SMark F. Adams PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */ 376d3d6bff4SMark F. Adams for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk; 37711e60469SMark F. Adams } 378038e3b61SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 379d3d6bff4SMark F. Adams ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 38011e60469SMark F. Adams CHKERRQ(ierr); 38192a756f0SMark F. Adams ierr = PetscFree( tidx ); CHKERRQ(ierr); 38211e60469SMark F. Adams /* 38311e60469SMark F. Adams Create a vector to contain the original vertex information for each element 38411e60469SMark F. Adams */ 385d3d6bff4SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 386d3d6bff4SMark F. Adams for( jj=0; jj<a_ndata_cols ; jj++ ) { 387d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs0 ; ii++) { 388d3d6bff4SMark F. Adams for( kk=0; kk<a_ndata_rows ; kk++ ) { 389d3d6bff4SMark F. Adams PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj; 390676e1743SMark F. Adams PetscScalar tt = (PetscScalar)data[ix]; 391676e1743SMark F. Adams ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES ); CHKERRQ(ierr); 392d3d6bff4SMark F. Adams } 393038e3b61SMark F. Adams } 394eb07cef2SMark F. Adams } 395eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 396eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 39711e60469SMark F. Adams /* 39811e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 39911e60469SMark F. Adams to the correct processor 40011e60469SMark F. Adams */ 40111e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 40211e60469SMark F. Adams CHKERRQ(ierr); 40311e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 40411e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40511e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40611e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 40711e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 40811e60469SMark F. Adams /* 40911e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 41011e60469SMark F. Adams */ 411eb07cef2SMark F. Adams ierr = PetscFree( *a_coarse_data ); CHKERRQ(ierr); 412d3d6bff4SMark F. Adams ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data ); CHKERRQ(ierr); 413eb07cef2SMark F. Adams 41411e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 415eb07cef2SMark F. Adams data = *a_coarse_data; 416d3d6bff4SMark F. Adams for( jj=0; jj<a_ndata_cols ; jj++ ) { 417d3d6bff4SMark F. Adams for( ii=0 ; ii<ncrs_new ; ii++) { 418d3d6bff4SMark F. Adams for( kk=0; kk<a_ndata_rows ; kk++ ) { 419d3d6bff4SMark F. Adams PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj; 420d3d6bff4SMark F. Adams data[ix] = PetscRealPart(array[jx]); 421d3d6bff4SMark F. Adams array[jx] = 1.e300; 422d3d6bff4SMark F. Adams } 423038e3b61SMark F. Adams } 424038e3b61SMark F. Adams } 42511e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 42611e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 42711e60469SMark F. Adams /* 42811e60469SMark F. Adams Invert for MatGetSubMatrix 42911e60469SMark F. Adams */ 430038e3b61SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr); 43111e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 43211e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 43311e60469SMark F. Adams /* A_crs output */ 434038e3b61SMark F. Adams ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 43511e60469SMark F. Adams CHKERRQ(ierr); 436eb07cef2SMark F. Adams 437038e3b61SMark F. Adams ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 438e33ef3b1SMark F. Adams Cmat = *a_Amat_crs; /* output */ 439038e3b61SMark F. Adams ierr = MatSetBlockSize( Cmat, a_cbs ); CHKERRQ(ierr); 440eb07cef2SMark F. Adams 44111e60469SMark F. Adams /* prolongator */ 44211e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 44311e60469SMark F. Adams { 44411e60469SMark F. Adams IS findices; 44511e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 44696434bc1SMark F. Adams #ifdef USE_R 44796434bc1SMark F. Adams ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew ); 44896434bc1SMark F. Adams #else 44911e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 45096434bc1SMark F. Adams #endif 45111e60469SMark F. Adams CHKERRQ(ierr); 452d61a3a7dSMark F. Adams 45311e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 45411e60469SMark F. Adams } 455d61a3a7dSMark F. Adams 4563530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 4573530afc2SMark F. Adams *a_P_inout = Pnew; /* output */ 45892a756f0SMark F. Adams 45911e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 46092a756f0SMark F. Adams ierr = PetscFree( counts ); CHKERRQ(ierr); 46192a756f0SMark F. Adams ierr = PetscFree( ranks ); CHKERRQ(ierr); 462e33ef3b1SMark F. Adams } 4635b89ad90SMark F. Adams 4645b89ad90SMark F. Adams PetscFunctionReturn(0); 4655b89ad90SMark F. Adams } 4665b89ad90SMark F. Adams 4675b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4685b89ad90SMark F. Adams /* 4695b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4705b89ad90SMark F. Adams by setting data structures and options. 4715b89ad90SMark F. Adams 4725b89ad90SMark F. Adams Input Parameter: 4735b89ad90SMark F. Adams . pc - the preconditioner context 4745b89ad90SMark F. Adams 4755b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 4765b89ad90SMark F. Adams 4775b89ad90SMark F. Adams Notes: 4785b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 4795b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 4805b89ad90SMark F. Adams */ 4815b89ad90SMark F. Adams #undef __FUNCT__ 4825b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 483eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc ) 4845b89ad90SMark F. Adams { 4855b89ad90SMark F. Adams PetscErrorCode ierr; 486eb07cef2SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 4875b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 48858471d46SMark F. Adams PC_MG_Levels **mglevels = mg->levels; 489eb07cef2SMark F. Adams Mat Amat = a_pc->mat, Pmat = a_pc->pmat; 490d3d6bff4SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend; 491eb07cef2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_pc)->comm; 4923530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 493038e3b61SMark F. Adams PetscBool isOK; 494587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 495587fa25dSMark F. Adams PetscReal *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS]; 496737a81a9SMark F. Adams MatInfo info; 4975ef31b24SMark F. Adams 4985b89ad90SMark F. Adams PetscFunctionBegin; 49903a628feSMark F. Adams pc_gamg->m_count++; 50058471d46SMark F. Adams if( a_pc->setupcalled > 0 ) { 50103a628feSMark F. Adams /* just do Galerkin grids */ 50258471d46SMark F. Adams Mat B,dA,dB; 50358471d46SMark F. Adams 50458471d46SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 50558471d46SMark F. Adams ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 50658471d46SMark F. Adams 50758471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 50858471d46SMark F. Adams ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 50958471d46SMark F. Adams /* (re)set to get dirty flag */ 51058471d46SMark F. Adams ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 51158471d46SMark F. Adams ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr); 51258471d46SMark F. Adams 51358471d46SMark F. Adams for (level=pc_gamg->m_Nlevels-2; level>-1; level--) { 51458471d46SMark F. Adams ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 51503a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 51603a628feSMark F. Adams if( pc_gamg->m_count == 2 ) { 51703a628feSMark F. Adams ierr = MatDestroy( &B ); CHKERRQ(ierr); 51803a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 51903a628feSMark F. Adams mglevels[level]->A = B; 52003a628feSMark F. Adams } 52103a628feSMark F. Adams else { 52258471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 52303a628feSMark F. Adams } 52458471d46SMark F. Adams ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 52558471d46SMark F. Adams dB = B; 52658471d46SMark F. Adams /* setup KSP/PC */ 52758471d46SMark F. Adams ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr); 52858471d46SMark F. Adams } 52958471d46SMark F. Adams 53058471d46SMark F. Adams #define PRINT_MATS !PETSC_TRUE 53158471d46SMark F. Adams /* plot levels - A */ 53258471d46SMark F. Adams if( PRINT_MATS ) { 53358471d46SMark F. Adams for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){ 53458471d46SMark F. Adams PetscViewer viewer; 53558471d46SMark F. Adams char fname[32]; KSP smoother; Mat Tmat, TTm; 53658471d46SMark F. Adams ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 53758471d46SMark F. Adams ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr); 53803a628feSMark F. Adams sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level); 53958471d46SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 54058471d46SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 54158471d46SMark F. Adams ierr = MatView( Tmat, viewer ); CHKERRQ(ierr); 54258471d46SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 54358471d46SMark F. Adams } 54458471d46SMark F. Adams } 54558471d46SMark F. Adams 54603a628feSMark F. Adams a_pc->setupcalled = 2; 54703a628feSMark F. Adams 54858471d46SMark F. Adams PetscFunctionReturn(0); 549eb07cef2SMark F. Adams } 550baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 551baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 5525b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 5535b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 5543530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 555eb07cef2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 556eb07cef2SMark F. Adams nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 557eb07cef2SMark F. Adams 558038e3b61SMark F. Adams /* get data of not around */ 559038e3b61SMark F. Adams if( pc_gamg->m_data == 0 && nloc > 0 ) { 560038e3b61SMark F. Adams ierr = PCSetCoordinates_GAMG( a_pc, -1, 0 ); CHKERRQ( ierr ); 561eb07cef2SMark F. Adams } 562eb07cef2SMark F. Adams data = pc_gamg->m_data; 563038e3b61SMark F. Adams 564eb07cef2SMark F. Adams /* Get A_i and R_i */ 565737a81a9SMark F. Adams ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 566bed7c2b9SMark 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", 5672c047e2dSMark F. Adams mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols, 5682c047e2dSMark F. Adams (int)(info.nz_used/(PetscReal)N),npe); 5698f4b7eb5SMark F. Adams for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 570c20e4228SMark F. Adams level < (GAMG_MAXLEVELS-1) && (level==0 || M>2*pc_gamg->m_min_eq_proc); /* && (npe==1 || nactivepe>1); */ 5710205a208SMark F. Adams level++ ){ 5725b89ad90SMark F. Adams level1 = level + 1; 573b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 574b4fbaa2aSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 575b4fbaa2aSMark F. Adams #endif 576b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 577b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 578b4fbaa2aSMark F. Adams #endif 579470e26f8SMark F. Adams ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_method, 580c20e4228SMark F. Adams level, pc_gamg->m_threshold, &bs, &Parr[level1], &coarse_data, &isOK, 581bed7c2b9SMark F. Adams &emaxs[level], pc_gamg->m_verbose ); 5825b89ad90SMark F. Adams CHKERRQ(ierr); 583d3d6bff4SMark F. Adams ierr = PetscFree( data ); CHKERRQ( ierr ); 584b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 585b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 586b4fbaa2aSMark F. Adams #endif 587baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 588baa4e9faSMark F. Adams if( isOK ) { 589b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 590b4fbaa2aSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 591b4fbaa2aSMark F. Adams #endif 592486a8d0bSJed Brown ierr = PCGAMGPartitionLevel(a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs, 593c20e4228SMark F. Adams &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1], 594c20e4228SMark F. Adams pc_gamg->m_avoid_repart, pc_gamg->m_min_eq_proc ); 5953530afc2SMark F. Adams CHKERRQ(ierr); 596b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 597b4fbaa2aSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 598b4fbaa2aSMark F. Adams #endif 5993530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 600737a81a9SMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 601bed7c2b9SMark 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", 6022c047e2dSMark F. Adams mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols, 6032c047e2dSMark F. Adams (int)(info.nz_used/(PetscReal)N),nactivepe); 604e33ef3b1SMark F. Adams /* coarse grids with SA can have zero row/cols from singleton aggregates */ 60558471d46SMark F. Adams /* aggregation method should gaurrentee this does not happen! */ 606737a81a9SMark F. Adams 607bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) { 608785cba28SMark F. Adams Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; 609785cba28SMark F. Adams v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */ 610e33ef3b1SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); 611785cba28SMark F. Adams nloceq = Iend-Istart; 612e33ef3b1SMark F. Adams ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); 613e33ef3b1SMark F. Adams ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); 614e33ef3b1SMark F. Adams ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); 615785cba28SMark F. Adams for(kk=0;kk<nloceq;kk++){ 616e33ef3b1SMark F. Adams if(data_arr[kk]==0.0) { 617785cba28SMark F. Adams id = kk + Istart; 618e33ef3b1SMark F. Adams ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); 619e33ef3b1SMark F. Adams CHKERRQ(ierr); 620ecd57171SMark F. Adams PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); 621e33ef3b1SMark F. Adams } 622e33ef3b1SMark F. Adams } 623e33ef3b1SMark F. Adams ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); 624e33ef3b1SMark F. Adams ierr = VecDestroy( &diag ); CHKERRQ(ierr); 625e33ef3b1SMark F. Adams ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 626e33ef3b1SMark F. Adams ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 627e33ef3b1SMark F. Adams } 628486a8d0bSJed Brown } else { 629be544d3cSMark F. Adams coarse_data = 0; 630baa4e9faSMark F. Adams break; 631baa4e9faSMark F. Adams } 632eb07cef2SMark F. Adams data = coarse_data; 633b4fbaa2aSMark F. Adams 634b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 635b4fbaa2aSMark F. Adams ierr = PetscLogStagePop(); CHKERRQ( ierr ); 636b4fbaa2aSMark F. Adams #endif 6375b89ad90SMark F. Adams } 638be544d3cSMark F. Adams if( coarse_data ) { 639eb07cef2SMark F. Adams ierr = PetscFree( coarse_data ); CHKERRQ( ierr ); 640be544d3cSMark F. Adams } 641bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 6425b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 6435b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 6445b89ad90SMark F. Adams fine_level = level; 645eb07cef2SMark F. Adams ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 6465b89ad90SMark F. Adams 647fc4362bfSMark F. Adams /* set default smoothers */ 648587fa25dSMark F. Adams for ( lidx=1, level = pc_gamg->m_Nlevels-2; 649587fa25dSMark F. Adams lidx <= fine_level; 650587fa25dSMark F. Adams lidx++, level--) { 6515745e0f5SMark F. Adams PetscReal emax, emin; 6525b89ad90SMark F. Adams KSP smoother; PC subpc; 653587fa25dSMark F. Adams ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 6545b89ad90SMark F. Adams ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 655038e3b61SMark F. Adams if( emaxs[level] > 0.0 ) emax=emaxs[level]; 656587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 657587fa25dSMark F. Adams KSP eksp; Mat Lmat = Aarr[level]; 658fc4362bfSMark F. Adams Vec bb, xx; PC pc; 659038e3b61SMark F. Adams 6605745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 6615745e0f5SMark F. Adams ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 662fc4362bfSMark F. Adams { 663fc4362bfSMark F. Adams PetscRandom rctx; 664fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 665fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 666fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 667fc4362bfSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 6685b89ad90SMark F. Adams } 669fc4362bfSMark F. Adams ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 670bed7c2b9SMark F. Adams ierr = KSPAppendOptionsPrefix( eksp, "eigen_estimate_"); CHKERRQ(ierr); 671e33ef3b1SMark F. Adams ierr = KSPSetType( eksp, KSPCG ); CHKERRQ(ierr); 672bed7c2b9SMark F. Adams ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 673fc4362bfSMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 67458471d46SMark F. Adams ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 675fc4362bfSMark F. Adams ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 67658471d46SMark F. Adams ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as eigen estimates op. */ 677038e3b61SMark F. Adams ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 678fc4362bfSMark F. Adams CHKERRQ(ierr); 679fc4362bfSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 680fc4362bfSMark F. Adams 681fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 682fc4362bfSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 683fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 684fc4362bfSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 685fc4362bfSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 686fc4362bfSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 687bed7c2b9SMark F. Adams if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e PC=%s\n",__FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER); 688fc4362bfSMark F. Adams } 689038e3b61SMark F. Adams { 690038e3b61SMark F. Adams PetscInt N1, N0, tt; 691038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 692038e3b61SMark F. Adams ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 693785cba28SMark F. Adams emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */ 694038e3b61SMark F. Adams emax *= 1.05; 695038e3b61SMark F. Adams } 696038e3b61SMark F. Adams 69758471d46SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); 698fc4362bfSMark F. Adams ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 6990e1b4bd6SMark F. Adams /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */ 7005745e0f5SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 701e33ef3b1SMark F. Adams ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); 7025745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 7035745e0f5SMark F. Adams } 7045745e0f5SMark F. Adams { 705ecd57171SMark F. Adams /* coarse grid */ 706ecd57171SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 7075745e0f5SMark F. Adams Mat Lmat = Aarr[pc_gamg->m_Nlevels-1]; 708eb07cef2SMark F. Adams ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 70958471d46SMark F. Adams ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 7105745e0f5SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 711ecd57171SMark F. Adams ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 712ecd57171SMark F. Adams ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 713ecd57171SMark F. Adams ierr = PCSetUp( subpc ); CHKERRQ(ierr); 714ecd57171SMark F. Adams ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 715ecd57171SMark F. Adams assert(ii==1); 716ecd57171SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 717ecd57171SMark F. Adams ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 718fc4362bfSMark F. Adams } 719737a81a9SMark F. Adams 720fc4362bfSMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 721eb07cef2SMark F. Adams ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr); 7225b89ad90SMark F. Adams { 7235b89ad90SMark F. Adams PetscBool galerkin; 724eb07cef2SMark F. Adams ierr = PCMGGetGalerkin( a_pc, &galerkin); CHKERRQ(ierr); 7255b89ad90SMark F. Adams if(galerkin){ 7265b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 7275b89ad90SMark F. Adams } 7285b89ad90SMark F. Adams } 7295745e0f5SMark F. Adams 73058471d46SMark F. Adams /* plot levels - R/P */ 73158471d46SMark F. Adams if( PRINT_MATS ) { 73258471d46SMark F. Adams for (level=pc_gamg->m_Nlevels-1;level>0;level--){ 73358471d46SMark F. Adams PetscViewer viewer; 73458471d46SMark F. Adams char fname[32]; 73503a628feSMark F. Adams sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level); 73611e60469SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 7375b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 738038e3b61SMark F. Adams ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr); 7395b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 74003a628feSMark F. Adams sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level); 741e33ef3b1SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 742e33ef3b1SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 743e33ef3b1SMark F. Adams ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr); 744e33ef3b1SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 7455b89ad90SMark F. Adams } 74658471d46SMark F. Adams } 74758471d46SMark F. Adams 74858471d46SMark F. Adams /* set interpolation between the levels, clean up */ 74958471d46SMark F. Adams for (lidx=0,level=pc_gamg->m_Nlevels-1; 75058471d46SMark F. Adams lidx<fine_level; 75158471d46SMark F. Adams lidx++, level--){ 75258471d46SMark F. Adams ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr); 753587fa25dSMark F. Adams ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 754587fa25dSMark F. Adams ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 7555b89ad90SMark F. Adams } 7565745e0f5SMark F. Adams 7575b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 758eb07cef2SMark F. Adams a_pc->setupcalled = 0; 759eb07cef2SMark F. Adams ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr ); 76003a628feSMark F. Adams a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ 76158471d46SMark F. Adams 76258471d46SMark F. Adams { 76358471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 76458471d46SMark F. Adams ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 76558471d46SMark F. Adams ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 76658471d46SMark F. Adams ierr = KSPSetUp( smoother ); CHKERRQ(ierr); 76758471d46SMark F. Adams } 768e33ef3b1SMark F. Adams 7695b89ad90SMark F. Adams PetscFunctionReturn(0); 7705b89ad90SMark F. Adams } 7715b89ad90SMark F. Adams 772eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 7735b89ad90SMark F. Adams /* 7745b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 7755b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 7765b89ad90SMark F. Adams 7775b89ad90SMark F. Adams Input Parameter: 7785b89ad90SMark F. Adams . pc - the preconditioner context 7795b89ad90SMark F. Adams 7805b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 7815b89ad90SMark F. Adams */ 7825b89ad90SMark F. Adams #undef __FUNCT__ 7835b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 7845b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 7855b89ad90SMark F. Adams { 7865b89ad90SMark F. Adams PetscErrorCode ierr; 7875b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7885b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 7895b89ad90SMark F. Adams 7905b89ad90SMark F. Adams PetscFunctionBegin; 7915b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 7925b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 7935b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 7945b89ad90SMark F. Adams PetscFunctionReturn(0); 7955b89ad90SMark F. Adams } 7965b89ad90SMark F. Adams 797676e1743SMark F. Adams 798676e1743SMark F. Adams #undef __FUNCT__ 799676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 800676e1743SMark F. Adams /*@ 801676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 802676e1743SMark F. Adams processor reduction. 803676e1743SMark F. Adams 804676e1743SMark F. Adams Not Collective on PC 805676e1743SMark F. Adams 806676e1743SMark F. Adams Input Parameters: 807676e1743SMark F. Adams . pc - the preconditioner context 808676e1743SMark F. Adams 809676e1743SMark F. Adams 810676e1743SMark F. Adams Options Database Key: 811676e1743SMark F. Adams . -pc_gamg_process_eq_limit 812676e1743SMark F. Adams 813676e1743SMark F. Adams Level: intermediate 814676e1743SMark F. Adams 815676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 816676e1743SMark F. Adams 817676e1743SMark F. Adams .seealso: () 818676e1743SMark F. Adams @*/ 819676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 820676e1743SMark F. Adams { 821676e1743SMark F. Adams PetscErrorCode ierr; 822676e1743SMark F. Adams 823676e1743SMark F. Adams PetscFunctionBegin; 824676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 825676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 826676e1743SMark F. Adams PetscFunctionReturn(0); 827676e1743SMark F. Adams } 828676e1743SMark F. Adams 829676e1743SMark F. Adams EXTERN_C_BEGIN 830676e1743SMark F. Adams #undef __FUNCT__ 831676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 832676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 833676e1743SMark F. Adams { 834c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 835c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 836676e1743SMark F. Adams 837676e1743SMark F. Adams PetscFunctionBegin; 838c20e4228SMark F. Adams if(n>0) pc_gamg->m_min_eq_proc = n; 839676e1743SMark F. Adams PetscFunctionReturn(0); 840676e1743SMark F. Adams } 841676e1743SMark F. Adams EXTERN_C_END 842676e1743SMark F. Adams 843676e1743SMark F. Adams #undef __FUNCT__ 844676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning" 845676e1743SMark F. Adams /*@ 846676e1743SMark F. Adams PCGAMGAvoidRepartitioning - Do not repartition the coarse grids 847676e1743SMark F. Adams 848676e1743SMark F. Adams Collective on PC 849676e1743SMark F. Adams 850676e1743SMark F. Adams Input Parameters: 851676e1743SMark F. Adams . pc - the preconditioner context 852676e1743SMark F. Adams 853676e1743SMark F. Adams 854676e1743SMark F. Adams Options Database Key: 855676e1743SMark F. Adams . -pc_gamg_avoid_repartitioning 856676e1743SMark F. Adams 857676e1743SMark F. Adams Level: intermediate 858676e1743SMark F. Adams 859676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 860676e1743SMark F. Adams 861676e1743SMark F. Adams .seealso: () 862676e1743SMark F. Adams @*/ 863676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n) 864676e1743SMark F. Adams { 865676e1743SMark F. Adams PetscErrorCode ierr; 866676e1743SMark F. Adams 867676e1743SMark F. Adams PetscFunctionBegin; 868676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 869676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 870676e1743SMark F. Adams PetscFunctionReturn(0); 871676e1743SMark F. Adams } 872676e1743SMark F. Adams 873676e1743SMark F. Adams EXTERN_C_BEGIN 874676e1743SMark F. Adams #undef __FUNCT__ 875676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG" 876676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n) 877676e1743SMark F. Adams { 878c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 879c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 880676e1743SMark F. Adams 881676e1743SMark F. Adams PetscFunctionBegin; 882c20e4228SMark F. Adams pc_gamg->m_avoid_repart = n; 883676e1743SMark F. Adams PetscFunctionReturn(0); 884676e1743SMark F. Adams } 885676e1743SMark F. Adams EXTERN_C_END 886676e1743SMark F. Adams 887676e1743SMark F. Adams #undef __FUNCT__ 8883542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 8893542efc5SMark F. Adams /*@ 8903542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 8913542efc5SMark F. Adams 8923542efc5SMark F. Adams Not collective on PC 8933542efc5SMark F. Adams 8943542efc5SMark F. Adams Input Parameters: 8953542efc5SMark F. Adams . pc - the preconditioner context 8963542efc5SMark F. Adams 8973542efc5SMark F. Adams 8983542efc5SMark F. Adams Options Database Key: 8993542efc5SMark F. Adams . -pc_gamg_threshold 9003542efc5SMark F. Adams 9013542efc5SMark F. Adams Level: intermediate 9023542efc5SMark F. Adams 9033542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 9043542efc5SMark F. Adams 9053542efc5SMark F. Adams .seealso: () 9063542efc5SMark F. Adams @*/ 9073542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 9083542efc5SMark F. Adams { 9093542efc5SMark F. Adams PetscErrorCode ierr; 9103542efc5SMark F. Adams 9113542efc5SMark F. Adams PetscFunctionBegin; 9123542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9133542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 9143542efc5SMark F. Adams PetscFunctionReturn(0); 9153542efc5SMark F. Adams } 9163542efc5SMark F. Adams 9173542efc5SMark F. Adams EXTERN_C_BEGIN 9183542efc5SMark F. Adams #undef __FUNCT__ 9193542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 9203542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 9213542efc5SMark F. Adams { 922c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 923c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9243542efc5SMark F. Adams 9253542efc5SMark F. Adams PetscFunctionBegin; 926c20e4228SMark F. Adams pc_gamg->m_threshold = n; 9273542efc5SMark F. Adams PetscFunctionReturn(0); 9283542efc5SMark F. Adams } 9293542efc5SMark F. Adams EXTERN_C_END 9303542efc5SMark F. Adams 9313542efc5SMark F. Adams #undef __FUNCT__ 932676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType" 933676e1743SMark F. Adams /*@ 934676e1743SMark F. Adams PCGAMGSetSolverType - Set solution method. 935676e1743SMark F. Adams 936676e1743SMark F. Adams Collective on PC 937676e1743SMark F. Adams 938676e1743SMark F. Adams Input Parameters: 939676e1743SMark F. Adams . pc - the preconditioner context 940676e1743SMark F. Adams 941676e1743SMark F. Adams 942676e1743SMark F. Adams Options Database Key: 9433542efc5SMark F. Adams . -pc_gamg_type 944676e1743SMark F. Adams 945676e1743SMark F. Adams Level: intermediate 946676e1743SMark F. Adams 947676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 948676e1743SMark F. Adams 949676e1743SMark F. Adams .seealso: () 950676e1743SMark F. Adams @*/ 951676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz ) 952676e1743SMark F. Adams { 953676e1743SMark F. Adams PetscErrorCode ierr; 954676e1743SMark F. Adams 955676e1743SMark F. Adams PetscFunctionBegin; 956676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 957676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz)); 958676e1743SMark F. Adams CHKERRQ(ierr); 959676e1743SMark F. Adams PetscFunctionReturn(0); 960676e1743SMark F. Adams } 961676e1743SMark F. Adams 962676e1743SMark F. Adams EXTERN_C_BEGIN 963676e1743SMark F. Adams #undef __FUNCT__ 964676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG" 965676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz ) 966676e1743SMark F. Adams { 967676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 968676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 969676e1743SMark F. Adams 970676e1743SMark F. Adams PetscFunctionBegin; 971676e1743SMark F. Adams if(sz < 64) strcpy(pc_gamg->m_type,str); 972676e1743SMark F. Adams PetscFunctionReturn(0); 973676e1743SMark F. Adams } 974676e1743SMark F. Adams EXTERN_C_END 975676e1743SMark F. Adams 9765b89ad90SMark F. Adams #undef __FUNCT__ 9775b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 9785b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 9795b89ad90SMark F. Adams { 980676e1743SMark F. Adams PetscErrorCode ierr; 981676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 982676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 983676e1743SMark F. Adams PetscBool flag; 9845b89ad90SMark F. Adams 9855b89ad90SMark F. Adams PetscFunctionBegin; 986676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 987676e1743SMark F. Adams { 988676e1743SMark F. Adams ierr = PetscOptionsString("-pc_gamg_type", 989470e26f8SMark F. Adams "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)", 990676e1743SMark F. Adams "PCGAMGSetSolverType", 991676e1743SMark F. Adams pc_gamg->m_type, 992676e1743SMark F. Adams pc_gamg->m_type, 993676e1743SMark F. Adams 64, 994676e1743SMark F. Adams &flag ); 995676e1743SMark F. Adams CHKERRQ(ierr); 996*d8c859f0SMark F. Adams if( flag && pc_gamg->m_data != 0 ) { 997*d8c859f0SMark F. Adams if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) || 998*d8c859f0SMark F. Adams (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) || 999*d8c859f0SMark F. Adams (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) { 1000*d8c859f0SMark F. Adams SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)"); 1001*d8c859f0SMark F. Adams } 1002*d8c859f0SMark F. Adams } 1003bed7c2b9SMark F. Adams 1004bed7c2b9SMark F. Adams /* -pc_gamg_verbose */ 1005bed7c2b9SMark 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); 10062a44bfbaSMark F. Adams 1007470e26f8SMark F. Adams if (flag && strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2; 1008470e26f8SMark F. Adams else if (flag && strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1; 1009470e26f8SMark F. Adams else pc_gamg->m_method = 0; 10102a44bfbaSMark F. Adams 1011c20e4228SMark F. Adams /* -pc_gamg_avoid_repartitioning */ 1012c20e4228SMark F. Adams pc_gamg->m_avoid_repart = PETSC_FALSE; 1013676e1743SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning", 1014676e1743SMark F. Adams "Do not repartion coarse grids (false)", 1015676e1743SMark F. Adams "PCGAMGAvoidRepartitioning", 1016c20e4228SMark F. Adams pc_gamg->m_avoid_repart, 1017c20e4228SMark F. Adams &pc_gamg->m_avoid_repart, 1018676e1743SMark F. Adams &flag); 1019676e1743SMark F. Adams CHKERRQ(ierr); 1020676e1743SMark F. Adams 1021c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 1022c20e4228SMark F. Adams pc_gamg->m_min_eq_proc = 600; 1023676e1743SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1024676e1743SMark F. Adams "Limit (goal) on number of equations per process on coarse grids", 1025676e1743SMark F. Adams "PCGAMGSetProcEqLim", 1026c20e4228SMark F. Adams pc_gamg->m_min_eq_proc, 1027c20e4228SMark F. Adams &pc_gamg->m_min_eq_proc, 1028676e1743SMark F. Adams &flag ); 1029676e1743SMark F. Adams CHKERRQ(ierr); 10303542efc5SMark F. Adams 1031c20e4228SMark F. Adams /* -pc_gamg_threshold */ 1032c20e4228SMark F. Adams pc_gamg->m_threshold = 0.05; 10333542efc5SMark F. Adams ierr = PetscOptionsReal("-pc_gamg_threshold", 10343542efc5SMark F. Adams "Relative threshold to use for dropping edges in aggregation graph", 10353542efc5SMark F. Adams "PCGAMGSetThreshold", 1036c20e4228SMark F. Adams pc_gamg->m_threshold, 1037c20e4228SMark F. Adams &pc_gamg->m_threshold, 10383542efc5SMark F. Adams &flag ); 10393542efc5SMark F. Adams CHKERRQ(ierr); 1040676e1743SMark F. Adams } 1041676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1042676e1743SMark F. Adams 10435b89ad90SMark F. Adams PetscFunctionReturn(0); 10445b89ad90SMark F. Adams } 10455b89ad90SMark F. Adams 10465b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 10475b89ad90SMark F. Adams /*MC 1048dc76db92SMark F. Adams PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two 1049dc76db92SMark F. Adams AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each 1050dc76db92SMark F. Adams vertex, and 2) smoothed aggregation. Smoothed aggregation (SA) can work without coordinates but it 1051dc76db92SMark F. Adams will generate some common non-trivial null spaces if coordinates are provided. The input fine grid matrix 1052dc76db92SMark F. Adams must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly. 1053dc76db92SMark F. Adams SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational 1054dc76db92SMark F. Adams modes, when coordinates are provide in 2D and 3D. 10555b89ad90SMark F. Adams 1056280d9858SJed Brown Options Database Keys: 10575b89ad90SMark F. Adams Multigrid options(inherited) 1058280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1059280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1060280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1061280d9858SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 10625b89ad90SMark F. Adams 10635b89ad90SMark F. Adams Level: intermediate 1064280d9858SJed Brown 10655b89ad90SMark F. Adams Concepts: multigrid 10665b89ad90SMark F. Adams 10675b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1068280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 10695b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 10705b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 10715b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 10725b89ad90SMark F. Adams M*/ 10735b89ad90SMark F. Adams 10745b89ad90SMark F. Adams EXTERN_C_BEGIN 10755b89ad90SMark F. Adams #undef __FUNCT__ 10765b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 10775b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 10785b89ad90SMark F. Adams { 10795b89ad90SMark F. Adams PetscErrorCode ierr; 10805b89ad90SMark F. Adams PC_GAMG *pc_gamg; 10815b89ad90SMark F. Adams PC_MG *mg; 10825ef31b24SMark F. Adams PetscClassId cookie; 10835ee9c036SSatish Balay #if defined PETSC_USE_LOG 10845ee9c036SSatish Balay static int count = 0; 10855ee9c036SSatish Balay #endif 10865b89ad90SMark F. Adams 10875b89ad90SMark F. Adams PetscFunctionBegin; 10885b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 10895b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 10905b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 10915b89ad90SMark F. Adams 10925b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 10935b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 109403a628feSMark F. Adams pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0; 10955b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 10965b89ad90SMark F. Adams mg->innerctx = pc_gamg; 10975b89ad90SMark F. Adams 10985b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 10995b89ad90SMark F. Adams 11005b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 11015b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 11025b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 11035b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 11045b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 11055b89ad90SMark F. Adams 11065b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 11075b89ad90SMark F. Adams "PCSetCoordinates_C", 11085b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 1109c97e1df0SMark F. Adams PCSetCoordinates_GAMG); 1110c97e1df0SMark F. Adams CHKERRQ(ierr); 1111c97e1df0SMark F. Adams 1112676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1113676e1743SMark F. Adams "PCGAMGSetProcEqLim_C", 1114676e1743SMark F. Adams "PCGAMGSetProcEqLim_GAMG", 1115676e1743SMark F. Adams PCGAMGSetProcEqLim_GAMG); 1116676e1743SMark F. Adams CHKERRQ(ierr); 1117676e1743SMark F. Adams 1118676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1119676e1743SMark F. Adams "PCGAMGAvoidRepartitioning_C", 1120676e1743SMark F. Adams "PCGAMGAvoidRepartitioning_GAMG", 1121676e1743SMark F. Adams PCGAMGAvoidRepartitioning_GAMG); 1122676e1743SMark F. Adams CHKERRQ(ierr); 1123676e1743SMark F. Adams 1124676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 11253542efc5SMark F. Adams "PCGAMGSetThreshold_C", 11263542efc5SMark F. Adams "PCGAMGSetThreshold_GAMG", 11273542efc5SMark F. Adams PCGAMGSetThreshold_GAMG); 11283542efc5SMark F. Adams CHKERRQ(ierr); 11293542efc5SMark F. Adams 11303542efc5SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1131676e1743SMark F. Adams "PCGAMGSetSolverType_C", 1132676e1743SMark F. Adams "PCGAMGSetSolverType_GAMG", 1133676e1743SMark F. Adams PCGAMGSetSolverType_GAMG); 1134676e1743SMark F. Adams CHKERRQ(ierr); 1135c97e1df0SMark F. Adams 1136b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG 1137785cba28SMark F. Adams if( count++ == 0 ) { 11385ef31b24SMark F. Adams PetscClassIdRegister("GAMG Setup",&cookie); 1139b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); 11402a44bfbaSMark F. Adams PetscLogEventRegister(" Graph", cookie, &gamg_setup_events[GRAPH]); 11412a44bfbaSMark F. Adams PetscLogEventRegister(" G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); 11422a44bfbaSMark F. Adams PetscLogEventRegister(" G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); 11432a44bfbaSMark F. Adams PetscLogEventRegister(" G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); 1144b4fbaa2aSMark F. Adams PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); 1145b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); 1146b4fbaa2aSMark F. Adams PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); 1147b4fbaa2aSMark F. Adams PetscLogEventRegister(" search&set", cookie, &gamg_setup_events[FIND_V]); 1148b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: init", cookie, &gamg_setup_events[SET7]); 1149b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); */ 1150b4fbaa2aSMark F. Adams PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); 1151b4fbaa2aSMark F. Adams PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); 1152b4fbaa2aSMark F. Adams PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_events[SET12]); 1153b4fbaa2aSMark F. Adams /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ 1154b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ 1155b4fbaa2aSMark F. Adams /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ 11565ef31b24SMark F. Adams 1157b4fbaa2aSMark F. Adams /* create timer stages */ 1158b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1159b4fbaa2aSMark F. Adams { 1160b4fbaa2aSMark F. Adams char str[32]; 1161b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d (finest)",0); 1162b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[0]); 1163b4fbaa2aSMark F. Adams PetscInt lidx; 1164b4fbaa2aSMark F. Adams for (lidx=1;lidx<9;lidx++){ 1165b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1166b4fbaa2aSMark F. Adams PetscLogStageRegister(str, &gamg_stages[lidx]); 1167b4fbaa2aSMark F. Adams } 1168b4fbaa2aSMark F. Adams } 1169b4fbaa2aSMark F. Adams #endif 1170b4fbaa2aSMark F. Adams } 1171b4fbaa2aSMark F. Adams #endif 11725b89ad90SMark F. Adams PetscFunctionReturn(0); 11735b89ad90SMark F. Adams } 11745b89ad90SMark F. Adams EXTERN_C_END 1175