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