Actual source code: ml.c

  1: #define PETSCKSP_DLL

  3: /* 
  4:    Provides an interface to the ML 4.0 smoothed Aggregation 
  5: */
 6:  #include private/pcimpl.h
 7:  #include src/ksp/pc/impls/mg/mgimpl.h
 8:  #include src/mat/impls/aij/seq/aij.h
 9:  #include src/mat/impls/aij/mpi/mpiaij.h

 12: #include <math.h> 
 13: #include "ml_include.h"

 16: /* The context (data structure) at each grid level */
 17: typedef struct {
 18:   Vec        x,b,r;           /* global vectors */
 19:   Mat        A,P,R;
 20:   KSP        ksp;
 21: } GridCtx;

 23: /* The context used to input PETSc matrix into ML at fine grid */
 24: typedef struct {
 25:   Mat          A,Aloc;
 26:   Vec          x,y;
 27:   ML_Operator  *mlmat;
 28:   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
 29: } FineGridCtx;

 31: /* The context associates a ML matrix with a PETSc shell matrix */
 32: typedef struct {
 33:   Mat          A;       /* PETSc shell matrix associated with mlmat */
 34:   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
 35:   Vec          y;
 36: } Mat_MLShell;

 38: /* Private context for the ML preconditioner */
 39: typedef struct {
 40:   ML           *ml_object;
 41:   ML_Aggregate *agg_object;
 42:   GridCtx      *gridctx;
 43:   FineGridCtx  *PetscMLdata;
 44:   PetscInt     fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme;
 45:   PetscReal    Threshold,DampingFactor;
 46:   PetscTruth   SpectralNormScheme_Anorm;
 47:   PetscMPIInt  size;
 48:   PetscErrorCode (*PCSetUp)(PC);
 49:   PetscErrorCode (*PCDestroy)(PC);
 50: } PC_ML;

 53:    int allocated_space,int columns[],double values[],int row_lengths[]);

 64: /* -------------------------------------------------------------------------- */
 65: /*
 66:    PCSetUp_ML - Prepares for the use of the ML preconditioner
 67:                     by setting data structures and options.   

 69:    Input Parameter:
 70: .  pc - the preconditioner context

 72:    Application Interface Routine: PCSetUp()

 74:    Notes:
 75:    The interface routine PCSetUp() is not usually called directly by
 76:    the user, but instead is called by PCApply() if necessary.
 77: */
 81: PetscErrorCode PCSetUp_ML(PC pc)
 82: {
 83:   PetscErrorCode  ierr;
 84:   PetscMPIInt     size;
 85:   FineGridCtx     *PetscMLdata;
 86:   ML              *ml_object;
 87:   ML_Aggregate    *agg_object;
 88:   ML_Operator     *mlmat;
 89:   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level;
 90:   Mat             A,Aloc;
 91:   GridCtx         *gridctx;
 92:   PC_ML           *pc_ml=PETSC_NULL;
 93:   PetscContainer  container;

 96:   PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);
 97:   if (container) {
 98:     PetscContainerGetPointer(container,(void **)&pc_ml);
 99:   } else {
100:     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
101:   }
102: 
103:   /* setup special features of PCML */
104:   /*--------------------------------*/
105:   /* covert A to Aloc to be used by ML at fine grid */
106:   A = pc->pmat;
107:   MPI_Comm_size(A->comm,&size);
108:   pc_ml->size = size;
109:   if (size > 1){
110:     MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);
111:   } else {
112:     Aloc = A;
113:   }

115:   /* create and initialize struct 'PetscMLdata' */
116:   PetscNew(FineGridCtx,&PetscMLdata);
117:   PetscMLdata->A    = A;
118:   PetscMLdata->Aloc = Aloc;
119:   PetscMalloc((Aloc->cmap.n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);
120:   pc_ml->PetscMLdata = PetscMLdata;

122:   VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);
123:   if (size == 1){
124:     VecSetSizes(PetscMLdata->x,A->cmap.n,A->cmap.n);
125:   } else {
126:     VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);
127:   }
128:   VecSetType(PetscMLdata->x,VECSEQ);

130:   VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);
131:   VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);
132:   VecSetType(PetscMLdata->y,VECSEQ);
133: 
134:   /* create ML discretization matrix at fine grid */
135:   MatGetSize(Aloc,&m,&nlocal_allcols);
136:   ML_Create(&ml_object,pc_ml->MaxNlevels);
137:   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
138:   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
139:   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);

