xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision fb6a8e6dfce6eb064cf7517c8c640342ca88b023)
1ab718edeSHong Zhang 
25582bec1SHong Zhang /*
32dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
47ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
57ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
65582bec1SHong Zhang */
7b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>   /*I "petscpc.h" I*/
8c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
9c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
10c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
11cb5d8e9eSHong Zhang 
125582bec1SHong Zhang #include <math.h>
132cf39c26SSatish Balay EXTERN_C_BEGIN
1468210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1568210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1668210224SSatish Balay #define HAVE_CONFIG_H
1768210224SSatish Balay #endif
18c6db04a5SJed Brown #include <ml_include.h>
195582bec1SHong Zhang EXTERN_C_END
205582bec1SHong Zhang 
21*fb6a8e6dSJed Brown typedef enum {PCML_NULLSPACE_AUTO,PCML_NULLSPACE_USER,PCML_NULLSPACE_BLOCK,PCML_NULLSPACE_SCALAR} PCMLNullSpaceType;
22*fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO","USER","BLOCK","SCALAR","PCMLNullSpaceType","PCML_NULLSPACE_",0};
23*fb6a8e6dSJed Brown 
245582bec1SHong Zhang /* The context (data structure) at each grid level */
255582bec1SHong Zhang typedef struct {
265582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
275582bec1SHong Zhang   Mat        A,P,R;
285582bec1SHong Zhang   KSP        ksp;
295582bec1SHong Zhang } GridCtx;
305582bec1SHong Zhang 
315582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
325582bec1SHong Zhang typedef struct {
33573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
34573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3524a42b14SHong Zhang   Vec          x,y;
365582bec1SHong Zhang   ML_Operator  *mlmat;
375582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
385582bec1SHong Zhang } FineGridCtx;
395582bec1SHong Zhang 
405582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
415582bec1SHong Zhang typedef struct {
425582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
435582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
4467d6f150SMatthew G Knepley   Vec          y, work;
455582bec1SHong Zhang } Mat_MLShell;
465582bec1SHong Zhang 
475582bec1SHong Zhang /* Private context for the ML preconditioner */
485582bec1SHong Zhang typedef struct {
495582bec1SHong Zhang   ML             *ml_object;
505582bec1SHong Zhang   ML_Aggregate   *agg_object;
515582bec1SHong Zhang   GridCtx        *gridctx;
525582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
53b5c8bdf8SJed Brown   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization;
54b5c8bdf8SJed Brown   PetscReal      Threshold,DampingFactor,EnergyMinimizationDropTol;
55ace3abfcSBarry Smith   PetscBool      SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable;
5648268eb4SJed Brown   PetscBool      reuse_interpolation;
57*fb6a8e6dSJed Brown   PCMLNullSpaceType nulltype;
58573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
595582bec1SHong Zhang } PC_ML;
6041ca0015SHong Zhang 
616562c4e1SBarry Smith #undef __FUNCT__
626562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow"
636562c4e1SBarry Smith static int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],int allocated_space, int columns[], double values[], int row_lengths[])
646562c4e1SBarry Smith {
656562c4e1SBarry Smith   PetscErrorCode ierr;
666562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
676562c4e1SBarry Smith   PetscScalar    *aa;
686562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
696562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
705582bec1SHong Zhang 
716562c4e1SBarry Smith 
726562c4e1SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,PETSC_NULL); if (ierr) return(0);
736562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
746562c4e1SBarry Smith     row   = requested_rows[i];
756562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
766562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
776562c4e1SBarry Smith     if ( (row >= 0) || (row <= (m-1)) ) {
786562c4e1SBarry Smith       aj = a->j + a->i[row];
796562c4e1SBarry Smith       aa = a->a + a->i[row];
806562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++){
816562c4e1SBarry Smith         columns[k]  = aj[j];
826562c4e1SBarry Smith         values[k++] = aa[j];
836562c4e1SBarry Smith       }
846562c4e1SBarry Smith     }
856562c4e1SBarry Smith   }
866562c4e1SBarry Smith   return(1);
876562c4e1SBarry Smith }
886562c4e1SBarry Smith 
896562c4e1SBarry Smith #undef __FUNCT__
906562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm"
916562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
926562c4e1SBarry Smith {
936562c4e1SBarry Smith   PetscErrorCode ierr;
946562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_data;
956562c4e1SBarry Smith   Mat            A=ml->A;
966562c4e1SBarry Smith   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
976562c4e1SBarry Smith   PetscMPIInt    size;
986562c4e1SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
996562c4e1SBarry Smith   PetscScalar    *array;
1006562c4e1SBarry Smith 
1016562c4e1SBarry Smith   PetscFunctionBegin;
1026562c4e1SBarry Smith   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1036562c4e1SBarry Smith   if (size == 1) return 0;
1046562c4e1SBarry Smith 
1056562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
1066562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1076562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1086562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1096562c4e1SBarry Smith   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
1106562c4e1SBarry Smith   for (i=in_length; i<out_length; i++){
1116562c4e1SBarry Smith     p[i] = array[i-in_length];
1126562c4e1SBarry Smith   }
1136562c4e1SBarry Smith   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
1146562c4e1SBarry Smith   PetscFunctionReturn(0);
1156562c4e1SBarry Smith }
1166562c4e1SBarry Smith 
1176562c4e1SBarry Smith #undef __FUNCT__
1186562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec"
1196562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
1206562c4e1SBarry Smith {
1216562c4e1SBarry Smith   PetscErrorCode ierr;
1226562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
1236562c4e1SBarry Smith   Mat            A=ml->A, Aloc=ml->Aloc;
1246562c4e1SBarry Smith   PetscMPIInt    size;
1256562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
1266562c4e1SBarry Smith   PetscInt       i;
1276562c4e1SBarry Smith 
1286562c4e1SBarry Smith   PetscFunctionBegin;
1296562c4e1SBarry Smith   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1306562c4e1SBarry Smith   if (size == 1){
1316562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
1326562c4e1SBarry Smith   } else {
1336562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
1346562c4e1SBarry Smith     PetscML_comm(pwork,ml);
1356562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
1366562c4e1SBarry Smith   }
1376562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
1386562c4e1SBarry Smith   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
1396562c4e1SBarry Smith   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
1406562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1416562c4e1SBarry Smith   PetscFunctionReturn(0);
1426562c4e1SBarry Smith }
1436562c4e1SBarry Smith 
1446562c4e1SBarry Smith #undef __FUNCT__
1456562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML"
1466562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1476562c4e1SBarry Smith {
1486562c4e1SBarry Smith   PetscErrorCode   ierr;
1496562c4e1SBarry Smith   Mat_MLShell      *shell;
1506562c4e1SBarry Smith   PetscScalar      *xarray,*yarray;
1516562c4e1SBarry Smith   PetscInt         x_length,y_length;
1526562c4e1SBarry Smith 
1536562c4e1SBarry Smith   PetscFunctionBegin;
1546562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
1556562c4e1SBarry Smith   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
1566562c4e1SBarry Smith   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
1576562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1586562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
1596562c4e1SBarry Smith   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
1606562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1616562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1626562c4e1SBarry Smith   PetscFunctionReturn(0);
1636562c4e1SBarry Smith }
1646562c4e1SBarry Smith 
1656562c4e1SBarry Smith #undef __FUNCT__
1666562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML"
16767d6f150SMatthew G Knepley /* Computes y = w + A * x
16867d6f150SMatthew G Knepley    It is possible that w == y, but not x == y
16967d6f150SMatthew G Knepley */
1706562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
1716562c4e1SBarry Smith {
1726562c4e1SBarry Smith   Mat_MLShell   *shell;
1736562c4e1SBarry Smith   PetscScalar   *xarray,*yarray;
1746562c4e1SBarry Smith   PetscInt       x_length,y_length;
17567d6f150SMatthew G Knepley   PetscErrorCode ierr;
1766562c4e1SBarry Smith 
1776562c4e1SBarry Smith   PetscFunctionBegin;
1786562c4e1SBarry Smith   ierr = MatShellGetContext(A, (void **) &shell);CHKERRQ(ierr);
17967d6f150SMatthew G Knepley   if (y == w) {
18067d6f150SMatthew G Knepley     if (!shell->work) {
18167d6f150SMatthew G Knepley       ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr);
18267d6f150SMatthew G Knepley     }
18367d6f150SMatthew G Knepley     ierr = VecGetArray(x,           &xarray);CHKERRQ(ierr);
18467d6f150SMatthew G Knepley     ierr = VecGetArray(shell->work, &yarray);CHKERRQ(ierr);
18567d6f150SMatthew G Knepley     x_length = shell->mlmat->invec_leng;
18667d6f150SMatthew G Knepley     y_length = shell->mlmat->outvec_leng;
18767d6f150SMatthew G Knepley     ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray);
18867d6f150SMatthew G Knepley     ierr = VecRestoreArray(x,           &xarray);CHKERRQ(ierr);
18967d6f150SMatthew G Knepley     ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr);
1903ba3408dSMatthew G Knepley     ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr);
19167d6f150SMatthew G Knepley   } else {
1926562c4e1SBarry Smith     ierr = VecGetArray(x, &xarray);CHKERRQ(ierr);
1936562c4e1SBarry Smith     ierr = VecGetArray(y, &yarray);CHKERRQ(ierr);
1946562c4e1SBarry Smith     x_length = shell->mlmat->invec_leng;
1956562c4e1SBarry Smith     y_length = shell->mlmat->outvec_leng;
1966562c4e1SBarry Smith     ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray);
1976562c4e1SBarry Smith     ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr);
1986562c4e1SBarry Smith     ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr);
1996562c4e1SBarry Smith     ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr);
20067d6f150SMatthew G Knepley   }
2016562c4e1SBarry Smith   PetscFunctionReturn(0);
2026562c4e1SBarry Smith }
2036562c4e1SBarry Smith 
2046562c4e1SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */
2056562c4e1SBarry Smith #undef __FUNCT__
2066562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML"
2076562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
2086562c4e1SBarry Smith {
2096562c4e1SBarry Smith   PetscErrorCode  ierr;
2106562c4e1SBarry Smith   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
2116562c4e1SBarry Smith   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
2126562c4e1SBarry Smith   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2136562c4e1SBarry Smith   PetscScalar     *aa=a->a,*ba=b->a,*ca;
2146562c4e1SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
2156562c4e1SBarry Smith   PetscInt        *ci,*cj,ncols;
2166562c4e1SBarry Smith 
2176562c4e1SBarry Smith   PetscFunctionBegin;
218e32f2f54SBarry Smith   if (am != an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
2196562c4e1SBarry Smith 
2206562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX){
2216562c4e1SBarry Smith     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2226562c4e1SBarry Smith     ci[0] = 0;
2236562c4e1SBarry Smith     for (i=0; i<am; i++){
2246562c4e1SBarry Smith       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
2256562c4e1SBarry Smith     }
2266562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
2276562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
2286562c4e1SBarry Smith 
2296562c4e1SBarry Smith     k = 0;
2306562c4e1SBarry Smith     for (i=0; i<am; i++){
2316562c4e1SBarry Smith       /* diagonal portion of A */
2326562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2336562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2346562c4e1SBarry Smith         cj[k]   = *aj++;
2356562c4e1SBarry Smith         ca[k++] = *aa++;
2366562c4e1SBarry Smith       }
2376562c4e1SBarry Smith       /* off-diagonal portion of A */
2386562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2396562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2406562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
2416562c4e1SBarry Smith         ca[k++] = *ba++;
2426562c4e1SBarry Smith       }
2436562c4e1SBarry Smith     }
244e32f2f54SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
2456562c4e1SBarry Smith 
2466562c4e1SBarry Smith     /* put together the new matrix */
2476562c4e1SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
2486562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
2496562c4e1SBarry Smith 
2506562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2516562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2526562c4e1SBarry Smith     mat = (Mat_SeqAIJ*)(*Aloc)->data;
2536562c4e1SBarry Smith     mat->free_a       = PETSC_TRUE;
2546562c4e1SBarry Smith     mat->free_ij      = PETSC_TRUE;
2556562c4e1SBarry Smith 
2566562c4e1SBarry Smith     mat->nonew    = 0;
2576562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX){
2586562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2596562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2606562c4e1SBarry Smith     for (i=0; i<am; i++) {
2616562c4e1SBarry Smith       /* diagonal portion of A */
2626562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2636562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2646562c4e1SBarry Smith       /* off-diagonal portion of A */
2656562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2666562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2676562c4e1SBarry Smith     }
2686562c4e1SBarry Smith   } else {
2690005702aSJed Brown     SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2706562c4e1SBarry Smith   }
2716562c4e1SBarry Smith   PetscFunctionReturn(0);
2726562c4e1SBarry Smith }
2736562c4e1SBarry Smith 
2746562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat);
2756562c4e1SBarry Smith #undef __FUNCT__
2766562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML"
2776562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2786562c4e1SBarry Smith {
2796562c4e1SBarry Smith   PetscErrorCode ierr;
2806562c4e1SBarry Smith   Mat_MLShell    *shell;
2816562c4e1SBarry Smith 
2826562c4e1SBarry Smith   PetscFunctionBegin;
2836562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
284601cad40SBrad Aagaard   ierr = VecDestroy(&shell->y);CHKERRQ(ierr);
285601cad40SBrad Aagaard   if (shell->work) {ierr = VecDestroy(&shell->work);CHKERRQ(ierr);}
2866562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2876562c4e1SBarry Smith   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
2886562c4e1SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
2896562c4e1SBarry Smith   PetscFunctionReturn(0);
2906562c4e1SBarry Smith }
2916562c4e1SBarry Smith 
2926562c4e1SBarry Smith #undef __FUNCT__
2936562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ"
2946562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2956562c4e1SBarry Smith {
2966562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
2976562c4e1SBarry Smith   PetscErrorCode        ierr;
298ae7fe62dSJed Brown   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = PETSC_NULL,nz_max;
2996562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
3006562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
3016562c4e1SBarry Smith 
3026562c4e1SBarry Smith   PetscFunctionBegin;
303e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3046562c4e1SBarry Smith   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
3056562c4e1SBarry Smith     if (reuse){
3066562c4e1SBarry Smith       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
3076562c4e1SBarry Smith       aij->i = ml_rowptr;
3086562c4e1SBarry Smith       aij->j = ml_cols;
3096562c4e1SBarry Smith       aij->a = ml_vals;
3106562c4e1SBarry Smith     } else {
3116562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
3126562c4e1SBarry Smith       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3136562c4e1SBarry Smith       for (i=0; i<m; i++) {
3146562c4e1SBarry Smith         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
3156562c4e1SBarry Smith       }
3166562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
3176562c4e1SBarry Smith       for (i=0; i<m; i++){
3186562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
3196562c4e1SBarry Smith         aj += nnz[i]; aa += nnz[i];
3206562c4e1SBarry Smith       }
3216562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
3226562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
3236562c4e1SBarry Smith     }
3246562c4e1SBarry Smith     PetscFunctionReturn(0);
3256562c4e1SBarry Smith   }
3266562c4e1SBarry Smith 
3276562c4e1SBarry Smith   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
328ae7fe62dSJed Brown   if (reuse) {
329ae7fe62dSJed Brown     for (nz_max=0,i=0; i<m; i++) nz_max = PetscMax(nz_max,ml_cols[i+1] - ml_cols[i] + 1);
330ae7fe62dSJed Brown   } else {
3316562c4e1SBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
3326562c4e1SBarry Smith     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3336562c4e1SBarry Smith     ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
3346562c4e1SBarry Smith 
3356562c4e1SBarry Smith     ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3366562c4e1SBarry Smith     nz_max = 1;
3376562c4e1SBarry Smith     for (i=0; i<m; i++) {
3386562c4e1SBarry Smith       nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
339ae7fe62dSJed Brown       if (nnz[i] > nz_max) nz_max = nnz[i];
3406562c4e1SBarry Smith     }
3416562c4e1SBarry Smith     ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
342ae7fe62dSJed Brown   }
3436562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
3446562c4e1SBarry Smith   for (i=0; i<m; i++) {
345ae7fe62dSJed Brown     PetscInt ncols;
3466562c4e1SBarry Smith     k = 0;
3476562c4e1SBarry Smith     /* diagonal entry */
3486562c4e1SBarry Smith     aj[k] = i; aa[k++] = ml_vals[i];
3496562c4e1SBarry Smith     /* off diagonal entries */
3506562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
3516562c4e1SBarry Smith       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
3526562c4e1SBarry Smith     }
353ae7fe62dSJed Brown     ncols = ml_cols[i+1] - ml_cols[i] + 1;
3546562c4e1SBarry Smith     /* sort aj and aa */
355ae7fe62dSJed Brown     ierr = PetscSortIntWithScalarArray(ncols,aj,aa);CHKERRQ(ierr);
356ae7fe62dSJed Brown     ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3576562c4e1SBarry Smith   }
3586562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3596562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3606562c4e1SBarry Smith 
3616562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3626562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
3636562c4e1SBarry Smith   PetscFunctionReturn(0);
3646562c4e1SBarry Smith }
3656562c4e1SBarry Smith 
3666562c4e1SBarry Smith #undef __FUNCT__
3676562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL"
3686562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3696562c4e1SBarry Smith {
3706562c4e1SBarry Smith   PetscErrorCode ierr;
3716562c4e1SBarry Smith   PetscInt       m,n;
3726562c4e1SBarry Smith   ML_Comm        *MLcomm;
3736562c4e1SBarry Smith   Mat_MLShell    *shellctx;
3746562c4e1SBarry Smith 
3756562c4e1SBarry Smith   PetscFunctionBegin;
3766562c4e1SBarry Smith   m = mlmat->outvec_leng;
3776562c4e1SBarry Smith   n = mlmat->invec_leng;
3786562c4e1SBarry Smith   if (!m || !n){
3796562c4e1SBarry Smith     newmat = PETSC_NULL;
3806562c4e1SBarry Smith     PetscFunctionReturn(0);
3816562c4e1SBarry Smith   }
3826562c4e1SBarry Smith 
3836562c4e1SBarry Smith   if (reuse){
3846562c4e1SBarry Smith     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
3856562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3866562c4e1SBarry Smith     PetscFunctionReturn(0);
3876562c4e1SBarry Smith   }
3886562c4e1SBarry Smith 
3896562c4e1SBarry Smith   MLcomm = mlmat->comm;
3906562c4e1SBarry Smith   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
3916562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3926562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
3936562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
3946562c4e1SBarry Smith   shellctx->A         = *newmat;
3956562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
39667d6f150SMatthew G Knepley   shellctx->work      = PETSC_NULL;
3979bb5392cSJed Brown   ierr = VecCreate(MLcomm->USR_comm,&shellctx->y);CHKERRQ(ierr);
3986562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
3996562c4e1SBarry Smith   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
4006562c4e1SBarry Smith   (*newmat)->ops->destroy = MatDestroy_ML;
4016562c4e1SBarry Smith   PetscFunctionReturn(0);
4026562c4e1SBarry Smith }
4036562c4e1SBarry Smith 
4046562c4e1SBarry Smith #undef __FUNCT__
4056562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ"
406ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
4076562c4e1SBarry Smith {
4086562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
4096562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*aj;
4106562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
4116562c4e1SBarry Smith   PetscErrorCode        ierr;
4126562c4e1SBarry Smith   PetscInt              i,j,k,*gordering;
413ae7fe62dSJed Brown   PetscInt              m=mlmat->outvec_leng,n,nz_max,row;
4146562c4e1SBarry Smith   Mat                   A;
4156562c4e1SBarry Smith 
4166562c4e1SBarry Smith   PetscFunctionBegin;
417e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
4186562c4e1SBarry Smith   n = mlmat->invec_leng;
419e32f2f54SBarry Smith   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
4206562c4e1SBarry Smith 
421ae7fe62dSJed Brown   if (reuse) {
422ae7fe62dSJed Brown     A = *newmat;
423ae7fe62dSJed Brown     for (nz_max=0,i=0; i<m; i++) nz_max = PetscMax(nz_max,ml_cols[i+1] - ml_cols[i] + 1);
424ae7fe62dSJed Brown   } else {
425ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
4266562c4e1SBarry Smith     ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
4276562c4e1SBarry Smith     ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
4286562c4e1SBarry Smith     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
4296562c4e1SBarry Smith     ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
4306562c4e1SBarry Smith 
4316562c4e1SBarry Smith     nz_max = 0;
4326562c4e1SBarry Smith     for (i=0; i<m; i++){
4336562c4e1SBarry Smith       nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
4346562c4e1SBarry Smith       if (nz_max < nnz[i]) nz_max = nnz[i];
4356562c4e1SBarry Smith       nnzA[i] = 1; /* diag */
4366562c4e1SBarry Smith       for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
4376562c4e1SBarry Smith         if (ml_cols[j] < m) nnzA[i]++;
4386562c4e1SBarry Smith       }
4396562c4e1SBarry Smith       nnzB[i] = nnz[i] - nnzA[i];
4406562c4e1SBarry Smith     }
4416562c4e1SBarry Smith     ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
442ae7fe62dSJed Brown     ierr = PetscFree3(nnzA,nnzB,nnz);
443ae7fe62dSJed Brown   }
4446562c4e1SBarry Smith 
4456562c4e1SBarry Smith   /* insert mat values -- remap row and column indices */
4466562c4e1SBarry Smith   nz_max++;
4476562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
4486562c4e1SBarry Smith   /* create global row numbering for a ML_Operator */
4496562c4e1SBarry Smith   ML_build_global_numbering(mlmat,&gordering,"rows");
4506562c4e1SBarry Smith   for (i=0; i<m; i++) {
451ae7fe62dSJed Brown     PetscInt ncols;
4526562c4e1SBarry Smith     row = gordering[i];
4536562c4e1SBarry Smith     k = 0;
4546562c4e1SBarry Smith     /* diagonal entry */
4556562c4e1SBarry Smith     aj[k] = row; aa[k++] = ml_vals[i];
4566562c4e1SBarry Smith     /* off diagonal entries */
4576562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
4586562c4e1SBarry Smith       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
4596562c4e1SBarry Smith     }
460ae7fe62dSJed Brown     ncols = ml_cols[i+1] - ml_cols[i] + 1;
461ae7fe62dSJed Brown     ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
4626562c4e1SBarry Smith   }
463eb4736aaSBarry Smith   ML_free(gordering);
4646562c4e1SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4656562c4e1SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4666562c4e1SBarry Smith   *newmat = A;
4676562c4e1SBarry Smith 
4686562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
4696562c4e1SBarry Smith   PetscFunctionReturn(0);
4706562c4e1SBarry Smith }
4716562c4e1SBarry Smith 
4726562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
47301da6913SBarry Smith #undef __FUNCT__
474a06653b4SBarry Smith #define __FUNCT__ "PCReset_ML"
47516336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
47601da6913SBarry Smith {
47701da6913SBarry Smith   PetscErrorCode  ierr;
478e0262f48SMatthew G Knepley   PC_MG           *mg = (PC_MG*)pc->data;
479e0262f48SMatthew G Knepley   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
48001da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
48101da6913SBarry Smith 
48201da6913SBarry Smith   PetscFunctionBegin;
48301da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
48401da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
48501da6913SBarry Smith 
48601da6913SBarry Smith   if (pc_ml->PetscMLdata) {
48701da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
488ae7fe62dSJed Brown     ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);
489ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr);
490ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr);
49101da6913SBarry Smith   }
49201da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
49301da6913SBarry Smith 
494f5a5dd59SJed Brown   if (pc_ml->gridctx) {
49501da6913SBarry Smith     for (level=0; level<fine_level; level++){
496601cad40SBrad Aagaard       if (pc_ml->gridctx[level].A){ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);}
497601cad40SBrad Aagaard       if (pc_ml->gridctx[level].P){ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);}
498601cad40SBrad Aagaard       if (pc_ml->gridctx[level].R){ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);}
499601cad40SBrad Aagaard       if (pc_ml->gridctx[level].x){ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);}
500601cad40SBrad Aagaard       if (pc_ml->gridctx[level].b){ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);}
501601cad40SBrad Aagaard       if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
50201da6913SBarry Smith     }
503f5a5dd59SJed Brown   }
50401da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
50501da6913SBarry Smith   PetscFunctionReturn(0);
50601da6913SBarry Smith }
5075582bec1SHong Zhang /* -------------------------------------------------------------------------- */
5085582bec1SHong Zhang /*
5095582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
5105582bec1SHong Zhang                     by setting data structures and options.
5115582bec1SHong Zhang 
5125582bec1SHong Zhang    Input Parameter:
5135582bec1SHong Zhang .  pc - the preconditioner context
5145582bec1SHong Zhang 
5155582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5165582bec1SHong Zhang 
5175582bec1SHong Zhang    Notes:
5185582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5195582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5205582bec1SHong Zhang */
5216ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
522a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
523c07bf074SBarry Smith 
5245582bec1SHong Zhang #undef __FUNCT__
5255582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
5266ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5275582bec1SHong Zhang {
5285582bec1SHong Zhang   PetscErrorCode  ierr;
529eef31507SHong Zhang   PetscMPIInt     size;
5305582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
5315582bec1SHong Zhang   ML              *ml_object;
5325582bec1SHong Zhang   ML_Aggregate    *agg_object;
5335582bec1SHong Zhang   ML_Operator     *mlmat;
5344f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5355582bec1SHong Zhang   Mat             A,Aloc;
5365582bec1SHong Zhang   GridCtx         *gridctx;
53701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
53801da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
539ace3abfcSBarry Smith   PetscBool       isSeq, isMPI;
540c07bf074SBarry Smith   KSP             smoother;
541c07bf074SBarry Smith   PC              subpc;
54248268eb4SJed Brown   PetscInt        mesh_level, old_mesh_level;
54348268eb4SJed Brown 
5445582bec1SHong Zhang 
5455582bec1SHong Zhang   PetscFunctionBegin;
54648268eb4SJed Brown   A = pc->pmat;
54748268eb4SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
54848268eb4SJed Brown 
549573998d7SHong Zhang   if (pc->setupcalled) {
55048268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
55148268eb4SJed Brown       /*
55248268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
55348268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
55448268eb4SJed Brown        */
55548268eb4SJed Brown       ml_object = pc_ml->ml_object;
55648268eb4SJed Brown       gridctx = pc_ml->gridctx;
55748268eb4SJed Brown       Nlevels = pc_ml->Nlevels;
55848268eb4SJed Brown       fine_level = Nlevels - 1;
55948268eb4SJed Brown       gridctx[fine_level].A = A;
56048268eb4SJed Brown 
56148268eb4SJed Brown       ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
56248268eb4SJed Brown       ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
56348268eb4SJed Brown       if (isMPI){
56448268eb4SJed Brown         ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
56548268eb4SJed Brown       } else if (isSeq) {
56648268eb4SJed Brown         Aloc = A;
567ae7fe62dSJed Brown         ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
56848268eb4SJed Brown       } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name);
56948268eb4SJed Brown 
57048268eb4SJed Brown       ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
57148268eb4SJed Brown       PetscMLdata = pc_ml->PetscMLdata;
572ae7fe62dSJed Brown       ierr = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr);
57348268eb4SJed Brown       PetscMLdata->A    = A;
57448268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
57548268eb4SJed Brown       ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
57648268eb4SJed Brown       ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
57748268eb4SJed Brown 
57848268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
57948268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
58048268eb4SJed Brown         old_mesh_level = mesh_level;
58148268eb4SJed Brown         mesh_level = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
58248268eb4SJed Brown 
58348268eb4SJed Brown         /* clean and regenerate A */
58448268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
58548268eb4SJed Brown         ML_Operator_Clean(mlmat);
58648268eb4SJed Brown         ML_Operator_Init(mlmat,ml_object->comm);
58748268eb4SJed Brown         ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level);
58848268eb4SJed Brown       }
58948268eb4SJed Brown 
59048268eb4SJed Brown       level = fine_level - 1;
59148268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
59248268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++){
59348268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
594ae7fe62dSJed Brown           ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
59548268eb4SJed Brown           level--;
59648268eb4SJed Brown         }
59748268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
59848268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++){
59948268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
600ae7fe62dSJed Brown           ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
60148268eb4SJed Brown           level--;
60248268eb4SJed Brown         }
60348268eb4SJed Brown       }
60448268eb4SJed Brown 
60548268eb4SJed Brown       for (level=0; level<fine_level; level++) {
60648268eb4SJed Brown         if (level > 0){
60748268eb4SJed Brown           ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
60848268eb4SJed Brown         }
6097721a10bSJed Brown         ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
61048268eb4SJed Brown       }
61148268eb4SJed Brown       ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
6127721a10bSJed Brown       ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
61348268eb4SJed Brown 
61448268eb4SJed Brown       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
61548268eb4SJed Brown       PetscFunctionReturn(0);
61648268eb4SJed Brown     } else {
617c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
61816336fedSMatthew G Knepley       ierr = PCReset_ML(pc);CHKERRQ(ierr);
619a06653b4SBarry Smith       ierr = PCReset_MG(pc);CHKERRQ(ierr);
620573998d7SHong Zhang     }
62148268eb4SJed Brown   }
622573998d7SHong Zhang 
6235582bec1SHong Zhang   /* setup special features of PCML */
6245582bec1SHong Zhang   /*--------------------------------*/
6255582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6265582bec1SHong Zhang   pc_ml->size = size;
627864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
628864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
629864b637dSMatthew Knepley   if (isMPI){
630db571536SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
631864b637dSMatthew Knepley   } else if (isSeq) {
6325582bec1SHong Zhang     Aloc = A;
633ae7fe62dSJed Brown     ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
634c5b2ea42SJed Brown   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name);
6355582bec1SHong Zhang 
6365582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
63738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
6385582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
639d0f46423SBarry Smith   ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
6405582bec1SHong Zhang 
64124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
642d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
64324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
64424a42b14SHong Zhang 
64524a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
646d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
64724a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
648573998d7SHong Zhang   PetscMLdata->A    = A;
649573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
65024a42b14SHong Zhang 
6515582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
65245cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
6535582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
6544f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
6555582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
656ead7dcbeSHong Zhang   ML_Comm_Set_UsrComm(ml_object->comm,((PetscObject)A)->comm);
657573998d7SHong Zhang   pc_ml->ml_object = ml_object;
6585582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
6595582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
6605582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
6615582bec1SHong Zhang 
662b5c8bdf8SJed Brown   ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO);
663b5c8bdf8SJed Brown 
6645582bec1SHong Zhang   /* aggregation */
6655582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
666573998d7SHong Zhang   pc_ml->agg_object = agg_object;
667573998d7SHong Zhang 
668*fb6a8e6dSJed Brown   {
669*fb6a8e6dSJed Brown     MatNullSpace mnull;
670*fb6a8e6dSJed Brown     ierr = MatGetNearNullSpace(A,&mnull);CHKERRQ(ierr);
671*fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
672*fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
673*fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
674*fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
675*fb6a8e6dSJed Brown     }
676*fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
677*fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
678*fb6a8e6dSJed Brown       PetscScalar *nullvec;
679*fb6a8e6dSJed Brown       const PetscScalar *v;
680*fb6a8e6dSJed Brown       PetscBool has_const;
681*fb6a8e6dSJed Brown       PetscInt i,j,mlocal,nvec;
682*fb6a8e6dSJed Brown       const Vec *vecs;
683*fb6a8e6dSJed Brown       if (!mnull) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
684*fb6a8e6dSJed Brown       ierr = MatGetLocalSize(Aloc,&mlocal,PETSC_NULL);CHKERRQ(ierr);
685*fb6a8e6dSJed Brown       ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr);
686*fb6a8e6dSJed Brown       ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof *nullvec,&nullvec);CHKERRQ(ierr);
687*fb6a8e6dSJed Brown       if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/m;
688*fb6a8e6dSJed Brown       for (i=0; i<nvec; i++) {
689*fb6a8e6dSJed Brown         ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
690*fb6a8e6dSJed Brown         for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j];
691*fb6a8e6dSJed Brown         ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
692*fb6a8e6dSJed Brown       }
693*fb6a8e6dSJed Brown       ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal);CHKERRQ(ierr);
694*fb6a8e6dSJed Brown       ierr = PetscFree(nullvec);CHKERRQ(ierr);
695*fb6a8e6dSJed Brown     } break;
696*fb6a8e6dSJed Brown     case PCML_NULLSPACE_BLOCK:
697*fb6a8e6dSJed Brown       ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
698*fb6a8e6dSJed Brown       break;
699*fb6a8e6dSJed Brown     case PCML_NULLSPACE_SCALAR:
700*fb6a8e6dSJed Brown       break;
701*fb6a8e6dSJed Brown     default: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unknown null space type");
702*fb6a8e6dSJed Brown     }
703*fb6a8e6dSJed Brown   }
7045582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
7055582bec1SHong Zhang   /* set options */
7065582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
7075582bec1SHong Zhang   case 1:
7085582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
7095582bec1SHong Zhang   case 2:
7105582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
7115582bec1SHong Zhang   case 3:
7125582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
7135582bec1SHong Zhang   }
7145582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
7155582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
7165582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
7177ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
7185582bec1SHong Zhang   }
719b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
720b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
721b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
722b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
723b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
724b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7255582bec1SHong Zhang 
726b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
7275582bec1SHong Zhang     Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
728b5c8bdf8SJed Brown   } else {
729b5c8bdf8SJed Brown     Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
730b5c8bdf8SJed Brown   }
73165e19b50SBarry Smith   if (Nlevels<=0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
732573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
733aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
734c07bf074SBarry Smith 
73597177400SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
736aa85bbbfSHong Zhang   /* set default smoothers */
737aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++){
738aa85bbbfSHong Zhang     if (size == 1){
739aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
740aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
741aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
742aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
743aa85bbbfSHong Zhang     } else {
744aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
745aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
746aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
747aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
748aa85bbbfSHong Zhang     }
749aa85bbbfSHong Zhang   }
75097177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
7515582bec1SHong Zhang 
7525582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
7535582bec1SHong Zhang   pc_ml->gridctx = gridctx;
7545582bec1SHong Zhang 
7555582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
7565582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
757e14861a4SHong Zhang   gridctx[fine_level].A = A;
758573998d7SHong Zhang 
759e14861a4SHong Zhang   level = fine_level - 1;
760ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
7615582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
762e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
763db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
764e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
765db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
766573998d7SHong Zhang 
767573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
768573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
7695582bec1SHong Zhang       level--;
7705582bec1SHong Zhang     }
771ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
7725582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
7735582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
774db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
775ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
776db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
777573998d7SHong Zhang 
7785582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
779ae7fe62dSJed Brown       ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
7805582bec1SHong Zhang       level--;
7815582bec1SHong Zhang     }
7825582bec1SHong Zhang   }
7835582bec1SHong Zhang 
784573998d7SHong Zhang   /* create vectors and ksp at all levels */
785ac346b81SHong Zhang   for (level=0; level<fine_level; level++){
786573998d7SHong Zhang     level1 = level + 1;
787e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
788d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
7895582bec1SHong Zhang     ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
79097177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
7915582bec1SHong Zhang 
792e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
793d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
7945582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
79597177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
796ac346b81SHong Zhang 
797e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
798d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
799ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
80097177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
801ac346b81SHong Zhang 
8025582bec1SHong Zhang     if (level == 0){
80397177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
8045582bec1SHong Zhang     } else {
80597177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
806573998d7SHong Zhang     }
807573998d7SHong Zhang   }
808573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
809573998d7SHong Zhang 
810573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
811573998d7SHong Zhang   for (level=0; level<fine_level; level++){
812573998d7SHong Zhang     level1 = level + 1;
813aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
814573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
815573998d7SHong Zhang     if (level > 0){
81697177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
8175582bec1SHong Zhang     }
8185582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
8195582bec1SHong Zhang   }
82097177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
821ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
8225582bec1SHong Zhang 
823c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
824c07bf074SBarry Smith   pc->setupcalled = 0;
8253751b4bdSBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8265582bec1SHong Zhang   PetscFunctionReturn(0);
8275582bec1SHong Zhang }
8285582bec1SHong Zhang 
8295582bec1SHong Zhang /* -------------------------------------------------------------------------- */
8305582bec1SHong Zhang /*
8315582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
8325582bec1SHong Zhang    that was created with PCCreate_ML().
8335582bec1SHong Zhang 
8345582bec1SHong Zhang    Input Parameter:
8355582bec1SHong Zhang .  pc - the preconditioner context
8365582bec1SHong Zhang 
8375582bec1SHong Zhang    Application Interface Routine: PCDestroy()
8385582bec1SHong Zhang */
8395582bec1SHong Zhang #undef __FUNCT__
8405582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
8416ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
8425582bec1SHong Zhang {
8435582bec1SHong Zhang   PetscErrorCode  ierr;
84401da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
84501da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
8465582bec1SHong Zhang 
8475582bec1SHong Zhang   PetscFunctionBegin;
84816336fedSMatthew G Knepley   ierr = PCReset_ML(pc);CHKERRQ(ierr);
84901da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
85001da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8515582bec1SHong Zhang   PetscFunctionReturn(0);
8525582bec1SHong Zhang }
8535582bec1SHong Zhang 
8545582bec1SHong Zhang #undef __FUNCT__
8555582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
8566ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
8575582bec1SHong Zhang {
8585582bec1SHong Zhang   PetscErrorCode  ierr;
8593751b4bdSBarry Smith   PetscInt        indx,PrintLevel;
8605582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
86101da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
86201da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
863b5c8bdf8SJed Brown   PetscMPIInt     size;
86488ff4cc7SJed Brown   MPI_Comm        comm = ((PetscObject)pc)->comm;
8655582bec1SHong Zhang 
8665582bec1SHong Zhang   PetscFunctionBegin;
86788ff4cc7SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
8685582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
8695582bec1SHong Zhang   PrintLevel    = 0;
8705582bec1SHong Zhang   indx          = 0;
8715582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
8725582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
873573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
874573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
8753751b4bdSBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
8765582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
877573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
878573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
879acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);CHKERRQ(ierr);
880acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,PETSC_NULL);CHKERRQ(ierr);
881acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,PETSC_NULL);CHKERRQ(ierr);
882*fb6a8e6dSJed Brown   ierr = PetscOptionsEnum("-pc_ml_nullspace","Which type of null space information to use","None",PCMLNullSpaceTypes,(PetscEnum)pc_ml->nulltype,(PetscEnum*)&pc_ml->nulltype,PETSC_NULL);CHKERRQ(ierr);
883b5c8bdf8SJed Brown   ierr = PetscOptionsInt("-pc_ml_EnergyMinimization","Energy minimization norm type (0=no minimization; see ML manual for 1,2,3; -1 and 4 undocumented)","None",pc_ml->EnergyMinimization,&pc_ml->EnergyMinimization,PETSC_NULL);CHKERRQ(ierr);
88448268eb4SJed Brown   ierr = PetscOptionsBool("-pc_ml_reuse_interpolation","Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)","None",pc_ml->reuse_interpolation,&pc_ml->reuse_interpolation,PETSC_NULL);CHKERRQ(ierr);
885b5c8bdf8SJed Brown   /*
886b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
887b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
888b5c8bdf8SJed Brown 
889b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
890b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
891b5c8bdf8SJed Brown   */
89288ff4cc7SJed Brown   if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
89388ff4cc7SJed Brown   if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
894b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2");CHKERRQ(ierr);}
895b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
896b5c8bdf8SJed Brown     ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,PETSC_NULL);CHKERRQ(ierr);
897b5c8bdf8SJed Brown   }
898b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
899b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
900acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,PETSC_NULL);CHKERRQ(ierr);
901b5c8bdf8SJed Brown   }
902b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
903b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
904acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_KeepAggInfo","Allows the preconditioner to be reused, or auxilliary matrices to be generated","None",pc_ml->KeepAggInfo,&pc_ml->KeepAggInfo,PETSC_NULL);CHKERRQ(ierr);
905b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
906b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
907acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_Reusable","Store intermedaiate data structures so that the multilevel hierarchy is reusable","None",pc_ml->Reusable,&pc_ml->Reusable,PETSC_NULL);CHKERRQ(ierr);
908b5c8bdf8SJed Brown   /*
909b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
910b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
911b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
912b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
913b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
914b5c8bdf8SJed Brown    */
915acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,PETSC_NULL);CHKERRQ(ierr);
9165582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
9175582bec1SHong Zhang   PetscFunctionReturn(0);
9185582bec1SHong Zhang }
9195582bec1SHong Zhang 
9205582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9215582bec1SHong Zhang /*
9225582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
9235582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
9245582bec1SHong Zhang    context, PC, that was created within PCCreate().
9255582bec1SHong Zhang 
9265582bec1SHong Zhang    Input Parameter:
9275582bec1SHong Zhang .  pc - the preconditioner context
9285582bec1SHong Zhang 
9295582bec1SHong Zhang    Application Interface Routine: PCCreate()
9305582bec1SHong Zhang */
9315582bec1SHong Zhang 
9325582bec1SHong Zhang /*MC
9331e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
9345582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
9356ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
9366ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
9375582bec1SHong Zhang 
9386ca4d86aSHong Zhang    Options Database Key:
9396ca4d86aSHong Zhang    Multigrid options(inherited)
9406ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
9416ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
9426ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
943ccd75124SBarry Smith    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
94451f519a2SBarry Smith    ML options:
945ccd75124SBarry Smith .  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
9466ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
9476ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
948f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
9496ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
9506ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
9517ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
9525582bec1SHong Zhang 
9535582bec1SHong Zhang    Level: intermediate
9545582bec1SHong Zhang 
9555582bec1SHong Zhang   Concepts: multigrid
9565582bec1SHong Zhang 
9575582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
95897177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
95997177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
96097177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
96197177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9625582bec1SHong Zhang M*/
9635582bec1SHong Zhang 
9645582bec1SHong Zhang EXTERN_C_BEGIN
9655582bec1SHong Zhang #undef __FUNCT__
9665582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
9677087cfbeSBarry Smith PetscErrorCode  PCCreate_ML(PC pc)
9685582bec1SHong Zhang {
9695582bec1SHong Zhang   PetscErrorCode  ierr;
9705582bec1SHong Zhang   PC_ML           *pc_ml;
97101da6913SBarry Smith   PC_MG           *mg;
9725582bec1SHong Zhang 
9735582bec1SHong Zhang   PetscFunctionBegin;
974573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
9755582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
97603bfa161SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
977e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
978e0f5d30fSBarry Smith   ierr = DMDestroy(&pc->dm);CHKERRQ(ierr);
979e0f5d30fSBarry Smith   mg = (PC_MG*)pc->data;
980c91913d3SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
9815582bec1SHong Zhang 
9825582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
98338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
98401da6913SBarry Smith   mg->innerctx = pc_ml;
9855582bec1SHong Zhang 
986573998d7SHong Zhang   pc_ml->ml_object     = 0;
987573998d7SHong Zhang   pc_ml->agg_object    = 0;
988573998d7SHong Zhang   pc_ml->gridctx       = 0;
989573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
990573998d7SHong Zhang   pc_ml->Nlevels       = -1;
991573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
992573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
9933751b4bdSBarry Smith   pc_ml->CoarsenScheme = 1;
994573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
995573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
996573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
997573998d7SHong Zhang   pc_ml->size          = 0;
998573998d7SHong Zhang 
9995582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
10005582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
10015582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1002a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
10035582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
10045582bec1SHong Zhang   PetscFunctionReturn(0);
10055582bec1SHong Zhang }
10065582bec1SHong Zhang EXTERN_C_END
1007