xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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*/
8b4876bcbSBarry Smith #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscksp.h" I*/
9c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
10c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
111e25c274SJed Brown #include <petscdm.h>            /* for DMDestroy(&pc->mg) hack */
12cb5d8e9eSHong Zhang 
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>
1939381ba2SJed Brown #include <ml_viz_stats.h>
205582bec1SHong Zhang EXTERN_C_END
215582bec1SHong Zhang 
22fb6a8e6dSJed Brown typedef enum {PCML_NULLSPACE_AUTO,PCML_NULLSPACE_USER,PCML_NULLSPACE_BLOCK,PCML_NULLSPACE_SCALAR} PCMLNullSpaceType;
23fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO","USER","BLOCK","SCALAR","PCMLNullSpaceType","PCML_NULLSPACE_",0};
24fb6a8e6dSJed Brown 
255582bec1SHong Zhang /* The context (data structure) at each grid level */
265582bec1SHong Zhang typedef struct {
275582bec1SHong Zhang   Vec x,b,r;                  /* global vectors */
285582bec1SHong Zhang   Mat A,P,R;
295582bec1SHong Zhang   KSP ksp;
3039381ba2SJed Brown   Vec coords;                 /* projected by ML, if PCSetCoordinates is called; values packed by node */
315582bec1SHong Zhang } GridCtx;
325582bec1SHong Zhang 
335582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
345582bec1SHong Zhang typedef struct {
35573998d7SHong Zhang   Mat         A;       /* Petsc matrix in aij format */
36573998d7SHong Zhang   Mat         Aloc;    /* local portion of A to be used by ML */
3724a42b14SHong Zhang   Vec         x,y;
385582bec1SHong Zhang   ML_Operator *mlmat;
395582bec1SHong Zhang   PetscScalar *pwork;  /* tmp array used by PetscML_comm() */
405582bec1SHong Zhang } FineGridCtx;
415582bec1SHong Zhang 
425582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
435582bec1SHong Zhang typedef struct {
445582bec1SHong Zhang   Mat         A;        /* PETSc shell matrix associated with mlmat */
455582bec1SHong Zhang   ML_Operator *mlmat;   /* ML matrix assorciated with A */
4667d6f150SMatthew G Knepley   Vec         y, work;
475582bec1SHong Zhang } Mat_MLShell;
485582bec1SHong Zhang 
495582bec1SHong Zhang /* Private context for the ML preconditioner */
505582bec1SHong Zhang typedef struct {
515582bec1SHong Zhang   ML                *ml_object;
525582bec1SHong Zhang   ML_Aggregate      *agg_object;
535582bec1SHong Zhang   GridCtx           *gridctx;
545582bec1SHong Zhang   FineGridCtx       *PetscMLdata;
5539381ba2SJed Brown   PetscInt          Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization,MinPerProc,PutOnSingleProc,RepartitionType,ZoltanScheme;
5639381ba2SJed Brown   PetscReal         Threshold,DampingFactor,EnergyMinimizationDropTol,MaxMinRatio,AuxThreshold;
5739381ba2SJed Brown   PetscBool         SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable,Repartition,Aux;
5848268eb4SJed Brown   PetscBool         reuse_interpolation;
59fb6a8e6dSJed Brown   PCMLNullSpaceType nulltype;
60573998d7SHong Zhang   PetscMPIInt       size; /* size of communicator for pc->pmat */
6139381ba2SJed Brown   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
6239381ba2SJed Brown   PetscInt          nloc;
6339381ba2SJed Brown   PetscReal         *coords; /* ML has a grid object for each level: the finest grid will point into coords */
645582bec1SHong Zhang } PC_ML;
6541ca0015SHong Zhang 
666562c4e1SBarry Smith #undef __FUNCT__
676562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow"
686562c4e1SBarry 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[])
696562c4e1SBarry Smith {
706562c4e1SBarry Smith   PetscErrorCode ierr;
716562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
726562c4e1SBarry Smith   PetscScalar    *aa;
736562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
746562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
755582bec1SHong Zhang 
760298fd71SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,NULL); if (ierr) return(0);
776562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
786562c4e1SBarry Smith     row            = requested_rows[i];
796562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
806562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
816562c4e1SBarry Smith     if ((row >= 0) || (row <= (m-1))) {
826562c4e1SBarry Smith       aj = a->j + a->i[row];
836562c4e1SBarry Smith       aa = a->a + a->i[row];
846562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++) {
856562c4e1SBarry Smith         columns[k]  = aj[j];
866562c4e1SBarry Smith         values[k++] = aa[j];
876562c4e1SBarry Smith       }
886562c4e1SBarry Smith     }
896562c4e1SBarry Smith   }
906562c4e1SBarry Smith   return(1);
916562c4e1SBarry Smith }
926562c4e1SBarry Smith 
936562c4e1SBarry Smith #undef __FUNCT__
946562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm"
956562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
966562c4e1SBarry Smith {
976562c4e1SBarry Smith   PetscErrorCode ierr;
986562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_data;
996562c4e1SBarry Smith   Mat            A   = ml->A;
1006562c4e1SBarry Smith   Mat_MPIAIJ     *a  = (Mat_MPIAIJ*)A->data;
1016562c4e1SBarry Smith   PetscMPIInt    size;
1026562c4e1SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
1036562c4e1SBarry Smith   PetscScalar    *array;
1046562c4e1SBarry Smith 
1056562c4e1SBarry Smith   PetscFunctionBegin;
106ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1076562c4e1SBarry Smith   if (size == 1) return 0;
1086562c4e1SBarry Smith 
1096562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
1106562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1116562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1126562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1136562c4e1SBarry Smith   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
1142fa5cd67SKarl Rupp   for (i=in_length; i<out_length; i++) p[i] = array[i-in_length];
1156562c4e1SBarry Smith   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
1166562c4e1SBarry Smith   PetscFunctionReturn(0);
1176562c4e1SBarry Smith }
1186562c4e1SBarry Smith 
1196562c4e1SBarry Smith #undef __FUNCT__
1206562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec"
1216562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
1226562c4e1SBarry Smith {
1236562c4e1SBarry Smith   PetscErrorCode ierr;
1246562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_Get_MyMatvecData(ML_data);
1256562c4e1SBarry Smith   Mat            A   = ml->A, Aloc=ml->Aloc;
1266562c4e1SBarry Smith   PetscMPIInt    size;
1276562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
1286562c4e1SBarry Smith   PetscInt       i;
1296562c4e1SBarry Smith 
1306562c4e1SBarry Smith   PetscFunctionBegin;
131ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1326562c4e1SBarry Smith   if (size == 1) {
1336562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
1346562c4e1SBarry Smith   } else {
1356562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
136b0250c70SBarry Smith     ierr = PetscML_comm(pwork,ml);CHKERRQ(ierr);
1376562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
1386562c4e1SBarry Smith   }
1396562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
1406562c4e1SBarry Smith   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
1416562c4e1SBarry Smith   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
1426562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1436562c4e1SBarry Smith   PetscFunctionReturn(0);
1446562c4e1SBarry Smith }
1456562c4e1SBarry Smith 
1466562c4e1SBarry Smith #undef __FUNCT__
1476562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML"
1486562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1496562c4e1SBarry Smith {
1506562c4e1SBarry Smith   PetscErrorCode ierr;
1516562c4e1SBarry Smith   Mat_MLShell    *shell;
1526562c4e1SBarry Smith   PetscScalar    *xarray,*yarray;
1536562c4e1SBarry Smith   PetscInt       x_length,y_length;
1546562c4e1SBarry Smith 
1556562c4e1SBarry Smith   PetscFunctionBegin;
1566562c4e1SBarry Smith   ierr     = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
1576562c4e1SBarry Smith   ierr     = VecGetArray(x,&xarray);CHKERRQ(ierr);
1586562c4e1SBarry Smith   ierr     = VecGetArray(y,&yarray);CHKERRQ(ierr);
1596562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1606562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
161815d23e5SBarry Smith   PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray));
1626562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1636562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1646562c4e1SBarry Smith   PetscFunctionReturn(0);
1656562c4e1SBarry Smith }
1666562c4e1SBarry Smith 
1676562c4e1SBarry Smith #undef __FUNCT__
1686562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML"
16967d6f150SMatthew G Knepley /* Computes y = w + A * x
17067d6f150SMatthew G Knepley    It is possible that w == y, but not x == y
17167d6f150SMatthew G Knepley */
1726562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
1736562c4e1SBarry Smith {
1746562c4e1SBarry Smith   Mat_MLShell    *shell;
1756562c4e1SBarry Smith   PetscScalar    *xarray,*yarray;
1766562c4e1SBarry Smith   PetscInt       x_length,y_length;
17767d6f150SMatthew G Knepley   PetscErrorCode ierr;
1786562c4e1SBarry Smith 
1796562c4e1SBarry Smith   PetscFunctionBegin;
1806562c4e1SBarry Smith   ierr = MatShellGetContext(A, (void**) &shell);CHKERRQ(ierr);
18167d6f150SMatthew G Knepley   if (y == w) {
18267d6f150SMatthew G Knepley     if (!shell->work) {
18367d6f150SMatthew G Knepley       ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr);
18467d6f150SMatthew G Knepley     }
18567d6f150SMatthew G Knepley     ierr     = VecGetArray(x,           &xarray);CHKERRQ(ierr);
18667d6f150SMatthew G Knepley     ierr     = VecGetArray(shell->work, &yarray);CHKERRQ(ierr);
18767d6f150SMatthew G Knepley     x_length = shell->mlmat->invec_leng;
18867d6f150SMatthew G Knepley     y_length = shell->mlmat->outvec_leng;
189815d23e5SBarry Smith     PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray));
19067d6f150SMatthew G Knepley     ierr = VecRestoreArray(x,           &xarray);CHKERRQ(ierr);
19167d6f150SMatthew G Knepley     ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr);
1923ba3408dSMatthew G Knepley     ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr);
19367d6f150SMatthew G Knepley   } else {
1946562c4e1SBarry Smith     ierr     = VecGetArray(x, &xarray);CHKERRQ(ierr);
1956562c4e1SBarry Smith     ierr     = VecGetArray(y, &yarray);CHKERRQ(ierr);
1966562c4e1SBarry Smith     x_length = shell->mlmat->invec_leng;
1976562c4e1SBarry Smith     y_length = shell->mlmat->outvec_leng;
198815d23e5SBarry Smith     PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray));
1996562c4e1SBarry Smith     ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr);
2006562c4e1SBarry Smith     ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr);
2016562c4e1SBarry Smith     ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr);
20267d6f150SMatthew G Knepley   }
2036562c4e1SBarry Smith   PetscFunctionReturn(0);
2046562c4e1SBarry Smith }
2056562c4e1SBarry Smith 
20679d04de1SBarry Smith /* newtype is ignored since only handles one case */
2076562c4e1SBarry Smith #undef __FUNCT__
2086562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML"
2096562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
2106562c4e1SBarry Smith {
2116562c4e1SBarry Smith   PetscErrorCode ierr;
2126562c4e1SBarry Smith   Mat_MPIAIJ     *mpimat=(Mat_MPIAIJ*)A->data;
2136562c4e1SBarry Smith   Mat_SeqAIJ     *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
2146562c4e1SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2156562c4e1SBarry Smith   PetscScalar    *aa=a->a,*ba=b->a,*ca;
2166562c4e1SBarry Smith   PetscInt       am =A->rmap->n,an=A->cmap->n,i,j,k;
2176562c4e1SBarry Smith   PetscInt       *ci,*cj,ncols;
2186562c4e1SBarry Smith 
2196562c4e1SBarry Smith   PetscFunctionBegin;
220e32f2f54SBarry 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);
2216562c4e1SBarry Smith 
2226562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
2236562c4e1SBarry Smith     ierr  = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2246562c4e1SBarry Smith     ci[0] = 0;
2252fa5cd67SKarl Rupp     for (i=0; i<am; i++) ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
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     }
268ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2696562c4e1SBarry Smith   PetscFunctionReturn(0);
2706562c4e1SBarry Smith }
2716562c4e1SBarry Smith 
2726562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat);
2736562c4e1SBarry Smith #undef __FUNCT__
2746562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML"
2756562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2766562c4e1SBarry Smith {
2776562c4e1SBarry Smith   PetscErrorCode ierr;
2786562c4e1SBarry Smith   Mat_MLShell    *shell;
2796562c4e1SBarry Smith 
2806562c4e1SBarry Smith   PetscFunctionBegin;
2816562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
282601cad40SBrad Aagaard   ierr = VecDestroy(&shell->y);CHKERRQ(ierr);
283601cad40SBrad Aagaard   if (shell->work) {ierr = VecDestroy(&shell->work);CHKERRQ(ierr);}
2846562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2856562c4e1SBarry Smith   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
2866562c4e1SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
2876562c4e1SBarry Smith   PetscFunctionReturn(0);
2886562c4e1SBarry Smith }
2896562c4e1SBarry Smith 
2906562c4e1SBarry Smith #undef __FUNCT__
2916562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ"
2926562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2936562c4e1SBarry Smith {
2946562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata*)mlmat->data;
2956562c4e1SBarry Smith   PetscErrorCode        ierr;
2960298fd71SBarry Smith   PetscInt              m       =mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = NULL,nz_max;
29739381ba2SJed Brown   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i;
2986562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2996562c4e1SBarry Smith 
3006562c4e1SBarry Smith   PetscFunctionBegin;
301e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3026562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
3036562c4e1SBarry Smith     if (reuse) {
3046562c4e1SBarry Smith       Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data;
3056562c4e1SBarry Smith       aij->i = ml_rowptr;
3066562c4e1SBarry Smith       aij->j = ml_cols;
3076562c4e1SBarry Smith       aij->a = ml_vals;
3086562c4e1SBarry Smith     } else {
3096562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
3106562c4e1SBarry Smith       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3112fa5cd67SKarl Rupp       for (i=0; i<m; i++) nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
3126562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
3136562c4e1SBarry Smith       for (i=0; i<m; i++) {
3146562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
3156562c4e1SBarry Smith         aj  += nnz[i]; aa += nnz[i];
3166562c4e1SBarry Smith       }
3176562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
3186562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
3196562c4e1SBarry Smith     }
3206562c4e1SBarry Smith     PetscFunctionReturn(0);
3216562c4e1SBarry Smith   }
3226562c4e1SBarry Smith 
32339381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
32439381ba2SJed Brown   ierr   = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
32539381ba2SJed Brown   if (!reuse) {
3266562c4e1SBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
3276562c4e1SBarry Smith     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3286562c4e1SBarry Smith     ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
32939381ba2SJed Brown     /* keep track of block size for A matrices */
33039381ba2SJed Brown     ierr = MatSetBlockSize (*newmat, mlmat->num_PDEs);CHKERRQ(ierr);
3316562c4e1SBarry Smith 
33239381ba2SJed Brown     ierr = PetscMalloc(m*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
3336562c4e1SBarry Smith     for (i=0; i<m; i++) {
334815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
3356562c4e1SBarry Smith     }
3366562c4e1SBarry Smith     ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
337ae7fe62dSJed Brown   }
3386562c4e1SBarry Smith   for (i=0; i<m; i++) {
339ae7fe62dSJed Brown     PetscInt ncols;
34039381ba2SJed Brown 
341815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
342ae7fe62dSJed Brown     ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3436562c4e1SBarry Smith   }
3446562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3456562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3466562c4e1SBarry Smith 
3476562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3486562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
3496562c4e1SBarry Smith   PetscFunctionReturn(0);
3506562c4e1SBarry Smith }
3516562c4e1SBarry Smith 
3526562c4e1SBarry Smith #undef __FUNCT__
3536562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL"
3546562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3556562c4e1SBarry Smith {
3566562c4e1SBarry Smith   PetscErrorCode ierr;
3576562c4e1SBarry Smith   PetscInt       m,n;
3586562c4e1SBarry Smith   ML_Comm        *MLcomm;
3596562c4e1SBarry Smith   Mat_MLShell    *shellctx;
3606562c4e1SBarry Smith 
3616562c4e1SBarry Smith   PetscFunctionBegin;
3626562c4e1SBarry Smith   m = mlmat->outvec_leng;
3636562c4e1SBarry Smith   n = mlmat->invec_leng;
3646562c4e1SBarry Smith 
3656562c4e1SBarry Smith   if (reuse) {
3666562c4e1SBarry Smith     ierr            = MatShellGetContext(*newmat,(void**)&shellctx);CHKERRQ(ierr);
3676562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3686562c4e1SBarry Smith     PetscFunctionReturn(0);
3696562c4e1SBarry Smith   }
3706562c4e1SBarry Smith 
3716562c4e1SBarry Smith   MLcomm = mlmat->comm;
3722fa5cd67SKarl Rupp 
3736562c4e1SBarry Smith   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
3746562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3756562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
3766562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
3772fa5cd67SKarl Rupp 
3786562c4e1SBarry Smith   shellctx->A         = *newmat;
3796562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
3800298fd71SBarry Smith   shellctx->work      = NULL;
3812fa5cd67SKarl Rupp 
3829bb5392cSJed Brown   ierr = VecCreate(MLcomm->USR_comm,&shellctx->y);CHKERRQ(ierr);
3836562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
3846562c4e1SBarry Smith   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
3852fa5cd67SKarl Rupp 
3866562c4e1SBarry Smith   (*newmat)->ops->destroy = MatDestroy_ML;
3876562c4e1SBarry Smith   PetscFunctionReturn(0);
3886562c4e1SBarry Smith }
3896562c4e1SBarry Smith 
3906562c4e1SBarry Smith #undef __FUNCT__
3916562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ"
392ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3936562c4e1SBarry Smith {
39439381ba2SJed Brown   PetscInt       *aj;
39539381ba2SJed Brown   PetscScalar    *aa;
3966562c4e1SBarry Smith   PetscErrorCode ierr;
39739381ba2SJed Brown   PetscInt       i,j,*gordering;
398ae7fe62dSJed Brown   PetscInt       m=mlmat->outvec_leng,n,nz_max,row;
3996562c4e1SBarry Smith   Mat            A;
4006562c4e1SBarry Smith 
4016562c4e1SBarry Smith   PetscFunctionBegin;
402e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
4036562c4e1SBarry Smith   n = mlmat->invec_leng;
404e32f2f54SBarry Smith   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
4056562c4e1SBarry Smith 
40639381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
40739381ba2SJed Brown   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
4082fa5cd67SKarl Rupp   if (reuse) A = *newmat;
4092fa5cd67SKarl Rupp   else {
410ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
4116562c4e1SBarry Smith     ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
4126562c4e1SBarry Smith     ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
4136562c4e1SBarry Smith     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
41439381ba2SJed Brown     /* keep track of block size for A matrices */
41539381ba2SJed Brown     ierr = MatSetBlockSize (A,mlmat->num_PDEs);CHKERRQ(ierr);
4166562c4e1SBarry Smith     ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
4176562c4e1SBarry Smith 
4186562c4e1SBarry Smith     for (i=0; i<m; i++) {
419815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
42039381ba2SJed Brown       nnzA[i] = 0;
42139381ba2SJed Brown       for (j=0; j<nnz[i]; j++) {
42239381ba2SJed Brown         if (aj[j] < m) nnzA[i]++;
4236562c4e1SBarry Smith       }
4246562c4e1SBarry Smith       nnzB[i] = nnz[i] - nnzA[i];
4256562c4e1SBarry Smith     }
4266562c4e1SBarry Smith     ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
427ae7fe62dSJed Brown     ierr = PetscFree3(nnzA,nnzB,nnz);
428ae7fe62dSJed Brown   }
4296562c4e1SBarry Smith   /* create global row numbering for a ML_Operator */
430815d23e5SBarry Smith   PetscStackCall("ML_build_global_numbering",ML_build_global_numbering(mlmat,&gordering,"rows"));
4316562c4e1SBarry Smith   for (i=0; i<m; i++) {
432ae7fe62dSJed Brown     PetscInt ncols;
4336562c4e1SBarry Smith     row = gordering[i];
43439381ba2SJed Brown 
435815d23e5SBarry Smith     PetscStackCall(",ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
4362fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
437ae7fe62dSJed Brown     ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
4386562c4e1SBarry Smith   }
439815d23e5SBarry Smith   PetscStackCall("ML_Operator_Getrow",ML_free(gordering));
4406562c4e1SBarry Smith   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4416562c4e1SBarry Smith   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4426562c4e1SBarry Smith   *newmat = A;
4436562c4e1SBarry Smith 
4446562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
4456562c4e1SBarry Smith   PetscFunctionReturn(0);
4466562c4e1SBarry Smith }
4476562c4e1SBarry Smith 
44839381ba2SJed Brown /* -------------------------------------------------------------------------- */
44939381ba2SJed Brown /*
45039381ba2SJed Brown    PCSetCoordinates_ML
45139381ba2SJed Brown 
45239381ba2SJed Brown    Input Parameter:
45339381ba2SJed Brown    .  pc - the preconditioner context
45439381ba2SJed Brown */
45539381ba2SJed Brown #undef __FUNCT__
45639381ba2SJed Brown #define __FUNCT__ "PCSetCoordinates_ML"
457f7a08781SBarry Smith static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
45839381ba2SJed Brown {
45939381ba2SJed Brown   PC_MG          *mg    = (PC_MG*)pc->data;
46039381ba2SJed Brown   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
46139381ba2SJed Brown   PetscErrorCode ierr;
46239381ba2SJed Brown   PetscInt       arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend;
46339381ba2SJed Brown   Mat            Amat = pc->pmat;
46439381ba2SJed Brown 
46539381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
46639381ba2SJed Brown   PetscFunctionBegin;
46739381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
46839381ba2SJed Brown   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
46939381ba2SJed Brown 
47039381ba2SJed Brown   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
47139381ba2SJed Brown   nloc = (Iend-my0)/bs;
47239381ba2SJed Brown 
473ce94432eSBarry Smith   if (nloc!=a_nloc) SETERRQ2(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Number of local blocks must locations = %d %d.",a_nloc,nloc);
474ce94432eSBarry Smith   if ((Iend-my0)%bs!=0) SETERRQ1(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
47539381ba2SJed Brown 
47639381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
47739381ba2SJed Brown   pc_ml->dim  = ndm;
47839381ba2SJed Brown   pc_ml->nloc = a_nloc;
47939381ba2SJed Brown   arrsz       = ndm * a_nloc;
48039381ba2SJed Brown 
48139381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
48239381ba2SJed Brown   if (pc_ml->coords==0 || (oldarrsz != arrsz)) {
48339381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
48439381ba2SJed Brown     ierr = PetscMalloc((arrsz)*sizeof(PetscReal), &pc_ml->coords);CHKERRQ(ierr);
48539381ba2SJed Brown   }
48639381ba2SJed Brown   for (kk=0; kk<arrsz; kk++) pc_ml->coords[kk] = -999.;
48739381ba2SJed Brown   /* copy data in - column oriented */
48839381ba2SJed Brown   for (kk = 0; kk < nloc; kk++) {
48939381ba2SJed Brown     for (ii = 0; ii < ndm; ii++) {
49039381ba2SJed Brown       pc_ml->coords[ii*nloc + kk] =  coords[kk*ndm + ii];
49139381ba2SJed Brown     }
49239381ba2SJed Brown   }
49339381ba2SJed Brown   PetscFunctionReturn(0);
49439381ba2SJed Brown }
49539381ba2SJed Brown 
4966562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
49701da6913SBarry Smith #undef __FUNCT__
498a06653b4SBarry Smith #define __FUNCT__ "PCReset_ML"
49916336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
50001da6913SBarry Smith {
50101da6913SBarry Smith   PetscErrorCode ierr;
502e0262f48SMatthew G Knepley   PC_MG          *mg    = (PC_MG*)pc->data;
503e0262f48SMatthew G Knepley   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
50439381ba2SJed Brown   PetscInt       level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim;
50501da6913SBarry Smith 
50601da6913SBarry Smith   PetscFunctionBegin;
50739381ba2SJed Brown   if (dim) {
50839381ba2SJed Brown     ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats*) pc_ml->ml_object->Grid[0].Grid;
50939381ba2SJed Brown 
51039381ba2SJed Brown     for (level=0; level<=fine_level; level++) {
51139381ba2SJed Brown       ierr = VecDestroy(&pc_ml->gridctx[level].coords);CHKERRQ(ierr);
51239381ba2SJed Brown     }
51339381ba2SJed Brown 
51439381ba2SJed Brown     grid_info->x = 0; /* do this so ML doesn't try to free coordinates */
51539381ba2SJed Brown     grid_info->y = 0;
51639381ba2SJed Brown     grid_info->z = 0;
51739381ba2SJed Brown 
518815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
51939381ba2SJed Brown   }
520815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Aggregate_Destroy(&pc_ml->agg_object));
521815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Destroy(&pc_ml->ml_object));
52201da6913SBarry Smith 
52301da6913SBarry Smith   if (pc_ml->PetscMLdata) {
52401da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
525ae7fe62dSJed Brown     ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);
526ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr);
527ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr);
52801da6913SBarry Smith   }
52901da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
53001da6913SBarry Smith 
531f5a5dd59SJed Brown   if (pc_ml->gridctx) {
53201da6913SBarry Smith     for (level=0; level<fine_level; level++) {
533601cad40SBrad Aagaard       if (pc_ml->gridctx[level].A) {ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);}
534601cad40SBrad Aagaard       if (pc_ml->gridctx[level].P) {ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);}
535601cad40SBrad Aagaard       if (pc_ml->gridctx[level].R) {ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);}
536601cad40SBrad Aagaard       if (pc_ml->gridctx[level].x) {ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);}
537601cad40SBrad Aagaard       if (pc_ml->gridctx[level].b) {ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);}
538601cad40SBrad Aagaard       if (pc_ml->gridctx[level+1].r) {ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
53901da6913SBarry Smith     }
540f5a5dd59SJed Brown   }
54101da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
54239381ba2SJed Brown   ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
5432fa5cd67SKarl Rupp 
54439381ba2SJed Brown   pc_ml->dim  = 0;
54539381ba2SJed Brown   pc_ml->nloc = 0;
54601da6913SBarry Smith   PetscFunctionReturn(0);
54701da6913SBarry Smith }
5485582bec1SHong Zhang /* -------------------------------------------------------------------------- */
5495582bec1SHong Zhang /*
5505582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
5515582bec1SHong Zhang                     by setting data structures and options.
5525582bec1SHong Zhang 
5535582bec1SHong Zhang    Input Parameter:
5545582bec1SHong Zhang .  pc - the preconditioner context
5555582bec1SHong Zhang 
5565582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5575582bec1SHong Zhang 
5585582bec1SHong Zhang    Notes:
5595582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5605582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5615582bec1SHong Zhang */
5626ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
563a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
564c07bf074SBarry Smith 
5655582bec1SHong Zhang #undef __FUNCT__
5665582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
5676ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5685582bec1SHong Zhang {
5695582bec1SHong Zhang   PetscErrorCode ierr;
570eef31507SHong Zhang   PetscMPIInt    size;
5715582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
5725582bec1SHong Zhang   ML             *ml_object;
5735582bec1SHong Zhang   ML_Aggregate   *agg_object;
5745582bec1SHong Zhang   ML_Operator    *mlmat;
5754f8eab3cSJed Brown   PetscInt       nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5765582bec1SHong Zhang   Mat            A,Aloc;
5775582bec1SHong Zhang   GridCtx        *gridctx;
57801da6913SBarry Smith   PC_MG          *mg    = (PC_MG*)pc->data;
57901da6913SBarry Smith   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
580ace3abfcSBarry Smith   PetscBool      isSeq, isMPI;
581c07bf074SBarry Smith   KSP            smoother;
582c07bf074SBarry Smith   PC             subpc;
58348268eb4SJed Brown   PetscInt       mesh_level, old_mesh_level;
58448268eb4SJed Brown 
5855582bec1SHong Zhang   PetscFunctionBegin;
58648268eb4SJed Brown   A    = pc->pmat;
587ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
58848268eb4SJed Brown 
589573998d7SHong Zhang   if (pc->setupcalled) {
59048268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
59148268eb4SJed Brown       /*
59248268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
59348268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
59448268eb4SJed Brown        */
59548268eb4SJed Brown       ml_object             = pc_ml->ml_object;
59648268eb4SJed Brown       gridctx               = pc_ml->gridctx;
59748268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
59848268eb4SJed Brown       fine_level            = Nlevels - 1;
59948268eb4SJed Brown       gridctx[fine_level].A = A;
60048268eb4SJed Brown 
601251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
602251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
60348268eb4SJed Brown       if (isMPI) {
6040298fd71SBarry Smith         ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
60548268eb4SJed Brown       } else if (isSeq) {
60648268eb4SJed Brown         Aloc = A;
607ae7fe62dSJed Brown         ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
608ce94432eSBarry Smith       } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name);
60948268eb4SJed Brown 
61048268eb4SJed Brown       ierr              = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
61148268eb4SJed Brown       PetscMLdata       = pc_ml->PetscMLdata;
612ae7fe62dSJed Brown       ierr              = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr);
61348268eb4SJed Brown       PetscMLdata->A    = A;
61448268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
615815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Destroy",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
616815d23e5SBarry Smith       PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
61748268eb4SJed Brown 
61848268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
61948268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
62048268eb4SJed Brown         old_mesh_level = mesh_level;
62148268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
62248268eb4SJed Brown 
62348268eb4SJed Brown         /* clean and regenerate A */
62448268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
625815d23e5SBarry Smith         PetscStackCall("ML_Operator_Clean",ML_Operator_Clean(mlmat));
626815d23e5SBarry Smith         PetscStackCall("ML_Operator_Init",ML_Operator_Init(mlmat,ml_object->comm));
627815d23e5SBarry Smith         PetscStackCall("ML_Gen_AmatrixRAP",ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
62848268eb4SJed Brown       }
62948268eb4SJed Brown 
63048268eb4SJed Brown       level = fine_level - 1;
63148268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
63248268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
63348268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
634ae7fe62dSJed Brown           ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
63548268eb4SJed Brown           level--;
63648268eb4SJed Brown         }
63748268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
63848268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
63948268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
640ae7fe62dSJed Brown           ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
64148268eb4SJed Brown           level--;
64248268eb4SJed Brown         }
64348268eb4SJed Brown       }
64448268eb4SJed Brown 
64548268eb4SJed Brown       for (level=0; level<fine_level; level++) {
64648268eb4SJed Brown         if (level > 0) {
64748268eb4SJed Brown           ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
64848268eb4SJed Brown         }
6497721a10bSJed Brown         ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
65048268eb4SJed Brown       }
65148268eb4SJed Brown       ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
6527721a10bSJed Brown       ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
65348268eb4SJed Brown 
65448268eb4SJed Brown       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
65548268eb4SJed Brown       PetscFunctionReturn(0);
65648268eb4SJed Brown     } else {
657c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
65816336fedSMatthew G Knepley       ierr = PCReset_ML(pc);CHKERRQ(ierr);
659a06653b4SBarry Smith       ierr = PCReset_MG(pc);CHKERRQ(ierr);
660573998d7SHong Zhang     }
66148268eb4SJed Brown   }
662573998d7SHong Zhang 
6635582bec1SHong Zhang   /* setup special features of PCML */
6645582bec1SHong Zhang   /*--------------------------------*/
6655582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6665582bec1SHong Zhang   pc_ml->size = size;
667251f4c67SDmitry Karpeev   ierr        = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
668251f4c67SDmitry Karpeev   ierr        = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
669864b637dSMatthew Knepley   if (isMPI) {
6700298fd71SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
671864b637dSMatthew Knepley   } else if (isSeq) {
6725582bec1SHong Zhang     Aloc = A;
673ae7fe62dSJed Brown     ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
674ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name);
6755582bec1SHong Zhang 
6765582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
67738f2d2fdSLisandro Dalcin   ierr               = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
6785582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
679d0f46423SBarry Smith   ierr               = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
6805582bec1SHong Zhang 
68124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
682d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
68324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
68424a42b14SHong Zhang 
68524a42b14SHong Zhang   ierr              = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
686d0f46423SBarry Smith   ierr              = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
68724a42b14SHong Zhang   ierr              = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
688573998d7SHong Zhang   PetscMLdata->A    = A;
689573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
69039381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
69139381ba2SJed Brown     PetscInt  i,j,dim=pc_ml->dim;
69239381ba2SJed Brown     PetscInt  nloc = pc_ml->nloc,nlocghost;
69339381ba2SJed Brown     PetscReal *ghostedcoords;
69439381ba2SJed Brown 
69539381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
69639381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
69739381ba2SJed Brown     ierr      = PetscMalloc(dim*nlocghost*sizeof(PetscReal),&ghostedcoords);CHKERRQ(ierr);
69839381ba2SJed Brown     for (i = 0; i < dim; i++) {
69939381ba2SJed Brown       /* copy coordinate values into first component of pwork */
70039381ba2SJed Brown       for (j = 0; j < nloc; j++) {
70139381ba2SJed Brown         PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
70239381ba2SJed Brown       }
70339381ba2SJed Brown       /* get the ghost values */
70439381ba2SJed Brown       ierr = PetscML_comm(PetscMLdata->pwork,PetscMLdata);CHKERRQ(ierr);
70539381ba2SJed Brown       /* write into the vector */
70639381ba2SJed Brown       for (j = 0; j < nlocghost; j++) {
70739381ba2SJed Brown         ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
70839381ba2SJed Brown       }
70939381ba2SJed Brown     }
71039381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
71139381ba2SJed Brown      * what ML needs */
71239381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
71339381ba2SJed Brown     pc_ml->coords = ghostedcoords;
71439381ba2SJed Brown   }
71524a42b14SHong Zhang 
7165582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
71745cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
7185582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
7194f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
720815d23e5SBarry Smith   PetscStackCall("ML_Create",ML_Create(&ml_object,pc_ml->MaxNlevels));
721ce94432eSBarry Smith   PetscStackCall("ML_Comm_Set_UsrComm",ML_Comm_Set_UsrComm(ml_object->comm,PetscObjectComm((PetscObject)A)));
722573998d7SHong Zhang   pc_ml->ml_object = ml_object;
723815d23e5SBarry Smith   PetscStackCall("ML_Init_Amatrix",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
724815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Getrow",ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols));
725815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
7265582bec1SHong Zhang 
727815d23e5SBarry Smith   PetscStackCall("ML_Set_Symmetrize",ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO));
728b5c8bdf8SJed Brown 
7295582bec1SHong Zhang   /* aggregation */
730815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Create",ML_Aggregate_Create(&agg_object));
731573998d7SHong Zhang   pc_ml->agg_object = agg_object;
732573998d7SHong Zhang 
733fb6a8e6dSJed Brown   {
734fb6a8e6dSJed Brown     MatNullSpace mnull;
735fb6a8e6dSJed Brown     ierr = MatGetNearNullSpace(A,&mnull);CHKERRQ(ierr);
736fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
737fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
738fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
739fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
740fb6a8e6dSJed Brown     }
741fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
742fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
743fb6a8e6dSJed Brown       PetscScalar       *nullvec;
744fb6a8e6dSJed Brown       const PetscScalar *v;
745fb6a8e6dSJed Brown       PetscBool         has_const;
7461c547e14SJed Brown       PetscInt          i,j,mlocal,nvec,M;
747fb6a8e6dSJed Brown       const Vec         *vecs;
7482fa5cd67SKarl Rupp 
749ce94432eSBarry Smith       if (!mnull) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
7500298fd71SBarry Smith       ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
7510298fd71SBarry Smith       ierr = MatGetLocalSize(Aloc,&mlocal,NULL);CHKERRQ(ierr);
752fb6a8e6dSJed Brown       ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr);
7538caf3d72SBarry Smith       ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr);
7541c547e14SJed Brown       if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M;
755fb6a8e6dSJed Brown       for (i=0; i<nvec; i++) {
756fb6a8e6dSJed Brown         ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
757fb6a8e6dSJed Brown         for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j];
758fb6a8e6dSJed Brown         ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
759fb6a8e6dSJed Brown       }
760815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Create",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal);CHKERRQ(ierr));
761fb6a8e6dSJed Brown       ierr = PetscFree(nullvec);CHKERRQ(ierr);
762fb6a8e6dSJed Brown     } break;
763fb6a8e6dSJed Brown     case PCML_NULLSPACE_BLOCK:
764815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_NullSpace",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr));
765fb6a8e6dSJed Brown       break;
766fb6a8e6dSJed Brown     case PCML_NULLSPACE_SCALAR:
767fb6a8e6dSJed Brown       break;
768ce94432eSBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unknown null space type");
769fb6a8e6dSJed Brown     }
770fb6a8e6dSJed Brown   }
771815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_MaxCoarseSize",ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize));
7725582bec1SHong Zhang   /* set options */
7735582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
7745582bec1SHong Zhang   case 1:
775815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_Coupled",ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));break;
7765582bec1SHong Zhang   case 2:
777815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_MIS",ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));break;
7785582bec1SHong Zhang   case 3:
779815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_METIS",ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));break;
7805582bec1SHong Zhang   }
781815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_Threshold",ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold));
782815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_DampingFactor",ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor));
7835582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm) {
784815d23e5SBarry Smith     PetscStackCall("ML_Set_SpectralNormScheme_Anorm",ML_Set_SpectralNormScheme_Anorm(ml_object));
7855582bec1SHong Zhang   }
786b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
787b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
788b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
789b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
790b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
791b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7925582bec1SHong Zhang 
79339381ba2SJed Brown   if (pc_ml->Aux) {
794ce94432eSBarry Smith     if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Auxiliary matrix requires coordinates");
79539381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
79639381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
79739381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
79839381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
79939381ba2SJed Brown   }
80039381ba2SJed Brown 
80139381ba2SJed Brown   if (pc_ml->dim) {
80239381ba2SJed Brown     PetscInt               i,dim = pc_ml->dim;
80339381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
80439381ba2SJed Brown     PetscInt               nlocghost;
80539381ba2SJed Brown 
80639381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
80739381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
80839381ba2SJed Brown 
809815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_VizAndStats_Setup(",ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
81039381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats*) ml_object->Grid[0].Grid;
81139381ba2SJed Brown     for (i = 0; i < dim; i++) {
81239381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
81339381ba2SJed Brown        * in pc_ml */
81439381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
81539381ba2SJed Brown       switch (i) {
81639381ba2SJed Brown       case 0: grid_info->x = pc_ml->coords + nlocghost * i; break;
81739381ba2SJed Brown       case 1: grid_info->y = pc_ml->coords + nlocghost * i; break;
81839381ba2SJed Brown       case 2: grid_info->z = pc_ml->coords + nlocghost * i; break;
819ce94432eSBarry Smith       default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
82039381ba2SJed Brown       }
82139381ba2SJed Brown     }
82239381ba2SJed Brown     grid_info->Ndim = dim;
82339381ba2SJed Brown   }
82439381ba2SJed Brown 
82539381ba2SJed Brown   /* repartitioning */
82639381ba2SJed Brown   if (pc_ml->Repartition) {
827815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Activate",ML_Repartition_Activate(ml_object));
828815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_LargestMinMaxRatio",ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio));
829815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_MinPerProc",ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc));
830815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_PutOnSingleProc",ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc));
83139381ba2SJed Brown #if 0                           /* Function not yet defined in ml-6.2 */
83239381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
83339381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
83439381ba2SJed Brown      * finest level (reflection default behavior in
83539381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
836815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
83739381ba2SJed Brown #endif
83839381ba2SJed Brown 
83939381ba2SJed Brown     if (!pc_ml->RepartitionType) {
84039381ba2SJed Brown       PetscInt i;
84139381ba2SJed Brown 
842ce94432eSBarry Smith       if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates");
843815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN));
844815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_Dimensions",ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
84539381ba2SJed Brown 
84639381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
84739381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Grid[i].Grid;
84839381ba2SJed Brown         grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
84939381ba2SJed Brown         /* defaults from ml_agg_info.c */
85039381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
85139381ba2SJed Brown         grid_info->zoltan_timers        = 0;
85239381ba2SJed Brown         grid_info->smoothing_steps      = 4;  /* only relevant to hypergraph / fast hypergraph */
85339381ba2SJed Brown       }
8542fa5cd67SKarl Rupp     } else {
855815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS));
85639381ba2SJed Brown     }
85739381ba2SJed Brown   }
85839381ba2SJed Brown 
859b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
860815d23e5SBarry Smith     PetscStackCall("ML_Gen_MGHierarchy_UsingAggregation",Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
861b5c8bdf8SJed Brown   } else {
862815d23e5SBarry Smith     PetscStackCall("ML_Gen_MultiLevelHierarchy_UsingAggregation",Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
863b5c8bdf8SJed Brown   }
864ce94432eSBarry Smith   if (Nlevels<=0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
865573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
866aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
867c07bf074SBarry Smith 
8680298fd71SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,NULL);CHKERRQ(ierr);
869aa85bbbfSHong Zhang   /* set default smoothers */
870aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++) {
871aa85bbbfSHong Zhang     ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
872aa85bbbfSHong Zhang     ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
873aa85bbbfSHong Zhang     ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
874aa85bbbfSHong Zhang     ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
875aa85bbbfSHong Zhang   }
876f2e59741SMatthew G Knepley   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
87797177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
878f2e59741SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8795582bec1SHong Zhang 
8805582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
8812fa5cd67SKarl Rupp 
8825582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8835582bec1SHong Zhang 
8845582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8855582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
886e14861a4SHong Zhang   gridctx[fine_level].A = A;
887573998d7SHong Zhang 
888e14861a4SHong Zhang   level = fine_level - 1;
889ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8905582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
891e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
892db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
893e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
894db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
895573998d7SHong Zhang 
896573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
897573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
8985582bec1SHong Zhang       level--;
8995582bec1SHong Zhang     }
900ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
9015582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
9025582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
903db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
904ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
905db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
906573998d7SHong Zhang 
9075582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
908ae7fe62dSJed Brown       ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
9095582bec1SHong Zhang       level--;
9105582bec1SHong Zhang     }
9115582bec1SHong Zhang   }
9125582bec1SHong Zhang 
913573998d7SHong Zhang   /* create vectors and ksp at all levels */
914ac346b81SHong Zhang   for (level=0; level<fine_level; level++) {
915573998d7SHong Zhang     level1 = level + 1;
916e64afeacSLisandro Dalcin     ierr   = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
917d0f46423SBarry Smith     ierr   = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
9185582bec1SHong Zhang     ierr   = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
91997177400SBarry Smith     ierr   = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
9205582bec1SHong Zhang 
921e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
922d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
9235582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
92497177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
925ac346b81SHong Zhang 
926e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
927d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
928ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
92997177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
930ac346b81SHong Zhang 
9315582bec1SHong Zhang     if (level == 0) {
93297177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
9335582bec1SHong Zhang     } else {
93497177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
935573998d7SHong Zhang     }
936573998d7SHong Zhang   }
937573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
938573998d7SHong Zhang 
939573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
940573998d7SHong Zhang   for (level=0; level<fine_level; level++) {
941573998d7SHong Zhang     level1 = level + 1;
942aea2a34eSBarry Smith     ierr   = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
943573998d7SHong Zhang     ierr   = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
944573998d7SHong Zhang     if (level > 0) {
94597177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
9465582bec1SHong Zhang     }
9475582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
9485582bec1SHong Zhang   }
94997177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
950ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
9515582bec1SHong Zhang 
95239381ba2SJed Brown   /* put coordinate info in levels */
95339381ba2SJed Brown   if (pc_ml->dim) {
95439381ba2SJed Brown     PetscInt  i,j,dim = pc_ml->dim;
95539381ba2SJed Brown     PetscInt  bs, nloc;
95639381ba2SJed Brown     PC        subpc;
95739381ba2SJed Brown     PetscReal *array;
95839381ba2SJed Brown 
95939381ba2SJed Brown     level = fine_level;
96039381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
961ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Amat[mllevel].to->Grid->Grid;
96239381ba2SJed Brown       MPI_Comm               comm       = ((PetscObject)gridctx[level].A)->comm;
96339381ba2SJed Brown 
96439381ba2SJed Brown       ierr  = MatGetBlockSize (gridctx[level].A, &bs);CHKERRQ(ierr);
9650298fd71SBarry Smith       ierr  = MatGetLocalSize (gridctx[level].A, NULL, &nloc);CHKERRQ(ierr);
96639381ba2SJed Brown       nloc /= bs; /* number of local nodes */
96739381ba2SJed Brown 
96839381ba2SJed Brown       ierr = VecCreate(comm,&gridctx[level].coords);CHKERRQ(ierr);
96939381ba2SJed Brown       ierr = VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE);CHKERRQ(ierr);
97039381ba2SJed Brown       ierr = VecSetType(gridctx[level].coords,VECMPI);CHKERRQ(ierr);
97139381ba2SJed Brown       ierr = VecGetArray(gridctx[level].coords,&array);CHKERRQ(ierr);
97239381ba2SJed Brown       for (j = 0; j < nloc; j++) {
97339381ba2SJed Brown         for (i = 0; i < dim; i++) {
97439381ba2SJed Brown           switch (i) {
97539381ba2SJed Brown           case 0: array[dim * j + i] = grid_info->x[j]; break;
97639381ba2SJed Brown           case 1: array[dim * j + i] = grid_info->y[j]; break;
97739381ba2SJed Brown           case 2: array[dim * j + i] = grid_info->z[j]; break;
978ce94432eSBarry Smith           default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
97939381ba2SJed Brown           }
98039381ba2SJed Brown         }
98139381ba2SJed Brown       }
98239381ba2SJed Brown 
98339381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
98439381ba2SJed Brown       ierr = KSPGetPC(gridctx[level].ksp,&subpc);CHKERRQ(ierr);
98539381ba2SJed Brown       ierr = PCSetCoordinates(subpc,dim,nloc,array);CHKERRQ(ierr);
98639381ba2SJed Brown       ierr = VecRestoreArray(gridctx[level].coords,&array);CHKERRQ(ierr);
98739381ba2SJed Brown       level--;
98839381ba2SJed Brown     }
98939381ba2SJed Brown   }
99039381ba2SJed Brown 
991c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
992c07bf074SBarry Smith   pc->setupcalled = 0;
9933751b4bdSBarry Smith   ierr            = PCSetUp_MG(pc);CHKERRQ(ierr);
9945582bec1SHong Zhang   PetscFunctionReturn(0);
9955582bec1SHong Zhang }
9965582bec1SHong Zhang 
9975582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9985582bec1SHong Zhang /*
9995582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
10005582bec1SHong Zhang    that was created with PCCreate_ML().
10015582bec1SHong Zhang 
10025582bec1SHong Zhang    Input Parameter:
10035582bec1SHong Zhang .  pc - the preconditioner context
10045582bec1SHong Zhang 
10055582bec1SHong Zhang    Application Interface Routine: PCDestroy()
10065582bec1SHong Zhang */
10075582bec1SHong Zhang #undef __FUNCT__
10085582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
10096ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
10105582bec1SHong Zhang {
10115582bec1SHong Zhang   PetscErrorCode ierr;
101201da6913SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
101301da6913SBarry Smith   PC_ML          *pc_ml= (PC_ML*)mg->innerctx;
10145582bec1SHong Zhang 
10155582bec1SHong Zhang   PetscFunctionBegin;
101616336fedSMatthew G Knepley   ierr = PCReset_ML(pc);CHKERRQ(ierr);
101701da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
101801da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
101900de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C","",NULL);CHKERRQ(ierr);
10205582bec1SHong Zhang   PetscFunctionReturn(0);
10215582bec1SHong Zhang }
10225582bec1SHong Zhang 
10235582bec1SHong Zhang #undef __FUNCT__
10245582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
10256ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
10265582bec1SHong Zhang {
10275582bec1SHong Zhang   PetscErrorCode ierr;
102839381ba2SJed Brown   PetscInt       indx,PrintLevel,partindx;
10295582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
103039381ba2SJed Brown   const char     *part[]   = {"Zoltan","ParMETIS"};
103139381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
103239381ba2SJed Brown   PetscInt   zidx;
103339381ba2SJed Brown   const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"};
103439381ba2SJed Brown #endif
103501da6913SBarry Smith   PC_MG       *mg    = (PC_MG*)pc->data;
103601da6913SBarry Smith   PC_ML       *pc_ml = (PC_ML*)mg->innerctx;
1037b5c8bdf8SJed Brown   PetscMPIInt size;
1038ce94432eSBarry Smith   MPI_Comm    comm;
10395582bec1SHong Zhang 
10405582bec1SHong Zhang   PetscFunctionBegin;
1041ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
104288ff4cc7SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
10435582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
10442fa5cd67SKarl Rupp 
10455582bec1SHong Zhang   PrintLevel = 0;
10465582bec1SHong Zhang   indx       = 0;
104739381ba2SJed Brown   partindx   = 0;
10482fa5cd67SKarl Rupp 
10490298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,NULL);CHKERRQ(ierr);
1050815d23e5SBarry Smith   PetscStackCall("ML_Set_PrintLeve",ML_Set_PrintLevel(PrintLevel));
10510298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,NULL);CHKERRQ(ierr);
10520298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,NULL);CHKERRQ(ierr);
10530298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,NULL);CHKERRQ(ierr);
10542fa5cd67SKarl Rupp 
10555582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
10562fa5cd67SKarl Rupp 
10570298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,NULL);CHKERRQ(ierr);
10580298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,NULL);CHKERRQ(ierr);
10590298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,NULL);CHKERRQ(ierr);
10600298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,NULL);CHKERRQ(ierr);
10610298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,NULL);CHKERRQ(ierr);
10620298fd71SBarry Smith   ierr = PetscOptionsEnum("-pc_ml_nullspace","Which type of null space information to use","None",PCMLNullSpaceTypes,(PetscEnum)pc_ml->nulltype,(PetscEnum*)&pc_ml->nulltype,NULL);CHKERRQ(ierr);
10630298fd71SBarry Smith   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,NULL);CHKERRQ(ierr);
10640298fd71SBarry Smith   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,NULL);CHKERRQ(ierr);
1065b5c8bdf8SJed Brown   /*
1066b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1067b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
1068b5c8bdf8SJed Brown 
1069b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1070b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
1071b5c8bdf8SJed Brown   */
107288ff4cc7SJed Brown   if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
107388ff4cc7SJed Brown   if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
1074b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2");CHKERRQ(ierr);}
1075b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
10760298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,NULL);CHKERRQ(ierr);
1077b5c8bdf8SJed Brown   }
1078b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
1079b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10800298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,NULL);CHKERRQ(ierr);
1081b5c8bdf8SJed Brown   }
1082b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1083b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
10840298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_KeepAggInfo","Allows the preconditioner to be reused, or auxilliary matrices to be generated","None",pc_ml->KeepAggInfo,&pc_ml->KeepAggInfo,NULL);CHKERRQ(ierr);
1085b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1086b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10870298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_Reusable","Store intermedaiate data structures so that the multilevel hierarchy is reusable","None",pc_ml->Reusable,&pc_ml->Reusable,NULL);CHKERRQ(ierr);
1088b5c8bdf8SJed Brown   /*
1089b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1090b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1091b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1092b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1093b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1094b5c8bdf8SJed Brown    */
10950298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,NULL);CHKERRQ(ierr);
10960298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the heirarchy","ML_Repartition_Activate",pc_ml->Repartition,&pc_ml->Repartition,NULL);CHKERRQ(ierr);
109739381ba2SJed Brown   if (pc_ml->Repartition) {
10980298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,NULL);CHKERRQ(ierr);
10990298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,NULL);CHKERRQ(ierr);
11000298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor","ML_Repartition_Set_PutOnSingleProc",pc_ml->PutOnSingleProc,&pc_ml->PutOnSingleProc,NULL);CHKERRQ(ierr);
110139381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
110239381ba2SJed Brown     partindx = 0;
11030298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,NULL);CHKERRQ(ierr);
11042fa5cd67SKarl Rupp 
110539381ba2SJed Brown     pc_ml->RepartitionType = partindx;
110639381ba2SJed Brown     if (!partindx) {
11075572b5bbSJed Brown       PetscInt zindx = 0;
11082fa5cd67SKarl Rupp 
11090298fd71SBarry Smith       ierr = PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,NULL);CHKERRQ(ierr);
11102fa5cd67SKarl Rupp 
111139381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
111239381ba2SJed Brown     }
111339381ba2SJed Brown #else
111439381ba2SJed Brown     partindx = 1;
11150298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,NULL);CHKERRQ(ierr);
1116ce94432eSBarry Smith     if (!partindx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan");
111739381ba2SJed Brown #endif
11180298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,NULL);CHKERRQ(ierr);
11190298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,NULL);CHKERRQ(ierr);
112039381ba2SJed Brown   }
11215582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
11225582bec1SHong Zhang   PetscFunctionReturn(0);
11235582bec1SHong Zhang }
11245582bec1SHong Zhang 
11255582bec1SHong Zhang /* -------------------------------------------------------------------------- */
11265582bec1SHong Zhang /*
11275582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
11285582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
11295582bec1SHong Zhang    context, PC, that was created within PCCreate().
11305582bec1SHong Zhang 
11315582bec1SHong Zhang    Input Parameter:
11325582bec1SHong Zhang .  pc - the preconditioner context
11335582bec1SHong Zhang 
11345582bec1SHong Zhang    Application Interface Routine: PCCreate()
11355582bec1SHong Zhang */
11365582bec1SHong Zhang 
11375582bec1SHong Zhang /*MC
11381e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
11395582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
11406ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
11416ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
11425582bec1SHong Zhang 
11436ca4d86aSHong Zhang    Options Database Key:
11446ca4d86aSHong Zhang    Multigrid options(inherited)
11456ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
11466ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
11476ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
11488c1c2452SJed Brown    -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
114951f519a2SBarry Smith    ML options:
1150ccd75124SBarry Smith .  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
11516ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
11526ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
1153f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
11546ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
11556ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
115639381ba2SJed Brown .  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
115739381ba2SJed Brown .  -pc_ml_repartition <false>: Allow ML to repartition levels of the heirarchy (ML_Repartition_Activate)
115839381ba2SJed Brown .  -pc_ml_repartitionMaxMinRatio <1.3>: Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio)
115939381ba2SJed Brown .  -pc_ml_repartitionMinPerProc <512>: Smallest repartitioned size (ML_Repartition_Set_MinPerProc)
116039381ba2SJed Brown .  -pc_ml_repartitionPutOnSingleProc <5000>: Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc)
116139381ba2SJed Brown .  -pc_ml_repartitionType <Zoltan>: Repartitioning library to use (ML_Repartition_Set_Partitioner)
116239381ba2SJed Brown .  -pc_ml_repartitionZoltanScheme <RCB>: Repartitioning scheme to use (None)
116339381ba2SJed Brown .  -pc_ml_Aux <false>: Aggregate using auxiliary coordinate-based laplacian (None)
116439381ba2SJed Brown -  -pc_ml_AuxThreshold <0.0>: Auxiliary smoother drop tol (None)
11655582bec1SHong Zhang 
11665582bec1SHong Zhang    Level: intermediate
11675582bec1SHong Zhang 
11685582bec1SHong Zhang   Concepts: multigrid
11695582bec1SHong Zhang 
11705582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
117197177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
117297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
117397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
117497177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11755582bec1SHong Zhang M*/
11765582bec1SHong Zhang 
11775582bec1SHong Zhang #undef __FUNCT__
11785582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
1179*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
11805582bec1SHong Zhang {
11815582bec1SHong Zhang   PetscErrorCode ierr;
11825582bec1SHong Zhang   PC_ML          *pc_ml;
118301da6913SBarry Smith   PC_MG          *mg;
11845582bec1SHong Zhang 
11855582bec1SHong Zhang   PetscFunctionBegin;
1186573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11875582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
118803bfa161SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
1189e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
1190e0f5d30fSBarry Smith   ierr         = DMDestroy(&pc->dm);CHKERRQ(ierr);
1191e0f5d30fSBarry Smith   mg           = (PC_MG*)pc->data;
1192c91913d3SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
11935582bec1SHong Zhang 
11945582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
119538f2d2fdSLisandro Dalcin   ierr         = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
119601da6913SBarry Smith   mg->innerctx = pc_ml;
11975582bec1SHong Zhang 
1198573998d7SHong Zhang   pc_ml->ml_object                = 0;
1199573998d7SHong Zhang   pc_ml->agg_object               = 0;
1200573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1201573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1202573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1203573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1204573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
12053751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1206573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1207573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0/3.0;
1208573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1209573998d7SHong Zhang   pc_ml->size                     = 0;
121039381ba2SJed Brown   pc_ml->dim                      = 0;
121139381ba2SJed Brown   pc_ml->nloc                     = 0;
121239381ba2SJed Brown   pc_ml->coords                   = 0;
121339381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
121439381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
121539381ba2SJed Brown   pc_ml->MinPerProc               = 512;
121639381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
121739381ba2SJed Brown   pc_ml->RepartitionType          = 0;
121839381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
121939381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
122039381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
122139381ba2SJed Brown 
122239381ba2SJed Brown   /* allow for coordinates to be passed */
122300de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C","PCSetCoordinates_ML",PCSetCoordinates_ML);CHKERRQ(ierr);
1224573998d7SHong Zhang 
12255582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
12265582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
12275582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1228a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
12295582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
12305582bec1SHong Zhang   PetscFunctionReturn(0);
12315582bec1SHong Zhang }
1232