xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision e14861a4f12c1b9107747fda3b05cd179800a721)
15582bec1SHong Zhang /*
25582bec1SHong Zhang    Provides an interface to the ML 3.0 smoothed Aggregation
35582bec1SHong Zhang */
45582bec1SHong Zhang #include "src/ksp/pc/pcimpl.h"   /*I "petscpc.h" I*/
55582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
65582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h"
75582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
85582bec1SHong Zhang EXTERN_C_BEGIN
95582bec1SHong Zhang #include <math.h>
105582bec1SHong Zhang #include "ml_include.h"
115582bec1SHong Zhang EXTERN_C_END
125582bec1SHong Zhang 
135582bec1SHong Zhang /* The context (data structure) at each grid level */
145582bec1SHong Zhang typedef struct {
155582bec1SHong Zhang   Vec        x,b,r;            /* global vectors */
165582bec1SHong Zhang   Mat        A,P,R;
175582bec1SHong Zhang   KSP        ksp;
185582bec1SHong Zhang } GridCtx;
195582bec1SHong Zhang 
205582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
215582bec1SHong Zhang typedef struct {
225582bec1SHong Zhang   Mat          A,Aloc;
2324a42b14SHong Zhang   Vec          x,y;
245582bec1SHong Zhang   ML_Operator  *mlmat;
255582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
265582bec1SHong Zhang   PetscInt     rlen_max,*cols; /* used by MatConvert_ML_SeqAIJ() */
275582bec1SHong Zhang   PetscScalar  *vals;          /* used by MatConvert_ML_SeqAIJ() */
285582bec1SHong Zhang } FineGridCtx;
295582bec1SHong Zhang 
305582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
315582bec1SHong Zhang typedef struct {
325582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
335582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
345582bec1SHong Zhang   Vec          y;
355582bec1SHong Zhang } Mat_MLShell;
365582bec1SHong Zhang 
375582bec1SHong Zhang /* Private context for the ML preconditioner */
385582bec1SHong Zhang typedef struct {
395582bec1SHong Zhang   ML           *ml_object;
405582bec1SHong Zhang   ML_Aggregate *agg_object;
415582bec1SHong Zhang   GridCtx      *gridctx;
425582bec1SHong Zhang   FineGridCtx  *PetscMLdata;
435582bec1SHong Zhang   PetscInt     fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme;
445582bec1SHong Zhang   PetscReal    Threshold,DampingFactor;
455582bec1SHong Zhang   PetscTruth   SpectralNormScheme_Anorm;
465582bec1SHong Zhang   PetscMPIInt  size;
475582bec1SHong Zhang 
485582bec1SHong Zhang   PetscErrorCode (*PCSetUp)(PC);
495582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
505582bec1SHong Zhang 
515582bec1SHong Zhang } PC_ML;
525582bec1SHong Zhang extern int PetscML_getrow(void *ML_data,int N_requested_rows,int requested_rows[],
535582bec1SHong Zhang             int allocated_space,int columns[],double values[],int row_lengths[]);
545582bec1SHong Zhang extern int PetscML_matvec(void *ML_data, int in_length, double p[], int out_length,double ap[]);
555582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
565582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
575582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
585582bec1SHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,const MatType,Mat*);
595582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
60*e14861a4SHong Zhang extern PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator*,Mat*);
61*e14861a4SHong Zhang extern PetscErrorCode MatConvert_ML_MPIAIJ(FineGridCtx*,Mat*);
625582bec1SHong Zhang extern PetscErrorCode MatConvert_ML_SHELL(ML_Operator*,Mat*);
635582bec1SHong Zhang 
645582bec1SHong Zhang /* -------------------------------------------------------------------------- */
655582bec1SHong Zhang /*
665582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
675582bec1SHong Zhang                     by setting data structures and options.
685582bec1SHong Zhang 
695582bec1SHong Zhang    Input Parameter:
705582bec1SHong Zhang .  pc - the preconditioner context
715582bec1SHong Zhang 
725582bec1SHong Zhang    Application Interface Routine: PCSetUp()
735582bec1SHong Zhang 
745582bec1SHong Zhang    Notes:
755582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
765582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
775582bec1SHong Zhang */
785582bec1SHong Zhang 
795582bec1SHong Zhang #undef __FUNCT__
805582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
815582bec1SHong Zhang static PetscErrorCode PCSetUp_ML(PC pc)
825582bec1SHong Zhang {
835582bec1SHong Zhang   PetscErrorCode       ierr;
845582bec1SHong Zhang   PetscMPIInt          size,rank;
855582bec1SHong Zhang   FineGridCtx          *PetscMLdata;
865582bec1SHong Zhang   ML                   *ml_object;
875582bec1SHong Zhang   ML_Aggregate         *agg_object;
885582bec1SHong Zhang   ML_Operator          *mlmat;
895582bec1SHong Zhang   PetscInt             nlocal_allcols,Nlevels,mllevel,level,m,fine_level;
905582bec1SHong Zhang   Mat                  A,Aloc;
915582bec1SHong Zhang   GridCtx              *gridctx;
925582bec1SHong Zhang   PC                   pc_coarse;
935582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
945582bec1SHong Zhang   PetscObjectContainer container;
955582bec1SHong Zhang 
965582bec1SHong Zhang   PetscFunctionBegin;
975582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
985582bec1SHong Zhang   if (container) {
995582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
1005582bec1SHong Zhang   } else {
1015582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
1025582bec1SHong Zhang   }
1035582bec1SHong Zhang 
1045582bec1SHong Zhang   /* setup special features of PCML */
1055582bec1SHong Zhang   /*--------------------------------*/
1065582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1075582bec1SHong Zhang   A = pc->pmat;
1085582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
109*e14861a4SHong Zhang   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); /* rm! */
1105582bec1SHong Zhang   pc_ml->size = size;
1115582bec1SHong Zhang   if (size > 1){
1125582bec1SHong Zhang     Aloc = PETSC_NULL;
1135582bec1SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,0,&Aloc);CHKERRQ(ierr);
1145582bec1SHong Zhang   } else {
1155582bec1SHong Zhang     Aloc = A;
1165582bec1SHong Zhang   }
1175582bec1SHong Zhang 
1185582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
1195582bec1SHong Zhang   ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1205582bec1SHong Zhang   ierr = PetscMemzero(PetscMLdata,sizeof(FineGridCtx));CHKERRQ(ierr);
1215582bec1SHong Zhang   PetscMLdata->A    = A;
1225582bec1SHong Zhang   PetscMLdata->Aloc = Aloc;
1235582bec1SHong Zhang   ierr = PetscMalloc((Aloc->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1245582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
1255582bec1SHong Zhang 
12624a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
12724a42b14SHong Zhang   if (size == 1){
12824a42b14SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,A->n,PETSC_DECIDE);CHKERRQ(ierr);
12924a42b14SHong Zhang   } else {
13024a42b14SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,Aloc->n,PETSC_DECIDE);CHKERRQ(ierr);
13124a42b14SHong Zhang   }
13224a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
13324a42b14SHong Zhang 
13424a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
13524a42b14SHong Zhang   ierr = VecSetSizes(PetscMLdata->y,A->m,PETSC_DECIDE);CHKERRQ(ierr);
13624a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
13724a42b14SHong Zhang 
1385582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
1395582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1405582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
1415582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1425582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1435582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1445582bec1SHong Zhang 
1455582bec1SHong Zhang   /* aggregation */
1465582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
1475582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1485582bec1SHong Zhang   /* set options */
1495582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1505582bec1SHong Zhang   case 1:
1515582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1525582bec1SHong Zhang   case 2:
1535582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1545582bec1SHong Zhang   case 3:
1555582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1565582bec1SHong Zhang   }
1575582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1585582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1595582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
1605582bec1SHong Zhang     ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object);
1615582bec1SHong Zhang   }
1625582bec1SHong Zhang 
1635582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
1645582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
1655582bec1SHong Zhang   ierr = MGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
1665582bec1SHong Zhang   pc_ml->ml_object  = ml_object;
1675582bec1SHong Zhang   pc_ml->agg_object = agg_object;
1685582bec1SHong Zhang 
1695582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
1705582bec1SHong Zhang   fine_level = Nlevels - 1;
1715582bec1SHong Zhang   pc_ml->gridctx = gridctx;
1725582bec1SHong Zhang   pc_ml->fine_level = fine_level;
1735582bec1SHong Zhang 
1745582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
1755582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
1765582bec1SHong Zhang   PetscMLdata->rlen_max = A->N;
1775582bec1SHong Zhang   ierr = PetscMalloc(PetscMLdata->rlen_max*(sizeof(PetscScalar)+sizeof(PetscInt)),&PetscMLdata->cols);CHKERRQ(ierr);
1785582bec1SHong Zhang   PetscMLdata->vals = (PetscScalar*)(PetscMLdata->cols + PetscMLdata->rlen_max);
179*e14861a4SHong Zhang 
180*e14861a4SHong Zhang   gridctx[fine_level].A = A;
181*e14861a4SHong Zhang   level = fine_level - 1;
182*e14861a4SHong Zhang   if (size == 1){ /* convert ML mat format into petsc seqaij format */
1835582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
184*e14861a4SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
185*e14861a4SHong Zhang       ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].P);CHKERRQ(ierr);
186*e14861a4SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
187*e14861a4SHong Zhang       ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
188*e14861a4SHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
189*e14861a4SHong Zhang       ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].R);CHKERRQ(ierr);
1905582bec1SHong Zhang       level--;
1915582bec1SHong Zhang     }
1925582bec1SHong Zhang   } else { /* convert ML mat format into petsc shell format */
1935582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
1945582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
1955582bec1SHong Zhang       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr);
1965582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
197*e14861a4SHong Zhang       /*
198*e14861a4SHong Zhang       if (mllevel==1) {
199*e14861a4SHong Zhang         ML_Operator_Print_UsingGlobalOrdering(mlmat,"Amat1",PETSC_NULL,PETSC_NULL);
200*e14861a4SHong Zhang       }
201*e14861a4SHong Zhang       */
2025582bec1SHong Zhang       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].A);CHKERRQ(ierr);
203*e14861a4SHong Zhang #ifndef TMP
204*e14861a4SHong Zhang       if (mllevel>0){
205*e14861a4SHong Zhang         PetscMLdata->mlmat  = &(ml_object->Amat[mllevel]);
206*e14861a4SHong Zhang         Mat A_tmp;
207*e14861a4SHong Zhang         ierr = MatConvert_ML_MPIAIJ(PetscMLdata,&A_tmp);CHKERRQ(ierr);
208*e14861a4SHong Zhang         /* ierr = MatView(A_tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
209*e14861a4SHong Zhang 
210*e14861a4SHong Zhang         Vec x,yy1,yy2;
211*e14861a4SHong Zhang         PetscInt am,an;
212*e14861a4SHong Zhang         ierr = MatGetLocalSize(A_tmp,&am,&an);CHKERRQ(ierr);
213*e14861a4SHong Zhang         ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
214*e14861a4SHong Zhang         ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
215*e14861a4SHong Zhang         ierr = VecSetFromOptions(x);CHKERRQ(ierr);
216*e14861a4SHong Zhang         ierr = VecCreate(PETSC_COMM_WORLD,&yy1);CHKERRQ(ierr);
217*e14861a4SHong Zhang         ierr = VecSetSizes(yy1,am,PETSC_DECIDE);CHKERRQ(ierr);
218*e14861a4SHong Zhang         ierr = VecSetFromOptions(yy1);CHKERRQ(ierr);
219*e14861a4SHong Zhang         ierr = VecDuplicate(yy1,&yy2);CHKERRQ(ierr);
220*e14861a4SHong Zhang 
221*e14861a4SHong Zhang         PetscRandom       rd;
222*e14861a4SHong Zhang         PetscReal         rnorm;
223*e14861a4SHong Zhang         PetscScalar       mone = -1.0;
224*e14861a4SHong Zhang         ierr = PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rd);CHKERRQ(ierr);
225*e14861a4SHong Zhang         ierr = VecSetRandom(rd,x);CHKERRQ(ierr);
226*e14861a4SHong Zhang         ierr = MatMult(gridctx[level].A,x,yy1);CHKERRQ(ierr);
227*e14861a4SHong Zhang         ierr = MatMult(A_tmp,x,yy2);CHKERRQ(ierr);
228*e14861a4SHong Zhang         ierr = VecAXPY(&mone,yy1,yy2);CHKERRQ(ierr);
229*e14861a4SHong Zhang         ierr = VecNorm(yy2,NORM_2,&rnorm);CHKERRQ(ierr);
230*e14861a4SHong Zhang         printf(" [%d] mllevel: %d rnorm %g\n",rank,mllevel,rnorm);
231*e14861a4SHong Zhang 
232*e14861a4SHong Zhang         ierr = MatDestroy(A_tmp);CHKERRQ(ierr);
233*e14861a4SHong Zhang         ierr = VecDestroy(x);CHKERRQ(ierr);
234*e14861a4SHong Zhang         ierr = VecDestroy(yy1);CHKERRQ(ierr);
235*e14861a4SHong Zhang         ierr = VecDestroy(yy2);CHKERRQ(ierr);
236*e14861a4SHong Zhang         ierr = PetscRandomDestroy(rd);CHKERRQ(ierr);
237*e14861a4SHong Zhang     }
238*e14861a4SHong Zhang #endif
2395582bec1SHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
2405582bec1SHong Zhang       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr);
2415582bec1SHong Zhang       level--;
2425582bec1SHong Zhang     }
2435582bec1SHong Zhang   }
2445582bec1SHong Zhang 
2455582bec1SHong Zhang   /* create coarse level and the interpolation between the levels */
2465582bec1SHong Zhang   level = fine_level;
2475582bec1SHong Zhang   while ( level >= 0 ){
2485582bec1SHong Zhang     if (level != fine_level){
2495582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr);
2505582bec1SHong Zhang       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->n,PETSC_DECIDE);CHKERRQ(ierr);
2515582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
2525582bec1SHong Zhang       ierr = MGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2535582bec1SHong Zhang 
2545582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr);
2555582bec1SHong Zhang       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr);
2565582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
2575582bec1SHong Zhang       ierr = MGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
2585582bec1SHong Zhang     }
2595582bec1SHong Zhang     ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].r);CHKERRQ(ierr);
2605582bec1SHong Zhang     ierr = VecSetSizes(gridctx[level].r,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr);
2615582bec1SHong Zhang     ierr = VecSetType(gridctx[level].r,VECMPI);CHKERRQ(ierr);
2625582bec1SHong Zhang     ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr);
2635582bec1SHong Zhang 
2645582bec1SHong Zhang     if (level == 0){
2655582bec1SHong Zhang       ierr = MGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2665582bec1SHong Zhang     } else {
2675582bec1SHong Zhang       ierr = MGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
2685582bec1SHong Zhang       ierr = MGSetResidual(pc,level,MGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2695582bec1SHong Zhang       if (level == fine_level){
2705582bec1SHong Zhang         ierr = KSPSetOptionsPrefix(gridctx[level].ksp,"mg_fine_");CHKERRQ(ierr);
2715582bec1SHong Zhang         ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr);
2725582bec1SHong Zhang       }
2735582bec1SHong Zhang     }
2745582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2755582bec1SHong Zhang 
2765582bec1SHong Zhang     if (level < fine_level){
2775582bec1SHong Zhang       if (size > 1){
2785582bec1SHong Zhang         ierr = KSPGetPC(gridctx[level].ksp,&pc_coarse);CHKERRQ(ierr);
2795582bec1SHong Zhang         ierr = PCSetType(pc_coarse,PCNONE);CHKERRQ(ierr);
2805582bec1SHong Zhang       }
2815582bec1SHong Zhang       ierr = MGSetInterpolate(pc,level+1,gridctx[level].P);CHKERRQ(ierr);
2825582bec1SHong Zhang       ierr = MGSetRestriction(pc,level+1,gridctx[level].R);CHKERRQ(ierr);
2835582bec1SHong Zhang     }
2845582bec1SHong Zhang     level--;
2855582bec1SHong Zhang   }
2865582bec1SHong Zhang 
2875582bec1SHong Zhang   /* now call PCSetUp_MG()         */
2885582bec1SHong Zhang   /*--------------------------------*/
2895582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
2905582bec1SHong Zhang   PetscFunctionReturn(0);
2915582bec1SHong Zhang }
2925582bec1SHong Zhang 
2935582bec1SHong Zhang #undef __FUNCT__
2945582bec1SHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML"
2955582bec1SHong Zhang PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr)
2965582bec1SHong Zhang {
2975582bec1SHong Zhang   PetscErrorCode       ierr;
2985582bec1SHong Zhang   PC_ML                *pc_ml = (PC_ML*)ptr;
2995582bec1SHong Zhang   PetscInt             level;
3005582bec1SHong Zhang 
3015582bec1SHong Zhang   PetscFunctionBegin;
3025582bec1SHong Zhang   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
3035582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
3045582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
3055582bec1SHong Zhang 
3065582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
30724a42b14SHong Zhang   ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);
30824a42b14SHong Zhang   ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);
309*e14861a4SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata->cols);CHKERRQ(ierr);
3105582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
3115582bec1SHong Zhang 
3125582bec1SHong Zhang   level = pc_ml->fine_level;
3135582bec1SHong Zhang   while ( level>= 0 ){
3145582bec1SHong Zhang     if (level != pc_ml->fine_level){
3155582bec1SHong Zhang       ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);
3165582bec1SHong Zhang       ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);
3175582bec1SHong Zhang       ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);
3185582bec1SHong Zhang       ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);
3195582bec1SHong Zhang       ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);
3205582bec1SHong Zhang     }
3215582bec1SHong Zhang     ierr = VecDestroy(pc_ml->gridctx[level].r);CHKERRQ(ierr);
3225582bec1SHong Zhang     level--;
3235582bec1SHong Zhang   }
3245582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
3255582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
3265582bec1SHong Zhang   PetscFunctionReturn(0);
3275582bec1SHong Zhang }
3285582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3295582bec1SHong Zhang /*
3305582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3315582bec1SHong Zhang    that was created with PCCreate_ML().
3325582bec1SHong Zhang 
3335582bec1SHong Zhang    Input Parameter:
3345582bec1SHong Zhang .  pc - the preconditioner context
3355582bec1SHong Zhang 
3365582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3375582bec1SHong Zhang */
3385582bec1SHong Zhang #undef __FUNCT__
3395582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3405582bec1SHong Zhang static PetscErrorCode PCDestroy_ML(PC pc)
3415582bec1SHong Zhang {
3425582bec1SHong Zhang   PetscErrorCode       ierr;
3435582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
3445582bec1SHong Zhang   PetscObjectContainer container;
3455582bec1SHong Zhang 
3465582bec1SHong Zhang   PetscFunctionBegin;
3475582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3485582bec1SHong Zhang   if (container) {
3495582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3505582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
3515582bec1SHong Zhang   } else {
3525582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3535582bec1SHong Zhang   }
3549cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
3555582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
3565582bec1SHong Zhang   ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
3575582bec1SHong Zhang 
3585582bec1SHong Zhang   ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
3595582bec1SHong Zhang   PetscFunctionReturn(0);
3605582bec1SHong Zhang }
3615582bec1SHong Zhang 
3625582bec1SHong Zhang #undef __FUNCT__
3635582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3645582bec1SHong Zhang static PetscErrorCode PCSetFromOptions_ML(PC pc)
3655582bec1SHong Zhang {
3665582bec1SHong Zhang   PetscErrorCode ierr;
3675582bec1SHong Zhang   PetscInt       indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
3685582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
3695582bec1SHong Zhang   PetscTruth     flg;
3705582bec1SHong Zhang   const char     *type[] = {"additive","multiplicative","full","cascade","kascade"};
3715582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3725582bec1SHong Zhang   PC_ML          *pc_ml=PETSC_NULL;
3735582bec1SHong Zhang   PetscObjectContainer container;
3745582bec1SHong Zhang 
3755582bec1SHong Zhang   PetscFunctionBegin;
3765582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3775582bec1SHong Zhang   if (container) {
3785582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3795582bec1SHong Zhang   } else {
3805582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3815582bec1SHong Zhang   }
3825582bec1SHong Zhang   ierr = PetscOptionsHead("MG options");CHKERRQ(ierr);
3835582bec1SHong Zhang   /* inherited MG options */
3845582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3855582bec1SHong Zhang   if (flg) {
3865582bec1SHong Zhang     ierr = MGSetCycles(pc,m);CHKERRQ(ierr);
3875582bec1SHong Zhang   }
3885582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3895582bec1SHong Zhang   if (flg) {
3905582bec1SHong Zhang     ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
3915582bec1SHong Zhang   }
3925582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3935582bec1SHong Zhang   if (flg) {
3945582bec1SHong Zhang     ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
3955582bec1SHong Zhang   }
3965582bec1SHong Zhang   ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr);
3975582bec1SHong Zhang   if (flg) {
3985582bec1SHong Zhang     MGType mg = (MGType) 0;
3995582bec1SHong Zhang     switch (indx) {
4005582bec1SHong Zhang     case 0:
4015582bec1SHong Zhang       mg = MGADDITIVE;
4025582bec1SHong Zhang       break;
4035582bec1SHong Zhang     case 1:
4045582bec1SHong Zhang       mg = MGMULTIPLICATIVE;
4055582bec1SHong Zhang       break;
4065582bec1SHong Zhang     case 2:
4075582bec1SHong Zhang       mg = MGFULL;
4085582bec1SHong Zhang       break;
4095582bec1SHong Zhang     case 3:
4105582bec1SHong Zhang       mg = MGKASKADE;
4115582bec1SHong Zhang       break;
4125582bec1SHong Zhang     case 4:
4135582bec1SHong Zhang       mg = MGKASKADE;
4145582bec1SHong Zhang       break;
4155582bec1SHong Zhang     }
4165582bec1SHong Zhang     ierr = MGSetType(pc,mg);CHKERRQ(ierr);
4175582bec1SHong Zhang   }
4185582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4195582bec1SHong Zhang 
4205582bec1SHong Zhang   /* ML options */
4215582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
4225582bec1SHong Zhang   /* set defaults */
4235582bec1SHong Zhang   PrintLevel    = 0;
4245582bec1SHong Zhang   MaxNlevels    = 10;
4255582bec1SHong Zhang   MaxCoarseSize = 1;
4265582bec1SHong Zhang   indx          = 0;
4275582bec1SHong Zhang   Threshold     = 0.0;
4285582bec1SHong Zhang   DampingFactor = 4.0/3.0;
4295582bec1SHong Zhang 
4305582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
4315582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
4325582bec1SHong Zhang 
4335582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
4345582bec1SHong Zhang   pc_ml->MaxNlevels = MaxNlevels;
4355582bec1SHong Zhang 
4365582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
4375582bec1SHong Zhang   pc_ml->MaxCoarseSize = MaxCoarseSize;
4385582bec1SHong Zhang 
4395582bec1SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
4405582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
4415582bec1SHong Zhang 
4425582bec1SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr);
4435582bec1SHong Zhang   pc_ml->DampingFactor = DampingFactor;
4445582bec1SHong Zhang 
4455582bec1SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr);
4465582bec1SHong Zhang   pc_ml->Threshold = Threshold;
4475582bec1SHong Zhang 
4485582bec1SHong Zhang   ierr = PetscOptionsLogical("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_FALSE);
4495582bec1SHong Zhang 
4505582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4515582bec1SHong Zhang   PetscFunctionReturn(0);
4525582bec1SHong Zhang }
4535582bec1SHong Zhang 
4545582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4555582bec1SHong Zhang /*
4565582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4575582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4585582bec1SHong Zhang    context, PC, that was created within PCCreate().
4595582bec1SHong Zhang 
4605582bec1SHong Zhang    Input Parameter:
4615582bec1SHong Zhang .  pc - the preconditioner context
4625582bec1SHong Zhang 
4635582bec1SHong Zhang    Application Interface Routine: PCCreate()
4645582bec1SHong Zhang */
4655582bec1SHong Zhang 
4665582bec1SHong Zhang /*MC
4675582bec1SHong Zhang      PCML - Use geometric multigrid preconditioning. This preconditioner requires you provide
4685582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4695582bec1SHong Zhang        operators are computed by ML and wrapped as PETSc shell matrices.
4705582bec1SHong Zhang 
4715582bec1SHong Zhang    Options Database Key: (not done yet!)
4725582bec1SHong Zhang +  -pc_mg_maxlevels <nlevels> - maximum number of levels including finest
4735582bec1SHong Zhang .  -pc_mg_cycles 1 or 2 - for V or W-cycle
4745582bec1SHong Zhang .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
4755582bec1SHong Zhang .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
4765582bec1SHong Zhang .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
4775582bec1SHong Zhang .  -pc_mg_monitor - print information on the multigrid convergence
4785582bec1SHong Zhang -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
4795582bec1SHong Zhang                         to the Socket viewer for reading from Matlab.
4805582bec1SHong Zhang 
4815582bec1SHong Zhang    Level: intermediate
4825582bec1SHong Zhang 
4835582bec1SHong Zhang   Concepts: multigrid
4845582bec1SHong Zhang 
4855582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
4865582bec1SHong Zhang            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
4875582bec1SHong Zhang            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
4885582bec1SHong Zhang            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
4895582bec1SHong Zhang            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()
4905582bec1SHong Zhang M*/
4915582bec1SHong Zhang 
4925582bec1SHong Zhang EXTERN_C_BEGIN
4935582bec1SHong Zhang #undef __FUNCT__
4945582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
4955582bec1SHong Zhang PetscErrorCode PCCreate_ML(PC pc)
4965582bec1SHong Zhang {
4975582bec1SHong Zhang   PetscErrorCode       ierr;
4985582bec1SHong Zhang   PC_ML                *pc_ml;
4995582bec1SHong Zhang   PetscObjectContainer container;
5005582bec1SHong Zhang 
5015582bec1SHong Zhang   PetscFunctionBegin;
5025582bec1SHong Zhang   /* initialize pc as PCMG */
5035582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
5045582bec1SHong Zhang 
5055582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
5065582bec1SHong Zhang   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
5079cb0ec93SHong Zhang   ierr = PetscMemzero(pc_ml,sizeof(PC_ML));CHKERRQ(ierr);
5085582bec1SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
5095582bec1SHong Zhang   ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
5109cb0ec93SHong Zhang   ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr);
5115582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
5125582bec1SHong Zhang 
5135582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
5145582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
5155582bec1SHong Zhang 
5165582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
5175582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
5185582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
5195582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
5205582bec1SHong Zhang   PetscFunctionReturn(0);
5215582bec1SHong Zhang }
5225582bec1SHong Zhang EXTERN_C_END
5235582bec1SHong Zhang 
5245582bec1SHong Zhang int PetscML_getrow(void *ML_data, int N_requested_rows, int requested_rows[],
5255582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
5265582bec1SHong Zhang {
5275582bec1SHong Zhang   PetscErrorCode ierr;
5285582bec1SHong Zhang   Mat            Aloc;
5295582bec1SHong Zhang   Mat_SeqAIJ     *a;
5305582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
5315582bec1SHong Zhang   PetscScalar    *aa;
5325582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5335582bec1SHong Zhang 
5345582bec1SHong Zhang   Aloc = ml->Aloc;
5355582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
5365582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
5375582bec1SHong Zhang 
5385582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
5395582bec1SHong Zhang     row   = requested_rows[i];
5405582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
5415582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
5425582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
5435582bec1SHong Zhang       aj = a->j + a->i[row];
5445582bec1SHong Zhang       aa = a->a + a->i[row];
5455582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
5465582bec1SHong Zhang         columns[k]  = aj[j];
5475582bec1SHong Zhang         values[k++] = aa[j];
5485582bec1SHong Zhang       }
5495582bec1SHong Zhang     }
5505582bec1SHong Zhang   }
5515582bec1SHong Zhang   return(1);
5525582bec1SHong Zhang }
5535582bec1SHong Zhang 
5545582bec1SHong Zhang int PetscML_matvec(void *ML_data,int in_length,double p[],int out_length,double ap[])
5555582bec1SHong Zhang {
5565582bec1SHong Zhang   PetscErrorCode ierr;
5575582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5585582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5595582bec1SHong Zhang   PetscMPIInt    size;
5605582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5615582bec1SHong Zhang   PetscInt       i;
5625582bec1SHong Zhang 
5635582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5645582bec1SHong Zhang   if (size == 1){
56524a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5665582bec1SHong Zhang   } else {
5675582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5685582bec1SHong Zhang     PetscML_comm(pwork,ml);
56924a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5705582bec1SHong Zhang   }
57124a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
57224a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
5735582bec1SHong Zhang   return 0;
5745582bec1SHong Zhang }
5755582bec1SHong Zhang 
5765582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5775582bec1SHong Zhang {
5785582bec1SHong Zhang   PetscErrorCode ierr;
5795582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5805582bec1SHong Zhang   Mat            A=ml->A;
5815582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5825582bec1SHong Zhang   PetscMPIInt    size;
5835582bec1SHong Zhang   PetscInt       i,in_length=A->m,out_length=ml->Aloc->n;
5845582bec1SHong Zhang   PetscScalar    *array;
5855582bec1SHong Zhang 
5865582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5875582bec1SHong Zhang   if (size == 1) return 0;
58824a42b14SHong Zhang 
58924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
59024a42b14SHong Zhang   ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
59124a42b14SHong Zhang   ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
5925582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5935582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5945582bec1SHong Zhang     p[i] = array[i-in_length];
5955582bec1SHong Zhang   }
5965582bec1SHong Zhang   return 0;
5975582bec1SHong Zhang }
5985582bec1SHong Zhang #undef __FUNCT__
5995582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
6005582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
6015582bec1SHong Zhang {
6025582bec1SHong Zhang   PetscErrorCode   ierr;
6035582bec1SHong Zhang   Mat_MLShell      *shell;
6045582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
6055582bec1SHong Zhang   PetscInt         x_length,y_length;
6065582bec1SHong Zhang 
6075582bec1SHong Zhang   PetscFunctionBegin;
6085582bec1SHong Zhang   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
6095582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6105582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6115582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6125582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6135582bec1SHong Zhang 
6145582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6155582bec1SHong Zhang 
6165582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6175582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
6185582bec1SHong Zhang   PetscFunctionReturn(0);
6195582bec1SHong Zhang }
6205582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
6215582bec1SHong Zhang #undef __FUNCT__
6225582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
6235582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
6245582bec1SHong Zhang {
6255582bec1SHong Zhang   PetscErrorCode    ierr;
6265582bec1SHong Zhang   Mat_MLShell       *shell;
6275582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
6285582bec1SHong Zhang   const PetscScalar one=1.0;
6295582bec1SHong Zhang   PetscInt          x_length,y_length;
6305582bec1SHong Zhang 
6315582bec1SHong Zhang   PetscFunctionBegin;
6325582bec1SHong Zhang   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
6335582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6345582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6355582bec1SHong Zhang 
6365582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6375582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6385582bec1SHong Zhang 
6395582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6405582bec1SHong Zhang 
6415582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6425582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
6435582bec1SHong Zhang   ierr = VecAXPY(&one,w,y);CHKERRQ(ierr);
6445582bec1SHong Zhang 
6455582bec1SHong Zhang   PetscFunctionReturn(0);
6465582bec1SHong Zhang }
6475582bec1SHong Zhang 
6485582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
6495582bec1SHong Zhang #undef __FUNCT__
6505582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
6515582bec1SHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,const MatType newtype,Mat *Aloc)
6525582bec1SHong Zhang {
6535582bec1SHong Zhang   PetscErrorCode  ierr;
6545582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6555582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6565582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6575582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
6585582bec1SHong Zhang   PetscInt        am=A->m,an=A->n,i,j,k;
6595582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6605582bec1SHong Zhang   MatReuse        scall=MAT_INITIAL_MATRIX;
6615582bec1SHong Zhang 
6625582bec1SHong Zhang   PetscFunctionBegin;
6635582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6645582bec1SHong Zhang 
6655582bec1SHong Zhang   if (*Aloc) scall = MAT_REUSE_MATRIX;
6665582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6675582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6685582bec1SHong Zhang     ci[0] = 0;
6695582bec1SHong Zhang     for (i=0; i<am; i++){
6705582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6715582bec1SHong Zhang     }
6725582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6735582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6745582bec1SHong Zhang 
6755582bec1SHong Zhang     k = 0;
6765582bec1SHong Zhang     for (i=0; i<am; i++){
6775582bec1SHong Zhang       /* diagonal portion of A */
6785582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6795582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6805582bec1SHong Zhang         cj[k]   = *aj++;
6815582bec1SHong Zhang         ca[k++] = *aa++;
6825582bec1SHong Zhang       }
6835582bec1SHong Zhang       /* off-diagonal portion of A */
6845582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6855582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6865582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6875582bec1SHong Zhang         ca[k++] = *ba++;
6885582bec1SHong Zhang       }
6895582bec1SHong Zhang     }
6905582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6915582bec1SHong Zhang 
6925582bec1SHong Zhang     /* put together the new matrix */
6935582bec1SHong Zhang     an = mpimat->A->n+mpimat->B->n;
6945582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6955582bec1SHong Zhang 
6965582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6975582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6985582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6995582bec1SHong Zhang     mat->freedata = PETSC_TRUE;
7005582bec1SHong Zhang     mat->nonew    = 0;
7015582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
7025582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
7035582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
7045582bec1SHong Zhang     for (i=0; i<am; i++) {
7055582bec1SHong Zhang       /* diagonal portion of A */
7065582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
7075582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
7085582bec1SHong Zhang       /* off-diagonal portion of A */
7095582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
7105582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
7115582bec1SHong Zhang     }
7125582bec1SHong Zhang   } else {
7135582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
7145582bec1SHong Zhang   }
7155582bec1SHong Zhang   PetscFunctionReturn(0);
7165582bec1SHong Zhang }
7175582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
7185582bec1SHong Zhang #undef __FUNCT__
7195582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
7205582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
7215582bec1SHong Zhang {
7225582bec1SHong Zhang   PetscErrorCode ierr;
7235582bec1SHong Zhang   Mat_MLShell    *shell;
7245582bec1SHong Zhang 
7255582bec1SHong Zhang   PetscFunctionBegin;
7265582bec1SHong Zhang   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
7275582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
7285582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
7295582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
7305582bec1SHong Zhang   PetscFunctionReturn(0);
7315582bec1SHong Zhang }
7325582bec1SHong Zhang 
733*e14861a4SHong Zhang extern PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt [],PetscScalar []);
7345582bec1SHong Zhang #undef __FUNCT__
7355582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SeqAIJ"
736*e14861a4SHong Zhang PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator *mlmat,Mat *newmat)
7375582bec1SHong Zhang {
738*e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
7395582bec1SHong Zhang   PetscErrorCode  ierr;
740*e14861a4SHong Zhang   PetscInt        m=mlmat->outvec_leng,n=mlmat->invec_leng,nnz[m],nz_max;
741*e14861a4SHong Zhang   PetscInt        *ml_cols=matdata->columns,*aj,i,j,k;
742*e14861a4SHong Zhang   PetscScalar     *ml_vals=matdata->values,*aa;
7435582bec1SHong Zhang 
7445582bec1SHong Zhang   PetscFunctionBegin;
745*e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
74624a42b14SHong Zhang   if (m != n){ /* pass array pointers if ml->mlmat is Pmat or Rmat */
747*e14861a4SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
74824a42b14SHong Zhang     PetscFunctionReturn(0);
74924a42b14SHong Zhang   }
7505582bec1SHong Zhang 
751*e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
7525582bec1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,m,n,PETSC_DECIDE,PETSC_DECIDE,newmat);CHKERRQ(ierr);
7535582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
754*e14861a4SHong Zhang 
755*e14861a4SHong Zhang   nz_max = 0;
7565582bec1SHong Zhang   for (i=0; i<m; i++) {
757*e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
758*e14861a4SHong Zhang     if (nnz[i] > nz_max) nz_max = nnz[i];
7595582bec1SHong Zhang   }
7605582bec1SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
7615582bec1SHong Zhang 
762*e14861a4SHong Zhang   nz_max++;
763*e14861a4SHong Zhang   ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
764*e14861a4SHong Zhang   aa = (PetscScalar*)(aj + nz_max);
76524a42b14SHong Zhang 
7665582bec1SHong Zhang   for (i=0; i<m; i++){
767*e14861a4SHong Zhang     k = 0;
768*e14861a4SHong Zhang     /* diagonal entry */
769*e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
770*e14861a4SHong Zhang     /* off diagonal entries */
771*e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
772*e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
77324a42b14SHong Zhang     }
774*e14861a4SHong Zhang     /* sort aj and aa ??? */
775*e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
776*e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7775582bec1SHong Zhang   }
7785582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7795582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
780*e14861a4SHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
7815582bec1SHong Zhang   PetscFunctionReturn(0);
7825582bec1SHong Zhang }
7835582bec1SHong Zhang 
7845582bec1SHong Zhang #undef __FUNCT__
7855582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SHELL"
7865582bec1SHong Zhang PetscErrorCode MatConvert_ML_SHELL(ML_Operator *mlmat,Mat *newmat)
7875582bec1SHong Zhang {
7885582bec1SHong Zhang   PetscErrorCode ierr;
7895582bec1SHong Zhang   PetscInt       m,n;
7905582bec1SHong Zhang   ML_Comm        *MLcomm;
7915582bec1SHong Zhang   Mat_MLShell    *shellctx;
7925582bec1SHong Zhang 
7935582bec1SHong Zhang   PetscFunctionBegin;
7945582bec1SHong Zhang   m = mlmat->outvec_leng;
7955582bec1SHong Zhang   n = mlmat->invec_leng;
7965582bec1SHong Zhang   if (!m || !n){
7975582bec1SHong Zhang     newmat = PETSC_NULL;
7985582bec1SHong Zhang   } else {
7995582bec1SHong Zhang     MLcomm = mlmat->comm;
8005582bec1SHong Zhang     ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
8015582bec1SHong Zhang     ierr = PetscMemzero(shellctx,sizeof(Mat_MLShell));CHKERRQ(ierr);
8025582bec1SHong Zhang     ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
8035582bec1SHong Zhang     ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void *)MatMult_ML);CHKERRQ(ierr);
8045582bec1SHong Zhang     ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void *)MatMultAdd_ML);CHKERRQ(ierr);
8055582bec1SHong Zhang     shellctx->A         = *newmat;
8065582bec1SHong Zhang     shellctx->mlmat     = mlmat;
8075582bec1SHong Zhang     ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
8085582bec1SHong Zhang     ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
8095582bec1SHong Zhang     ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
8105582bec1SHong Zhang     (*newmat)->ops->destroy = MatDestroy_ML;
8115582bec1SHong Zhang   }
8125582bec1SHong Zhang   PetscFunctionReturn(0);
8135582bec1SHong Zhang }
814*e14861a4SHong Zhang 
815*e14861a4SHong Zhang #undef __FUNCT__
816*e14861a4SHong Zhang #define __FUNCT__ "MatConvert_ML_MPIAIJ"
817*e14861a4SHong Zhang PetscErrorCode MatConvert_ML_MPIAIJ(FineGridCtx *ml,Mat *newmat)
818*e14861a4SHong Zhang {
819*e14861a4SHong Zhang   PetscErrorCode  ierr;
820*e14861a4SHong Zhang   ML_Operator     *mat=ml->mlmat;
821*e14861a4SHong Zhang   PetscInt        i,j,*cols=ml->cols,cstart,cend,rlen_max,*gordering;
822*e14861a4SHong Zhang   PetscInt        m=mat->outvec_leng,n,nnzA[m],nnzB[m],nnz[m],nz_max,row,col,*cols_tmp;
823*e14861a4SHong Zhang   PetscScalar     *vals=ml->vals;
824*e14861a4SHong Zhang   Mat             C;
825*e14861a4SHong Zhang   Mat_MPIAIJ      *c;
826*e14861a4SHong Zhang   PetscMPIInt     rank;
827*e14861a4SHong Zhang 
828*e14861a4SHong Zhang   PetscFunctionBegin;
829*e14861a4SHong Zhang   if ( mat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mat->getrow = NULL");
830*e14861a4SHong Zhang 
831*e14861a4SHong Zhang   ML_build_global_numbering(mat,mat->comm,&gordering);
832*e14861a4SHong Zhang   ierr = MPI_Comm_rank(ml->A->comm,&rank);CHKERRQ(ierr);
833*e14861a4SHong Zhang   n = mat->invec_leng;
834*e14861a4SHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] m: %d, %d; n: %d,\n",rank,m,mat->getrow->Nrows,n);CHKERRQ(ierr);*/
835*e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
836*e14861a4SHong Zhang 
837*e14861a4SHong Zhang   ierr = MatCreate(ml->A->comm,m,n,PETSC_DECIDE,PETSC_DECIDE,&C);CHKERRQ(ierr);
838*e14861a4SHong Zhang   ierr = MatSetType(C,MATMPIAIJ);CHKERRQ(ierr);
839*e14861a4SHong Zhang   c        = (Mat_MPIAIJ*)C->data;
840*e14861a4SHong Zhang   cstart   = c->cstart;
841*e14861a4SHong Zhang   cend     = c->cend;
842*e14861a4SHong Zhang   rlen_max = C->N;
843*e14861a4SHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]  cstart/end: %d %d, C->N: %d\n",rank,cstart,cend,C->N); */
844*e14861a4SHong Zhang 
845*e14861a4SHong Zhang   ierr = PetscMalloc(nz_max*sizeof(int),&cols_tmp);CHKERRQ(ierr);
846*e14861a4SHong Zhang 
847*e14861a4SHong Zhang   nz_max = 0;
848*e14861a4SHong Zhang   for (i=0; i<m; i++){
849*e14861a4SHong Zhang     ML_get_matrix_row(mat,1,&i,&rlen_max,&cols,&vals,&nnz[i],0);
850*e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
851*e14861a4SHong Zhang     nnzA[i] = 0;
852*e14861a4SHong Zhang     for (j=0; j<nnz[i]; j++){
853*e14861a4SHong Zhang       col = gordering[cols[j]];
854*e14861a4SHong Zhang       if (cstart <= col && col < cend ) nnzA[i]++;
855*e14861a4SHong Zhang     }
856*e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
857*e14861a4SHong Zhang   }
858*e14861a4SHong Zhang   ierr = MatMPIAIJSetPreallocation(C,0,nnzA,0,nnzB);CHKERRQ(ierr);
859*e14861a4SHong Zhang 
860*e14861a4SHong Zhang   /* insert values -- remap row and column indices */
861*e14861a4SHong Zhang   for (i=0; i<m; i++){
862*e14861a4SHong Zhang     ML_get_matrix_row(mat,1,&i,&rlen_max,&cols_tmp,&vals,&nnz[i],0);
863*e14861a4SHong Zhang     for (j=0; j<nnz[i]; j++){
864*e14861a4SHong Zhang       cols[j] = gordering[cols_tmp[j]];
865*e14861a4SHong Zhang     }
866*e14861a4SHong Zhang     row = gordering[i];
867*e14861a4SHong Zhang     ierr = MatSetValues(C,1,&row,nnz[i],cols,ml->vals,INSERT_VALUES);CHKERRQ(ierr);
868*e14861a4SHong Zhang   }
869*e14861a4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
870*e14861a4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
871*e14861a4SHong Zhang   *newmat = C;
872*e14861a4SHong Zhang 
873*e14861a4SHong Zhang   ierr = PetscFree(cols_tmp);CHKERRQ(ierr);
874*e14861a4SHong Zhang   PetscFunctionReturn(0);
875*e14861a4SHong Zhang }
876*e14861a4SHong Zhang 
877*e14861a4SHong Zhang /* -----------------------------------------------------------------------*/
878*e14861a4SHong Zhang #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;}
879*e14861a4SHong Zhang 
880*e14861a4SHong Zhang #undef __FUNCT__
881*e14861a4SHong Zhang #define __FUNCT__ "PetscSortIntWithScalarArray_Private"
882*e14861a4SHong Zhang /*
883*e14861a4SHong Zhang    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
884*e14861a4SHong Zhang    Assumes 0 origin for v, number of elements = right+1 (right is index of
885*e14861a4SHong Zhang    right-most member).
886*e14861a4SHong Zhang */
887*e14861a4SHong Zhang static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right)
888*e14861a4SHong Zhang {
889*e14861a4SHong Zhang   PetscErrorCode ierr;
890*e14861a4SHong Zhang   PetscInt       i,vl,last,tmp;
891*e14861a4SHong Zhang   PetscScalar    stmp;
892*e14861a4SHong Zhang 
893*e14861a4SHong Zhang   PetscFunctionBegin;
894*e14861a4SHong Zhang   if (right <= 1) {
895*e14861a4SHong Zhang     if (right == 1) {
896*e14861a4SHong Zhang       if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp);
897*e14861a4SHong Zhang     }
898*e14861a4SHong Zhang     PetscFunctionReturn(0);
899*e14861a4SHong Zhang   }
900*e14861a4SHong Zhang   SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp);
901*e14861a4SHong Zhang   vl   = v[0];
902*e14861a4SHong Zhang   last = 0;
903*e14861a4SHong Zhang   for (i=1; i<=right; i++) {
904*e14861a4SHong Zhang     if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);}
905*e14861a4SHong Zhang   }
906*e14861a4SHong Zhang   SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp);
907*e14861a4SHong Zhang   ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr);
908*e14861a4SHong Zhang   ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
909*e14861a4SHong Zhang   PetscFunctionReturn(0);
910*e14861a4SHong Zhang }
911*e14861a4SHong Zhang 
912*e14861a4SHong Zhang #undef __FUNCT__
913*e14861a4SHong Zhang #define __FUNCT__ "PetscSortIntWithScalarArray"
914*e14861a4SHong Zhang /*@
915*e14861a4SHong Zhang    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
916*e14861a4SHong Zhang        changes a second array to match the sorted first array.
917*e14861a4SHong Zhang 
918*e14861a4SHong Zhang    Not Collective
919*e14861a4SHong Zhang 
920*e14861a4SHong Zhang    Input Parameters:
921*e14861a4SHong Zhang +  n  - number of values
922*e14861a4SHong Zhang .  i  - array of integers
923*e14861a4SHong Zhang -  I - second array of integers
924*e14861a4SHong Zhang 
925*e14861a4SHong Zhang    Level: intermediate
926*e14861a4SHong Zhang 
927*e14861a4SHong Zhang    Concepts: sorting^ints with array
928*e14861a4SHong Zhang 
929*e14861a4SHong Zhang .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
930*e14861a4SHong Zhang @*/
931*e14861a4SHong Zhang PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar I[])
932*e14861a4SHong Zhang {
933*e14861a4SHong Zhang   PetscErrorCode ierr;
934*e14861a4SHong Zhang   PetscInt       j,k,tmp,ik;
935*e14861a4SHong Zhang   PetscScalar    stmp;
936*e14861a4SHong Zhang 
937*e14861a4SHong Zhang   PetscFunctionBegin;
938*e14861a4SHong Zhang   if (n<8) {
939*e14861a4SHong Zhang     for (k=0; k<n; k++) {
940*e14861a4SHong Zhang       ik = i[k];
941*e14861a4SHong Zhang       for (j=k+1; j<n; j++) {
942*e14861a4SHong Zhang 	if (ik > i[j]) {
943*e14861a4SHong Zhang 	  SWAP2IntScalar(i[k],i[j],I[k],I[j],tmp,stmp);
944*e14861a4SHong Zhang 	  ik = i[k];
945*e14861a4SHong Zhang 	}
946*e14861a4SHong Zhang       }
947*e14861a4SHong Zhang     }
948*e14861a4SHong Zhang   } else {
949*e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray_Private(i,I,n-1);CHKERRQ(ierr);
950*e14861a4SHong Zhang   }
951*e14861a4SHong Zhang   PetscFunctionReturn(0);
952*e14861a4SHong Zhang }
953