15b89ad90SMark F. Adams /* 25b89ad90SMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4*5ef31b24SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> 55b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */ 65b89ad90SMark F. Adams typedef struct gamg_TAG { 75b89ad90SMark F. Adams PetscInt m_dim; 85b89ad90SMark F. Adams PetscInt m_Nlevels; 95b89ad90SMark F. Adams PetscInt m_data_sz; 105b89ad90SMark F. Adams PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 115b89ad90SMark F. Adams } PC_GAMG; 125b89ad90SMark F. Adams 13*5ef31b24SMark F. Adams #define TOP_GRID_LIM 1000 14*5ef31b24SMark F. Adams 155b89ad90SMark F. Adams /* -----------------------------------------------------------------------------*/ 165b89ad90SMark F. Adams #undef __FUNCT__ 175b89ad90SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 185b89ad90SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 195b89ad90SMark F. Adams { 205b89ad90SMark F. Adams PetscErrorCode ierr; 215b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 225b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 235b89ad90SMark F. Adams 245b89ad90SMark F. Adams PetscFunctionBegin; 255b89ad90SMark F. Adams ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 265b89ad90SMark F. Adams PetscFunctionReturn(0); 275b89ad90SMark F. Adams } 285b89ad90SMark F. Adams 295b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 305b89ad90SMark F. Adams /* 315b89ad90SMark F. Adams PCSetCoordinates_GAMG 325b89ad90SMark F. Adams 335b89ad90SMark F. Adams Input Parameter: 345b89ad90SMark F. Adams . pc - the preconditioner context 355b89ad90SMark F. Adams */ 36a92563c5SMark F. Adams EXTERN_C_BEGIN 375b89ad90SMark F. Adams #undef __FUNCT__ 385b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 396c237d78SBarry Smith PetscErrorCode PCSetCoordinates_GAMG(PC pc, const int ndm,PetscReal *coords ) 405b89ad90SMark F. Adams { 415b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 425b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 436c237d78SBarry Smith PetscErrorCode ierr; 446c237d78SBarry Smith PetscInt bs, my0, tt; 456c237d78SBarry Smith Mat mat = pc->pmat; 466c237d78SBarry Smith PetscInt arrsz; 475b89ad90SMark F. Adams 485b89ad90SMark F. Adams PetscFunctionBegin; 495b89ad90SMark F. Adams ierr = MatGetBlockSize( mat, &bs ); CHKERRQ( ierr ); 505b89ad90SMark F. Adams ierr = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr); 516c237d78SBarry Smith arrsz = (tt-my0)/bs*ndm; 525b89ad90SMark F. Adams 535b89ad90SMark F. Adams // put coordinates 546c237d78SBarry Smith if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 555b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 566c237d78SBarry Smith ierr = PetscMalloc(arrsz*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 575b89ad90SMark F. Adams } 585b89ad90SMark F. Adams 595b89ad90SMark F. Adams /* copy data in */ 605b89ad90SMark F. Adams for(tt=0;tt<arrsz;tt++){ 615b89ad90SMark F. Adams pc_gamg->m_data[tt] = coords[tt]; 625b89ad90SMark F. Adams } 635b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 645b89ad90SMark F. Adams pc_gamg->m_dim = ndm; 655b89ad90SMark F. Adams 665b89ad90SMark F. Adams PetscFunctionReturn(0); 675b89ad90SMark F. Adams } 68a92563c5SMark F. Adams EXTERN_C_END 695b89ad90SMark F. Adams 705b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 715b89ad90SMark F. Adams /* 7211e60469SMark F. Adams partitionLevel 735b89ad90SMark F. Adams 745b89ad90SMark F. Adams Input Parameter: 753530afc2SMark F. Adams . a_Amat_fine - matrix on this fine (k) level 7611e60469SMark F. Adams . a_dime - 2 or 3 773530afc2SMark F. Adams In/Output Parameter: 783530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 793530afc2SMark F. Adams . a_coarse_crds - coordinates that need to be moved 803530afc2SMark F. Adams . a_active_proc - number of active procs 8111e60469SMark F. Adams Output Parameter: 823530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 835b89ad90SMark F. Adams */ 845cb416c2SMark F. Adams 855b89ad90SMark F. Adams #undef __FUNCT__ 8611e60469SMark F. Adams #define __FUNCT__ "partitionLevel" 873530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine, 8811e60469SMark F. Adams PetscInt a_dim, 893530afc2SMark F. Adams Mat *a_P_inout, 903530afc2SMark F. Adams PetscReal **a_coarse_crds, 913530afc2SMark F. Adams PetscMPIInt *a_active_proc, 923530afc2SMark F. Adams Mat *a_Amat_crs 9311e60469SMark F. Adams ) 945b89ad90SMark F. Adams { 955b89ad90SMark F. Adams PetscErrorCode ierr; 963530afc2SMark F. Adams Mat Amat, Pnew, Pold = *a_P_inout; 9711e60469SMark F. Adams IS new_indices,isnum; 983530afc2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 99*5ef31b24SMark F. Adams PetscMPIInt nactive_procs,mype,npe; 10011e60469SMark F. Adams PetscInt Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,bs=1; /* bs ??? */ 1015b89ad90SMark F. Adams 1025b89ad90SMark F. Adams PetscFunctionBegin; 10311e60469SMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 10411e60469SMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 10511e60469SMark F. Adams /* RAP */ 1063530afc2SMark F. Adams ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Amat ); CHKERRQ(ierr); 1073530afc2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart0, &Iend0 ); CHKERRQ(ierr); /* x2 size of 'a_coarse_crds' */ 10811e60469SMark F. Adams ncrs0 = (Iend0 - Istart0)/bs; 109acadaa72SMark F. Adams 11011e60469SMark F. Adams /* Repartition Amat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 11111e60469SMark F. Adams { 1123530afc2SMark F. Adams PetscInt neq,N,counts[npe]; 11311e60469SMark F. Adams IS isnewproc; 114*5ef31b24SMark F. Adams PetscMPIInt new_npe,targ_npe; 1153530afc2SMark F. Adams 1163530afc2SMark F. Adams ierr = MatGetSize( Amat, &neq, &N );CHKERRQ(ierr); 1173530afc2SMark F. Adams #define MIN_EQ_PROC 100 118*5ef31b24SMark F. Adams nactive_procs = *a_active_proc; 1193530afc2SMark F. Adams targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */ 120*5ef31b24SMark F. Adams if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */ 121*5ef31b24SMark F. Adams else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */ 1223530afc2SMark F. Adams else { 123*5ef31b24SMark F. Adams PetscMPIInt factstart,fact; 1243530afc2SMark F. Adams new_npe = -9999; 125*5ef31b24SMark F. Adams factstart = nactive_procs; 1263530afc2SMark F. Adams for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */ 127*5ef31b24SMark F. Adams if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) { 128*5ef31b24SMark F. Adams new_npe = nactive_procs/fact; 1293530afc2SMark F. Adams } 1303530afc2SMark F. Adams } 1313530afc2SMark F. Adams assert(new_npe != -9999); 1323530afc2SMark F. Adams } 1333530afc2SMark F. Adams *a_active_proc = new_npe; /* output for next time */ 134*5ef31b24SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t%s num active procs = %d (N loc = %d)\n",__FUNCT__,new_npe,ncrs0); 13511e60469SMark F. Adams { /* partition: get 'isnewproc' */ 13611e60469SMark F. Adams MatPartitioning mpart; 137*5ef31b24SMark F. Adams Mat adj; 138*5ef31b24SMark F. Adams const PetscInt *is_idx; 139*5ef31b24SMark F. Adams PetscInt is_sz,kk,jj,old_fact=(npe/nactive_procs),new_fact=(npe/new_npe),*isnewproc_idx; 140c9a0b8beSMark F. Adams /* create sub communicator */ 141c9a0b8beSMark F. Adams MPI_Comm cm,new_comm; 142*5ef31b24SMark F. Adams int membershipKey = mype % old_fact; 143c9a0b8beSMark F. Adams ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr); 144c9a0b8beSMark F. Adams ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr); 145c9a0b8beSMark F. Adams ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 146c9a0b8beSMark F. Adams 147*5ef31b24SMark F. Adams /* MatPartitioningApply call MatConvert, which is collective */ 148*5ef31b24SMark F. Adams ierr = MatConvert(Amat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);CHKERRQ(ierr); 149*5ef31b24SMark F. Adams if( membershipKey == 0 ) { 150*5ef31b24SMark F. Adams /* hack to fix global data that pmetis.c uses in 'adj' */ 151*5ef31b24SMark F. Adams for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) { 152*5ef31b24SMark F. Adams adj->rmap->range[jj] = adj->rmap->range[kk]; 153*5ef31b24SMark F. Adams } 154*5ef31b24SMark F. Adams ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr); 155*5ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 15611e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 157*5ef31b24SMark F. Adams ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr); 15811e60469SMark F. Adams ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 15911e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 160*5ef31b24SMark F. Adams /* collect IS info */ 161*5ef31b24SMark F. Adams ierr = ISGetLocalSize( isnewproc, &is_sz ); CHKERRQ(ierr); 162*5ef31b24SMark F. Adams ierr = PetscMalloc( is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr); 163*5ef31b24SMark F. Adams ierr = ISGetIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 164*5ef31b24SMark F. Adams for(kk=0;kk<is_sz;kk++){ 165c9a0b8beSMark F. Adams /* spread partitioning across machine - probably the right thing to do but machine specific */ 166*5ef31b24SMark F. Adams isnewproc_idx[kk] = is_idx[kk] * new_fact; 167*5ef31b24SMark F. Adams } 168*5ef31b24SMark F. Adams ierr = ISRestoreIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 169*5ef31b24SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 170*5ef31b24SMark F. Adams } 171*5ef31b24SMark F. Adams else { 172*5ef31b24SMark F. Adams isnewproc_idx = 0; 173*5ef31b24SMark F. Adams is_sz = 0; 174*5ef31b24SMark F. Adams } 175*5ef31b24SMark F. Adams ierr = MatDestroy( &adj ); CHKERRQ(ierr); 176*5ef31b24SMark F. Adams ierr = MPI_Comm_free( &new_comm ); CHKERRQ(ierr); 177*5ef31b24SMark F. Adams ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc ); 178*5ef31b24SMark F. Adams if( membershipKey == 0 ) { 179*5ef31b24SMark F. Adams ierr = PetscFree( isnewproc_idx ); CHKERRQ(ierr); 180*5ef31b24SMark F. Adams } 18111e60469SMark F. Adams } 182c9a0b8beSMark F. Adams 18311e60469SMark F. Adams /* 18411e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 18511e60469SMark F. Adams */ 18611e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 18711e60469SMark F. Adams /* 18811e60469SMark F. Adams Determine how many elements are assigned to each processor 18911e60469SMark F. Adams */ 19011e60469SMark F. Adams ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr); 19111e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 19211e60469SMark F. Adams ncrs_new = counts[mype]; 19311e60469SMark F. Adams } 194*5ef31b24SMark F. Adams 19511e60469SMark F. Adams { /* Create a vector to contain the newly ordered element information */ 19611e60469SMark F. Adams const PetscInt *idx; 19711e60469SMark F. Adams PetscInt i,j,k; 19811e60469SMark F. Adams IS isscat; 19911e60469SMark F. Adams PetscScalar *array; 20011e60469SMark F. Adams Vec src_crd, dest_crd; 2013530afc2SMark F. Adams PetscReal *coords = *a_coarse_crds; 20211e60469SMark F. Adams VecScatter vecscat; 203acadaa72SMark F. Adams 20411e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 20511e60469SMark F. Adams ierr = VecSetSizes( dest_crd, a_dim*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 20611e60469SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/ 20711e60469SMark F. Adams /* 20811e60469SMark F. Adams There are three data items per cell (element), the integer vertex numbers of its three 20911e60469SMark F. Adams coordinates (we convert to double to use the scatter) (one can think 21011e60469SMark F. Adams of the vectors of having a block size of 3, then there is one index in idx[] for each element) 21111e60469SMark F. Adams */ 21211e60469SMark F. Adams { 21311e60469SMark F. Adams PetscInt tidx[ncrs0*a_dim]; 21411e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 21511e60469SMark F. Adams for (i=0,j=0; i<ncrs0; i++) { 21611e60469SMark F. Adams for (k=0; k<a_dim; k++) tidx[j++] = idx[i]*a_dim + k; 21711e60469SMark F. Adams } 21811e60469SMark F. Adams ierr = ISCreateGeneral( wcomm, a_dim*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 21911e60469SMark F. Adams CHKERRQ(ierr); 22011e60469SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 22111e60469SMark F. Adams } 22211e60469SMark F. Adams /* 22311e60469SMark F. Adams Create a vector to contain the original vertex information for each element 22411e60469SMark F. Adams */ 22511e60469SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, a_dim*ncrs0, &src_crd ); CHKERRQ(ierr); 22611e60469SMark F. Adams ierr = VecGetArray( src_crd, &array ); CHKERRQ(ierr); 22711e60469SMark F. Adams for (i=0; i<a_dim*ncrs0; i++) array[i] = coords[i]; 22811e60469SMark F. Adams ierr = VecRestoreArray( src_crd, &array ); CHKERRQ(ierr); 22911e60469SMark F. Adams /* 23011e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 23111e60469SMark F. Adams to the correct processor 23211e60469SMark F. Adams */ 23311e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 23411e60469SMark F. Adams CHKERRQ(ierr); 23511e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 23611e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 23711e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 23811e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 23911e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 24011e60469SMark F. Adams /* 24111e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 24211e60469SMark F. Adams */ 2433530afc2SMark F. Adams ierr = PetscFree( *a_coarse_crds ); CHKERRQ(ierr); 2443530afc2SMark F. Adams ierr = PetscMalloc( a_dim*ncrs_new*sizeof(PetscReal), a_coarse_crds ); CHKERRQ(ierr); 2453530afc2SMark F. Adams coords = *a_coarse_crds; /* convient */ 24611e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 24711e60469SMark F. Adams for (i=0; i<a_dim*ncrs_new; i++) coords[i] = PetscRealPart(array[i]); 24811e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 24911e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 25011e60469SMark F. Adams } 25111e60469SMark F. Adams /* 25211e60469SMark F. Adams Invert for MatGetSubMatrix 25311e60469SMark F. Adams */ 25411e60469SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new, &new_indices ); CHKERRQ(ierr); 25511e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 25611e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 25711e60469SMark F. Adams /* A_crs output */ 2583530afc2SMark F. Adams ierr = MatGetSubMatrix( Amat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 25911e60469SMark F. Adams CHKERRQ(ierr); 26011e60469SMark F. Adams ierr = MatDestroy( &Amat ); CHKERRQ(ierr); 2613530afc2SMark F. Adams Amat = *a_Amat_crs; 26211e60469SMark F. Adams /* prolongator */ 26311e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 26411e60469SMark F. Adams { 26511e60469SMark F. Adams IS findices; 26611e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 26711e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 26811e60469SMark F. Adams CHKERRQ(ierr); 26911e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 27011e60469SMark F. Adams } 2713530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 2723530afc2SMark F. Adams *a_P_inout = Pnew; /* output */ 273acadaa72SMark F. Adams 27411e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 2755b89ad90SMark F. Adams 2765b89ad90SMark F. Adams PetscFunctionReturn(0); 2775b89ad90SMark F. Adams } 2785b89ad90SMark F. Adams 279*5ef31b24SMark F. Adams #define GAMG_MAXLEVELS 20 280*5ef31b24SMark F. Adams #if defined(PETSC_USE_LOG) 281*5ef31b24SMark F. Adams PetscLogStage gamg_stages[20]; 282*5ef31b24SMark F. Adams #endif 2835b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 2845b89ad90SMark F. Adams /* 2855b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 2865b89ad90SMark F. Adams by setting data structures and options. 2875b89ad90SMark F. Adams 2885b89ad90SMark F. Adams Input Parameter: 2895b89ad90SMark F. Adams . pc - the preconditioner context 2905b89ad90SMark F. Adams 2915b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 2925b89ad90SMark F. Adams 2935b89ad90SMark F. Adams Notes: 2945b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 2955b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 2965b89ad90SMark F. Adams */ 2975b89ad90SMark F. Adams #undef __FUNCT__ 2985b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 2995b89ad90SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc ) 3005b89ad90SMark F. Adams { 3015b89ad90SMark F. Adams PetscErrorCode ierr; 3025b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 3035b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 3045b89ad90SMark F. Adams Mat Amat = pc->mat, Pmat = pc->pmat; 3055b89ad90SMark F. Adams PetscBool isSeq, isMPI; 3065b89ad90SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, lidx; 3075b89ad90SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 3083530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 3090205a208SMark F. Adams PetscBool isOK; 3105b89ad90SMark F. Adams 311*5ef31b24SMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data; 312*5ef31b24SMark F. Adams 3135b89ad90SMark F. Adams PetscFunctionBegin; 314baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 315baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 3165b89ad90SMark F. Adams if (pc->setupcalled){ 3175b89ad90SMark F. Adams /* no state data in GAMG to destroy (now) */ 3185b89ad90SMark F. Adams ierr = PCReset_MG(pc); CHKERRQ(ierr); 3195b89ad90SMark F. Adams } 3206c237d78SBarry Smith if (!pc_gamg->m_data) SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates"); 3215b89ad90SMark F. Adams /* setup special features of PCGAMG */ 3225b89ad90SMark F. Adams ierr = PetscTypeCompare((PetscObject)Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 3235b89ad90SMark F. Adams ierr = PetscTypeCompare((PetscObject)Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 3245b89ad90SMark F. Adams if (isMPI) { 3255b89ad90SMark F. Adams } else if (isSeq) { 3265b89ad90SMark F. Adams } else SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with GAMG. GAMG can only handle AIJ matrices.",((PetscObject)Amat)->type_name); 3275b89ad90SMark F. Adams 3285b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 3295b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 3305b89ad90SMark F. Adams if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs); 3315b89ad90SMark F. Adams 3325b89ad90SMark F. Adams /* Get A_i and R_i */ 3333530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 3343530afc2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s level %d N=%d\n",0,__FUNCT__,0,N); 3350205a208SMark F. Adams for (level=0, Aarr[0] = Pmat, nactivepe = npe; 336*5ef31b24SMark F. Adams level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM); /* hard wired stopping logic */ 3370205a208SMark F. Adams level++ ) { 3385b89ad90SMark F. Adams level1 = level + 1; 339*5ef31b24SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr); 3400f9369f8SMark F. Adams ierr = createProlongation( Aarr[level], crds, pc_gamg->m_dim, 341*5ef31b24SMark F. Adams &Parr[level1], &coarse_crds, &isOK ); 3425b89ad90SMark F. Adams CHKERRQ(ierr); 343*5ef31b24SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr); 344*5ef31b24SMark F. Adams 3455b89ad90SMark F. Adams ierr = PetscFree( crds ); CHKERRQ( ierr ); 346baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 347baa4e9faSMark F. Adams if( isOK ) { 348*5ef31b24SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr); 349*5ef31b24SMark F. Adams ierr = partitionLevel( Aarr[level], pc_gamg->m_dim, &Parr[level1], &coarse_crds, &nactivepe, &Aarr[level1] ); 3503530afc2SMark F. Adams CHKERRQ(ierr); 351*5ef31b24SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr); 3523530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 3533530afc2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done level %d N=%d, active pe = %d\n",0,__FUNCT__,level1,N,nactivepe); 354baa4e9faSMark F. Adams } 355baa4e9faSMark F. Adams else{ 356baa4e9faSMark F. Adams break; 357baa4e9faSMark F. Adams } 35811e60469SMark F. Adams crds = coarse_crds; 3595b89ad90SMark F. Adams } 3605b89ad90SMark F. Adams ierr = PetscFree( coarse_crds ); CHKERRQ( ierr ); 361baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 3625b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 3635b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 3645b89ad90SMark F. Adams fine_level = level; 3655b89ad90SMark F. Adams ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 3665b89ad90SMark F. Adams 3675b89ad90SMark F. Adams /* set default smoothers */ 3685b89ad90SMark F. Adams PetscReal emax = 2.0, emin; 3695b89ad90SMark F. Adams for (level=1; level<=fine_level; level++){ 3705b89ad90SMark F. Adams KSP smoother; PC subpc; 3715b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 3725b89ad90SMark F. Adams ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr); 3730f9369f8SMark F. Adams emin = emax/10.0; /* fix!!! */ 3745b89ad90SMark F. Adams ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 3755b89ad90SMark F. Adams ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr); 3765b89ad90SMark F. Adams ierr = KSPChebychevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); /* need auto !!!!*/ 3775b89ad90SMark F. Adams } 3785b89ad90SMark F. Adams ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 3795b89ad90SMark F. Adams { 3805b89ad90SMark F. Adams PetscBool galerkin; 3815b89ad90SMark F. Adams ierr = PCMGGetGalerkin( pc, &galerkin); CHKERRQ(ierr); 3825b89ad90SMark F. Adams if(galerkin){ 3835b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 3845b89ad90SMark F. Adams } 3855b89ad90SMark F. Adams } 386*5ef31b24SMark F. Adams { 387*5ef31b24SMark F. Adams char str[32]; 388*5ef31b24SMark F. Adams sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1); 389*5ef31b24SMark F. Adams PetscLogStageRegister(str, &gamg_stages[fine_level]); 390*5ef31b24SMark F. Adams } 3915b89ad90SMark F. Adams /* create coarse level and the interpolation between the levels */ 3925b89ad90SMark F. Adams for (level=0,lidx=pc_gamg->m_Nlevels-1; level<fine_level; level++,lidx--){ 3935b89ad90SMark F. Adams level1 = level + 1; 394*5ef31b24SMark F. Adams ierr = PCMGSetInterpolation(pc,level1,Parr[lidx]);CHKERRQ(ierr); 3956c237d78SBarry Smith if( !PETSC_TRUE ) { 39611e60469SMark F. Adams PetscViewer viewer; char fname[32]; 39711e60469SMark F. Adams sprintf(fname,"Amat_%d.m",lidx); 39811e60469SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 3995b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 40011e60469SMark F. Adams ierr = MatView( Aarr[lidx], viewer ); CHKERRQ(ierr); 4015b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 4025b89ad90SMark F. Adams } 403*5ef31b24SMark F. Adams ierr = MatDestroy( &Parr[lidx] ); CHKERRQ(ierr); 404a92563c5SMark F. Adams { 4055b89ad90SMark F. Adams KSP smoother; 4065b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr); 4075b89ad90SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN ); 4085b89ad90SMark F. Adams CHKERRQ(ierr); 409*5ef31b24SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 410*5ef31b24SMark F. Adams } 411a92563c5SMark F. Adams ierr = MatDestroy( &Aarr[lidx] ); CHKERRQ(ierr); 412*5ef31b24SMark F. Adams { 413*5ef31b24SMark F. Adams char str[32]; 414*5ef31b24SMark F. Adams sprintf(str,"MG Level %d (%d)",level+1,lidx-1); 415*5ef31b24SMark F. Adams PetscLogStageRegister(str, &gamg_stages[lidx-1]); 416a92563c5SMark F. Adams } 4175b89ad90SMark F. Adams } 4185b89ad90SMark F. Adams { /* fine level (no P) */ 4195b89ad90SMark F. Adams KSP smoother; 4205b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr); 4215b89ad90SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN ); 4225b89ad90SMark F. Adams CHKERRQ(ierr); 423*5ef31b24SMark F. Adams ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 4245b89ad90SMark F. Adams } 4255b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 4265b89ad90SMark F. Adams pc->setupcalled = 0; 4275b89ad90SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 4285b89ad90SMark F. Adams PetscFunctionReturn(0); 4295b89ad90SMark F. Adams } 4305b89ad90SMark F. Adams 4315b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4325b89ad90SMark F. Adams /* 4335b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 4345b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 4355b89ad90SMark F. Adams 4365b89ad90SMark F. Adams Input Parameter: 4375b89ad90SMark F. Adams . pc - the preconditioner context 4385b89ad90SMark F. Adams 4395b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 4405b89ad90SMark F. Adams */ 4415b89ad90SMark F. Adams #undef __FUNCT__ 4425b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 4435b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 4445b89ad90SMark F. Adams { 4455b89ad90SMark F. Adams PetscErrorCode ierr; 4465b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 4475b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 4485b89ad90SMark F. Adams 4495b89ad90SMark F. Adams PetscFunctionBegin; 4505b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 4515b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 4525b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 4535b89ad90SMark F. Adams PetscFunctionReturn(0); 4545b89ad90SMark F. Adams } 4555b89ad90SMark F. Adams 4565b89ad90SMark F. Adams #undef __FUNCT__ 4575b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 4585b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 4595b89ad90SMark F. Adams { 4605b89ad90SMark F. Adams /* PetscErrorCode ierr; */ 4615b89ad90SMark F. Adams /* PC_MG *mg = (PC_MG*)pc->data; */ 4625b89ad90SMark F. Adams /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 4635b89ad90SMark F. Adams /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 4645b89ad90SMark F. Adams 4655b89ad90SMark F. Adams PetscFunctionBegin; 4665b89ad90SMark F. Adams PetscFunctionReturn(0); 4675b89ad90SMark F. Adams } 4685b89ad90SMark F. Adams 4695b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4705b89ad90SMark F. Adams /* 4715b89ad90SMark F. Adams PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 4725b89ad90SMark F. Adams 4735b89ad90SMark F. Adams Input Parameter: 4745b89ad90SMark F. Adams . pc - the preconditioner context 4755b89ad90SMark F. Adams 4765b89ad90SMark F. Adams Application Interface Routine: PCCreate() 4775b89ad90SMark F. Adams 4785b89ad90SMark F. Adams */ 4795b89ad90SMark F. Adams /* MC 4805b89ad90SMark F. Adams PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 4815b89ad90SMark F. Adams fine grid discretization matrix and coordinates on the fine grid. 4825b89ad90SMark F. Adams 4835b89ad90SMark F. Adams Options Database Key: 4845b89ad90SMark F. Adams Multigrid options(inherited) 4855b89ad90SMark F. Adams + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 4865b89ad90SMark F. Adams . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 4875b89ad90SMark F. Adams . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 4885b89ad90SMark F. Adams -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 4895b89ad90SMark F. Adams GAMG options: 4905b89ad90SMark F. Adams 4915b89ad90SMark F. Adams Level: intermediate 4925b89ad90SMark F. Adams Concepts: multigrid 4935b89ad90SMark F. Adams 4945b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 4955b89ad90SMark F. Adams PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 4965b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 4975b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 4985b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4995b89ad90SMark F. Adams M */ 5005b89ad90SMark F. Adams 5015b89ad90SMark F. Adams EXTERN_C_BEGIN 5025b89ad90SMark F. Adams #undef __FUNCT__ 5035b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 5045b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 5055b89ad90SMark F. Adams { 5065b89ad90SMark F. Adams PetscErrorCode ierr; 5075b89ad90SMark F. Adams PC_GAMG *pc_gamg; 5085b89ad90SMark F. Adams PC_MG *mg; 509*5ef31b24SMark F. Adams PetscClassId cookie; 5105b89ad90SMark F. Adams 5115b89ad90SMark F. Adams PetscFunctionBegin; 5125b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 5135b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 5145b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 5155b89ad90SMark F. Adams 5165b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 5175b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 5185b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 5195b89ad90SMark F. Adams mg->innerctx = pc_gamg; 5205b89ad90SMark F. Adams 5215b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 5225b89ad90SMark F. Adams 5235b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 5245b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 5255b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 5265b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 5275b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 5285b89ad90SMark F. Adams 5295b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 5305b89ad90SMark F. Adams "PCSetCoordinates_C", 5315b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 5325b89ad90SMark F. Adams PCSetCoordinates_GAMG);CHKERRQ(ierr); 533*5ef31b24SMark F. Adams 534*5ef31b24SMark F. Adams PetscClassIdRegister("GAMG Setup",&cookie); 535*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-createProl", cookie, &gamg_setup_stages[SET1]); 536*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-partLevel", cookie, &gamg_setup_stages[SET2]); 537*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-MIS Graph", cookie, &gamg_setup_stages[SET3]); 538*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-MIS-Agg", cookie, &gamg_setup_stages[SET4]); 539*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-growSupp", cookie, &gamg_setup_stages[SET5]); 540*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-tri-Prol", cookie, &gamg_setup_stages[SET6]); 541*5ef31b24SMark F. Adams PetscLogEventRegister("GAMG-find-prol", cookie, &gamg_setup_stages[FIND_V]); 542*5ef31b24SMark F. Adams 5435b89ad90SMark F. Adams PetscFunctionReturn(0); 5445b89ad90SMark F. Adams } 5455b89ad90SMark F. Adams EXTERN_C_END 546