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