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