141:   /* aggregation */
142:   ML_Aggregate_Create(&agg_object);
143:   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
144:   /* set options */
145:   switch (pc_ml->CoarsenScheme) {
146:   case 1:
147:     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
148:   case 2:
149:     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
150:   case 3:
151:     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
152:   }
153:   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
154:   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
155:   if (pc_ml->SpectralNormScheme_Anorm){
156:     ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object);
157:   }
158: 
159:   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
160:   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
161:   PCMGSetLevels(pc,Nlevels,PETSC_NULL);
162:   PCSetFromOptions_MG(pc); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
163:   pc_ml->ml_object  = ml_object;
164:   pc_ml->agg_object = agg_object;

166:   PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);
167:   fine_level = Nlevels - 1;
168:   pc_ml->gridctx = gridctx;
169:   pc_ml->fine_level = fine_level;

171:   /* wrap ML matrices by PETSc shell matrices at coarsened grids. 
172:      Level 0 is the finest grid for ML, but coarsest for PETSc! */
173:   gridctx[fine_level].A = A;
174:   level = fine_level - 1;
175:   if (size == 1){ /* convert ML P, R and A into seqaij format */
176:     for (mllevel=1; mllevel<Nlevels; mllevel++){
177:       mlmat  = &(ml_object->Pmat[mllevel]);
178:       MatWrapML_SeqAIJ(mlmat,&gridctx[level].P);
179:       mlmat  = &(ml_object->Amat[mllevel]);
180:       MatWrapML_SeqAIJ(mlmat,&gridctx[level].A);
181:       mlmat  = &(ml_object->Rmat[mllevel-1]);
182:       MatWrapML_SeqAIJ(mlmat,&gridctx[level].R);
183:       level--;
184:     }
185:   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
186:     for (mllevel=1; mllevel<Nlevels; mllevel++){
187:       mlmat  = &(ml_object->Pmat[mllevel]);
188:       MatWrapML_SHELL(mlmat,&gridctx[level].P);
189:       mlmat  = &(ml_object->Rmat[mllevel-1]);
190:       MatWrapML_SHELL(mlmat,&gridctx[level].R);
191:       mlmat  = &(ml_object->Amat[mllevel]);
192:       MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);
193:       level--;
194:     }
195:   }

197:   /* create coarse level and the interpolation between the levels */
198:   for (level=0; level<fine_level; level++){
199:     VecCreate(gridctx[level].A->comm,&gridctx[level].x);
200:     VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);
201:     VecSetType(gridctx[level].x,VECMPI);
202:     PCMGSetX(pc,level,gridctx[level].x);
203: 
204:     VecCreate(gridctx[level].A->comm,&gridctx[level].b);
205:     VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);
206:     VecSetType(gridctx[level].b,VECMPI);
207:     PCMGSetRhs(pc,level,gridctx[level].b);
208: 
209:     level1 = level + 1;
210:     VecCreate(gridctx[level1].A->comm,&gridctx[level1].r);
211:     VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);
212:     VecSetType(gridctx[level1].r,VECMPI);
213:     PCMGSetR(pc,level1,gridctx[level1].r);

215:     PCMGSetInterpolate(pc,level1,gridctx[level].P);
216:     PCMGSetRestriction(pc,level1,gridctx[level].R);

218:     if (level == 0){
219:       PCMGGetCoarseSolve(pc,&gridctx[level].ksp);
220:     } else {
221:       PCMGGetSmoother(pc,level,&gridctx[level].ksp);
222:       PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);
223:     }
224:     KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);
225:   }
226:   PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);
227:   PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);
228:   KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);
229: 
230:   /* now call PCSetUp_MG()         */
231:   /*--------------------------------*/
232:   (*pc_ml->PCSetUp)(pc);
233:   return(0);
234: }

