xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 2e68589b902172b31894c5911c32a0b8b89b52fd)
1*2e68589bSMark F. Adams /*
2*2e68589bSMark F. Adams  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3*2e68589bSMark F. Adams  */
4*2e68589bSMark F. Adams 
5*2e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
6*2e68589bSMark F. Adams #include <private/kspimpl.h>
7*2e68589bSMark F. Adams 
8*2e68589bSMark F. Adams #include <assert.h>
9*2e68589bSMark F. Adams #include <petscblaslapack.h>
10*2e68589bSMark F. Adams 
11*2e68589bSMark F. Adams typedef struct {
12*2e68589bSMark F. Adams   PetscInt smooths;
13*2e68589bSMark F. Adams } PC_GAMG_AGG;
14*2e68589bSMark F. Adams 
15*2e68589bSMark F. Adams #undef __FUNCT__
16*2e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths"
17*2e68589bSMark F. Adams /*@
18*2e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
19*2e68589bSMark F. Adams 
20*2e68589bSMark F. Adams    Not Collective on PC
21*2e68589bSMark F. Adams 
22*2e68589bSMark F. Adams    Input Parameters:
23*2e68589bSMark F. Adams .  pc - the preconditioner context
24*2e68589bSMark F. Adams 
25*2e68589bSMark F. Adams    Options Database Key:
26*2e68589bSMark F. Adams .  -pc_gamg_agg_nsmooths
27*2e68589bSMark F. Adams 
28*2e68589bSMark F. Adams    Level: intermediate
29*2e68589bSMark F. Adams 
30*2e68589bSMark F. Adams    Concepts: Aggregation AMG preconditioner
31*2e68589bSMark F. Adams 
32*2e68589bSMark F. Adams .seealso: ()
33*2e68589bSMark F. Adams @*/
34*2e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
35*2e68589bSMark F. Adams {
36*2e68589bSMark F. Adams   PetscErrorCode ierr;
37*2e68589bSMark F. Adams 
38*2e68589bSMark F. Adams   PetscFunctionBegin;
39*2e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
40*2e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
41*2e68589bSMark F. Adams   PetscFunctionReturn(0);
42*2e68589bSMark F. Adams }
43*2e68589bSMark F. Adams 
44*2e68589bSMark F. Adams EXTERN_C_BEGIN
45*2e68589bSMark F. Adams #undef __FUNCT__
46*2e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths_GAMG"
47*2e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n)
48*2e68589bSMark F. Adams {
49*2e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
50*2e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
51*2e68589bSMark F. Adams   PC_GAMG_AGG      *pc_gamg_sa = (PC_GAMG_AGG*)pc_gamg->subctx;
52*2e68589bSMark F. Adams 
53*2e68589bSMark F. Adams   PetscFunctionBegin;
54*2e68589bSMark F. Adams   assert(n>=0);
55*2e68589bSMark F. Adams   pc_gamg_sa->smooths = n;
56*2e68589bSMark F. Adams   PetscFunctionReturn(0);
57*2e68589bSMark F. Adams }
58*2e68589bSMark F. Adams EXTERN_C_END
59*2e68589bSMark F. Adams 
60*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
61*2e68589bSMark F. Adams /*
62*2e68589bSMark F. Adams    PCSetFromOptions_GAMG_AGG
63*2e68589bSMark F. Adams 
64*2e68589bSMark F. Adams   Input Parameter:
65*2e68589bSMark F. Adams    . pc -
66*2e68589bSMark F. Adams */
67*2e68589bSMark F. Adams #undef __FUNCT__
68*2e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
69*2e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GAMG_AGG( PC pc )
70*2e68589bSMark F. Adams {
71*2e68589bSMark F. Adams   PetscErrorCode  ierr;
72*2e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
73*2e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
74*2e68589bSMark F. Adams   PC_GAMG_AGG      *pc_gamg_sa = (PC_GAMG_AGG*)pc_gamg->subctx;
75*2e68589bSMark F. Adams   PetscBool        flag;
76*2e68589bSMark F. Adams 
77*2e68589bSMark F. Adams   PetscFunctionBegin;
78*2e68589bSMark F. Adams   /* call base class */
79*2e68589bSMark F. Adams   ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr);
80*2e68589bSMark F. Adams 
81*2e68589bSMark F. Adams   ierr = PetscOptionsHead("GAMG-AGG options"); CHKERRQ(ierr);
82*2e68589bSMark F. Adams   {
83*2e68589bSMark F. Adams     /* -pc_gamg_agg_nsmooths */
84*2e68589bSMark F. Adams     pc_gamg_sa->smooths = 0;
85*2e68589bSMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths",
86*2e68589bSMark F. Adams                            "smoothing steps for smoothed aggregation, usually 1 (0)",
87*2e68589bSMark F. Adams                            "PCGAMGSetNSmooths",
88*2e68589bSMark F. Adams                            pc_gamg_sa->smooths,
89*2e68589bSMark F. Adams                            &pc_gamg_sa->smooths,
90*2e68589bSMark F. Adams                            &flag);
91*2e68589bSMark F. Adams     CHKERRQ(ierr);
92*2e68589bSMark F. Adams   }
93*2e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
94*2e68589bSMark F. Adams 
95*2e68589bSMark F. Adams   if( pc_gamg->verbose ) {
96*2e68589bSMark F. Adams     PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done\n",0,__FUNCT__);
97*2e68589bSMark F. Adams   }
98*2e68589bSMark F. Adams 
99*2e68589bSMark F. Adams   PetscFunctionReturn(0);
100*2e68589bSMark F. Adams }
101*2e68589bSMark F. Adams 
102*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
103*2e68589bSMark F. Adams /*
104*2e68589bSMark F. Adams    PCDestroy_AGG
105*2e68589bSMark F. Adams 
106*2e68589bSMark F. Adams   Input Parameter:
107*2e68589bSMark F. Adams    . pc -
108*2e68589bSMark F. Adams */
109*2e68589bSMark F. Adams #undef __FUNCT__
110*2e68589bSMark F. Adams #define __FUNCT__ "PCDestroy_AGG"
111*2e68589bSMark F. Adams PetscErrorCode PCDestroy_AGG( PC pc )
112*2e68589bSMark F. Adams {
113*2e68589bSMark F. Adams   PetscErrorCode  ierr;
114*2e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
115*2e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
116*2e68589bSMark F. Adams   PC_GAMG_AGG      *pc_gamg_sa = (PC_GAMG_AGG*)pc_gamg->subctx;
117*2e68589bSMark F. Adams 
118*2e68589bSMark F. Adams   PetscFunctionBegin;
119*2e68589bSMark F. Adams   if( pc_gamg_sa ) {
120*2e68589bSMark F. Adams     ierr = PetscFree(pc_gamg_sa);CHKERRQ(ierr);
121*2e68589bSMark F. Adams     pc_gamg_sa = 0;
122*2e68589bSMark F. Adams   }
123*2e68589bSMark F. Adams 
124*2e68589bSMark F. Adams   /* call base class */
125*2e68589bSMark F. Adams   ierr = PCDestroy_GAMG( pc );CHKERRQ(ierr);
126*2e68589bSMark F. Adams 
127*2e68589bSMark F. Adams   PetscFunctionReturn(0);
128*2e68589bSMark F. Adams }
129*2e68589bSMark F. Adams 
130*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
131*2e68589bSMark F. Adams /*
132*2e68589bSMark F. Adams    PCSetCoordinates_AGG
133*2e68589bSMark F. Adams 
134*2e68589bSMark F. Adams    Input Parameter:
135*2e68589bSMark F. Adams    .  pc - the preconditioner context
136*2e68589bSMark F. Adams */
137*2e68589bSMark F. Adams EXTERN_C_BEGIN
138*2e68589bSMark F. Adams #undef __FUNCT__
139*2e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
140*2e68589bSMark F. Adams PetscErrorCode PCSetCoordinates_AGG( PC pc, PetscInt ndm, PetscReal *coords )
141*2e68589bSMark F. Adams {
142*2e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
143*2e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
144*2e68589bSMark F. Adams   PetscErrorCode ierr;
145*2e68589bSMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
146*2e68589bSMark F. Adams   Mat            Amat = pc->pmat;
147*2e68589bSMark F. Adams 
148*2e68589bSMark F. Adams   PetscFunctionBegin;
149*2e68589bSMark F. Adams   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
150*2e68589bSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
151*2e68589bSMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
152*2e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
153*2e68589bSMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
154*2e68589bSMark F. Adams 
155*2e68589bSMark F. Adams   /* SA: null space vectors */
156*2e68589bSMark F. Adams   if( coords && bs==1 ) pc_gamg->data_cols = 1; /* scalar w/ coords and SA (not needed) */
157*2e68589bSMark F. Adams   else if( coords ) pc_gamg->data_cols = (ndm==2 ? 3 : 6); /* elasticity */
158*2e68589bSMark F. Adams   else pc_gamg->data_cols = bs; /* no data, force SA with constant null space vectors */
159*2e68589bSMark F. Adams   pc_gamg->data_rows = bs;
160*2e68589bSMark F. Adams 
161*2e68589bSMark F. Adams   arrsz = nloc*pc_gamg->data_rows*pc_gamg->data_cols;
162*2e68589bSMark F. Adams 
163*2e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
164*2e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
165*2e68589bSMark F. Adams     ierr = PetscFree( pc_gamg->data );  CHKERRQ(ierr);
166*2e68589bSMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr);
167*2e68589bSMark F. Adams   }
168*2e68589bSMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->data[kk] = -999999999.;
169*2e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
170*2e68589bSMark F. Adams   /* copy data in - column oriented */
171*2e68589bSMark F. Adams   for(kk=0;kk<nloc;kk++){
172*2e68589bSMark F. Adams     const PetscInt M = Iend - my0;
173*2e68589bSMark F. Adams     PetscReal *data = &pc_gamg->data[kk*bs];
174*2e68589bSMark F. Adams     if( pc_gamg->data_cols==1 ) *data = 1.0;
175*2e68589bSMark F. Adams     else {
176*2e68589bSMark F. Adams       for(ii=0;ii<bs;ii++)
177*2e68589bSMark F. Adams         for(jj=0;jj<bs;jj++)
178*2e68589bSMark F. Adams           if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
179*2e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
180*2e68589bSMark F. Adams       if( coords ) {
181*2e68589bSMark F. Adams         if( ndm == 2 ){ /* rotational modes */
182*2e68589bSMark F. Adams           data += 2*M;
183*2e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
184*2e68589bSMark F. Adams           data[1] =  coords[2*kk];
185*2e68589bSMark F. Adams         }
186*2e68589bSMark F. Adams         else {
187*2e68589bSMark F. Adams           data += 3*M;
188*2e68589bSMark F. Adams           data[0] = 0.0;               data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
189*2e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  coords[3*kk];
190*2e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
191*2e68589bSMark F. Adams         }
192*2e68589bSMark F. Adams       }
193*2e68589bSMark F. Adams     }
194*2e68589bSMark F. Adams   }
195*2e68589bSMark F. Adams   assert(pc_gamg->data[arrsz] == -99.);
196*2e68589bSMark F. Adams 
197*2e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
198*2e68589bSMark F. Adams 
199*2e68589bSMark F. Adams   PetscFunctionReturn(0);
200*2e68589bSMark F. Adams }
201*2e68589bSMark F. Adams EXTERN_C_END
202*2e68589bSMark F. Adams 
203*2e68589bSMark F. Adams 
204*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
205*2e68589bSMark F. Adams /*
206*2e68589bSMark F. Adams    PCSetData_AGG
207*2e68589bSMark F. Adams 
208*2e68589bSMark F. Adams   Input Parameter:
209*2e68589bSMark F. Adams    . pc -
210*2e68589bSMark F. Adams */
211*2e68589bSMark F. Adams #undef __FUNCT__
212*2e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
213*2e68589bSMark F. Adams PetscErrorCode PCSetData_AGG( PC pc )
214*2e68589bSMark F. Adams {
215*2e68589bSMark F. Adams   PetscErrorCode  ierr;
216*2e68589bSMark F. Adams   PetscFunctionBegin;
217*2e68589bSMark F. Adams   ierr = PCSetCoordinates_AGG( pc, -1, PETSC_NULL ); CHKERRQ(ierr);
218*2e68589bSMark F. Adams   PetscFunctionReturn(0);
219*2e68589bSMark F. Adams }
220*2e68589bSMark F. Adams 
221*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
222*2e68589bSMark F. Adams /*
223*2e68589bSMark F. Adams    PCCreateGAMG_AGG
224*2e68589bSMark F. Adams 
225*2e68589bSMark F. Adams   Input Parameter:
226*2e68589bSMark F. Adams    . pc -
227*2e68589bSMark F. Adams */
228*2e68589bSMark F. Adams #undef __FUNCT__
229*2e68589bSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
230*2e68589bSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG( PC pc )
231*2e68589bSMark F. Adams {
232*2e68589bSMark F. Adams   PetscErrorCode  ierr;
233*2e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
234*2e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
235*2e68589bSMark F. Adams   PC_GAMG_AGG      *pc_gamg_sa;
236*2e68589bSMark F. Adams 
237*2e68589bSMark F. Adams   PetscFunctionBegin;
238*2e68589bSMark F. Adams   /* create sub context for SA */
239*2e68589bSMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_sa ); CHKERRQ(ierr);
240*2e68589bSMark F. Adams   assert(!pc_gamg->subctx);
241*2e68589bSMark F. Adams   pc_gamg->subctx = pc_gamg_sa;
242*2e68589bSMark F. Adams 
243*2e68589bSMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
244*2e68589bSMark F. Adams   pc->ops->destroy        = PCDestroy_AGG;
245*2e68589bSMark F. Adams   /* reset does not do anything; setup not virtual */
246*2e68589bSMark F. Adams 
247*2e68589bSMark F. Adams   /* set internal function pointers */
248*2e68589bSMark F. Adams   pc_gamg->createprolongator = PCGAMGcreateProl_AGG;
249*2e68589bSMark F. Adams   pc_gamg->createdefaultdata = PCSetData_AGG;
250*2e68589bSMark F. Adams 
251*2e68589bSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
252*2e68589bSMark F. Adams                                             "PCSetCoordinates_C",
253*2e68589bSMark F. Adams                                             "PCSetCoordinates_AGG",
254*2e68589bSMark F. Adams                                             PCSetCoordinates_AGG);
255*2e68589bSMark F. Adams   PetscFunctionReturn(0);
256*2e68589bSMark F. Adams }
257*2e68589bSMark F. Adams 
258*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
259*2e68589bSMark F. Adams /*
260*2e68589bSMark F. Adams  formProl0
261*2e68589bSMark F. Adams 
262*2e68589bSMark F. Adams    Input Parameter:
263*2e68589bSMark F. Adams    . selected - list of selected local ID, includes selected ghosts
264*2e68589bSMark F. Adams    . locals_llist - linked list with aggregates
265*2e68589bSMark F. Adams    . bs - block size
266*2e68589bSMark F. Adams    . nSAvec - num columns of new P
267*2e68589bSMark F. Adams    . my0crs - global index of locals
268*2e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
269*2e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
270*2e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
271*2e68589bSMark F. Adams   Output Parameter:
272*2e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
273*2e68589bSMark F. Adams    . a_Prol - prolongation operator
274*2e68589bSMark F. Adams */
275*2e68589bSMark F. Adams #undef __FUNCT__
276*2e68589bSMark F. Adams #define __FUNCT__ "formProl0"
277*2e68589bSMark F. Adams PetscErrorCode formProl0(IS selected, /* list of selected local ID, includes selected ghosts */
278*2e68589bSMark F. Adams                          IS locals_llist, /* linked list from selected vertices of aggregate unselected vertices */
279*2e68589bSMark F. Adams                          const PetscInt bs,
280*2e68589bSMark F. Adams                          const PetscInt nSAvec,
281*2e68589bSMark F. Adams                          const PetscInt my0crs,
282*2e68589bSMark F. Adams 			 const PetscInt data_stride,
283*2e68589bSMark F. Adams 			 PetscReal data_in[],
284*2e68589bSMark F. Adams                          const PetscInt flid_fgid[],
285*2e68589bSMark F. Adams                          PetscReal **a_data_out,
286*2e68589bSMark F. Adams                          Mat a_Prol /* prolongation operator (output)*/
287*2e68589bSMark F. Adams                          )
288*2e68589bSMark F. Adams {
289*2e68589bSMark F. Adams   PetscErrorCode ierr;
290*2e68589bSMark F. Adams   PetscInt  Istart,Iend,nFineLoc,clid,flid,aggID,kk,jj,ii,nLocalSelected,ndone,nSelected;
291*2e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)a_Prol)->comm;
292*2e68589bSMark F. Adams   PetscMPIInt    mype, npe;
293*2e68589bSMark F. Adams   const PetscInt *selected_idx,*llist_idx;
294*2e68589bSMark F. Adams   PetscReal      *out_data;
295*2e68589bSMark F. Adams 
296*2e68589bSMark F. Adams   PetscFunctionBegin;
297*2e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
298*2e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
299*2e68589bSMark F. Adams   ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend );    CHKERRQ(ierr);
300*2e68589bSMark F. Adams   nFineLoc = (Iend-Istart)/bs; assert((Iend-Istart)%bs==0);
301*2e68589bSMark F. Adams 
302*2e68589bSMark F. Adams   ierr = ISGetLocalSize( selected, &nSelected );        CHKERRQ(ierr);
303*2e68589bSMark F. Adams   ierr = ISGetIndices( selected, &selected_idx );       CHKERRQ(ierr);
304*2e68589bSMark F. Adams   for(kk=0,nLocalSelected=0;kk<nSelected;kk++){
305*2e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
306*2e68589bSMark F. Adams     if(lid<nFineLoc) nLocalSelected++;
307*2e68589bSMark F. Adams   }
308*2e68589bSMark F. Adams 
309*2e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
310*2e68589bSMark F. Adams #define DATA_OUT_STRIDE (nLocalSelected*nSAvec)
311*2e68589bSMark F. Adams   ierr = PetscMalloc( (DATA_OUT_STRIDE*nSAvec+1)*sizeof(PetscReal), &out_data ); CHKERRQ(ierr);
312*2e68589bSMark F. Adams   for(ii=0;ii<DATA_OUT_STRIDE*nSAvec+1;ii++) out_data[ii]=1.e300;
313*2e68589bSMark F. Adams   *a_data_out = out_data; /* output - stride nLocalSelected*nSAvec */
314*2e68589bSMark F. Adams 
315*2e68589bSMark F. Adams   /* find points and set prolongation */
316*2e68589bSMark F. Adams   ndone = 0;
317*2e68589bSMark F. Adams   ierr = ISGetIndices( locals_llist, &llist_idx );      CHKERRQ(ierr);
318*2e68589bSMark F. Adams   for( clid = 0 ; clid < nLocalSelected ; clid++ ){
319*2e68589bSMark F. Adams     PetscInt cgid = my0crs + clid, cids[100];
320*2e68589bSMark F. Adams 
321*2e68589bSMark F. Adams     /* count agg */
322*2e68589bSMark F. Adams     aggID = 0;
323*2e68589bSMark F. Adams     flid = selected_idx[clid]; assert(flid != -1);
324*2e68589bSMark F. Adams     do{
325*2e68589bSMark F. Adams       aggID++;
326*2e68589bSMark F. Adams     } while( (flid=llist_idx[flid]) != -1 );
327*2e68589bSMark F. Adams 
328*2e68589bSMark F. Adams     /* get block */
329*2e68589bSMark F. Adams     {
330*2e68589bSMark F. Adams       PetscBLASInt   asz=aggID,M=asz*bs,N=nSAvec,INFO;
331*2e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs;
332*2e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
333*2e68589bSMark F. Adams       PetscInt       *fids;
334*2e68589bSMark F. Adams 
335*2e68589bSMark F. Adams       ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr);
336*2e68589bSMark F. Adams       ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr);
337*2e68589bSMark F. Adams       ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr);
338*2e68589bSMark F. Adams       ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr);
339*2e68589bSMark F. Adams       ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr);
340*2e68589bSMark F. Adams 
341*2e68589bSMark F. Adams       flid = selected_idx[clid];
342*2e68589bSMark F. Adams       aggID = 0;
343*2e68589bSMark F. Adams       do{
344*2e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
345*2e68589bSMark F. Adams         PetscReal *data = &data_in[flid*bs];
346*2e68589bSMark F. Adams         for( kk = ii = 0; ii < bs ; ii++ ) {
347*2e68589bSMark F. Adams           for( jj = 0; jj < N ; jj++ ) {
348*2e68589bSMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = data[jj*data_stride + ii];
349*2e68589bSMark F. Adams           }
350*2e68589bSMark F. Adams         }
351*2e68589bSMark F. Adams 
352*2e68589bSMark F. Adams         /* set fine IDs */
353*2e68589bSMark F. Adams         for(kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
354*2e68589bSMark F. Adams 
355*2e68589bSMark F. Adams         aggID++;
356*2e68589bSMark F. Adams       }while( (flid=llist_idx[flid]) != -1 );
357*2e68589bSMark F. Adams 
358*2e68589bSMark F. Adams       /* pad with zeros */
359*2e68589bSMark F. Adams       for( ii = asz*bs; ii < Mdata ; ii++ ) {
360*2e68589bSMark F. Adams 	for( jj = 0; jj < N ; jj++, kk++ ) {
361*2e68589bSMark F. Adams 	  qqc[jj*Mdata + ii] = .0;
362*2e68589bSMark F. Adams 	}
363*2e68589bSMark F. Adams       }
364*2e68589bSMark F. Adams 
365*2e68589bSMark F. Adams       ndone += aggID;
366*2e68589bSMark F. Adams       /* QR */
367*2e68589bSMark F. Adams       LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
368*2e68589bSMark F. Adams       if( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRS error");
369*2e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
370*2e68589bSMark F. Adams       {
371*2e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
372*2e68589bSMark F. Adams         for( jj = 0; jj < nSAvec ; jj++ ) {
373*2e68589bSMark F. Adams           for( ii = 0; ii < nSAvec ; ii++ ) {
374*2e68589bSMark F. Adams             assert(data[jj*DATA_OUT_STRIDE + ii] == 1.e300);
375*2e68589bSMark F. Adams             if( ii <= jj ) data[jj*DATA_OUT_STRIDE + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
376*2e68589bSMark F. Adams 	    else data[jj*DATA_OUT_STRIDE + ii] = 0.;
377*2e68589bSMark F. Adams           }
378*2e68589bSMark F. Adams         }
379*2e68589bSMark F. Adams       }
380*2e68589bSMark F. Adams 
381*2e68589bSMark F. Adams       /* get Q - row oriented */
382*2e68589bSMark F. Adams       LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
383*2e68589bSMark F. Adams       if( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR error arg %d",-INFO);
384*2e68589bSMark F. Adams 
385*2e68589bSMark F. Adams       for( ii = 0 ; ii < M ; ii++ ){
386*2e68589bSMark F. Adams         for( jj = 0 ; jj < N ; jj++ ) {
387*2e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
388*2e68589bSMark F. Adams         }
389*2e68589bSMark F. Adams       }
390*2e68589bSMark F. Adams 
391*2e68589bSMark F. Adams       /* add diagonal block of P0 */
392*2e68589bSMark F. Adams       for(kk=0;kk<N;kk++) cids[kk] = N*cgid + kk; /* global col IDs in P0 */
393*2e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr);
394*2e68589bSMark F. Adams 
395*2e68589bSMark F. Adams       ierr = PetscFree( qqc );  CHKERRQ(ierr);
396*2e68589bSMark F. Adams       ierr = PetscFree( qqr );  CHKERRQ(ierr);
397*2e68589bSMark F. Adams       ierr = PetscFree( TAU );  CHKERRQ(ierr);
398*2e68589bSMark F. Adams       ierr = PetscFree( WORK );  CHKERRQ(ierr);
399*2e68589bSMark F. Adams       ierr = PetscFree( fids );  CHKERRQ(ierr);
400*2e68589bSMark F. Adams     } /* scoping */
401*2e68589bSMark F. Adams   } /* for all coarse nodes */
402*2e68589bSMark F. Adams   assert(out_data[nSAvec*DATA_OUT_STRIDE]==1.e300);
403*2e68589bSMark F. Adams 
404*2e68589bSMark F. Adams /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); /\* synchronous version *\/ */
405*2e68589bSMark F. Adams /* MatGetSize( a_Prol, &kk, &jj ); */
406*2e68589bSMark F. Adams /* PetscPrintf(PETSC_COMM_WORLD," **** [%d]%s %d total done, N=%d (%d local done)\n",mype,__FUNCT__,ii,kk/bs,ndone); */
407*2e68589bSMark F. Adams 
408*2e68589bSMark F. Adams   ierr = ISRestoreIndices( selected, &selected_idx );     CHKERRQ(ierr);
409*2e68589bSMark F. Adams   ierr = ISRestoreIndices( locals_llist, &llist_idx );     CHKERRQ(ierr);
410*2e68589bSMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411*2e68589bSMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
412*2e68589bSMark F. Adams 
413*2e68589bSMark F. Adams   PetscFunctionReturn(0);
414*2e68589bSMark F. Adams }
415*2e68589bSMark F. Adams 
416*2e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
417*2e68589bSMark F. Adams /*
418*2e68589bSMark F. Adams    PCGAMGcreateProl_AGG
419*2e68589bSMark F. Adams 
420*2e68589bSMark F. Adams   Input Parameter:
421*2e68589bSMark F. Adams    . pc - this
422*2e68589bSMark F. Adams    . Amat - matrix on this fine level
423*2e68589bSMark F. Adams    . data[nloc*data_sz(in)]
424*2e68589bSMark F. Adams   Output Parameter:
425*2e68589bSMark F. Adams    . a_P_out - prolongation operator to the next level
426*2e68589bSMark F. Adams    . a_data_out - data of coarse grid points (num local columns in 'a_P_out')
427*2e68589bSMark F. Adams */
428*2e68589bSMark F. Adams #undef __FUNCT__
429*2e68589bSMark F. Adams #define __FUNCT__ "PCGAMGcreateProl_AGG"
430*2e68589bSMark F. Adams PetscErrorCode PCGAMGcreateProl_AGG( PC pc,
431*2e68589bSMark F. Adams                                     const Mat Amat,
432*2e68589bSMark F. Adams                                     const PetscReal data[],
433*2e68589bSMark F. Adams                                     Mat *a_P_out,
434*2e68589bSMark F. Adams                                     PetscReal **a_data_out
435*2e68589bSMark F. Adams                                    )
436*2e68589bSMark F. Adams {
437*2e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
438*2e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
439*2e68589bSMark F. Adams   PC_GAMG_AGG    *pc_gamg_sa = (PC_GAMG_AGG*)pc_gamg->subctx;
440*2e68589bSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
441*2e68589bSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cols;
442*2e68589bSMark F. Adams   PetscErrorCode ierr;
443*2e68589bSMark F. Adams   PetscInt       Istart,Iend,nloc,jj,kk,my0,nLocalSelected,NN,MM,bs_in;
444*2e68589bSMark F. Adams   Mat            Prol, Gmat, AuxMat;
445*2e68589bSMark F. Adams   PetscMPIInt    mype, npe;
446*2e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
447*2e68589bSMark F. Adams   IS             permIS, llist_1, selected_1;
448*2e68589bSMark F. Adams   const PetscInt *selected_idx,col_bs=data_cols;
449*2e68589bSMark F. Adams 
450*2e68589bSMark F. Adams   PetscFunctionBegin;
451*2e68589bSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
452*2e68589bSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
453*2e68589bSMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
454*2e68589bSMark F. Adams   ierr = MatGetSize( Amat, &MM, &NN );  CHKERRQ(ierr);
455*2e68589bSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs_in ); CHKERRQ( ierr );
456*2e68589bSMark F. Adams   nloc = (Iend-Istart)/bs_in; my0 = Istart/bs_in; assert((Iend-Istart)%bs_in==0);
457*2e68589bSMark F. Adams 
458*2e68589bSMark F. Adams   ierr = createGraph( pc, Amat, &Gmat, &AuxMat, &permIS ); CHKERRQ(ierr);
459*2e68589bSMark F. Adams 
460*2e68589bSMark F. Adams   /* SELECT COARSE POINTS */
461*2e68589bSMark F. Adams #if defined PETSC_USE_LOG
462*2e68589bSMark F. Adams   ierr = PetscLogEventBegin(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
463*2e68589bSMark F. Adams #endif
464*2e68589bSMark F. Adams 
465*2e68589bSMark F. Adams   ierr = maxIndSetAgg( permIS, Gmat, AuxMat, PETSC_TRUE, &selected_1, &llist_1 );
466*2e68589bSMark F. Adams   CHKERRQ(ierr);
467*2e68589bSMark F. Adams   ierr = MatDestroy( &AuxMat );  CHKERRQ(ierr);
468*2e68589bSMark F. Adams 
469*2e68589bSMark F. Adams #if defined PETSC_USE_LOG
470*2e68589bSMark F. Adams   ierr = PetscLogEventEnd(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
471*2e68589bSMark F. Adams #endif
472*2e68589bSMark F. Adams   ierr = ISDestroy(&permIS); CHKERRQ(ierr);
473*2e68589bSMark F. Adams 
474*2e68589bSMark F. Adams   /* get 'nLocalSelected' */
475*2e68589bSMark F. Adams   ierr = ISGetLocalSize( selected_1, &NN );        CHKERRQ(ierr);
476*2e68589bSMark F. Adams   ierr = ISGetIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
477*2e68589bSMark F. Adams   for(kk=0,nLocalSelected=0;kk<NN;kk++){
478*2e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
479*2e68589bSMark F. Adams     if(lid<nloc) nLocalSelected++;
480*2e68589bSMark F. Adams   }
481*2e68589bSMark F. Adams   ierr = ISRestoreIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
482*2e68589bSMark F. Adams 
483*2e68589bSMark F. Adams   /* create prolongator, create P matrix */
484*2e68589bSMark F. Adams   ierr = MatCreateMPIAIJ(wcomm,
485*2e68589bSMark F. Adams                          nloc*bs_in, nLocalSelected*col_bs,
486*2e68589bSMark F. Adams                          PETSC_DETERMINE, PETSC_DETERMINE,
487*2e68589bSMark F. Adams                          data_cols, PETSC_NULL, data_cols, PETSC_NULL,
488*2e68589bSMark F. Adams                          &Prol );
489*2e68589bSMark F. Adams   CHKERRQ(ierr);
490*2e68589bSMark F. Adams 
491*2e68589bSMark F. Adams   /* can get all points "removed" */
492*2e68589bSMark F. Adams   ierr =  MatGetSize( Prol, &kk, &NN ); CHKERRQ(ierr);
493*2e68589bSMark F. Adams   if( NN==0 ) {
494*2e68589bSMark F. Adams     if( verbose ) {
495*2e68589bSMark F. Adams       PetscPrintf(PETSC_COMM_WORLD,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__);
496*2e68589bSMark F. Adams     }
497*2e68589bSMark F. Adams     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
498*2e68589bSMark F. Adams     ierr = ISDestroy( &llist_1 ); CHKERRQ(ierr);
499*2e68589bSMark F. Adams     ierr = ISDestroy( &selected_1 ); CHKERRQ(ierr);
500*2e68589bSMark F. Adams     ierr = MatDestroy( &Gmat );  CHKERRQ(ierr);
501*2e68589bSMark F. Adams     *a_P_out = PETSC_NULL;  /* out */
502*2e68589bSMark F. Adams     PetscFunctionReturn(0);
503*2e68589bSMark F. Adams   }
504*2e68589bSMark F. Adams 
505*2e68589bSMark F. Adams   { /* SA */
506*2e68589bSMark F. Adams     PetscReal *data_w_ghost;
507*2e68589bSMark F. Adams     PetscInt  myCrs0, nbnodes=0, *flid_fgid;
508*2e68589bSMark F. Adams 
509*2e68589bSMark F. Adams     /* create global vector of data in 'data_w_ghost' */
510*2e68589bSMark F. Adams #if defined PETSC_USE_LOG
511*2e68589bSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
512*2e68589bSMark F. Adams #endif
513*2e68589bSMark F. Adams     if (npe > 1) {
514*2e68589bSMark F. Adams       PetscReal *tmp_gdata,*tmp_ldata,*tp2;
515*2e68589bSMark F. Adams 
516*2e68589bSMark F. Adams       ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr);
517*2e68589bSMark F. Adams       for( jj = 0 ; jj < data_cols ; jj++ ){
518*2e68589bSMark F. Adams         for( kk = 0 ; kk < bs_in ; kk++) {
519*2e68589bSMark F. Adams           PetscInt ii,nnodes;
520*2e68589bSMark F. Adams           const PetscReal *tp = data + jj*bs_in*nloc + kk;
521*2e68589bSMark F. Adams           for( ii = 0 ; ii < nloc ; ii++, tp += bs_in ){
522*2e68589bSMark F. Adams             tmp_ldata[ii] = *tp;
523*2e68589bSMark F. Adams           }
524*2e68589bSMark F. Adams           ierr = getDataWithGhosts( Gmat, 1, tmp_ldata, &nnodes, &tmp_gdata );
525*2e68589bSMark F. Adams           CHKERRQ(ierr);
526*2e68589bSMark F. Adams           if(jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
527*2e68589bSMark F. Adams             ierr = PetscMalloc( nnodes*bs_in*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr);
528*2e68589bSMark F. Adams             nbnodes = bs_in*nnodes;
529*2e68589bSMark F. Adams           }
530*2e68589bSMark F. Adams           tp2 = data_w_ghost + jj*bs_in*nnodes + kk;
531*2e68589bSMark F. Adams           for( ii = 0 ; ii < nnodes ; ii++, tp2 += bs_in ){
532*2e68589bSMark F. Adams             *tp2 = tmp_gdata[ii];
533*2e68589bSMark F. Adams           }
534*2e68589bSMark F. Adams           ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr);
535*2e68589bSMark F. Adams         }
536*2e68589bSMark F. Adams       }
537*2e68589bSMark F. Adams       ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr);
538*2e68589bSMark F. Adams     }
539*2e68589bSMark F. Adams     else {
540*2e68589bSMark F. Adams       nbnodes = bs_in*nloc;
541*2e68589bSMark F. Adams       data_w_ghost = (PetscReal*)data;
542*2e68589bSMark F. Adams     }
543*2e68589bSMark F. Adams 
544*2e68589bSMark F. Adams     /* scan my coarse zero gid */
545*2e68589bSMark F. Adams     MPI_Scan( &nLocalSelected, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm );
546*2e68589bSMark F. Adams     myCrs0 -= nLocalSelected;
547*2e68589bSMark F. Adams 
548*2e68589bSMark F. Adams     /* get P0 */
549*2e68589bSMark F. Adams     if( npe > 1 ){
550*2e68589bSMark F. Adams       PetscReal *fid_glid_loc,*fiddata;
551*2e68589bSMark F. Adams       PetscInt nnodes;
552*2e68589bSMark F. Adams 
553*2e68589bSMark F. Adams       ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr);
554*2e68589bSMark F. Adams       for(kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
555*2e68589bSMark F. Adams       ierr = getDataWithGhosts(Gmat, 1, fid_glid_loc, &nnodes, &fiddata);
556*2e68589bSMark F. Adams       CHKERRQ(ierr);
557*2e68589bSMark F. Adams       ierr = PetscMalloc( nnodes*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
558*2e68589bSMark F. Adams       for(kk=0;kk<nnodes;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
559*2e68589bSMark F. Adams       ierr = PetscFree( fiddata ); CHKERRQ(ierr);
560*2e68589bSMark F. Adams       assert(nnodes==nbnodes/bs_in);
561*2e68589bSMark F. Adams       ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr);
562*2e68589bSMark F. Adams     }
563*2e68589bSMark F. Adams     else {
564*2e68589bSMark F. Adams       ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
565*2e68589bSMark F. Adams       for(kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk;
566*2e68589bSMark F. Adams     }
567*2e68589bSMark F. Adams     ierr = MatDestroy( &Gmat );  CHKERRQ(ierr);
568*2e68589bSMark F. Adams #if defined PETSC_USE_LOG
569*2e68589bSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
570*2e68589bSMark F. Adams #endif
571*2e68589bSMark F. Adams 
572*2e68589bSMark F. Adams     ierr = formProl0(selected_1,llist_1,bs_in,data_cols,myCrs0,nbnodes,
573*2e68589bSMark F. Adams                      data_w_ghost,flid_fgid,a_data_out,Prol);
574*2e68589bSMark F. Adams     CHKERRQ(ierr);
575*2e68589bSMark F. Adams 
576*2e68589bSMark F. Adams     if (npe > 1) ierr = PetscFree( data_w_ghost );      CHKERRQ(ierr);
577*2e68589bSMark F. Adams     ierr = PetscFree( flid_fgid ); CHKERRQ(ierr);
578*2e68589bSMark F. Adams 
579*2e68589bSMark F. Adams     /* smooth P0 */
580*2e68589bSMark F. Adams     for( jj = 0 ; jj < pc_gamg_sa->smooths ; jj++ ){
581*2e68589bSMark F. Adams       Mat tMat;
582*2e68589bSMark F. Adams       Vec diag;
583*2e68589bSMark F. Adams       PetscReal alpha, emax, emin;
584*2e68589bSMark F. Adams #if defined PETSC_USE_LOG
585*2e68589bSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
586*2e68589bSMark F. Adams #endif
587*2e68589bSMark F. Adams       if( jj == 0 ) {
588*2e68589bSMark F. Adams         KSP eksp;
589*2e68589bSMark F. Adams         Vec bb, xx;
590*2e68589bSMark F. Adams         PC pc;
591*2e68589bSMark F. Adams         ierr = MatGetVecs( Amat, &bb, 0 );         CHKERRQ(ierr);
592*2e68589bSMark F. Adams         ierr = MatGetVecs( Amat, &xx, 0 );         CHKERRQ(ierr);
593*2e68589bSMark F. Adams         {
594*2e68589bSMark F. Adams           PetscRandom    rctx;
595*2e68589bSMark F. Adams           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
596*2e68589bSMark F. Adams           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
597*2e68589bSMark F. Adams           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
598*2e68589bSMark F. Adams           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
599*2e68589bSMark F. Adams         }
600*2e68589bSMark F. Adams         ierr = KSPCreate(wcomm,&eksp);                            CHKERRQ(ierr);
601*2e68589bSMark F. Adams         /* ierr = KSPSetType( eksp, KSPCG );                         CHKERRQ(ierr); */
602*2e68589bSMark F. Adams         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
603*2e68589bSMark F. Adams         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
604*2e68589bSMark F. Adams         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE );    CHKERRQ(ierr);
605*2e68589bSMark F. Adams         ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN );
606*2e68589bSMark F. Adams         CHKERRQ( ierr );
607*2e68589bSMark F. Adams         ierr = KSPGetPC( eksp, &pc );                              CHKERRQ( ierr );
608*2e68589bSMark F. Adams         ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr);  /* smoother */
609*2e68589bSMark F. Adams         ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);
610*2e68589bSMark F. Adams         CHKERRQ(ierr);
611*2e68589bSMark F. Adams         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
612*2e68589bSMark F. Adams         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE );        CHKERRQ(ierr);
613*2e68589bSMark F. Adams 
614*2e68589bSMark F. Adams         /* solve - keep stuff out of logging */
615*2e68589bSMark F. Adams         ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
616*2e68589bSMark F. Adams         ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
617*2e68589bSMark F. Adams         ierr = KSPSolve( eksp, bb, xx );                              CHKERRQ(ierr);
618*2e68589bSMark F. Adams         ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
619*2e68589bSMark F. Adams         ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
620*2e68589bSMark F. Adams 
621*2e68589bSMark F. Adams         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
622*2e68589bSMark F. Adams         if( verbose ) {
623*2e68589bSMark F. Adams           PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
624*2e68589bSMark F. Adams                       __FUNCT__,emax,emin,PCJACOBI);
625*2e68589bSMark F. Adams         }
626*2e68589bSMark F. Adams         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
627*2e68589bSMark F. Adams         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
628*2e68589bSMark F. Adams         ierr = KSPDestroy( &eksp );     CHKERRQ(ierr);
629*2e68589bSMark F. Adams 
630*2e68589bSMark F. Adams         if( pc_gamg->emax_id == -1 ) {
631*2e68589bSMark F. Adams           ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id );
632*2e68589bSMark F. Adams           assert(pc_gamg->emax_id != -1 );
633*2e68589bSMark F. Adams         }
634*2e68589bSMark F. Adams         ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr);
635*2e68589bSMark F. Adams       }
636*2e68589bSMark F. Adams 
637*2e68589bSMark F. Adams       /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
638*2e68589bSMark F. Adams       ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat );   CHKERRQ(ierr);
639*2e68589bSMark F. Adams       ierr = MatGetVecs( Amat, &diag, 0 );    CHKERRQ(ierr);
640*2e68589bSMark F. Adams       ierr = MatGetDiagonal( Amat, diag );    CHKERRQ(ierr); /* effectively PCJACOBI */
641*2e68589bSMark F. Adams       ierr = VecReciprocal( diag );         CHKERRQ(ierr);
642*2e68589bSMark F. Adams       ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr);
643*2e68589bSMark F. Adams       ierr = VecDestroy( &diag );           CHKERRQ(ierr);
644*2e68589bSMark F. Adams       alpha = -1.5/emax;
645*2e68589bSMark F. Adams       ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN );           CHKERRQ(ierr);
646*2e68589bSMark F. Adams       ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
647*2e68589bSMark F. Adams       Prol = tMat;
648*2e68589bSMark F. Adams #if defined PETSC_USE_LOG
649*2e68589bSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
650*2e68589bSMark F. Adams #endif
651*2e68589bSMark F. Adams     }
652*2e68589bSMark F. Adams   } /* scoping - SA code */
653*2e68589bSMark F. Adams 
654*2e68589bSMark F. Adams   /* attach block size of columns */
655*2e68589bSMark F. Adams   if( pc_gamg->col_bs_id == -1 ) {
656*2e68589bSMark F. Adams     ierr = PetscObjectComposedDataRegister( &pc_gamg->col_bs_id );
657*2e68589bSMark F. Adams     assert(pc_gamg->col_bs_id != -1 );
658*2e68589bSMark F. Adams   }
659*2e68589bSMark F. Adams   ierr = PetscObjectComposedDataSetInt( (PetscObject)Prol, pc_gamg->col_bs_id, data_cols ); CHKERRQ(ierr);
660*2e68589bSMark F. Adams 
661*2e68589bSMark F. Adams   *a_P_out = Prol;  /* out */
662*2e68589bSMark F. Adams 
663*2e68589bSMark F. Adams   ierr = ISDestroy( &llist_1 ); CHKERRQ(ierr);
664*2e68589bSMark F. Adams   ierr = ISDestroy( &selected_1 ); CHKERRQ(ierr);
665*2e68589bSMark F. Adams 
666*2e68589bSMark F. Adams   PetscFunctionReturn(0);
667*2e68589bSMark F. Adams }
668