xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision fb967f244476fae796298bc2a6a7dcdee8e846be)
1 #define PETSCKSP_DLL
2 
3 /*
4    Provides an interface to the ML 4.0 smoothed Aggregation
5 */
6 #include "private/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   PetscContainer  container;
94 
95   PetscFunctionBegin;
96   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
97   if (container) {
98     ierr = PetscContainerGetPointer(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->cmap.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->cmap.n,A->cmap.n);CHKERRQ(ierr);
125   } else {
126     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.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->rmap.n,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->cmap.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->rmap.n,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->rmap.n,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 
230   /* now call PCSetUp_MG()         */
231   /*--------------------------------*/
232   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
233   PetscFunctionReturn(0);
234 }
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "PetscContainerDestroy_PC_ML"
238 PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr)
239 {
240   PetscErrorCode       ierr;
241   PC_ML                *pc_ml = (PC_ML*)ptr;
242   PetscInt             level;
243 
244   PetscFunctionBegin;
245   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
246   ML_Aggregate_Destroy(&pc_ml->agg_object);
247   ML_Destroy(&pc_ml->ml_object);
248 
249   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
250   if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
251   if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
252   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
253 
254   for (level=0; level<pc_ml->fine_level; level++){
255     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
256     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
257     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
258     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
259     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
260     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
261   }
262   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
263   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
264   PetscFunctionReturn(0);
265 }
266 /* -------------------------------------------------------------------------- */
267 /*
268    PCDestroy_ML - Destroys the private context for the ML preconditioner
269    that was created with PCCreate_ML().
270 
271    Input Parameter:
272 .  pc - the preconditioner context
273 
274    Application Interface Routine: PCDestroy()
275 */
276 #undef __FUNCT__
277 #define __FUNCT__ "PCDestroy_ML"
278 PetscErrorCode PCDestroy_ML(PC pc)
279 {
280   PetscErrorCode  ierr;
281   PC_ML           *pc_ml=PETSC_NULL;
282   PetscContainer  container;
283 
284   PetscFunctionBegin;
285   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
286   if (container) {
287     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
288     pc->ops->destroy = pc_ml->PCDestroy;
289   } else {
290     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
291   }
292   /* detach pc and PC_ML and dereference container */
293   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
294   ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
295 
296   ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "PCSetFromOptions_ML"
302 PetscErrorCode PCSetFromOptions_ML(PC pc)
303 {
304   PetscErrorCode  ierr;
305   PetscInt        indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
306   PetscReal       Threshold,DampingFactor;
307   PetscTruth      flg;
308   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
309   PC_ML           *pc_ml=PETSC_NULL;
310   PetscContainer  container;
311   PCMGType        mgtype;
312 
313   PetscFunctionBegin;
314   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
315   if (container) {
316     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
317   } else {
318     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
319   }
320 
321   /* inherited MG options */
322   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
323     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
324     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
325     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
326     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
327   ierr = PetscOptionsTail();CHKERRQ(ierr);
328 
329   /* ML options */
330   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
331   /* set defaults */
332   PrintLevel    = 0;
333   MaxNlevels    = 10;
334   MaxCoarseSize = 1;
335   indx          = 0;
336   Threshold     = 0.0;
337   DampingFactor = 4.0/3.0;
338 
339   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
340   ML_Set_PrintLevel(PrintLevel);
341 
342   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
343   pc_ml->MaxNlevels = MaxNlevels;
344 
345   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
346   pc_ml->MaxCoarseSize = MaxCoarseSize;
347 
348   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
349   pc_ml->CoarsenScheme = indx;
350 
351   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr);
352   pc_ml->DampingFactor = DampingFactor;
353 
354   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr);
355   pc_ml->Threshold = Threshold;
356 
357   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
358 
359   ierr = PetscOptionsTail();CHKERRQ(ierr);
360   PetscFunctionReturn(0);
361 }
362 
363 /* -------------------------------------------------------------------------- */
364 /*
365    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
366    and sets this as the private data within the generic preconditioning
367    context, PC, that was created within PCCreate().
368 
369    Input Parameter:
370 .  pc - the preconditioner context
371 
372    Application Interface Routine: PCCreate()
373 */
374 
375 /*MC
376      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
377        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
378        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
379        and the restriction/interpolation operators wrapped as PETSc shell matrices.
380 
381    Options Database Key:
382    Multigrid options(inherited)
383 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
384 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
385 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
386 -  -pc_mg_type <multiplicative> (one of) additive multiplicative full cascade kascade
387 
388    ML options:
389 +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
390 .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
391 .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
392 .  -pc_ml_CoarsenScheme <Uncoupled> (one of) Uncoupled Coupled MIS METIS
393 .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
394 .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
395 -  -pc_ml_SpectralNormScheme_Anorm: <false> Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm)
396 
397    Level: intermediate
398 
399   Concepts: multigrid
400 
401 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
402            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
403            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
404            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
405            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
406 M*/
407 
408 EXTERN_C_BEGIN
409 #undef __FUNCT__
410 #define __FUNCT__ "PCCreate_ML"
411 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
412 {
413   PetscErrorCode  ierr;
414   PC_ML           *pc_ml;
415   PetscContainer  container;
416 
417   PetscFunctionBegin;
418   /* initialize pc as PCMG */
419   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
420 
421   /* create a supporting struct and attach it to pc */
422   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
423   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
424   ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
425   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
426   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
427 
428   pc_ml->PCSetUp   = pc->ops->setup;
429   pc_ml->PCDestroy = pc->ops->destroy;
430 
431   /* overwrite the pointers of PCMG by the functions of PCML */
432   pc->ops->setfromoptions = PCSetFromOptions_ML;
433   pc->ops->setup          = PCSetUp_ML;
434   pc->ops->destroy        = PCDestroy_ML;
435   PetscFunctionReturn(0);
436 }
437 EXTERN_C_END
438 
439 int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
440    int allocated_space, int columns[], double values[], int row_lengths[])
441 {
442   PetscErrorCode ierr;
443   Mat            Aloc;
444   Mat_SeqAIJ     *a;
445   PetscInt       m,i,j,k=0,row,*aj;
446   PetscScalar    *aa;
447   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
448 
449   Aloc = ml->Aloc;
450   a    = (Mat_SeqAIJ*)Aloc->data;
451   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
452 
453   for (i = 0; i<N_requested_rows; i++) {
454     row   = requested_rows[i];
455     row_lengths[i] = a->ilen[row];
456     if (allocated_space < k+row_lengths[i]) return(0);
457     if ( (row >= 0) || (row <= (m-1)) ) {
458       aj = a->j + a->i[row];
459       aa = a->a + a->i[row];
460       for (j=0; j<row_lengths[i]; j++){
461         columns[k]  = aj[j];
462         values[k++] = aa[j];
463       }
464     }
465   }
466   return(1);
467 }
468 
469 int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
470 {
471   PetscErrorCode ierr;
472   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
473   Mat            A=ml->A, Aloc=ml->Aloc;
474   PetscMPIInt    size;
475   PetscScalar    *pwork=ml->pwork;
476   PetscInt       i;
477 
478   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
479   if (size == 1){
480     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
481   } else {
482     for (i=0; i<in_length; i++) pwork[i] = p[i];
483     PetscML_comm(pwork,ml);
484     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
485   }
486   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
487   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
488   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
489   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
490   return 0;
491 }
492 
493 int PetscML_comm(double p[],void *ML_data)
494 {
495   PetscErrorCode ierr;
496   FineGridCtx    *ml=(FineGridCtx*)ML_data;
497   Mat            A=ml->A;
498   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
499   PetscMPIInt    size;
500   PetscInt       i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n;
501   PetscScalar    *array;
502 
503   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
504   if (size == 1) return 0;
505 
506   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
507   ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
508   ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
509   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
510   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
511   for (i=in_length; i<out_length; i++){
512     p[i] = array[i-in_length];
513   }
514   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
515   return 0;
516 }
517 #undef __FUNCT__
518 #define __FUNCT__ "MatMult_ML"
519 PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
520 {
521   PetscErrorCode   ierr;
522   Mat_MLShell      *shell;
523   PetscScalar      *xarray,*yarray;
524   PetscInt         x_length,y_length;
525 
526   PetscFunctionBegin;
527   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
528   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
529   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
530   x_length = shell->mlmat->invec_leng;
531   y_length = shell->mlmat->outvec_leng;
532 
533   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
534 
535   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
536   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
537   PetscFunctionReturn(0);
538 }
539 /* MatMultAdd_ML -  Compute y = w + A*x */
540 #undef __FUNCT__
541 #define __FUNCT__ "MatMultAdd_ML"
542 PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
543 {
544   PetscErrorCode    ierr;
545   Mat_MLShell       *shell;
546   PetscScalar       *xarray,*yarray;
547   PetscInt          x_length,y_length;
548 
549   PetscFunctionBegin;
550   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
551   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
552   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
553 
554   x_length = shell->mlmat->invec_leng;
555   y_length = shell->mlmat->outvec_leng;
556 
557   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
558 
559   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
560   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
561   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
562 
563   PetscFunctionReturn(0);
564 }
565 
566 /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
567 #undef __FUNCT__
568 #define __FUNCT__ "MatConvert_MPIAIJ_ML"
569 PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
570 {
571   PetscErrorCode  ierr;
572   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
573   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
574   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
575   PetscScalar     *aa=a->a,*ba=b->a,*ca;
576   PetscInt        am=A->rmap.n,an=A->cmap.n,i,j,k;
577   PetscInt        *ci,*cj,ncols;
578 
579   PetscFunctionBegin;
580   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
581 
582   if (scall == MAT_INITIAL_MATRIX){
583     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
584     ci[0] = 0;
585     for (i=0; i<am; i++){
586       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
587     }
588     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
589     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
590 
591     k = 0;
592     for (i=0; i<am; i++){
593       /* diagonal portion of A */
594       ncols = ai[i+1] - ai[i];
595       for (j=0; j<ncols; j++) {
596         cj[k]   = *aj++;
597         ca[k++] = *aa++;
598       }
599       /* off-diagonal portion of A */
600       ncols = bi[i+1] - bi[i];
601       for (j=0; j<ncols; j++) {
602         cj[k]   = an + (*bj); bj++;
603         ca[k++] = *ba++;
604       }
605     }
606     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
607 
608     /* put together the new matrix */
609     an = mpimat->A->cmap.n+mpimat->B->cmap.n;
610     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
611 
612     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
613     /* Since these are PETSc arrays, change flags to free them as necessary. */
614     mat = (Mat_SeqAIJ*)(*Aloc)->data;
615     mat->free_a       = PETSC_TRUE;
616     mat->free_ij      = PETSC_TRUE;
617 
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   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
649   PetscFunctionReturn(0);
650 }
651 
652 #undef __FUNCT__
653 #define __FUNCT__ "MatWrapML_SeqAIJ"
654 PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,Mat *newmat)
655 {
656   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
657   PetscErrorCode        ierr;
658   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
659   PetscInt              *ml_cols=matdata->columns,*aj,i,j,k;
660   PetscScalar           *ml_vals=matdata->values,*aa;
661 
662   PetscFunctionBegin;
663   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
664   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
665     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
666     PetscFunctionReturn(0);
667   }
668 
669   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
670   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
671   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
672   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
673   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
674 
675   nz_max = 0;
676   for (i=0; i<m; i++) {
677     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
678     if (nnz[i] > nz_max) nz_max = nnz[i];
679   }
680   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
681   ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* check! */
682 
683   nz_max++;
684   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
685   aa = (PetscScalar*)(aj + nz_max);
686 
687   for (i=0; i<m; i++){
688     k = 0;
689     /* diagonal entry */
690     aj[k] = i; aa[k++] = ml_vals[i];
691     /* off diagonal entries */
692     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
693       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
694     }
695     /* sort aj and aa */
696     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
697     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
698   }
699   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
700   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
701   ierr = PetscFree(aj);CHKERRQ(ierr);
702   ierr = PetscFree(nnz);CHKERRQ(ierr);
703   PetscFunctionReturn(0);
704 }
705 
706 #undef __FUNCT__
707 #define __FUNCT__ "MatWrapML_SHELL"
708 PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,Mat *newmat)
709 {
710   PetscErrorCode ierr;
711   PetscInt       m,n;
712   ML_Comm        *MLcomm;
713   Mat_MLShell    *shellctx;
714 
715   PetscFunctionBegin;
716   m = mlmat->outvec_leng;
717   n = mlmat->invec_leng;
718   if (!m || !n){
719     newmat = PETSC_NULL;
720   } else {
721     MLcomm = mlmat->comm;
722     ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
723     ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
724     ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
725     ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
726     shellctx->A         = *newmat;
727     shellctx->mlmat     = mlmat;
728     ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
729     ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
730     ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
731     (*newmat)->ops->destroy = MatDestroy_ML;
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "MatWrapML_MPIAIJ"
738 PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
739 {
740   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
741   PetscInt              *ml_cols=matdata->columns,*aj;
742   PetscScalar           *ml_vals=matdata->values,*aa;
743   PetscErrorCode        ierr;
744   PetscInt              i,j,k,*gordering;
745   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
746   Mat                   A;
747 
748   PetscFunctionBegin;
749   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
750   n = mlmat->invec_leng;
751   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
752 
753   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
754   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
755   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
756   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
757 
758   nz_max = 0;
759   for (i=0; i<m; i++){
760     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
761     if (nz_max < nnz[i]) nz_max = nnz[i];
762     nnzA[i] = 1; /* diag */
763     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
764       if (ml_cols[j] < m) nnzA[i]++;
765     }
766     nnzB[i] = nnz[i] - nnzA[i];
767   }
768   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
769 
770   /* insert mat values -- remap row and column indices */
771   nz_max++;
772   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
773   aa = (PetscScalar*)(aj + nz_max);
774   ML_build_global_numbering(mlmat,&gordering);
775   for (i=0; i<m; i++){
776     row = gordering[i];
777     k = 0;
778     /* diagonal entry */
779     aj[k] = row; aa[k++] = ml_vals[i];
780     /* off diagonal entries */
781     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
782       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
783     }
784     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
785   }
786   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
787   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
788   *newmat = A;
789 
790   ierr = PetscFree3(nnzA,nnzB,nnz);
791   ierr = PetscFree(aj);CHKERRQ(ierr);
792   PetscFunctionReturn(0);
793 }
794