238: PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr)
239: {
240:   PetscErrorCode       ierr;
241:   PC_ML                *pc_ml = (PC_ML*)ptr;
242:   PetscInt             level;

245:   if (pc_ml->size > 1){MatDestroy(pc_ml->PetscMLdata->Aloc);}
246:   ML_Aggregate_Destroy(&pc_ml->agg_object);
247:   ML_Destroy(&pc_ml->ml_object);

249:   PetscFree(pc_ml->PetscMLdata->pwork);
250:   if (pc_ml->PetscMLdata->x){VecDestroy(pc_ml->PetscMLdata->x);}
251:   if (pc_ml->PetscMLdata->y){VecDestroy(pc_ml->PetscMLdata->y);}
252:   PetscFree(pc_ml->PetscMLdata);

254:   for (level=0; level<pc_ml->fine_level; level++){
255:     if (pc_ml->gridctx[level].A){MatDestroy(pc_ml->gridctx[level].A);}
256:     if (pc_ml->gridctx[level].P){MatDestroy(pc_ml->gridctx[level].P);}
257:     if (pc_ml->gridctx[level].R){MatDestroy(pc_ml->gridctx[level].R);}
258:     if (pc_ml->gridctx[level].x){VecDestroy(pc_ml->gridctx[level].x);}
259:     if (pc_ml->gridctx[level].b){VecDestroy(pc_ml->gridctx[level].b);}
260:     if (pc_ml->gridctx[level+1].r){VecDestroy(pc_ml->gridctx[level+1].r);}
261:   }
262:   PetscFree(pc_ml->gridctx);
263:   PetscFree(pc_ml);
264:   return(0);
265: }
266: /* -------------------------------------------------------------------------- */
267: /*
268:    PCDestroy_ML - Destroys the private context for the ML preconditioner
269:    that was created with PCCreate_ML().

271:    Input Parameter:
272: .  pc - the preconditioner context

274:    Application Interface Routine: PCDestroy()
275: */
278: PetscErrorCode PCDestroy_ML(PC pc)
279: {
280:   PetscErrorCode  ierr;
281:   PC_ML           *pc_ml=PETSC_NULL;
282:   PetscContainer  container;

285:   PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);
286:   if (container) {
287:     PetscContainerGetPointer(container,(void **)&pc_ml);
288:     pc->ops->destroy = pc_ml->PCDestroy;
289:   } else {
290:     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
291:   }
292:   /* detach pc and PC_ML and dereference container */
293:   PetscObjectCompose((PetscObject)pc,"PC_ML",0);
294:   (*pc->ops->destroy)(pc);

296:   PetscContainerDestroy(container);
297:   return(0);
298: }

302: PetscErrorCode PCSetFromOptions_ML(PC pc)
303: {
304:   PetscErrorCode  ierr;
305:   PetscInt        indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
306:   PetscReal       Threshold,DampingFactor;
307:   PetscTruth      flg;
308:   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
309:   PC_ML           *pc_ml=PETSC_NULL;
310:   PetscContainer  container;
311:   PCMGType        mgtype;

314:   PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);
315:   if (container) {
316:     PetscContainerGetPointer(container,(void **)&pc_ml);
317:   } else {
318:     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
319:   }

321:   /* inherited MG options */
322:   PetscOptionsHead("Multigrid options(inherited)");
323:     PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);
324:     PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);
325:     PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);
326:     PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);
327:   PetscOptionsTail();

329:   /* ML options */
330:   PetscOptionsHead("ML options");
331:   /* set defaults */
332:   PrintLevel    = 0;
333:   MaxNlevels    = 10;
334:   MaxCoarseSize = 1;
335:   indx          = 0;
336:   Threshold     = 0.0;
337:   DampingFactor = 4.0/3.0;
338: 
339:   PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);
340:   ML_Set_PrintLevel(PrintLevel);

342:   PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);
343:   pc_ml->MaxNlevels = MaxNlevels;

345:   PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);
346:   pc_ml->MaxCoarseSize = MaxCoarseSize;

348:   PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);
349:   pc_ml->CoarsenScheme = indx;

351:   PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);
352:   pc_ml->DampingFactor = DampingFactor;
353: 
354:   PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);
355:   pc_ml->Threshold = Threshold;

357:   PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
358: 
359:   PetscOptionsTail();
360:   return(0);
361: }

363: /* -------------------------------------------------------------------------- */
364: /*
365:    PCCreate_ML - Creates a ML preconditioner context, PC_ML, 
366:    and sets this as the private data within the generic preconditioning 
367:    context, PC, that was created within PCCreate().

369:    Input Parameter:
370: .  pc - the preconditioner context

372:    Application Interface Routine: PCCreate()
373: */

