xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 038e3b619a263fb3ff92c5f8a0dd3ca839869fd8)
15b89ad90SMark F. Adams /*
25b89ad90SMark F. Adams  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
45ef31b24SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>
5f96513f1SMatthew G Knepley 
6f96513f1SMatthew G Knepley PetscLogEvent gamg_setup_stages[NUM_SET];
7f96513f1SMatthew G Knepley 
85b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */
95b89ad90SMark F. Adams typedef struct gamg_TAG {
105b89ad90SMark F. Adams   PetscInt       m_dim;
115b89ad90SMark F. Adams   PetscInt       m_Nlevels;
125b89ad90SMark F. Adams   PetscInt       m_data_sz;
135b89ad90SMark F. Adams   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
145b89ad90SMark F. Adams } PC_GAMG;
155b89ad90SMark F. Adams 
16fc4362bfSMark F. Adams #define TOP_GRID_LIM 100
175ef31b24SMark 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);
29eb07cef2SMark F. Adams   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
305b89ad90SMark F. Adams   PetscFunctionReturn(0);
315b89ad90SMark F. Adams }
325b89ad90SMark F. Adams 
335b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
345b89ad90SMark F. Adams /*
355b89ad90SMark F. Adams    PCSetCoordinates_GAMG
365b89ad90SMark F. Adams 
375b89ad90SMark F. Adams    Input Parameter:
385b89ad90SMark F. Adams    .  pc - the preconditioner context
395b89ad90SMark F. Adams */
40a92563c5SMark F. Adams EXTERN_C_BEGIN
415b89ad90SMark F. Adams #undef __FUNCT__
425b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
43eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
445b89ad90SMark F. Adams {
45eb07cef2SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
465b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
476c237d78SBarry Smith   PetscErrorCode ierr;
48*038e3b61SMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,sz,M;
49*038e3b61SMark F. Adams   Mat            Amat = a_pc->pmat;
50eb07cef2SMark F. Adams   PetscBool      useSA = PETSC_FALSE, flag;
51eb07cef2SMark F. Adams   char           str[16];
525b89ad90SMark F. Adams 
535b89ad90SMark F. Adams   PetscFunctionBegin;
54*038e3b61SMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
55*038e3b61SMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &kk ); CHKERRQ(ierr);
56*038e3b61SMark F. Adams   nloc = (kk-my0)/bs; M = kk - my0;
57*038e3b61SMark F. Adams   if((kk-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
58*038e3b61SMark F. Adams 
59eb07cef2SMark F. Adams   ierr  = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,16,&flag);    CHKERRQ( ierr );
60eb07cef2SMark F. Adams   useSA = (PetscBool)(flag && strcmp(str,"sa") == 0);
61*038e3b61SMark F. Adams 
62*038e3b61SMark F. Adams   if(a_coords == 0) useSA = PETSC_TRUE; /* use SA if no data */
63*038e3b61SMark F. Adams   if( !useSA ) sz = a_ndm; /* coordinates */
64*038e3b61SMark F. Adams   else{ /* SA: null space vectors */
65*038e3b61SMark F. Adams     if(a_coords != 0) sz = (a_ndm==2 ? 3*bs : 6*bs); /* elasticity */
66*038e3b61SMark F. Adams     else sz = bs*bs; /* no data, force SA with constant null space vectors */
67*038e3b61SMark F. Adams     if(a_coords != 0 && bs==1 ) sz = 1; /* scalar problem with coords and SA (not needed) */
68*038e3b61SMark F. Adams   }
69eb07cef2SMark F. Adams   arrsz = nloc*sz;
705b89ad90SMark F. Adams 
71*038e3b61SMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
726c237d78SBarry Smith   if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) {
735b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
74eb07cef2SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr);
755b89ad90SMark F. Adams   }
76*038e3b61SMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
77eb07cef2SMark F. Adams   pc_gamg->m_data[arrsz] = -99.;
78*038e3b61SMark F. Adams   /* copy data in - column oriented */
79eb07cef2SMark F. Adams   if( useSA ) {
80*038e3b61SMark F. Adams     for(kk=0;kk<nloc;kk++){
81*038e3b61SMark F. Adams       PetscReal *data = &pc_gamg->m_data[kk*bs];
82*038e3b61SMark F. Adams       if( sz==1 ) *data = 1.0;
83*038e3b61SMark F. Adams       else {
84*038e3b61SMark F. Adams         for(ii=0;ii<bs;ii++)
85*038e3b61SMark F. Adams 	  for(jj=0;jj<bs;jj++)
86*038e3b61SMark F. Adams 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
87*038e3b61SMark F. Adams 	    else data[ii*M + jj] = 0.0;
88eb07cef2SMark F. Adams         if( a_coords != 0 ) {
89*038e3b61SMark F. Adams           if( a_ndm == 2 ){ /* rotational modes */
90*038e3b61SMark F. Adams             data += 2*M;
91*038e3b61SMark F. Adams             data[0] = -a_coords[2*kk+1];
92*038e3b61SMark F. Adams             data[1] =  a_coords[2*kk];
935b89ad90SMark F. Adams           }
94eb07cef2SMark F. Adams           else {
95*038e3b61SMark F. Adams             data += 3*M;
96*038e3b61SMark F. Adams             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
97*038e3b61SMark F. Adams             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
98*038e3b61SMark F. Adams             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
99*038e3b61SMark F. Adams           }
100eb07cef2SMark F. Adams         }
101eb07cef2SMark F. Adams       }
102eb07cef2SMark F. Adams     }
103eb07cef2SMark F. Adams   }
104eb07cef2SMark F. Adams   else {
105*038e3b61SMark F. Adams     for( kk = 0 ; kk < nloc ; kk++ ){
106*038e3b61SMark F. Adams       for( ii = 0 ; ii < a_ndm ; ii++ ) {
107*038e3b61SMark F. Adams         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
108eb07cef2SMark F. Adams       }
109eb07cef2SMark F. Adams     }
110*038e3b61SMark F. Adams   }
111*038e3b61SMark F. Adams 
112eb07cef2SMark F. Adams   assert(pc_gamg->m_data[arrsz] == -99.);
113*038e3b61SMark F. Adams   for(kk=0;kk<arrsz;kk++){
114*038e3b61SMark F. Adams     assert(pc_gamg->m_data[kk] != -999.); // debug
115*038e3b61SMark F. Adams   }
116*038e3b61SMark F. Adams 
1175b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
118eb07cef2SMark F. Adams   pc_gamg->m_dim = a_ndm;
1195b89ad90SMark F. Adams 
1205b89ad90SMark F. Adams   PetscFunctionReturn(0);
1215b89ad90SMark F. Adams }
122a92563c5SMark F. Adams EXTERN_C_END
1235b89ad90SMark F. Adams 
1245b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1255b89ad90SMark F. Adams /*
12611e60469SMark F. Adams    partitionLevel
1275b89ad90SMark F. Adams 
1285b89ad90SMark F. Adams    Input Parameter:
1293530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
130*038e3b61SMark F. Adams    . a_data_sz - size of data to move (coarse grid)
1313530afc2SMark F. Adams    In/Output Parameter:
1323530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
133eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
1343530afc2SMark F. Adams    . a_active_proc - number of active procs
13511e60469SMark F. Adams    Output Parameter:
1363530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1375b89ad90SMark F. Adams */
1385cb416c2SMark F. Adams 
1395b89ad90SMark F. Adams #undef __FUNCT__
14011e60469SMark F. Adams #define __FUNCT__ "partitionLevel"
1413530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine,
142eb07cef2SMark F. Adams                                PetscInt a_data_sz,
143*038e3b61SMark F. Adams 			       PetscInt a_cbs,
1443530afc2SMark F. Adams                                Mat *a_P_inout,
145eb07cef2SMark F. Adams                                PetscReal **a_coarse_data,
1463530afc2SMark F. Adams                                PetscMPIInt *a_active_proc,
1473530afc2SMark F. Adams                                Mat *a_Amat_crs
14811e60469SMark F. Adams                                )
1495b89ad90SMark F. Adams {
1505b89ad90SMark F. Adams   PetscErrorCode   ierr;
151*038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
15211e60469SMark F. Adams   IS               new_indices,isnum;
1533530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
1545ef31b24SMark F. Adams   PetscMPIInt      nactive_procs,mype,npe;
155*038e3b61SMark F. Adams   PetscInt         Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,fbs;
1565b89ad90SMark F. Adams 
1575b89ad90SMark F. Adams   PetscFunctionBegin;
15811e60469SMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
15911e60469SMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
160*038e3b61SMark F. Adams   ierr = MatGetBlockSize( a_Amat_fine, &fbs ); CHKERRQ(ierr);
16111e60469SMark F. Adams   /* RAP */
162*038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
163*038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
164acadaa72SMark F. Adams 
165*038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
166*038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert( (Iend0-Istart0)%a_cbs == 0 );
167*038e3b61SMark F. Adams 
168*038e3b61SMark F. Adams   /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
16911e60469SMark F. Adams   {
170*038e3b61SMark F. Adams     PetscInt        neq,NN,counts[npe];
17111e60469SMark F. Adams     IS              isnewproc;
1725ef31b24SMark F. Adams     PetscMPIInt     new_npe,targ_npe;
1733530afc2SMark F. Adams 
174*038e3b61SMark F. Adams     ierr = MatGetSize( Cmat, &neq, &NN );CHKERRQ(ierr);
1753530afc2SMark F. Adams #define MIN_EQ_PROC 100
1765ef31b24SMark F. Adams     nactive_procs = *a_active_proc;
1773530afc2SMark F. Adams     targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */
1785ef31b24SMark F. Adams     if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */
1795ef31b24SMark F. Adams     else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */
1803530afc2SMark F. Adams     else {
1815ef31b24SMark F. Adams       PetscMPIInt     factstart,fact;
1823530afc2SMark F. Adams       new_npe = -9999;
1835ef31b24SMark F. Adams       factstart = nactive_procs;
1843530afc2SMark F. Adams       for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */
1855ef31b24SMark F. Adams         if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) {
1865ef31b24SMark F. Adams           new_npe = nactive_procs/fact;
1873530afc2SMark F. Adams         }
1883530afc2SMark F. Adams       }
189*038e3b61SMark F. Adams       if(new_npe == -9999) new_npe = nactive_procs;
1903530afc2SMark F. Adams     }
1913530afc2SMark F. Adams     *a_active_proc = new_npe; /* output for next time */
19211e60469SMark F. Adams     { /* partition: get 'isnewproc' */
19311e60469SMark F. Adams       MatPartitioning  mpart;
1945ef31b24SMark F. Adams       Mat              adj;
1955ef31b24SMark F. Adams       const PetscInt  *is_idx;
196eb07cef2SMark F. Adams       PetscInt         is_sz,kk,jj,ii,old_fact=(npe/nactive_procs),new_fact=(npe/new_npe),*isnewproc_idx;
197c9a0b8beSMark F. Adams       /* create sub communicator  */
198c9a0b8beSMark F. Adams       MPI_Comm cm,new_comm;
1995ef31b24SMark F. Adams       int membershipKey = mype % old_fact;
200c9a0b8beSMark F. Adams       ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr);
201c9a0b8beSMark F. Adams       ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
202c9a0b8beSMark F. Adams       ierr = MPI_Comm_free( &cm );                            CHKERRQ(ierr);
203c9a0b8beSMark F. Adams 
2045ef31b24SMark F. Adams       /* MatPartitioningApply call MatConvert, which is collective */
205*038e3b61SMark F. Adams       if( a_cbs==1) {
206*038e3b61SMark F. Adams         ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
207eb07cef2SMark F. Adams       }
208eb07cef2SMark F. Adams       else {
209*038e3b61SMark F. Adams 	/* make a scalar matrix to partition */
210eb07cef2SMark F. Adams         Mat tMat;
211eb07cef2SMark F. Adams         PetscInt Ii,ncols; const PetscScalar *vals; const PetscInt *idx;
212eb07cef2SMark F. Adams         ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
213eb07cef2SMark F. Adams                                 PETSC_DETERMINE, PETSC_DETERMINE,
214eb07cef2SMark F. Adams                                 25, PETSC_NULL, 10, PETSC_NULL,
215eb07cef2SMark F. Adams                                 &tMat );
216eb07cef2SMark F. Adams 
217eb07cef2SMark F. Adams         for ( Ii = Istart0; Ii < Iend0; Ii++ ) {
218*038e3b61SMark F. Adams           PetscInt dest_row = Ii/a_cbs;
219*038e3b61SMark F. Adams           ierr = MatGetRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
220eb07cef2SMark F. Adams           for( jj = 0 ; jj < ncols ; jj++ ){
221*038e3b61SMark F. Adams             PetscInt dest_col = idx[jj]/a_cbs;
222eb07cef2SMark F. Adams             PetscScalar v = 1.0;
223eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
224eb07cef2SMark F. Adams           }
225*038e3b61SMark F. Adams           ierr = MatRestoreRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
226eb07cef2SMark F. Adams         }
227eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
228eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
229eb07cef2SMark F. Adams 
230eb07cef2SMark F. Adams         ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
231eb07cef2SMark F. Adams 
232eb07cef2SMark F. Adams         ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
233eb07cef2SMark F. Adams       }
2345ef31b24SMark F. Adams       if( membershipKey == 0 ) {
235*038e3b61SMark F. Adams 	if(ncrs0==0)SETERRQ(wcomm,PETSC_ERR_LIB,"zero local nodes -- increase min");
2365ef31b24SMark F. Adams         /* hack to fix global data that pmetis.c uses in 'adj' */
2375ef31b24SMark F. Adams         for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) {
2385ef31b24SMark F. Adams           adj->rmap->range[jj] = adj->rmap->range[kk];
2395ef31b24SMark F. Adams         }
2405ef31b24SMark F. Adams         ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr);
2415ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
24211e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr);
2435ef31b24SMark F. Adams         ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr);
24411e60469SMark F. Adams         ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
24511e60469SMark F. Adams         ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr);
2465ef31b24SMark F. Adams         /* collect IS info */
2475ef31b24SMark F. Adams         ierr = ISGetLocalSize( isnewproc, &is_sz );        CHKERRQ(ierr);
248*038e3b61SMark F. Adams         ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
2495ef31b24SMark F. Adams         ierr = ISGetIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
250eb07cef2SMark F. Adams         /* spread partitioning across machine - probably the right thing to do but machine spec. */
251eb07cef2SMark F. Adams         for(kk=0,jj=0;kk<is_sz;kk++){
252*038e3b61SMark F. Adams           for(ii=0 ; ii<a_cbs ; ii++, jj++ ) { /* expand for equation level by 'bs' */
253eb07cef2SMark F. Adams             isnewproc_idx[jj] = is_idx[kk] * new_fact;
254eb07cef2SMark F. Adams           }
2555ef31b24SMark F. Adams         }
2565ef31b24SMark F. Adams         ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
2575ef31b24SMark F. Adams         ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
258*038e3b61SMark F. Adams         is_sz *= a_cbs;
2595ef31b24SMark F. Adams       }
2605ef31b24SMark F. Adams       else {
2615ef31b24SMark F. Adams         isnewproc_idx = 0;
2625ef31b24SMark F. Adams         is_sz = 0;
2635ef31b24SMark F. Adams       }
2645ef31b24SMark F. Adams       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
2655ef31b24SMark F. Adams       ierr = MPI_Comm_free( &new_comm );    CHKERRQ(ierr);
2665ef31b24SMark F. Adams       ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
2675ef31b24SMark F. Adams       if( membershipKey == 0 ) {
2685ef31b24SMark F. Adams         ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
2695ef31b24SMark F. Adams       }
27011e60469SMark F. Adams     }
27111e60469SMark F. Adams     /*
27211e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
27311e60469SMark F. Adams      */
27411e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
27511e60469SMark F. Adams     /*
27611e60469SMark F. Adams      Determine how many elements are assigned to each processor
27711e60469SMark F. Adams      */
27811e60469SMark F. Adams     ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr);
27911e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
280*038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
28111e60469SMark F. Adams   }
2825ef31b24SMark F. Adams 
28311e60469SMark F. Adams   { /* Create a vector to contain the newly ordered element information */
28411e60469SMark F. Adams     const PetscInt *idx;
285*038e3b61SMark F. Adams     PetscInt        ii,jj,kk;
28611e60469SMark F. Adams     IS              isscat;
28711e60469SMark F. Adams     PetscScalar    *array;
28811e60469SMark F. Adams     Vec             src_crd, dest_crd;
289eb07cef2SMark F. Adams     PetscReal      *data = *a_coarse_data;
29011e60469SMark F. Adams     VecScatter      vecscat;
291eb07cef2SMark F. Adams     PetscInt        tidx[ncrs0*a_data_sz];
292acadaa72SMark F. Adams 
29311e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
294eb07cef2SMark F. Adams     ierr = VecSetSizes( dest_crd, a_data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
29511e60469SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/
29611e60469SMark F. Adams     /*
297eb07cef2SMark F. Adams      There are 'a_data_sz' data items per node, (one can think of the vectors of having a
298*038e3b61SMark F. Adams      block size of 'a_data_sz').  Note, ISs are expanded into equation space by 'a_cbs'.
29911e60469SMark F. Adams      */
30011e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
301*038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
302*038e3b61SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs;
303*038e3b61SMark F. Adams       for( kk=0; kk<a_data_sz ; kk++, jj++) tidx[jj] = id*a_data_sz + kk;
30411e60469SMark F. Adams     }
305*038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
306eb07cef2SMark F. Adams     ierr = ISCreateGeneral( wcomm, a_data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
30711e60469SMark F. Adams     CHKERRQ(ierr);
30811e60469SMark F. Adams     /*
30911e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
31011e60469SMark F. Adams      */
311eb07cef2SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, a_data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
312*038e3b61SMark F. Adams     for(ii=0, jj=0; ii<ncrs0 ; ii++) {
313*038e3b61SMark F. Adams       for( kk=0; kk < a_data_sz ; kk++, jj++) {
314*038e3b61SMark F. Adams         ierr = VecSetValues( src_crd, 1, &jj, &data[kk*ncrs0 + ii], INSERT_VALUES );  CHKERRQ(ierr);
315*038e3b61SMark F. Adams       }
316eb07cef2SMark F. Adams     }
317eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
318eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
31911e60469SMark F. Adams     /*
32011e60469SMark F. Adams      Scatter the element vertex information (still in the original vertex ordering)
32111e60469SMark F. Adams      to the correct processor
32211e60469SMark F. Adams      */
32311e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
32411e60469SMark F. Adams     CHKERRQ(ierr);
32511e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
32611e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32711e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32811e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
32911e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
33011e60469SMark F. Adams     /*
33111e60469SMark F. Adams      Put the element vertex data into a new allocation of the gdata->ele
33211e60469SMark F. Adams      */
333eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
334eb07cef2SMark F. Adams     ierr = PetscMalloc( a_data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
335*038e3b61SMark F. Adams     VecGetLocalSize( dest_crd, &ii ); assert(ii==a_data_sz*ncrs_new);
336eb07cef2SMark F. Adams 
33711e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
338eb07cef2SMark F. Adams     data = *a_coarse_data;
339*038e3b61SMark F. Adams     for(ii=0, jj=0; ii<ncrs_new ; ii++) {
340*038e3b61SMark F. Adams       for( kk=0; kk < a_data_sz ; kk++, jj++) {
341*038e3b61SMark F. Adams         data[kk*ncrs_new + ii] = PetscRealPart(array[jj]);
342*038e3b61SMark F. Adams       }
343*038e3b61SMark F. Adams     }
34411e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
34511e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
34611e60469SMark F. Adams   }
34711e60469SMark F. Adams   /*
34811e60469SMark F. Adams    Invert for MatGetSubMatrix
34911e60469SMark F. Adams    */
350*038e3b61SMark F. Adams   ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
35111e60469SMark F. Adams   ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
35211e60469SMark F. Adams   ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
35311e60469SMark F. Adams   /* A_crs output */
354*038e3b61SMark F. Adams   ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
35511e60469SMark F. Adams   CHKERRQ(ierr);
356eb07cef2SMark F. Adams 
357*038e3b61SMark F. Adams   ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
358*038e3b61SMark F. Adams   Cmat = *a_Amat_crs;
359*038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
360eb07cef2SMark F. Adams 
36111e60469SMark F. Adams   /* prolongator */
36211e60469SMark F. Adams   ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
36311e60469SMark F. Adams   {
36411e60469SMark F. Adams     IS findices;
36511e60469SMark F. Adams     ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
36611e60469SMark F. Adams     ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
36711e60469SMark F. Adams     CHKERRQ(ierr);
36811e60469SMark F. Adams     ierr = ISDestroy( &findices ); CHKERRQ(ierr);
36911e60469SMark F. Adams   }
3703530afc2SMark F. Adams   ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
3713530afc2SMark F. Adams   *a_P_inout = Pnew; /* output */
37211e60469SMark F. Adams   ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
3735b89ad90SMark F. Adams 
3745b89ad90SMark F. Adams   PetscFunctionReturn(0);
3755b89ad90SMark F. Adams }
3765b89ad90SMark F. Adams 
377*038e3b61SMark F. Adams #define GAMG_MAXLEVELS 2
3785ef31b24SMark F. Adams #if defined(PETSC_USE_LOG)
3795ef31b24SMark F. Adams PetscLogStage  gamg_stages[20];
3805ef31b24SMark F. Adams #endif
3815b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
3825b89ad90SMark F. Adams /*
3835b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
3845b89ad90SMark F. Adams                     by setting data structures and options.
3855b89ad90SMark F. Adams 
3865b89ad90SMark F. Adams    Input Parameter:
3875b89ad90SMark F. Adams .  pc - the preconditioner context
3885b89ad90SMark F. Adams 
3895b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
3905b89ad90SMark F. Adams 
3915b89ad90SMark F. Adams    Notes:
3925b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
3935b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
3945b89ad90SMark F. Adams */
3955b89ad90SMark F. Adams #undef __FUNCT__
3965b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
397eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
3985b89ad90SMark F. Adams {
3995b89ad90SMark F. Adams   PetscErrorCode  ierr;
400eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4015b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
402eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
403eb07cef2SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, data_sz, Istart, Iend;
404eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4053530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
406*038e3b61SMark F. Adams   PetscBool        isOK;
407587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
408587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
4095ef31b24SMark F. Adams 
4105b89ad90SMark F. Adams   PetscFunctionBegin;
411eb07cef2SMark F. Adams   if( a_pc->setupcalled ) {
412eb07cef2SMark F. Adams     /* no state data in GAMG to destroy */
413eb07cef2SMark F. Adams     ierr = PCReset_MG( a_pc ); CHKERRQ(ierr);
414eb07cef2SMark F. Adams   }
415baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
416baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
4175b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
4185b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
4193530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
420eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
421eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
422eb07cef2SMark F. Adams 
423*038e3b61SMark F. Adams   /* get data of not around */
424*038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
425*038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
426eb07cef2SMark F. Adams   }
427eb07cef2SMark F. Adams   data = pc_gamg->m_data;
428*038e3b61SMark F. Adams   data_sz = pc_gamg->m_data_sz/nloc; assert(data!=0 || data_sz==0);
429*038e3b61SMark F. Adams 
430eb07cef2SMark F. Adams   /* Get A_i and R_i */
431eb07cef2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, data size=%d m_data_sz=%d nloc=%d\n",0,__FUNCT__,0,N,data_sz,pc_gamg->m_data_sz,nloc);
4328f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
433eb07cef2SMark F. Adams         level < GAMG_MAXLEVELS-1 && (level==0 || M/bs>TOP_GRID_LIM) && (npe==1 || nactivepe>1);
4340205a208SMark F. Adams         level++ ) {
4355b89ad90SMark F. Adams     level1 = level + 1;
4365ef31b24SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
437*038e3b61SMark F. Adams     ierr = createProlongation( Aarr[level], data, pc_gamg->m_dim, &bs, &data_sz,
438587fa25dSMark F. Adams                                &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
4395b89ad90SMark F. Adams     CHKERRQ(ierr);
4405ef31b24SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
4415ef31b24SMark F. Adams 
442eb07cef2SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
443baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
444baa4e9faSMark F. Adams     if( isOK ) {
4455ef31b24SMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
446*038e3b61SMark F. Adams       ierr = partitionLevel( Aarr[level], data_sz, bs,
447eb07cef2SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
4483530afc2SMark F. Adams       CHKERRQ(ierr);
4495ef31b24SMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
4503530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
451eb07cef2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, %d data size, %d active pes\n",0,__FUNCT__,level1,N,data_sz,nactivepe);
452baa4e9faSMark F. Adams     }
453baa4e9faSMark F. Adams     else{
454baa4e9faSMark F. Adams       break;
455baa4e9faSMark F. Adams     }
456eb07cef2SMark F. Adams     data = coarse_data;
4575b89ad90SMark F. Adams   }
458eb07cef2SMark F. Adams   ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
459baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
4605b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
4615b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
4625b89ad90SMark F. Adams   fine_level = level;
463eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
4645b89ad90SMark F. Adams 
465fc4362bfSMark F. Adams   /* set default smoothers */
466587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
467587fa25dSMark F. Adams         lidx <= fine_level;
468587fa25dSMark F. Adams         lidx++, level--) {
4695745e0f5SMark F. Adams     PetscReal emax, emin;
4705b89ad90SMark F. Adams     KSP smoother; PC subpc;
471587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
4725b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
473*038e3b61SMark F. Adams     if( emaxs[level] > 0.0 ) emax = emaxs[level];
474587fa25dSMark F. Adams     else{ /* eigen estimate 'emax' */
475587fa25dSMark F. Adams       KSP eksp; Mat Lmat = Aarr[level];
476fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
477*038e3b61SMark F. Adams 
4785745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
4795745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
480fc4362bfSMark F. Adams       {
481fc4362bfSMark F. Adams 	PetscRandom    rctx;
482fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
483fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
484fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
485fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
4865b89ad90SMark F. Adams       }
487fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
488fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
4895745e0f5SMark F. Adams       ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr );
490fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
491fc4362bfSMark F. Adams       ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr); /* should be same as above */
492*038e3b61SMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
493fc4362bfSMark F. Adams       CHKERRQ(ierr);
494fc4362bfSMark F. Adams       ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr);
495fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
496fc4362bfSMark F. Adams 
497fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
498fc4362bfSMark F. Adams       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
499fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
500fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
501fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
502fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
503fc4362bfSMark F. Adams     }
504*038e3b61SMark F. Adams     {
505*038e3b61SMark F. Adams       PetscInt N1, N0, tt;
506*038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
507*038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
508*038e3b61SMark F. Adams       emin = 1.5*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
509*038e3b61SMark F. Adams       PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s max eigen = %e min = %e (N=%d)\n",__FUNCT__,emax,emin,N1/bs);
510*038e3b61SMark F. Adams       emax *= 1.05;
511*038e3b61SMark F. Adams     }
512*038e3b61SMark F. Adams 
513587fa25dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN );
514fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
5155745e0f5SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
5165745e0f5SMark F. Adams     ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
5175745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
5185745e0f5SMark F. Adams   }
5195745e0f5SMark F. Adams   {
5205745e0f5SMark F. Adams     KSP smoother; /* coarse grid */
5215745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
522eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
5235745e0f5SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN );
5245745e0f5SMark F. Adams     CHKERRQ(ierr);
5255745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
526fc4362bfSMark F. Adams   }
527fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
528eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
5295b89ad90SMark F. Adams   {
5305b89ad90SMark F. Adams     PetscBool galerkin;
531eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
5325b89ad90SMark F. Adams     if(galerkin){
5335b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
5345b89ad90SMark F. Adams     }
5355b89ad90SMark F. Adams   }
5365745e0f5SMark F. Adams 
5375745e0f5SMark F. Adams   /* set interpolation between the levels, create timer stages, clean up */
5388f4b7eb5SMark F. Adams   if( PETSC_FALSE ) {
5395ef31b24SMark F. Adams     char str[32];
5405ef31b24SMark F. Adams     sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1);
5415ef31b24SMark F. Adams     PetscLogStageRegister(str, &gamg_stages[fine_level]);
5425ef31b24SMark F. Adams   }
543587fa25dSMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
544587fa25dSMark F. Adams        lidx<fine_level;
545587fa25dSMark F. Adams        lidx++, level--){
546587fa25dSMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
5476c237d78SBarry Smith     if( !PETSC_TRUE ) {
54811e60469SMark F. Adams       PetscViewer viewer; char fname[32];
549*038e3b61SMark F. Adams       sprintf(fname,"Pmat_%d.m",level);
55011e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
5515b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
552*038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
5535b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
5545b89ad90SMark F. Adams     }
555587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
556587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
5578f4b7eb5SMark F. Adams     if( PETSC_FALSE ) {
5585ef31b24SMark F. Adams       char str[32];
559587fa25dSMark F. Adams       sprintf(str,"MG Level %d (%d)",lidx+1,level-1);
560587fa25dSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[level-1]);
561a92563c5SMark F. Adams     }
5625b89ad90SMark F. Adams   }
5635745e0f5SMark F. Adams 
5645b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
565eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
566eb07cef2SMark F. Adams   ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr);
5675b89ad90SMark F. Adams   PetscFunctionReturn(0);
5685b89ad90SMark F. Adams }
5695b89ad90SMark F. Adams 
570eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
5715b89ad90SMark F. Adams /*
5725b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
5735b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
5745b89ad90SMark F. Adams 
5755b89ad90SMark F. Adams    Input Parameter:
5765b89ad90SMark F. Adams .  pc - the preconditioner context
5775b89ad90SMark F. Adams 
5785b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
5795b89ad90SMark F. Adams */
5805b89ad90SMark F. Adams #undef __FUNCT__
5815b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
5825b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
5835b89ad90SMark F. Adams {
5845b89ad90SMark F. Adams   PetscErrorCode  ierr;
5855b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
5865b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
5875b89ad90SMark F. Adams 
5885b89ad90SMark F. Adams   PetscFunctionBegin;
5895b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
5905b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
5915b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
5925b89ad90SMark F. Adams   PetscFunctionReturn(0);
5935b89ad90SMark F. Adams }
5945b89ad90SMark F. Adams 
5955b89ad90SMark F. Adams #undef __FUNCT__
5965b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
5975b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
5985b89ad90SMark F. Adams {
5995b89ad90SMark F. Adams   /* PetscErrorCode  ierr; */
6005b89ad90SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
6015b89ad90SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
6025b89ad90SMark F. Adams   /* MPI_Comm        comm = ((PetscObject)pc)->comm; */
6035b89ad90SMark F. Adams 
6045b89ad90SMark F. Adams   PetscFunctionBegin;
6055b89ad90SMark F. Adams   PetscFunctionReturn(0);
6065b89ad90SMark F. Adams }
6075b89ad90SMark F. Adams 
6085b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
6095b89ad90SMark F. Adams /*
6105b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
6115b89ad90SMark F. Adams 
6125b89ad90SMark F. Adams    Input Parameter:
6135b89ad90SMark F. Adams .  pc - the preconditioner context
6145b89ad90SMark F. Adams 
6155b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
6165b89ad90SMark F. Adams 
6175b89ad90SMark F. Adams   */
6185b89ad90SMark F. Adams  /* MC
6195b89ad90SMark F. Adams      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
6205b89ad90SMark F. Adams        fine grid discretization matrix and coordinates on the fine grid.
6215b89ad90SMark F. Adams 
6225b89ad90SMark F. Adams    Options Database Key:
6235b89ad90SMark F. Adams    Multigrid options(inherited)
6245b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
6255b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
6265b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
6275b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
6285b89ad90SMark F. Adams    GAMG options:
6295b89ad90SMark F. Adams 
6305b89ad90SMark F. Adams    Level: intermediate
6315b89ad90SMark F. Adams   Concepts: multigrid
6325b89ad90SMark F. Adams 
6335b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
6345b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
6355b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
6365b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
6375b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
6385b89ad90SMark F. Adams M */
6395b89ad90SMark F. Adams 
6405b89ad90SMark F. Adams EXTERN_C_BEGIN
6415b89ad90SMark F. Adams #undef __FUNCT__
6425b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
6435b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
6445b89ad90SMark F. Adams {
6455b89ad90SMark F. Adams   PetscErrorCode  ierr;
6465b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
6475b89ad90SMark F. Adams   PC_MG           *mg;
6485ef31b24SMark F. Adams   PetscClassId     cookie;
6495b89ad90SMark F. Adams 
6505b89ad90SMark F. Adams   PetscFunctionBegin;
6515b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
6525b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
6535b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
6545b89ad90SMark F. Adams 
6555b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
6565b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
6575b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
6585b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
6595b89ad90SMark F. Adams 
6605b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
6615b89ad90SMark F. Adams 
6625b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
6635b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
6645b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
6655b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
6665b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
6675b89ad90SMark F. Adams 
6685b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
6695b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
6705b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
6715b89ad90SMark F. Adams 					    PCSetCoordinates_GAMG);CHKERRQ(ierr);
6725ef31b24SMark F. Adams 
6735ef31b24SMark F. Adams   PetscClassIdRegister("GAMG Setup",&cookie);
6745ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-createProl", cookie, &gamg_setup_stages[SET1]);
6755ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-partLevel", cookie, &gamg_setup_stages[SET2]);
6765ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-MIS Graph", cookie, &gamg_setup_stages[SET3]);
6775ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-MIS-Agg", cookie, &gamg_setup_stages[SET4]);
6785ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-growSupp", cookie, &gamg_setup_stages[SET5]);
6795ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-tri-Prol", cookie, &gamg_setup_stages[SET6]);
6805ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-find-prol", cookie, &gamg_setup_stages[FIND_V]);
6815ef31b24SMark F. Adams 
6825b89ad90SMark F. Adams   PetscFunctionReturn(0);
6835b89ad90SMark F. Adams }
6845b89ad90SMark F. Adams EXTERN_C_END
685