15b89ad90SMark F. Adams /* 25b89ad90SMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 45b89ad90SMark F. Adams #include <private/pcimpl.h> /*I "petscpc.h" I*/ 55b89ad90SMark F. Adams #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 65b89ad90SMark F. Adams #include <../src/mat/impls/aij/seq/aij.h> 75b89ad90SMark F. Adams #include <../src/mat/impls/aij/mpi/mpiaij.h> 811e60469SMark F. Adams #include <assert.h> 95b89ad90SMark F. Adams 105b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */ 115b89ad90SMark F. Adams typedef struct gamg_TAG { 125b89ad90SMark F. Adams PetscInt m_dim; 135b89ad90SMark F. Adams PetscInt m_Nlevels; 145b89ad90SMark F. Adams PetscInt m_data_sz; 155b89ad90SMark F. Adams PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 165b89ad90SMark F. Adams } PC_GAMG; 175b89ad90SMark F. Adams 185b89ad90SMark F. Adams /* -----------------------------------------------------------------------------*/ 195b89ad90SMark F. Adams #undef __FUNCT__ 205b89ad90SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 215b89ad90SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 225b89ad90SMark F. Adams { 235b89ad90SMark F. Adams PetscErrorCode ierr; 245b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 255b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 265b89ad90SMark F. Adams 275b89ad90SMark F. Adams PetscFunctionBegin; 285b89ad90SMark F. Adams ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 295b89ad90SMark F. Adams PetscFunctionReturn(0); 305b89ad90SMark F. Adams } 315b89ad90SMark F. Adams 325b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 335b89ad90SMark F. Adams /* 345b89ad90SMark F. Adams PCSetCoordinates_GAMG 355b89ad90SMark F. Adams 365b89ad90SMark F. Adams Input Parameter: 375b89ad90SMark F. Adams . pc - the preconditioner context 385b89ad90SMark F. Adams */ 39a92563c5SMark F. Adams EXTERN_C_BEGIN 405b89ad90SMark F. Adams #undef __FUNCT__ 415b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 426c237d78SBarry Smith PetscErrorCode PCSetCoordinates_GAMG(PC pc, const int ndm,PetscReal *coords ) 435b89ad90SMark F. Adams { 445b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 455b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 466c237d78SBarry Smith PetscErrorCode ierr; 476c237d78SBarry Smith PetscInt bs, my0, tt; 486c237d78SBarry Smith Mat mat = pc->pmat; 496c237d78SBarry Smith PetscInt arrsz; 505b89ad90SMark F. Adams 515b89ad90SMark F. Adams PetscFunctionBegin; 525b89ad90SMark F. Adams ierr = MatGetBlockSize( mat, &bs ); CHKERRQ( ierr ); 535b89ad90SMark F. Adams ierr = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr); 546c237d78SBarry Smith arrsz = (tt-my0)/bs*ndm; 555b89ad90SMark F. Adams 565b89ad90SMark F. Adams // put coordinates 576c237d78SBarry Smith if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 585b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 596c237d78SBarry Smith ierr = PetscMalloc(arrsz*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 605b89ad90SMark F. Adams } 615b89ad90SMark F. Adams 625b89ad90SMark F. Adams /* copy data in */ 635b89ad90SMark F. Adams for(tt=0;tt<arrsz;tt++){ 645b89ad90SMark F. Adams pc_gamg->m_data[tt] = coords[tt]; 655b89ad90SMark F. Adams } 665b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 675b89ad90SMark F. Adams pc_gamg->m_dim = ndm; 685b89ad90SMark F. Adams 695b89ad90SMark F. Adams PetscFunctionReturn(0); 705b89ad90SMark F. Adams } 71a92563c5SMark F. Adams EXTERN_C_END 725b89ad90SMark F. Adams 735b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 745b89ad90SMark F. Adams /* 7511e60469SMark F. Adams partitionLevel 765b89ad90SMark F. Adams 775b89ad90SMark F. Adams Input Parameter: 78*3530afc2SMark F. Adams . a_Amat_fine - matrix on this fine (k) level 7911e60469SMark F. Adams . a_dime - 2 or 3 80*3530afc2SMark F. Adams In/Output Parameter: 81*3530afc2SMark F. Adams . a_P_inout - prolongation operator to the next level (k-1) 82*3530afc2SMark F. Adams . a_coarse_crds - coordinates that need to be moved 83*3530afc2SMark F. Adams . a_active_proc - number of active procs 8411e60469SMark F. Adams Output Parameter: 85*3530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 865b89ad90SMark F. Adams */ 875b89ad90SMark F. Adams #undef __FUNCT__ 8811e60469SMark F. Adams #define __FUNCT__ "partitionLevel" 89*3530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine, 9011e60469SMark F. Adams PetscInt a_dim, 91*3530afc2SMark F. Adams Mat *a_P_inout, 92*3530afc2SMark F. Adams PetscReal **a_coarse_crds, 93*3530afc2SMark F. Adams PetscMPIInt *a_active_proc, 94*3530afc2SMark F. Adams Mat *a_Amat_crs 9511e60469SMark F. Adams ) 965b89ad90SMark F. Adams { 975b89ad90SMark F. Adams PetscErrorCode ierr; 98*3530afc2SMark F. Adams Mat Amat, Pnew, Pold = *a_P_inout; 9911e60469SMark F. Adams IS new_indices,isnum; 100*3530afc2SMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 101*3530afc2SMark F. Adams PetscMPIInt nactive,mype,npe; 10211e60469SMark F. Adams PetscInt Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,bs=1; /* bs ??? */ 1035b89ad90SMark F. Adams 1045b89ad90SMark F. Adams PetscFunctionBegin; 10511e60469SMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 10611e60469SMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 10711e60469SMark F. Adams /* RAP */ 108*3530afc2SMark F. Adams ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Amat ); CHKERRQ(ierr); 109*3530afc2SMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart0, &Iend0 ); CHKERRQ(ierr); /* x2 size of 'a_coarse_crds' */ 11011e60469SMark F. Adams ncrs0 = (Iend0 - Istart0)/bs; 111acadaa72SMark F. Adams 11211e60469SMark F. Adams /* Repartition Amat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 11311e60469SMark F. Adams { 114*3530afc2SMark F. Adams PetscInt neq,N,counts[npe]; 11511e60469SMark F. Adams IS isnewproc; 116*3530afc2SMark F. Adams PetscMPIInt factstart,fact,new_npe,targ_npe; 117*3530afc2SMark F. Adams 118*3530afc2SMark F. Adams ierr = MatGetSize( Amat, &neq, &N );CHKERRQ(ierr); 119*3530afc2SMark F. Adams #define MIN_EQ_PROC 100 120*3530afc2SMark F. Adams nactive = *a_active_proc; 121*3530afc2SMark F. Adams targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */ 122*3530afc2SMark F. Adams if( targ_npe == 0 || neq < 1000 ) new_npe = 1; /* chop coarsest grid */ 123*3530afc2SMark F. Adams else if (targ_npe > nactive ) new_npe = nactive; /* no change */ 124*3530afc2SMark F. Adams else { 125*3530afc2SMark F. Adams new_npe = -9999; 126*3530afc2SMark F. Adams factstart = nactive; 127*3530afc2SMark F. Adams for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */ 128*3530afc2SMark F. Adams if( nactive%fact==0 && neq/(nactive/fact) > MIN_EQ_PROC ) { 129*3530afc2SMark F. Adams new_npe = nactive/fact; 130*3530afc2SMark F. Adams } 131*3530afc2SMark F. Adams } 132*3530afc2SMark F. Adams assert(new_npe != -9999); 133*3530afc2SMark F. Adams } 134*3530afc2SMark F. Adams *a_active_proc = new_npe; /* output for next time */ 135*3530afc2SMark F. Adams 13611e60469SMark F. Adams { /* partition: get 'isnewproc' */ 13711e60469SMark F. Adams MatPartitioning mpart; 13811e60469SMark F. Adams ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr); 13911e60469SMark F. Adams ierr = MatPartitioningSetAdjacency( mpart, Amat ); CHKERRQ(ierr); 14011e60469SMark F. Adams ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 14111e60469SMark F. Adams ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 14211e60469SMark F. Adams ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 14311e60469SMark F. Adams } 14411e60469SMark F. Adams /* 14511e60469SMark F. Adams Create an index set from the isnewproc index set to indicate the mapping TO 14611e60469SMark F. Adams */ 14711e60469SMark F. Adams ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 14811e60469SMark F. Adams /* 14911e60469SMark F. Adams Determine how many elements are assigned to each processor 15011e60469SMark F. Adams */ 15111e60469SMark F. Adams ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr); 15211e60469SMark F. Adams ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 15311e60469SMark F. Adams ncrs_new = counts[mype]; 15411e60469SMark F. Adams } 15511e60469SMark F. Adams { /* Create a vector to contain the newly ordered element information */ 15611e60469SMark F. Adams const PetscInt *idx; 15711e60469SMark F. Adams PetscInt i,j,k; 15811e60469SMark F. Adams IS isscat; 15911e60469SMark F. Adams PetscScalar *array; 16011e60469SMark F. Adams Vec src_crd, dest_crd; 161*3530afc2SMark F. Adams PetscReal *coords = *a_coarse_crds; 16211e60469SMark F. Adams VecScatter vecscat; 163acadaa72SMark F. Adams 16411e60469SMark F. Adams ierr = VecCreate( wcomm, &dest_crd ); 16511e60469SMark F. Adams ierr = VecSetSizes( dest_crd, a_dim*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 16611e60469SMark F. Adams ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/ 16711e60469SMark F. Adams /* 16811e60469SMark F. Adams There are three data items per cell (element), the integer vertex numbers of its three 16911e60469SMark F. Adams coordinates (we convert to double to use the scatter) (one can think 17011e60469SMark F. Adams of the vectors of having a block size of 3, then there is one index in idx[] for each element) 17111e60469SMark F. Adams */ 17211e60469SMark F. Adams { 17311e60469SMark F. Adams PetscInt tidx[ncrs0*a_dim]; 17411e60469SMark F. Adams ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 17511e60469SMark F. Adams for (i=0,j=0; i<ncrs0; i++) { 17611e60469SMark F. Adams for (k=0; k<a_dim; k++) tidx[j++] = idx[i]*a_dim + k; 17711e60469SMark F. Adams } 17811e60469SMark F. Adams ierr = ISCreateGeneral( wcomm, a_dim*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 17911e60469SMark F. Adams CHKERRQ(ierr); 18011e60469SMark F. Adams ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 18111e60469SMark F. Adams } 18211e60469SMark F. Adams /* 18311e60469SMark F. Adams Create a vector to contain the original vertex information for each element 18411e60469SMark F. Adams */ 18511e60469SMark F. Adams ierr = VecCreateSeq( PETSC_COMM_SELF, a_dim*ncrs0, &src_crd ); CHKERRQ(ierr); 18611e60469SMark F. Adams ierr = VecGetArray( src_crd, &array ); CHKERRQ(ierr); 18711e60469SMark F. Adams for (i=0; i<a_dim*ncrs0; i++) array[i] = coords[i]; 18811e60469SMark F. Adams ierr = VecRestoreArray( src_crd, &array ); CHKERRQ(ierr); 18911e60469SMark F. Adams /* 19011e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 19111e60469SMark F. Adams to the correct processor 19211e60469SMark F. Adams */ 19311e60469SMark F. Adams ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 19411e60469SMark F. Adams CHKERRQ(ierr); 19511e60469SMark F. Adams ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 19611e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 19711e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 19811e60469SMark F. Adams ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 19911e60469SMark F. Adams ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 20011e60469SMark F. Adams /* 20111e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 20211e60469SMark F. Adams */ 203*3530afc2SMark F. Adams ierr = PetscFree( *a_coarse_crds ); CHKERRQ(ierr); 204*3530afc2SMark F. Adams ierr = PetscMalloc( a_dim*ncrs_new*sizeof(PetscReal), a_coarse_crds ); CHKERRQ(ierr); 205*3530afc2SMark F. Adams coords = *a_coarse_crds; /* convient */ 20611e60469SMark F. Adams ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 20711e60469SMark F. Adams for (i=0; i<a_dim*ncrs_new; i++) coords[i] = PetscRealPart(array[i]); 20811e60469SMark F. Adams ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 20911e60469SMark F. Adams ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 21011e60469SMark F. Adams } 21111e60469SMark F. Adams /* 21211e60469SMark F. Adams Invert for MatGetSubMatrix 21311e60469SMark F. Adams */ 21411e60469SMark F. Adams ierr = ISInvertPermutation( isnum, ncrs_new, &new_indices ); CHKERRQ(ierr); 21511e60469SMark F. Adams ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 21611e60469SMark F. Adams ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 21711e60469SMark F. Adams /* A_crs output */ 218*3530afc2SMark F. Adams ierr = MatGetSubMatrix( Amat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 21911e60469SMark F. Adams CHKERRQ(ierr); 22011e60469SMark F. Adams ierr = MatDestroy( &Amat ); CHKERRQ(ierr); 221*3530afc2SMark F. Adams Amat = *a_Amat_crs; 22211e60469SMark F. Adams /* prolongator */ 22311e60469SMark F. Adams ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 22411e60469SMark F. Adams { 22511e60469SMark F. Adams IS findices; 22611e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 22711e60469SMark F. Adams ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 22811e60469SMark F. Adams CHKERRQ(ierr); 22911e60469SMark F. Adams ierr = ISDestroy( &findices ); CHKERRQ(ierr); 23011e60469SMark F. Adams } 231*3530afc2SMark F. Adams ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 232*3530afc2SMark F. Adams *a_P_inout = Pnew; /* output */ 233acadaa72SMark F. Adams 23411e60469SMark F. Adams ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 2355b89ad90SMark F. Adams 2365b89ad90SMark F. Adams PetscFunctionReturn(0); 2375b89ad90SMark F. Adams } 2385b89ad90SMark F. Adams 2395b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 2405b89ad90SMark F. Adams /* 2415b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 2425b89ad90SMark F. Adams by setting data structures and options. 2435b89ad90SMark F. Adams 2445b89ad90SMark F. Adams Input Parameter: 2455b89ad90SMark F. Adams . pc - the preconditioner context 2465b89ad90SMark F. Adams 2475b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 2485b89ad90SMark F. Adams 2495b89ad90SMark F. Adams Notes: 2505b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 2515b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 2525b89ad90SMark F. Adams */ 2535b89ad90SMark F. Adams extern PetscErrorCode PCSetFromOptions_MG(PC); 2545b89ad90SMark F. Adams extern PetscErrorCode PCReset_MG(PC); 2550f9369f8SMark F. Adams extern PetscErrorCode createProlongation( Mat, PetscReal [], const PetscInt, 256baa4e9faSMark F. Adams Mat *, PetscReal **, PetscBool *a_isOK ); 2575b89ad90SMark F. Adams #undef __FUNCT__ 2585b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 2595b89ad90SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc ) 2605b89ad90SMark F. Adams { 2615b89ad90SMark F. Adams PetscErrorCode ierr; 2625b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 2635b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 2645b89ad90SMark F. Adams Mat Amat = pc->mat, Pmat = pc->pmat; 2655b89ad90SMark F. Adams PetscBool isSeq, isMPI; 2665b89ad90SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, lidx; 2675b89ad90SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 268*3530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 2695b89ad90SMark F. Adams 2705b89ad90SMark F. Adams PetscFunctionBegin; 271baa4e9faSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 272baa4e9faSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 2735b89ad90SMark F. Adams if (pc->setupcalled){ 2745b89ad90SMark F. Adams /* no state data in GAMG to destroy (now) */ 2755b89ad90SMark F. Adams ierr = PCReset_MG(pc); CHKERRQ(ierr); 2765b89ad90SMark F. Adams } 2776c237d78SBarry Smith if (!pc_gamg->m_data) SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates"); 2785b89ad90SMark F. Adams /* setup special features of PCGAMG */ 2795b89ad90SMark F. Adams ierr = PetscTypeCompare((PetscObject)Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 2805b89ad90SMark F. Adams ierr = PetscTypeCompare((PetscObject)Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 2815b89ad90SMark F. Adams if (isMPI) { 2825b89ad90SMark F. Adams } else if (isSeq) { 2835b89ad90SMark 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); 2845b89ad90SMark F. Adams 2855b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 2865b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 2875b89ad90SMark F. Adams if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs); 2885b89ad90SMark F. Adams 2895b89ad90SMark F. Adams /* Get A_i and R_i */ 290*3530afc2SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 291*3530afc2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s level %d N=%d\n",0,__FUNCT__,0,N); 29211e60469SMark F. Adams #define GAMG_MAXLEVELS 20 2935b89ad90SMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Rarr[GAMG_MAXLEVELS]; PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data; 294baa4e9faSMark F. Adams PetscBool isOK; 295*3530afc2SMark F. Adams for (level=0, Aarr[0] = Pmat, nactivepe = npe; level < GAMG_MAXLEVELS-1; level++ ){ 296*3530afc2SMark F. Adams if( nactivepe == 1 ) { 2975b89ad90SMark F. Adams break; 2985b89ad90SMark F. Adams } 2995b89ad90SMark F. Adams level1 = level + 1; 3000f9369f8SMark F. Adams ierr = createProlongation( Aarr[level], crds, pc_gamg->m_dim, 301baa4e9faSMark F. Adams &Rarr[level1], &coarse_crds, &isOK ); 3025b89ad90SMark F. Adams CHKERRQ(ierr); 3035b89ad90SMark F. Adams ierr = PetscFree( crds ); CHKERRQ( ierr ); 304baa4e9faSMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 305baa4e9faSMark F. Adams if( isOK ) { 306*3530afc2SMark F. Adams ierr = partitionLevel( Aarr[level], pc_gamg->m_dim, &Rarr[level1], &coarse_crds, &nactivepe, &Aarr[level1] ); 307*3530afc2SMark F. Adams CHKERRQ(ierr); 308*3530afc2SMark F. Adams ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 309*3530afc2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done level %d N=%d, active pe = %d\n",0,__FUNCT__,level1,N,nactivepe); 310baa4e9faSMark F. Adams } 311baa4e9faSMark F. Adams else{ 312baa4e9faSMark F. Adams break; 313baa4e9faSMark F. Adams } 31411e60469SMark F. Adams crds = coarse_crds; 3155b89ad90SMark F. Adams } 3165b89ad90SMark F. Adams ierr = PetscFree( coarse_crds ); CHKERRQ( ierr ); 317baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 3185b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 3195b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 3205b89ad90SMark F. Adams fine_level = level; 3215b89ad90SMark F. Adams ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 3225b89ad90SMark F. Adams 3235b89ad90SMark F. Adams /* set default smoothers */ 3245b89ad90SMark F. Adams PetscReal emax = 2.0, emin; 3255b89ad90SMark F. Adams for (level=1; level<=fine_level; level++){ 3265b89ad90SMark F. Adams KSP smoother; PC subpc; 3275b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 3285b89ad90SMark F. Adams ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr); 3290f9369f8SMark F. Adams emin = emax/10.0; /* fix!!! */ 3305b89ad90SMark F. Adams ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 3315b89ad90SMark F. Adams ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr); 3325b89ad90SMark F. Adams ierr = KSPChebychevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); /* need auto !!!!*/ 3335b89ad90SMark F. Adams } 3345b89ad90SMark F. Adams ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 3355b89ad90SMark F. Adams { 3365b89ad90SMark F. Adams PetscBool galerkin; 3375b89ad90SMark F. Adams ierr = PCMGGetGalerkin( pc, &galerkin); CHKERRQ(ierr); 3385b89ad90SMark F. Adams if(galerkin){ 3395b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 3405b89ad90SMark F. Adams } 3415b89ad90SMark F. Adams } 3425b89ad90SMark F. Adams /* create coarse level and the interpolation between the levels */ 3435b89ad90SMark F. Adams for (level=0,lidx=pc_gamg->m_Nlevels-1; level<fine_level; level++,lidx--){ 3445b89ad90SMark F. Adams level1 = level + 1; 3455b89ad90SMark F. Adams ierr = PCMGSetInterpolation(pc,level1,Rarr[lidx]);CHKERRQ(ierr); 3466c237d78SBarry Smith if(!PETSC_TRUE) { 34711e60469SMark F. Adams PetscViewer viewer; char fname[32]; 34811e60469SMark F. Adams sprintf(fname,"Amat_%d.m",lidx); 34911e60469SMark F. Adams ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 3505b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 35111e60469SMark F. Adams ierr = MatView( Aarr[lidx], viewer ); CHKERRQ(ierr); 3525b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 3535b89ad90SMark F. Adams } 354a92563c5SMark F. Adams ierr = MatDestroy( &Rarr[lidx] ); CHKERRQ(ierr); 355a92563c5SMark F. Adams { 3565b89ad90SMark F. Adams KSP smoother; 3575b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr); 3585b89ad90SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN ); 3595b89ad90SMark F. Adams CHKERRQ(ierr); 360a92563c5SMark F. Adams ierr = MatDestroy( &Aarr[lidx] ); CHKERRQ(ierr); 361a92563c5SMark F. Adams } 3625b89ad90SMark F. Adams } 3635b89ad90SMark F. Adams { /* fine level (no P) */ 3645b89ad90SMark F. Adams KSP smoother; 3655b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr); 3665b89ad90SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN ); 3675b89ad90SMark F. Adams CHKERRQ(ierr); 3685b89ad90SMark F. Adams } 3695b89ad90SMark F. Adams 3705b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 3715b89ad90SMark F. Adams pc->setupcalled = 0; 3725b89ad90SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 3735b89ad90SMark F. Adams PetscFunctionReturn(0); 3745b89ad90SMark F. Adams } 3755b89ad90SMark F. Adams 3765b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 3775b89ad90SMark F. Adams /* 3785b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 3795b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 3805b89ad90SMark F. Adams 3815b89ad90SMark F. Adams Input Parameter: 3825b89ad90SMark F. Adams . pc - the preconditioner context 3835b89ad90SMark F. Adams 3845b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 3855b89ad90SMark F. Adams */ 3865b89ad90SMark F. Adams #undef __FUNCT__ 3875b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 3885b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 3895b89ad90SMark F. Adams { 3905b89ad90SMark F. Adams PetscErrorCode ierr; 3915b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 3925b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 3935b89ad90SMark F. Adams 3945b89ad90SMark F. Adams PetscFunctionBegin; 3955b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 3965b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 3975b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 3985b89ad90SMark F. Adams PetscFunctionReturn(0); 3995b89ad90SMark F. Adams } 4005b89ad90SMark F. Adams 4015b89ad90SMark F. Adams #undef __FUNCT__ 4025b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 4035b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 4045b89ad90SMark F. Adams { 4055b89ad90SMark F. Adams /* PetscErrorCode ierr; */ 4065b89ad90SMark F. Adams /* PC_MG *mg = (PC_MG*)pc->data; */ 4075b89ad90SMark F. Adams /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 4085b89ad90SMark F. Adams /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 4095b89ad90SMark F. Adams 4105b89ad90SMark F. Adams PetscFunctionBegin; 4115b89ad90SMark F. Adams PetscFunctionReturn(0); 4125b89ad90SMark F. Adams } 4135b89ad90SMark F. Adams 4145b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4155b89ad90SMark F. Adams /* 4165b89ad90SMark F. Adams PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 4175b89ad90SMark F. Adams 4185b89ad90SMark F. Adams Input Parameter: 4195b89ad90SMark F. Adams . pc - the preconditioner context 4205b89ad90SMark F. Adams 4215b89ad90SMark F. Adams Application Interface Routine: PCCreate() 4225b89ad90SMark F. Adams 4235b89ad90SMark F. Adams */ 4245b89ad90SMark F. Adams /* MC 4255b89ad90SMark F. Adams PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 4265b89ad90SMark F. Adams fine grid discretization matrix and coordinates on the fine grid. 4275b89ad90SMark F. Adams 4285b89ad90SMark F. Adams Options Database Key: 4295b89ad90SMark F. Adams Multigrid options(inherited) 4305b89ad90SMark F. Adams + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 4315b89ad90SMark F. Adams . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 4325b89ad90SMark F. Adams . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 4335b89ad90SMark F. Adams -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 4345b89ad90SMark F. Adams GAMG options: 4355b89ad90SMark F. Adams 4365b89ad90SMark F. Adams Level: intermediate 4375b89ad90SMark F. Adams Concepts: multigrid 4385b89ad90SMark F. Adams 4395b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 4405b89ad90SMark F. Adams PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 4415b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 4425b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 4435b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4445b89ad90SMark F. Adams M */ 4455b89ad90SMark F. Adams 4465b89ad90SMark F. Adams EXTERN_C_BEGIN 4475b89ad90SMark F. Adams #undef __FUNCT__ 4485b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 4495b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 4505b89ad90SMark F. Adams { 4515b89ad90SMark F. Adams PetscErrorCode ierr; 4525b89ad90SMark F. Adams PC_GAMG *pc_gamg; 4535b89ad90SMark F. Adams PC_MG *mg; 4545b89ad90SMark F. Adams 4555b89ad90SMark F. Adams PetscFunctionBegin; 4565b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 4575b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4585b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 4595b89ad90SMark F. Adams 4605b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 4615b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 4625b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 4635b89ad90SMark F. Adams mg->innerctx = pc_gamg; 4645b89ad90SMark F. Adams 4655b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 4665b89ad90SMark F. Adams 4675b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 4685b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 4695b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 4705b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 4715b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 4725b89ad90SMark F. Adams 4735b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 4745b89ad90SMark F. Adams "PCSetCoordinates_C", 4755b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 4765b89ad90SMark F. Adams PCSetCoordinates_GAMG);CHKERRQ(ierr); 4775b89ad90SMark F. Adams PetscFunctionReturn(0); 4785b89ad90SMark F. Adams } 4795b89ad90SMark F. Adams EXTERN_C_END 480