375: /*MC
376:      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 
377:        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 
378:        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
379:        and the restriction/interpolation operators wrapped as PETSc shell matrices.

381:    Options Database Key: 
382:    Multigrid options(inherited)
383: +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
384: .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
385: .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
386: -  -pc_mg_type <multiplicative> (one of) additive multiplicative full cascade kascade
387:    
388:    ML options:
389: +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
390: .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
391: .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
392: .  -pc_ml_CoarsenScheme <Uncoupled> (one of) Uncoupled Coupled MIS METIS
393: .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
394: .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
395: -  -pc_ml_SpectralNormScheme_Anorm: <false> Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm)

397:    Level: intermediate

399:   Concepts: multigrid
400:  
401: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 
402:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
403:            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
404:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
405:            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()      
406: M*/

411: PetscErrorCode  PCCreate_ML(PC pc)
412: {
413:   PetscErrorCode  ierr;
414:   PC_ML           *pc_ml;
415:   PetscContainer  container;

418:   /* initialize pc as PCMG */
419:   PCSetType(pc,PCMG); /* calls PCCreate_MG() and MGCreate_Private() */

421:   /* create a supporting struct and attach it to pc */
422:   PetscNew(PC_ML,&pc_ml);
423:   PetscContainerCreate(PETSC_COMM_SELF,&container);
424:   PetscContainerSetPointer(container,pc_ml);
425:   PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);
426:   PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);
427: 
428:   pc_ml->PCSetUp   = pc->ops->setup;
429:   pc_ml->PCDestroy = pc->ops->destroy;

431:   /* overwrite the pointers of PCMG by the functions of PCML */
432:   pc->ops->setfromoptions = PCSetFromOptions_ML;
433:   pc->ops->setup          = PCSetUp_ML;
434:   pc->ops->destroy        = PCDestroy_ML;
435:   return(0);
436: }

439: int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
440:    int allocated_space, int columns[], double values[], int row_lengths[])
441: {
443:   Mat            Aloc;
444:   Mat_SeqAIJ     *a;
445:   PetscInt       m,i,j,k=0,row,*aj;
446:   PetscScalar    *aa;
447:   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);

449:   Aloc = ml->Aloc;
450:   a    = (Mat_SeqAIJ*)Aloc->data;
451:   MatGetSize(Aloc,&m,PETSC_NULL);

453:   for (i = 0; i<N_requested_rows; i++) {
454:     row   = requested_rows[i];
455:     row_lengths[i] = a->ilen[row];
456:     if (allocated_space < k+row_lengths[i]) return(0);
457:     if ( (row >= 0) || (row <= (m-1)) ) {
458:       aj = a->j + a->i[row];
459:       aa = a->a + a->i[row];
460:       for (j=0; j<row_lengths[i]; j++){
461:         columns[k]  = aj[j];
462:         values[k++] = aa[j];
463:       }
464:     }
465:   }
466:   return(1);
467: }

469: int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
470: {
472:   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
473:   Mat            A=ml->A, Aloc=ml->Aloc;
474:   PetscMPIInt    size;
475:   PetscScalar    *pwork=ml->pwork;
476:   PetscInt       i;

478:   MPI_Comm_size(A->comm,&size);
479:   if (size == 1){
480:     VecPlaceArray(ml->x,p);
481:   } else {
482:     for (i=0; i<in_length; i++) pwork[i] = p[i];
483:     PetscML_comm(pwork,ml);
484:     VecPlaceArray(ml->x,pwork);
485:   }
486:   VecPlaceArray(ml->y,ap);
487:   MatMult(Aloc,ml->x,ml->y);
488:   VecResetArray(ml->x);
489:   VecResetArray(ml->y);
490:   return 0;
491: }

