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