xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 5582bec17b042ee1c19decb9d092fe7ad142b625)
1 /*
2    Provides an interface to the ML 3.0 smoothed Aggregation
3 */
4 #include "src/ksp/pc/pcimpl.h"   /*I "petscpc.h" I*/
5 #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
6 #include "src/mat/impls/aij/seq/aij.h"
7 #include "src/mat/impls/aij/mpi/mpiaij.h"
8 EXTERN_C_BEGIN
9 #include <math.h>
10 #include "ml_include.h"
11 EXTERN_C_END
12 
13 /* The context (data structure) at each grid level */
14 typedef struct {
15   Vec        x,b,r;            /* global vectors */
16   Mat        A,P,R;
17   KSP        ksp;
18 } GridCtx;
19 
20 /* The context used to input PETSc matrix into ML at fine grid */
21 typedef struct {
22   Mat          A,Aloc;
23   ML_Operator  *mlmat;
24   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
25   PetscInt     rlen_max,*cols; /* used by MatConvert_ML_SeqAIJ() */
26   PetscScalar  *vals;          /* used by MatConvert_ML_SeqAIJ() */
27 } FineGridCtx;
28 
29 /* The context associates a ML matrix with a PETSc shell matrix */
30 typedef struct {
31   Mat          A;       /* PETSc shell matrix associated with mlmat */
32   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
33   Vec          y;
34 } Mat_MLShell;
35 
36 /* Private context for the ML preconditioner */
37 typedef struct {
38   ML           *ml_object;
39   ML_Aggregate *agg_object;
40   GridCtx      *gridctx;
41   FineGridCtx  *PetscMLdata;
42   PetscInt     fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme;
43   PetscReal    Threshold,DampingFactor;
44   PetscTruth   SpectralNormScheme_Anorm;
45   PetscMPIInt  size;
46 
47   PetscErrorCode (*PCSetUp)(PC);
48   PetscErrorCode (*PCDestroy)(PC);
49 
50 } PC_ML;
51 extern int PetscML_getrow(void *ML_data,int N_requested_rows,int requested_rows[],
52             int allocated_space,int columns[],double values[],int row_lengths[]);
53 extern int PetscML_matvec(void *ML_data, int in_length, double p[], int out_length,double ap[]);
54 extern int PetscML_comm(double x[], void *ML_data);
55 extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
56 extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
57 extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,const MatType,Mat*);
58 extern PetscErrorCode MatDestroy_ML(Mat);
59 extern PetscErrorCode MatConvert_ML_SeqAIJ(FineGridCtx*,Mat*);
60 extern PetscErrorCode MatConvert_ML_SHELL(ML_Operator*,Mat*);
61 
62 /* -------------------------------------------------------------------------- */
63 /*
64    PCSetUp_ML - Prepares for the use of the ML preconditioner
65                     by setting data structures and options.
66 
67    Input Parameter:
68 .  pc - the preconditioner context
69 
70    Application Interface Routine: PCSetUp()
71 
72    Notes:
73    The interface routine PCSetUp() is not usually called directly by
74    the user, but instead is called by PCApply() if necessary.
75 */
76 
77 #undef __FUNCT__
78 #define __FUNCT__ "PCSetUp_ML"
79 static PetscErrorCode PCSetUp_ML(PC pc)
80 {
81   PetscErrorCode       ierr;
82   PetscMPIInt          size,rank;
83   FineGridCtx          *PetscMLdata;
84   ML                   *ml_object;
85   ML_Aggregate         *agg_object;
86   ML_Operator          *mlmat;
87   PetscInt             nlocal_allcols,Nlevels,mllevel,level,m,fine_level;
88   Mat                  A,Aloc;
89   GridCtx              *gridctx;
90   PC                   pc_coarse;
91   PC_ML                *pc_ml=PETSC_NULL;
92   PetscObjectContainer container;
93 
94   PetscFunctionBegin;
95   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
96   if (container) {
97     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
98   } else {
99     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
100   }
101 
102   /* setup special features of PCML */
103   /*--------------------------------*/
104   /* covert A to Aloc to be used by ML at fine grid */
105   A = pc->pmat;
106   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
107   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr);
108   pc_ml->size = size;
109   if (size > 1){
110     Aloc = PETSC_NULL;
111     ierr = MatConvert_MPIAIJ_ML(A,0,&Aloc);CHKERRQ(ierr);
112   } else {
113     Aloc = A;
114   }
115 
116   /* create and initialize struct 'PetscMLdata' */
117   ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
118   ierr = PetscMemzero(PetscMLdata,sizeof(FineGridCtx));CHKERRQ(ierr);
119   PetscMLdata->A    = A;
120   PetscMLdata->Aloc = Aloc;
121   ierr = PetscMalloc((Aloc->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
122   pc_ml->PetscMLdata = PetscMLdata;
123 
124   /* create ML discretization matrix at fine grid */
125   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
126   ML_Create(&ml_object,pc_ml->MaxNlevels);
127   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
128   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
129   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
130 
131   /* aggregation */
132   ML_Aggregate_Create(&agg_object);
133   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
134   /* set options */
135   switch (pc_ml->CoarsenScheme) {
136   case 1:
137     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
138   case 2:
139     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
140   case 3:
141     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
142   }
143   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
144   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
145   if (pc_ml->SpectralNormScheme_Anorm){
146     ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object);
147   }
148 
149   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
150   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
151   ierr = MGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
152   pc_ml->ml_object  = ml_object;
153   pc_ml->agg_object = agg_object;
154 
155   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
156   fine_level = Nlevels - 1;
157   pc_ml->gridctx = gridctx;
158   pc_ml->fine_level = fine_level;
159 
160   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
161      Level 0 is the finest grid for ML, but coarsest for PETSc! */
162   gridctx[fine_level].A = A;
163   level = fine_level - 1;
164   if (size == 1){ /* convert ML mat format into petsc seqaij format */
165     PetscMLdata->rlen_max = A->N;
166     ierr = PetscMalloc(PetscMLdata->rlen_max*(sizeof(PetscScalar)+sizeof(PetscInt)),&PetscMLdata->cols);CHKERRQ(ierr);
167     PetscMLdata->vals = (PetscScalar*)(PetscMLdata->cols + PetscMLdata->rlen_max);
168     for (mllevel=1; mllevel<Nlevels; mllevel++){
169       PetscMLdata->mlmat  = &(ml_object->Pmat[mllevel]);
170       ierr = MatConvert_ML_SeqAIJ(PetscMLdata,&gridctx[level].P);CHKERRQ(ierr);
171       PetscMLdata->mlmat  = &(ml_object->Amat[mllevel]);
172       ierr = MatConvert_ML_SeqAIJ(PetscMLdata,&gridctx[level].A);CHKERRQ(ierr);
173       PetscMLdata->mlmat  = &(ml_object->Rmat[mllevel-1]);
174       ierr = MatConvert_ML_SeqAIJ(PetscMLdata,&gridctx[level].R);CHKERRQ(ierr);
175       level--;
176     }
177   } else { /* convert ML mat format into petsc shell format */
178     for (mllevel=1; mllevel<Nlevels; mllevel++){
179       mlmat  = &(ml_object->Pmat[mllevel]);
180       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr);
181       mlmat  = &(ml_object->Amat[mllevel]);
182       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].A);CHKERRQ(ierr);
183       mlmat  = &(ml_object->Rmat[mllevel-1]);
184       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr);
185       level--;
186     }
187   }
188 
189   /* create coarse level and the interpolation between the levels */
190   level = fine_level;
191   while ( level >= 0 ){
192     if (level != fine_level){
193       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr);
194       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->n,PETSC_DECIDE);CHKERRQ(ierr);
195       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
196       ierr = MGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
197 
198       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr);
199       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr);
200       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
201       ierr = MGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
202     }
203     ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].r);CHKERRQ(ierr);
204     ierr = VecSetSizes(gridctx[level].r,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr);
205     ierr = VecSetType(gridctx[level].r,VECMPI);CHKERRQ(ierr);
206     ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr);
207 
208     if (level == 0){
209       ierr = MGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
210     } else {
211       ierr = MGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
212       ierr = MGSetResidual(pc,level,MGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
213       if (level == fine_level){
214         ierr = KSPSetOptionsPrefix(gridctx[level].ksp,"mg_fine_");CHKERRQ(ierr);
215         ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr);
216       }
217     }
218     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
219 
220     if (level < fine_level){
221       if (size > 1){
222         ierr = KSPGetPC(gridctx[level].ksp,&pc_coarse);CHKERRQ(ierr);
223         ierr = PCSetType(pc_coarse,PCNONE);CHKERRQ(ierr);
224       }
225       ierr = MGSetInterpolate(pc,level+1,gridctx[level].P);CHKERRQ(ierr);
226       ierr = MGSetRestriction(pc,level+1,gridctx[level].R);CHKERRQ(ierr);
227     }
228     level--;
229   }
230 
231   /* now call PCSetUp_MG()         */
232   /*--------------------------------*/
233   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML"
239 PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr)
240 {
241   PetscErrorCode       ierr;
242   PC_ML                *pc_ml = (PC_ML*)ptr;
243   PetscInt             level;
244 
245   PetscFunctionBegin;
246   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
247   ML_Aggregate_Destroy(&pc_ml->agg_object);
248   ML_Destroy(&pc_ml->ml_object);
249 
250   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
251   if (pc_ml->size == 1){ierr = PetscFree(pc_ml->PetscMLdata->cols);CHKERRQ(ierr);}
252   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
253 
254   level = pc_ml->fine_level;
255   while ( level>= 0 ){
256     if (level != pc_ml->fine_level){
257       ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);
258       ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);
259       ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);
260       ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);
261       ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);
262     }
263     ierr = VecDestroy(pc_ml->gridctx[level].r);CHKERRQ(ierr);
264     level--;
265   }
266   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
267   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 /* -------------------------------------------------------------------------- */
271 extern PetscErrorCode PetscObjectContainerSetUserDestroy(PetscObjectContainer,PetscErrorCode (*)(void*));
272 /*
273    PCDestroy_ML - Destroys the private context for the ML preconditioner
274    that was created with PCCreate_ML().
275 
276    Input Parameter:
277 .  pc - the preconditioner context
278 
279    Application Interface Routine: PCDestroy()
280 */
281 #undef __FUNCT__
282 #define __FUNCT__ "PCDestroy_ML"
283 static PetscErrorCode PCDestroy_ML(PC pc)
284 {
285   PetscErrorCode       ierr;
286   PC_ML                *pc_ml=PETSC_NULL;
287   PetscObjectContainer container;
288 
289   PetscFunctionBegin;
290   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
291   if (container) {
292     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
293     pc->ops->destroy = pc_ml->PCDestroy;
294   } else {
295     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
296   }
297   /* detach pc and PC_ML and dereference container ??? */
298   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
299   ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
300 
301   ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr);
302   ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "PCSetFromOptions_ML"
308 static PetscErrorCode PCSetFromOptions_ML(PC pc)
309 {
310   PetscErrorCode ierr;
311   PetscInt       indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
312   PetscReal      Threshold,DampingFactor;
313   PetscTruth     flg;
314   const char     *type[] = {"additive","multiplicative","full","cascade","kascade"};
315   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
316   PC_ML          *pc_ml=PETSC_NULL;
317   PetscObjectContainer container;
318 
319   PetscFunctionBegin;
320   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
321   if (container) {
322     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
323   } else {
324     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
325   }
326   ierr = PetscOptionsHead("MG options");CHKERRQ(ierr);
327   /* inherited MG options */
328   ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
329   if (flg) {
330     ierr = MGSetCycles(pc,m);CHKERRQ(ierr);
331   }
332   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
333   if (flg) {
334     ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
335   }
336   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
337   if (flg) {
338     ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
339   }
340   ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr);
341   if (flg) {
342     MGType mg = (MGType) 0;
343     switch (indx) {
344     case 0:
345       mg = MGADDITIVE;
346       break;
347     case 1:
348       mg = MGMULTIPLICATIVE;
349       break;
350     case 2:
351       mg = MGFULL;
352       break;
353     case 3:
354       mg = MGKASKADE;
355       break;
356     case 4:
357       mg = MGKASKADE;
358       break;
359     }
360     ierr = MGSetType(pc,mg);CHKERRQ(ierr);
361   }
362   ierr = PetscOptionsTail();CHKERRQ(ierr);
363 
364   /* ML options */
365   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
366   /* set defaults */
367   PrintLevel    = 0;
368   MaxNlevels    = 10;
369   MaxCoarseSize = 1;
370   indx          = 0;
371   Threshold     = 0.0;
372   DampingFactor = 4.0/3.0;
373 
374   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
375   ML_Set_PrintLevel(PrintLevel);
376 
377   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
378   pc_ml->MaxNlevels = MaxNlevels;
379 
380   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
381   pc_ml->MaxCoarseSize = MaxCoarseSize;
382 
383   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
384   pc_ml->CoarsenScheme = indx;
385 
386   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr);
387   pc_ml->DampingFactor = DampingFactor;
388 
389   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr);
390   pc_ml->Threshold = Threshold;
391 
392   ierr = PetscOptionsLogical("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_FALSE);
393 
394   ierr = PetscOptionsTail();CHKERRQ(ierr);
395   PetscFunctionReturn(0);
396 }
397 
398 /* -------------------------------------------------------------------------- */
399 /*
400    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
401    and sets this as the private data within the generic preconditioning
402    context, PC, that was created within PCCreate().
403 
404    Input Parameter:
405 .  pc - the preconditioner context
406 
407    Application Interface Routine: PCCreate()
408 */
409 
410 /*MC
411      PCML - Use geometric multigrid preconditioning. This preconditioner requires you provide
412        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
413        operators are computed by ML and wrapped as PETSc shell matrices.
414 
415    Options Database Key: (not done yet!)
416 +  -pc_mg_maxlevels <nlevels> - maximum number of levels including finest
417 .  -pc_mg_cycles 1 or 2 - for V or W-cycle
418 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
419 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
420 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
421 .  -pc_mg_monitor - print information on the multigrid convergence
422 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
423                         to the Socket viewer for reading from Matlab.
424 
425    Level: intermediate
426 
427   Concepts: multigrid
428 
429 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
430            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
431            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
432            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
433            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()
434 M*/
435 
436 EXTERN_C_BEGIN
437 #undef __FUNCT__
438 #define __FUNCT__ "PCCreate_ML"
439 PetscErrorCode PCCreate_ML(PC pc)
440 {
441   PetscErrorCode       ierr;
442   PC_ML                *pc_ml;
443   PetscObjectContainer container;
444 
445   PetscFunctionBegin;
446   /* initialize pc as PCMG */
447   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
448 
449   /* create a supporting struct and attach it to pc */
450   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
451   /* ierr = PetscMemzero(pc_ml,sizeof(PC_ML));CHKERRQ(ierr); */
452   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
453   ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
454   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
455 
456   pc_ml->PCSetUp   = pc->ops->setup;
457   pc_ml->PCDestroy = pc->ops->destroy;
458 
459   /* overwrite the pointers of PCMG by the functions of PCML */
460   pc->ops->setfromoptions = PCSetFromOptions_ML;
461   pc->ops->setup          = PCSetUp_ML;
462   pc->ops->destroy        = PCDestroy_ML;
463   PetscFunctionReturn(0);
464 }
465 EXTERN_C_END
466 
467 int PetscML_getrow(void *ML_data, int N_requested_rows, int requested_rows[],
468    int allocated_space, int columns[], double values[], int row_lengths[])
469 {
470   PetscErrorCode ierr;
471   Mat            Aloc;
472   Mat_SeqAIJ     *a;
473   PetscInt       m,i,j,k=0,row,*aj;
474   PetscScalar    *aa;
475   FineGridCtx    *ml=(FineGridCtx*)ML_data;
476 
477   Aloc = ml->Aloc;
478   a    = (Mat_SeqAIJ*)Aloc->data;
479   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
480 
481   for (i = 0; i<N_requested_rows; i++) {
482     row   = requested_rows[i];
483     row_lengths[i] = a->ilen[row];
484     if (allocated_space < k+row_lengths[i]) return(0);
485     if ( (row >= 0) || (row <= (m-1)) ) {
486       aj = a->j + a->i[row];
487       aa = a->a + a->i[row];
488       for (j=0; j<row_lengths[i]; j++){
489         columns[k]  = aj[j];
490         values[k++] = aa[j];
491       }
492     }
493   }
494   return(1);
495 }
496 
497 int PetscML_matvec(void *ML_data,int in_length,double p[],int out_length,double ap[])
498 {
499   PetscErrorCode ierr;
500   FineGridCtx    *ml=(FineGridCtx*)ML_data;
501   Mat            A=ml->A, Aloc=ml->Aloc;
502   Vec            x,y;
503   PetscMPIInt    size;
504   PetscScalar    *pwork=ml->pwork;
505   PetscInt       i;
506 
507   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
508   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,out_length,ap,&y);CHKERRQ(ierr);
509   if (size == 1){
510     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,in_length,p,&x);CHKERRQ(ierr);
511   } else {
512     for (i=0; i<in_length; i++) pwork[i] = p[i];
513     PetscML_comm(pwork,ml);
514     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Aloc->n,pwork,&x);CHKERRQ(ierr);
515   }
516   ierr = MatMult(Aloc,x,y);CHKERRQ(ierr);
517   ierr = VecDestroy(x);CHKERRQ(ierr);
518   ierr = VecDestroy(y);CHKERRQ(ierr);
519   return 0;
520 }
521 
522 int PetscML_comm(double p[],void *ML_data)
523 {
524   PetscErrorCode ierr;
525   FineGridCtx    *ml=(FineGridCtx*)ML_data;
526   Mat            A=ml->A;
527   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
528   Vec            x;
529   PetscMPIInt    size;
530   PetscInt       i,in_length=A->m,out_length=ml->Aloc->n;
531   PetscScalar    *array;
532 
533   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
534   if (size == 1) return 0;
535   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,in_length,p,&x);CHKERRQ(ierr);
536   ierr = VecScatterBegin(x,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
537   ierr = VecScatterEnd(x,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
538   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
539   for (i=in_length; i<out_length; i++){
540     p[i] = array[i-in_length];
541   }
542   ierr = VecDestroy(x);CHKERRQ(ierr);
543   return 0;
544 }
545 #undef __FUNCT__
546 #define __FUNCT__ "MatMult_ML"
547 PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
548 {
549   PetscErrorCode   ierr;
550   Mat_MLShell      *shell;
551   PetscScalar      *xarray,*yarray;
552   PetscInt         x_length,y_length;
553 
554   PetscFunctionBegin;
555   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
556   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
557   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
558   x_length = shell->mlmat->invec_leng;
559   y_length = shell->mlmat->outvec_leng;
560 
561   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
562 
563   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
564   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
565   PetscFunctionReturn(0);
566 }
567 /* MatMultAdd_ML -  Compute y = w + A*x */
568 #undef __FUNCT__
569 #define __FUNCT__ "MatMultAdd_ML"
570 PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
571 {
572   PetscErrorCode    ierr;
573   Mat_MLShell       *shell;
574   PetscScalar       *xarray,*yarray;
575   const PetscScalar one=1.0;
576   PetscInt          x_length,y_length;
577 
578   PetscFunctionBegin;
579   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
580   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
581   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
582 
583   x_length = shell->mlmat->invec_leng;
584   y_length = shell->mlmat->outvec_leng;
585 
586   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
587 
588   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
589   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
590   ierr = VecAXPY(&one,w,y);CHKERRQ(ierr);
591 
592   PetscFunctionReturn(0);
593 }
594 
595 /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
596 #undef __FUNCT__
597 #define __FUNCT__ "MatConvert_MPIAIJ_ML"
598 PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,const MatType newtype,Mat *Aloc)
599 {
600   PetscErrorCode  ierr;
601   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
602   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
603   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
604   PetscScalar     *aa=a->a,*ba=b->a,*ca;
605   PetscInt        am=A->m,an=A->n,i,j,k;
606   PetscInt        *ci,*cj,ncols;
607   MatReuse        scall=MAT_INITIAL_MATRIX;
608 
609   PetscFunctionBegin;
610   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
611 
612   if (*Aloc) scall = MAT_REUSE_MATRIX;
613   if (scall == MAT_INITIAL_MATRIX){
614     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
615     ci[0] = 0;
616     for (i=0; i<am; i++){
617       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
618     }
619     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
620     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
621 
622     k = 0;
623     for (i=0; i<am; i++){
624       /* diagonal portion of A */
625       ncols = ai[i+1] - ai[i];
626       for (j=0; j<ncols; j++) {
627         cj[k]   = *aj++;
628         ca[k++] = *aa++;
629       }
630       /* off-diagonal portion of A */
631       ncols = bi[i+1] - bi[i];
632       for (j=0; j<ncols; j++) {
633         cj[k]   = an + (*bj); bj++;
634         ca[k++] = *ba++;
635       }
636     }
637     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
638 
639     /* put together the new matrix */
640     an = mpimat->A->n+mpimat->B->n;
641     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
642 
643     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
644     /* Since these are PETSc arrays, change flags to free them as necessary. */
645     mat = (Mat_SeqAIJ*)(*Aloc)->data;
646     mat->freedata = PETSC_TRUE;
647     mat->nonew    = 0;
648   } else if (scall == MAT_REUSE_MATRIX){
649     mat=(Mat_SeqAIJ*)(*Aloc)->data;
650     ci = mat->i; cj = mat->j; ca = mat->a;
651     for (i=0; i<am; i++) {
652       /* diagonal portion of A */
653       ncols = ai[i+1] - ai[i];
654       for (j=0; j<ncols; j++) *ca++ = *aa++;
655       /* off-diagonal portion of A */
656       ncols = bi[i+1] - bi[i];
657       for (j=0; j<ncols; j++) *ca++ = *ba++;
658     }
659   } else {
660     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
661   }
662   PetscFunctionReturn(0);
663 }
664 extern PetscErrorCode MatDestroy_Shell(Mat);
665 #undef __FUNCT__
666 #define __FUNCT__ "MatDestroy_ML"
667 PetscErrorCode MatDestroy_ML(Mat A)
668 {
669   PetscErrorCode ierr;
670   Mat_MLShell    *shell;
671 
672   PetscFunctionBegin;
673   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
674   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
675   ierr = PetscFree(shell);CHKERRQ(ierr);
676   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "MatConvert_ML_SeqAIJ"
682 PetscErrorCode MatConvert_ML_SeqAIJ(FineGridCtx *ml,Mat *newmat)
683 {
684   PetscErrorCode  ierr;
685   PetscInt        i;
686   ML_Operator     *mat=ml->mlmat;
687   PetscInt        m=mat->outvec_leng,n= mat->invec_leng,nnz[m];
688 
689   PetscFunctionBegin;
690   if ( mat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mat->getrow = NULL");
691 
692   ierr = MatCreate(PETSC_COMM_SELF,m,n,PETSC_DECIDE,PETSC_DECIDE,newmat);CHKERRQ(ierr);
693   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
694   for (i=0; i<m; i++){
695     ML_get_matrix_row(mat,1,&i,&ml->rlen_max,&ml->cols,&ml->vals,&nnz[i],0);
696   }
697   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
698 
699   for (i=0; i<m; i++){
700     ML_get_matrix_row(mat,1,&i,&ml->rlen_max,&ml->cols,&ml->vals,&nnz[i],0);
701     ierr = MatSetValues(*newmat,1,&i,nnz[i],ml->cols,ml->vals,INSERT_VALUES);CHKERRQ(ierr);
702   }
703   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
704   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
705   PetscFunctionReturn(0);
706 }
707 
708 #undef __FUNCT__
709 #define __FUNCT__ "MatConvert_ML_SHELL"
710 PetscErrorCode MatConvert_ML_SHELL(ML_Operator *mlmat,Mat *newmat)
711 {
712   PetscErrorCode ierr;
713   PetscInt       m,n;
714   ML_Comm        *MLcomm;
715   Mat_MLShell    *shellctx;
716 
717   PetscFunctionBegin;
718   m = mlmat->outvec_leng;
719   n = mlmat->invec_leng;
720   if (!m || !n){
721     newmat = PETSC_NULL;
722   } else {
723     MLcomm = mlmat->comm;
724     ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
725     ierr = PetscMemzero(shellctx,sizeof(Mat_MLShell));CHKERRQ(ierr);
726     ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
727     ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void *)MatMult_ML);CHKERRQ(ierr);
728     ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void *)MatMultAdd_ML);CHKERRQ(ierr);
729     shellctx->A         = *newmat;
730     shellctx->mlmat     = mlmat;
731     ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
732     ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
733     ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
734     (*newmat)->ops->destroy = MatDestroy_ML;
735   }
736   PetscFunctionReturn(0);
737 }
738