493: int PetscML_comm(double p[],void *ML_data)
494: {
496:   FineGridCtx    *ml=(FineGridCtx*)ML_data;
497:   Mat            A=ml->A;
498:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
499:   PetscMPIInt    size;
500:   PetscInt       i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n;
501:   PetscScalar    *array;

503:   MPI_Comm_size(A->comm,&size);
504:   if (size == 1) return 0;
505: 
506:   VecPlaceArray(ml->y,p);
507:   VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);
508:   VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);
509:   VecResetArray(ml->y);
510:   VecGetArray(a->lvec,&array);
511:   for (i=in_length; i<out_length; i++){
512:     p[i] = array[i-in_length];
513:   }
514:   VecRestoreArray(a->lvec,&array);
515:   return 0;
516: }
519: PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
520: {
521:   PetscErrorCode   ierr;
522:   Mat_MLShell      *shell;
523:   PetscScalar      *xarray,*yarray;
524:   PetscInt         x_length,y_length;
525: 
527:   MatShellGetContext(A,(void **)&shell);
528:   VecGetArray(x,&xarray);
529:   VecGetArray(y,&yarray);
530:   x_length = shell->mlmat->invec_leng;
531:   y_length = shell->mlmat->outvec_leng;

533:   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);

535:   VecRestoreArray(x,&xarray);
536:   VecRestoreArray(y,&yarray);
537:   return(0);
538: }
539: /* MatMultAdd_ML -  Compute y = w + A*x */
542: PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
543: {
544:   PetscErrorCode    ierr;
545:   Mat_MLShell       *shell;
546:   PetscScalar       *xarray,*yarray;
547:   PetscInt          x_length,y_length;
548: 
550:   MatShellGetContext(A,(void **)&shell);
551:   VecGetArray(x,&xarray);
552:   VecGetArray(y,&yarray);

554:   x_length = shell->mlmat->invec_leng;
555:   y_length = shell->mlmat->outvec_leng;

557:   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);

559:   VecRestoreArray(x,&xarray);
560:   VecRestoreArray(y,&yarray);
561:   VecAXPY(y,1.0,w);

563:   return(0);
564: }

566: /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
569: PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
570: {
571:   PetscErrorCode  ierr;
572:   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
573:   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
574:   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
575:   PetscScalar     *aa=a->a,*ba=b->a,*ca;
576:   PetscInt        am=A->rmap.n,an=A->cmap.n,i,j,k;
577:   PetscInt        *ci,*cj,ncols;

580:   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);

582:   if (scall == MAT_INITIAL_MATRIX){
583:     PetscMalloc((1+am)*sizeof(PetscInt),&ci);
584:     ci[0] = 0;
585:     for (i=0; i<am; i++){
586:       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
587:     }
588:     PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);
589:     PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);

591:     k = 0;
592:     for (i=0; i<am; i++){
593:       /* diagonal portion of A */
594:       ncols = ai[i+1] - ai[i];
595:       for (j=0; j<ncols; j++) {
596:         cj[k]   = *aj++;
597:         ca[k++] = *aa++;
598:       }
599:       /* off-diagonal portion of A */
600:       ncols = bi[i+1] - bi[i];
601:       for (j=0; j<ncols; j++) {
602:         cj[k]   = an + (*bj); bj++;
603:         ca[k++] = *ba++;
604:       }
605:     }
606:     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);

608:     /* put together the new matrix */
609:     an = mpimat->A->cmap.n+mpimat->B->cmap.n;
610:     MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);

612:     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
613:     /* Since these are PETSc arrays, change flags to free them as necessary. */
614:     mat = (Mat_SeqAIJ*)(*Aloc)->data;
615:     mat->free_a       = PETSC_TRUE;
616:     mat->free_ij      = PETSC_TRUE;

618:     mat->nonew    = 0;
619:   } else if (scall == MAT_REUSE_MATRIX){
620:     mat=(Mat_SeqAIJ*)(*Aloc)->data;
621:     ci = mat->i; cj = mat->j; ca = mat->a;
622:     for (i=0; i<am; i++) {
623:       /* diagonal portion of A */
624:       ncols = ai[i+1] - ai[i];
625:       for (j=0; j<ncols; j++) *ca++ = *aa++;
626:       /* off-diagonal portion of A */
627:       ncols = bi[i+1] - bi[i];
628:       for (j=0; j<ncols; j++) *ca++ = *ba++;
629:     }
630:   } else {
631:     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
632:   }
633:   return(0);
634: }
638: PetscErrorCode MatDestroy_ML(Mat A)
639: {
641:   Mat_MLShell    *shell;

644:   MatShellGetContext(A,(void **)&shell);
645:   VecDestroy(shell->y);
646:   PetscFree(shell);
647:   MatDestroy_Shell(A);
648:   PetscObjectChangeTypeName((PetscObject)A,0);
649:   return(0);
650: }

