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