654: PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,Mat *newmat)
655: {
656:   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
657:   PetscErrorCode        ierr;
658:   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
659:   PetscInt              *ml_cols=matdata->columns,*aj,i,j,k;
660:   PetscScalar           *ml_vals=matdata->values,*aa;
661: 
663:   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
664:   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
665:     MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);
666:     return(0);
667:   }

669:   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
670:   MatCreate(PETSC_COMM_SELF,newmat);
671:   MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);
672:   MatSetType(*newmat,MATSEQAIJ);
673:   PetscMalloc((m+1)*sizeof(PetscInt),&nnz);

675:   nz_max = 0;
676:   for (i=0; i<m; i++) {
677:     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
678:     if (nnz[i] > nz_max) nz_max = nnz[i];
679:   }
680:   MatSeqAIJSetPreallocation(*newmat,0,nnz);
681:   MatSetOption(*newmat,MAT_COLUMNS_SORTED); /* check! */

683:   nz_max++;
684:   PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);
685:   aa = (PetscScalar*)(aj + nz_max);

687:   for (i=0; i<m; i++){
688:     k = 0;
689:     /* diagonal entry */
690:     aj[k] = i; aa[k++] = ml_vals[i];
691:     /* off diagonal entries */
692:     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
693:       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
694:     }
695:     /* sort aj and aa */
696:     PetscSortIntWithScalarArray(nnz[i],aj,aa);
697:     MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);
698:   }
699:   MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
700:   MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
701:   PetscFree(aj);
702:   PetscFree(nnz);
703:   return(0);
704: }

708: PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,Mat *newmat)
709: {
711:   PetscInt       m,n;
712:   ML_Comm        *MLcomm;
713:   Mat_MLShell    *shellctx;

716:   m = mlmat->outvec_leng;
717:   n = mlmat->invec_leng;
718:   if (!m || !n){
719:     newmat = PETSC_NULL;
720:   } else {
721:     MLcomm = mlmat->comm;
722:     PetscNew(Mat_MLShell,&shellctx);
723:     MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);
724:     MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);
725:     MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);
726:     shellctx->A         = *newmat;
727:     shellctx->mlmat     = mlmat;
728:     VecCreate(PETSC_COMM_WORLD,&shellctx->y);
729:     VecSetSizes(shellctx->y,m,PETSC_DECIDE);
730:     VecSetFromOptions(shellctx->y);
731:     (*newmat)->ops->destroy = MatDestroy_ML;
732:   }
733:   return(0);
734: }

738: PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
739: {
740:   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
741:   PetscInt              *ml_cols=matdata->columns,*aj;
742:   PetscScalar           *ml_vals=matdata->values,*aa;
743:   PetscErrorCode        ierr;
744:   PetscInt              i,j,k,*gordering;
745:   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
746:   Mat                   A;

749:   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
750:   n = mlmat->invec_leng;
751:   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);

753:   MatCreate(mlmat->comm->USR_comm,&A);
754:   MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);
755:   MatSetType(A,MATMPIAIJ);
756:   PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);
757: 
758:   nz_max = 0;
759:   for (i=0; i<m; i++){
760:     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
761:     if (nz_max < nnz[i]) nz_max = nnz[i];
762:     nnzA[i] = 1; /* diag */
763:     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
764:       if (ml_cols[j] < m) nnzA[i]++;
765:     }
766:     nnzB[i] = nnz[i] - nnzA[i];
767:   }
768:   MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);

770:   /* insert mat values -- remap row and column indices */
771:   nz_max++;
772:   PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);
773:   aa = (PetscScalar*)(aj + nz_max);
774:   ML_build_global_numbering(mlmat,&gordering);
775:   for (i=0; i<m; i++){
776:     row = gordering[i];
777:     k = 0;
778:     /* diagonal entry */
779:     aj[k] = row; aa[k++] = ml_vals[i];
780:     /* off diagonal entries */
781:     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
782:       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
783:     }
784:     MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);
785:   }
786:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
787:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
788:   *newmat = A;

790:   PetscFree3(nnzA,nnzB,nnz);
791:   PetscFree(aj);
792:   return(0);
793: }