xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 6562c4e16f1ba21ea29f57503f81105a84770009)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
42dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
57ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
67ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
75582bec1SHong Zhang */
86356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
97c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
107c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
117c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
12cb5d8e9eSHong Zhang 
135582bec1SHong Zhang #include <math.h>
142cf39c26SSatish Balay EXTERN_C_BEGIN
1568210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1668210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1768210224SSatish Balay #define HAVE_CONFIG_H
1868210224SSatish Balay #endif
195582bec1SHong Zhang #include "ml_include.h"
205582bec1SHong Zhang EXTERN_C_END
215582bec1SHong Zhang 
225582bec1SHong Zhang /* The context (data structure) at each grid level */
235582bec1SHong Zhang typedef struct {
245582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
255582bec1SHong Zhang   Mat        A,P,R;
265582bec1SHong Zhang   KSP        ksp;
275582bec1SHong Zhang } GridCtx;
285582bec1SHong Zhang 
295582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
305582bec1SHong Zhang typedef struct {
31573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
32573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3324a42b14SHong Zhang   Vec          x,y;
345582bec1SHong Zhang   ML_Operator  *mlmat;
355582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
365582bec1SHong Zhang } FineGridCtx;
375582bec1SHong Zhang 
385582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
395582bec1SHong Zhang typedef struct {
405582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
415582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
425582bec1SHong Zhang   Vec          y;
435582bec1SHong Zhang } Mat_MLShell;
445582bec1SHong Zhang 
455582bec1SHong Zhang /* Private context for the ML preconditioner */
465582bec1SHong Zhang typedef struct {
475582bec1SHong Zhang   ML             *ml_object;
485582bec1SHong Zhang   ML_Aggregate   *agg_object;
495582bec1SHong Zhang   GridCtx        *gridctx;
505582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
51573998d7SHong Zhang   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme;
525582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
535582bec1SHong Zhang   PetscTruth     SpectralNormScheme_Anorm;
54573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
555582bec1SHong Zhang } PC_ML;
5641ca0015SHong Zhang 
57*6562c4e1SBarry Smith #undef __FUNCT__
58*6562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow"
59*6562c4e1SBarry 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[])
60*6562c4e1SBarry Smith {
61*6562c4e1SBarry Smith   PetscErrorCode ierr;
62*6562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
63*6562c4e1SBarry Smith   PetscScalar    *aa;
64*6562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
65*6562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
665582bec1SHong Zhang 
67*6562c4e1SBarry Smith 
68*6562c4e1SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,PETSC_NULL); if (ierr) return(0);
69*6562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
70*6562c4e1SBarry Smith     row   = requested_rows[i];
71*6562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
72*6562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
73*6562c4e1SBarry Smith     if ( (row >= 0) || (row <= (m-1)) ) {
74*6562c4e1SBarry Smith       aj = a->j + a->i[row];
75*6562c4e1SBarry Smith       aa = a->a + a->i[row];
76*6562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++){
77*6562c4e1SBarry Smith         columns[k]  = aj[j];
78*6562c4e1SBarry Smith         values[k++] = aa[j];
79*6562c4e1SBarry Smith       }
80*6562c4e1SBarry Smith     }
81*6562c4e1SBarry Smith   }
82*6562c4e1SBarry Smith   return(1);
83*6562c4e1SBarry Smith }
84*6562c4e1SBarry Smith 
85*6562c4e1SBarry Smith #undef __FUNCT__
86*6562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm"
87*6562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
88*6562c4e1SBarry Smith {
89*6562c4e1SBarry Smith   PetscErrorCode ierr;
90*6562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_data;
91*6562c4e1SBarry Smith   Mat            A=ml->A;
92*6562c4e1SBarry Smith   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
93*6562c4e1SBarry Smith   PetscMPIInt    size;
94*6562c4e1SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
95*6562c4e1SBarry Smith   PetscScalar    *array;
96*6562c4e1SBarry Smith 
97*6562c4e1SBarry Smith   PetscFunctionBegin;
98*6562c4e1SBarry Smith   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
99*6562c4e1SBarry Smith   if (size == 1) return 0;
100*6562c4e1SBarry Smith 
101*6562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
102*6562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
103*6562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
104*6562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
105*6562c4e1SBarry Smith   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
106*6562c4e1SBarry Smith   for (i=in_length; i<out_length; i++){
107*6562c4e1SBarry Smith     p[i] = array[i-in_length];
108*6562c4e1SBarry Smith   }
109*6562c4e1SBarry Smith   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
110*6562c4e1SBarry Smith   PetscFunctionReturn(0);
111*6562c4e1SBarry Smith }
112*6562c4e1SBarry Smith 
113*6562c4e1SBarry Smith #undef __FUNCT__
114*6562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec"
115*6562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
116*6562c4e1SBarry Smith {
117*6562c4e1SBarry Smith   PetscErrorCode ierr;
118*6562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
119*6562c4e1SBarry Smith   Mat            A=ml->A, Aloc=ml->Aloc;
120*6562c4e1SBarry Smith   PetscMPIInt    size;
121*6562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
122*6562c4e1SBarry Smith   PetscInt       i;
123*6562c4e1SBarry Smith 
124*6562c4e1SBarry Smith   PetscFunctionBegin;
125*6562c4e1SBarry Smith   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
126*6562c4e1SBarry Smith   if (size == 1){
127*6562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
128*6562c4e1SBarry Smith   } else {
129*6562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
130*6562c4e1SBarry Smith     PetscML_comm(pwork,ml);
131*6562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
132*6562c4e1SBarry Smith   }
133*6562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
134*6562c4e1SBarry Smith   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
135*6562c4e1SBarry Smith   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
136*6562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
137*6562c4e1SBarry Smith   PetscFunctionReturn(0);
138*6562c4e1SBarry Smith }
139*6562c4e1SBarry Smith 
140*6562c4e1SBarry Smith #undef __FUNCT__
141*6562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML"
142*6562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
143*6562c4e1SBarry Smith {
144*6562c4e1SBarry Smith   PetscErrorCode   ierr;
145*6562c4e1SBarry Smith   Mat_MLShell      *shell;
146*6562c4e1SBarry Smith   PetscScalar      *xarray,*yarray;
147*6562c4e1SBarry Smith   PetscInt         x_length,y_length;
148*6562c4e1SBarry Smith 
149*6562c4e1SBarry Smith   PetscFunctionBegin;
150*6562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
151*6562c4e1SBarry Smith   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
152*6562c4e1SBarry Smith   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
153*6562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
154*6562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
155*6562c4e1SBarry Smith   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
156*6562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
157*6562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
158*6562c4e1SBarry Smith   PetscFunctionReturn(0);
159*6562c4e1SBarry Smith }
160*6562c4e1SBarry Smith 
161*6562c4e1SBarry Smith #undef __FUNCT__
162*6562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML"
163*6562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
164*6562c4e1SBarry Smith {
165*6562c4e1SBarry Smith   PetscErrorCode    ierr;
166*6562c4e1SBarry Smith   Mat_MLShell       *shell;
167*6562c4e1SBarry Smith   PetscScalar       *xarray,*yarray;
168*6562c4e1SBarry Smith   PetscInt          x_length,y_length;
169*6562c4e1SBarry Smith 
170*6562c4e1SBarry Smith   PetscFunctionBegin;
171*6562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
172*6562c4e1SBarry Smith   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
173*6562c4e1SBarry Smith   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
174*6562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
175*6562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
176*6562c4e1SBarry Smith   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
177*6562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
178*6562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
179*6562c4e1SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
180*6562c4e1SBarry Smith   PetscFunctionReturn(0);
181*6562c4e1SBarry Smith }
182*6562c4e1SBarry Smith 
183*6562c4e1SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */
184*6562c4e1SBarry Smith #undef __FUNCT__
185*6562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML"
186*6562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
187*6562c4e1SBarry Smith {
188*6562c4e1SBarry Smith   PetscErrorCode  ierr;
189*6562c4e1SBarry Smith   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
190*6562c4e1SBarry Smith   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
191*6562c4e1SBarry Smith   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
192*6562c4e1SBarry Smith   PetscScalar     *aa=a->a,*ba=b->a,*ca;
193*6562c4e1SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
194*6562c4e1SBarry Smith   PetscInt        *ci,*cj,ncols;
195*6562c4e1SBarry Smith 
196*6562c4e1SBarry Smith   PetscFunctionBegin;
197*6562c4e1SBarry Smith   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
198*6562c4e1SBarry Smith 
199*6562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX){
200*6562c4e1SBarry Smith     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
201*6562c4e1SBarry Smith     ci[0] = 0;
202*6562c4e1SBarry Smith     for (i=0; i<am; i++){
203*6562c4e1SBarry Smith       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
204*6562c4e1SBarry Smith     }
205*6562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
206*6562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
207*6562c4e1SBarry Smith 
208*6562c4e1SBarry Smith     k = 0;
209*6562c4e1SBarry Smith     for (i=0; i<am; i++){
210*6562c4e1SBarry Smith       /* diagonal portion of A */
211*6562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
212*6562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
213*6562c4e1SBarry Smith         cj[k]   = *aj++;
214*6562c4e1SBarry Smith         ca[k++] = *aa++;
215*6562c4e1SBarry Smith       }
216*6562c4e1SBarry Smith       /* off-diagonal portion of A */
217*6562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
218*6562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
219*6562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
220*6562c4e1SBarry Smith         ca[k++] = *ba++;
221*6562c4e1SBarry Smith       }
222*6562c4e1SBarry Smith     }
223*6562c4e1SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
224*6562c4e1SBarry Smith 
225*6562c4e1SBarry Smith     /* put together the new matrix */
226*6562c4e1SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
227*6562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
228*6562c4e1SBarry Smith 
229*6562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
230*6562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
231*6562c4e1SBarry Smith     mat = (Mat_SeqAIJ*)(*Aloc)->data;
232*6562c4e1SBarry Smith     mat->free_a       = PETSC_TRUE;
233*6562c4e1SBarry Smith     mat->free_ij      = PETSC_TRUE;
234*6562c4e1SBarry Smith 
235*6562c4e1SBarry Smith     mat->nonew    = 0;
236*6562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX){
237*6562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
238*6562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
239*6562c4e1SBarry Smith     for (i=0; i<am; i++) {
240*6562c4e1SBarry Smith       /* diagonal portion of A */
241*6562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
242*6562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
243*6562c4e1SBarry Smith       /* off-diagonal portion of A */
244*6562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
245*6562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
246*6562c4e1SBarry Smith     }
247*6562c4e1SBarry Smith   } else {
248*6562c4e1SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
249*6562c4e1SBarry Smith   }
250*6562c4e1SBarry Smith   PetscFunctionReturn(0);
251*6562c4e1SBarry Smith }
252*6562c4e1SBarry Smith 
253*6562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat);
254*6562c4e1SBarry Smith #undef __FUNCT__
255*6562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML"
256*6562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
257*6562c4e1SBarry Smith {
258*6562c4e1SBarry Smith   PetscErrorCode ierr;
259*6562c4e1SBarry Smith   Mat_MLShell    *shell;
260*6562c4e1SBarry Smith 
261*6562c4e1SBarry Smith   PetscFunctionBegin;
262*6562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
263*6562c4e1SBarry Smith   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
264*6562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
265*6562c4e1SBarry Smith   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
266*6562c4e1SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
267*6562c4e1SBarry Smith   PetscFunctionReturn(0);
268*6562c4e1SBarry Smith }
269*6562c4e1SBarry Smith 
270*6562c4e1SBarry Smith #undef __FUNCT__
271*6562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ"
272*6562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
273*6562c4e1SBarry Smith {
274*6562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
275*6562c4e1SBarry Smith   PetscErrorCode        ierr;
276*6562c4e1SBarry Smith   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
277*6562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
278*6562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
279*6562c4e1SBarry Smith 
280*6562c4e1SBarry Smith   PetscFunctionBegin;
281*6562c4e1SBarry Smith   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
282*6562c4e1SBarry Smith   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
283*6562c4e1SBarry Smith     if (reuse){
284*6562c4e1SBarry Smith       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
285*6562c4e1SBarry Smith       aij->i = ml_rowptr;
286*6562c4e1SBarry Smith       aij->j = ml_cols;
287*6562c4e1SBarry Smith       aij->a = ml_vals;
288*6562c4e1SBarry Smith     } else {
289*6562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
290*6562c4e1SBarry Smith       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
291*6562c4e1SBarry Smith       for (i=0; i<m; i++) {
292*6562c4e1SBarry Smith         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
293*6562c4e1SBarry Smith       }
294*6562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
295*6562c4e1SBarry Smith       for (i=0; i<m; i++){
296*6562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
297*6562c4e1SBarry Smith         aj += nnz[i]; aa += nnz[i];
298*6562c4e1SBarry Smith       }
299*6562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
300*6562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
301*6562c4e1SBarry Smith     }
302*6562c4e1SBarry Smith     PetscFunctionReturn(0);
303*6562c4e1SBarry Smith   }
304*6562c4e1SBarry Smith 
305*6562c4e1SBarry Smith   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
306*6562c4e1SBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
307*6562c4e1SBarry Smith   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
308*6562c4e1SBarry Smith   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
309*6562c4e1SBarry Smith 
310*6562c4e1SBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
311*6562c4e1SBarry Smith   nz_max = 1;
312*6562c4e1SBarry Smith   for (i=0; i<m; i++) {
313*6562c4e1SBarry Smith     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
314*6562c4e1SBarry Smith     if (nnz[i] > nz_max) nz_max += nnz[i];
315*6562c4e1SBarry Smith   }
316*6562c4e1SBarry Smith 
317*6562c4e1SBarry Smith   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
318*6562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
319*6562c4e1SBarry Smith   for (i=0; i<m; i++){
320*6562c4e1SBarry Smith     k = 0;
321*6562c4e1SBarry Smith     /* diagonal entry */
322*6562c4e1SBarry Smith     aj[k] = i; aa[k++] = ml_vals[i];
323*6562c4e1SBarry Smith     /* off diagonal entries */
324*6562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
325*6562c4e1SBarry Smith       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
326*6562c4e1SBarry Smith     }
327*6562c4e1SBarry Smith     /* sort aj and aa */
328*6562c4e1SBarry Smith     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
329*6562c4e1SBarry Smith     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
330*6562c4e1SBarry Smith   }
331*6562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
332*6562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
333*6562c4e1SBarry Smith 
334*6562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
335*6562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
336*6562c4e1SBarry Smith   PetscFunctionReturn(0);
337*6562c4e1SBarry Smith }
338*6562c4e1SBarry Smith 
339*6562c4e1SBarry Smith #undef __FUNCT__
340*6562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL"
341*6562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
342*6562c4e1SBarry Smith {
343*6562c4e1SBarry Smith   PetscErrorCode ierr;
344*6562c4e1SBarry Smith   PetscInt       m,n;
345*6562c4e1SBarry Smith   ML_Comm        *MLcomm;
346*6562c4e1SBarry Smith   Mat_MLShell    *shellctx;
347*6562c4e1SBarry Smith 
348*6562c4e1SBarry Smith   PetscFunctionBegin;
349*6562c4e1SBarry Smith   m = mlmat->outvec_leng;
350*6562c4e1SBarry Smith   n = mlmat->invec_leng;
351*6562c4e1SBarry Smith   if (!m || !n){
352*6562c4e1SBarry Smith     newmat = PETSC_NULL;
353*6562c4e1SBarry Smith     PetscFunctionReturn(0);
354*6562c4e1SBarry Smith   }
355*6562c4e1SBarry Smith 
356*6562c4e1SBarry Smith   if (reuse){
357*6562c4e1SBarry Smith     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
358*6562c4e1SBarry Smith     shellctx->mlmat = mlmat;
359*6562c4e1SBarry Smith     PetscFunctionReturn(0);
360*6562c4e1SBarry Smith   }
361*6562c4e1SBarry Smith 
362*6562c4e1SBarry Smith   MLcomm = mlmat->comm;
363*6562c4e1SBarry Smith   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
364*6562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
365*6562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
366*6562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
367*6562c4e1SBarry Smith   shellctx->A         = *newmat;
368*6562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
369*6562c4e1SBarry Smith   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
370*6562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
371*6562c4e1SBarry Smith   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
372*6562c4e1SBarry Smith   (*newmat)->ops->destroy = MatDestroy_ML;
373*6562c4e1SBarry Smith   PetscFunctionReturn(0);
374*6562c4e1SBarry Smith }
375*6562c4e1SBarry Smith 
376*6562c4e1SBarry Smith #undef __FUNCT__
377*6562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ"
378*6562c4e1SBarry Smith static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
379*6562c4e1SBarry Smith {
380*6562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
381*6562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*aj;
382*6562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
383*6562c4e1SBarry Smith   PetscErrorCode        ierr;
384*6562c4e1SBarry Smith   PetscInt              i,j,k,*gordering;
385*6562c4e1SBarry Smith   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
386*6562c4e1SBarry Smith   Mat                   A;
387*6562c4e1SBarry Smith 
388*6562c4e1SBarry Smith   PetscFunctionBegin;
389*6562c4e1SBarry Smith   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
390*6562c4e1SBarry Smith   n = mlmat->invec_leng;
391*6562c4e1SBarry Smith   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
392*6562c4e1SBarry Smith 
393*6562c4e1SBarry Smith   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
394*6562c4e1SBarry Smith   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
395*6562c4e1SBarry Smith   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
396*6562c4e1SBarry Smith   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
397*6562c4e1SBarry Smith 
398*6562c4e1SBarry Smith   nz_max = 0;
399*6562c4e1SBarry Smith   for (i=0; i<m; i++){
400*6562c4e1SBarry Smith     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
401*6562c4e1SBarry Smith     if (nz_max < nnz[i]) nz_max = nnz[i];
402*6562c4e1SBarry Smith     nnzA[i] = 1; /* diag */
403*6562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
404*6562c4e1SBarry Smith       if (ml_cols[j] < m) nnzA[i]++;
405*6562c4e1SBarry Smith     }
406*6562c4e1SBarry Smith     nnzB[i] = nnz[i] - nnzA[i];
407*6562c4e1SBarry Smith   }
408*6562c4e1SBarry Smith   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
409*6562c4e1SBarry Smith 
410*6562c4e1SBarry Smith   /* insert mat values -- remap row and column indices */
411*6562c4e1SBarry Smith   nz_max++;
412*6562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
413*6562c4e1SBarry Smith   /* create global row numbering for a ML_Operator */
414*6562c4e1SBarry Smith   ML_build_global_numbering(mlmat,&gordering,"rows");
415*6562c4e1SBarry Smith   for (i=0; i<m; i++){
416*6562c4e1SBarry Smith     row = gordering[i];
417*6562c4e1SBarry Smith     k = 0;
418*6562c4e1SBarry Smith     /* diagonal entry */
419*6562c4e1SBarry Smith     aj[k] = row; aa[k++] = ml_vals[i];
420*6562c4e1SBarry Smith     /* off diagonal entries */
421*6562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
422*6562c4e1SBarry Smith       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
423*6562c4e1SBarry Smith     }
424*6562c4e1SBarry Smith     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
425*6562c4e1SBarry Smith   }
426*6562c4e1SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
427*6562c4e1SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
428*6562c4e1SBarry Smith   *newmat = A;
429*6562c4e1SBarry Smith 
430*6562c4e1SBarry Smith   ierr = PetscFree3(nnzA,nnzB,nnz);
431*6562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
432*6562c4e1SBarry Smith   PetscFunctionReturn(0);
433*6562c4e1SBarry Smith }
434*6562c4e1SBarry Smith 
435*6562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
43601da6913SBarry Smith #undef __FUNCT__
4373751b4bdSBarry Smith #define __FUNCT__ "PCDestroy_ML_Private"
4383751b4bdSBarry Smith PetscErrorCode PCDestroy_ML_Private(void *ptr)
43901da6913SBarry Smith {
44001da6913SBarry Smith   PetscErrorCode  ierr;
44101da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)ptr;
44201da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
44301da6913SBarry Smith 
44401da6913SBarry Smith   PetscFunctionBegin;
44501da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
44601da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
44701da6913SBarry Smith 
44801da6913SBarry Smith   if (pc_ml->PetscMLdata) {
44901da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
45001da6913SBarry Smith     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
45101da6913SBarry Smith     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
45201da6913SBarry Smith     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
45301da6913SBarry Smith   }
45401da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
45501da6913SBarry Smith 
45601da6913SBarry Smith   for (level=0; level<fine_level; level++){
45701da6913SBarry Smith     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
45801da6913SBarry Smith     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
45901da6913SBarry Smith     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
46001da6913SBarry Smith     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
46101da6913SBarry Smith     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
46201da6913SBarry Smith     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
46301da6913SBarry Smith   }
46401da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
46501da6913SBarry Smith   PetscFunctionReturn(0);
46601da6913SBarry Smith }
4675582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4685582bec1SHong Zhang /*
4695582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
4705582bec1SHong Zhang                     by setting data structures and options.
4715582bec1SHong Zhang 
4725582bec1SHong Zhang    Input Parameter:
4735582bec1SHong Zhang .  pc - the preconditioner context
4745582bec1SHong Zhang 
4755582bec1SHong Zhang    Application Interface Routine: PCSetUp()
4765582bec1SHong Zhang 
4775582bec1SHong Zhang    Notes:
4785582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
4795582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
4805582bec1SHong Zhang */
4816ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
482c07bf074SBarry Smith extern PetscErrorCode PCDestroy_MG_Private(PC);
483c07bf074SBarry Smith 
4845582bec1SHong Zhang #undef __FUNCT__
4855582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
4866ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
4875582bec1SHong Zhang {
4885582bec1SHong Zhang   PetscErrorCode  ierr;
489eef31507SHong Zhang   PetscMPIInt     size;
4905582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
4915582bec1SHong Zhang   ML              *ml_object;
4925582bec1SHong Zhang   ML_Aggregate    *agg_object;
4935582bec1SHong Zhang   ML_Operator     *mlmat;
4944f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
4955582bec1SHong Zhang   Mat             A,Aloc;
4965582bec1SHong Zhang   GridCtx         *gridctx;
49701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
49801da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
499864b637dSMatthew Knepley   PetscTruth      isSeq, isMPI;
500c07bf074SBarry Smith   KSP             smoother;
501c07bf074SBarry Smith   PC              subpc;
5025582bec1SHong Zhang 
5035582bec1SHong Zhang   PetscFunctionBegin;
504573998d7SHong Zhang   if (pc->setupcalled){
505c07bf074SBarry Smith     /* since ML can change the size of vectors/matrices at any level we must destroy everything */
5063751b4bdSBarry Smith     ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
507c07bf074SBarry Smith     ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr);
508573998d7SHong Zhang   }
509573998d7SHong Zhang 
5105582bec1SHong Zhang   /* setup special features of PCML */
5115582bec1SHong Zhang   /*--------------------------------*/
5125582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
5135582bec1SHong Zhang   A = pc->pmat;
5147adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5155582bec1SHong Zhang   pc_ml->size = size;
516864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
517864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
518864b637dSMatthew Knepley   if (isMPI){
519db571536SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
520864b637dSMatthew Knepley   } else if (isSeq) {
5215582bec1SHong Zhang     Aloc = A;
522864b637dSMatthew Knepley   } else {
523864b637dSMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices.");
5245582bec1SHong Zhang   }
5255582bec1SHong Zhang 
5265582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
52738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
5285582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
529d0f46423SBarry Smith   ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
5305582bec1SHong Zhang 
53124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
532d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
53324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
53424a42b14SHong Zhang 
53524a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
536d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
53724a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
538573998d7SHong Zhang   PetscMLdata->A    = A;
539573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
54024a42b14SHong Zhang 
5415582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
54245cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
5435582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
5444f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
5455582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
546573998d7SHong Zhang   pc_ml->ml_object = ml_object;
5475582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
5485582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
5495582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
5505582bec1SHong Zhang 
5515582bec1SHong Zhang   /* aggregation */
5525582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
553573998d7SHong Zhang   pc_ml->agg_object = agg_object;
554573998d7SHong Zhang 
5554f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
5565582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
5575582bec1SHong Zhang   /* set options */
5585582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
5595582bec1SHong Zhang   case 1:
5605582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
5615582bec1SHong Zhang   case 2:
5625582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
5635582bec1SHong Zhang   case 3:
5645582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
5655582bec1SHong Zhang   }
5665582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
5675582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
5685582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
5697ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
5705582bec1SHong Zhang   }
5715582bec1SHong Zhang 
5725582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
5735582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
574573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
575aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
576c07bf074SBarry Smith 
57797177400SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
578aa85bbbfSHong Zhang   /* set default smoothers */
579aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++){
580aa85bbbfSHong Zhang     if (size == 1){
581aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
582aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
583aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
584aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
585aa85bbbfSHong Zhang     } else {
586aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
587aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
588aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
589aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
590aa85bbbfSHong Zhang     }
591aa85bbbfSHong Zhang   }
59297177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
5935582bec1SHong Zhang 
5945582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
5955582bec1SHong Zhang   pc_ml->gridctx = gridctx;
5965582bec1SHong Zhang 
5975582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
5985582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
599e14861a4SHong Zhang   gridctx[fine_level].A = A;
600573998d7SHong Zhang 
601e14861a4SHong Zhang   level = fine_level - 1;
602ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
6035582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
604e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
605db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
606e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
607db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
608573998d7SHong Zhang 
609573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
610573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
6115582bec1SHong Zhang       level--;
6125582bec1SHong Zhang     }
613ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
6145582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
6155582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
616db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
617ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
618db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
619573998d7SHong Zhang 
6205582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
621eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
6225582bec1SHong Zhang       level--;
6235582bec1SHong Zhang     }
6245582bec1SHong Zhang   }
6255582bec1SHong Zhang 
626573998d7SHong Zhang   /* create vectors and ksp at all levels */
627ac346b81SHong Zhang   for (level=0; level<fine_level; level++){
628573998d7SHong Zhang     level1 = level + 1;
629e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
630d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
6315582bec1SHong Zhang     ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
63297177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
6335582bec1SHong Zhang 
634e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
635d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
6365582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
63797177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
638ac346b81SHong Zhang 
639e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
640d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
641ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
64297177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
643ac346b81SHong Zhang 
6445582bec1SHong Zhang     if (level == 0){
64597177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
6465582bec1SHong Zhang     } else {
64797177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
648573998d7SHong Zhang     }
649573998d7SHong Zhang   }
650573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
651573998d7SHong Zhang 
652573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
653573998d7SHong Zhang   for (level=0; level<fine_level; level++){
654573998d7SHong Zhang     level1 = level + 1;
655aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
656573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
657573998d7SHong Zhang     if (level > 0){
65897177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
6595582bec1SHong Zhang     }
6605582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
6615582bec1SHong Zhang   }
66297177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
663ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
6645582bec1SHong Zhang 
665c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
666c07bf074SBarry Smith   pc->setupcalled = 0;
6673751b4bdSBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6685582bec1SHong Zhang   PetscFunctionReturn(0);
6695582bec1SHong Zhang }
6705582bec1SHong Zhang 
6715582bec1SHong Zhang /* -------------------------------------------------------------------------- */
6725582bec1SHong Zhang /*
6735582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
6745582bec1SHong Zhang    that was created with PCCreate_ML().
6755582bec1SHong Zhang 
6765582bec1SHong Zhang    Input Parameter:
6775582bec1SHong Zhang .  pc - the preconditioner context
6785582bec1SHong Zhang 
6795582bec1SHong Zhang    Application Interface Routine: PCDestroy()
6805582bec1SHong Zhang */
6815582bec1SHong Zhang #undef __FUNCT__
6825582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
6836ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
6845582bec1SHong Zhang {
6855582bec1SHong Zhang   PetscErrorCode  ierr;
68601da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
68701da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
6885582bec1SHong Zhang 
6895582bec1SHong Zhang   PetscFunctionBegin;
6903751b4bdSBarry Smith   ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
69101da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
69201da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
6935582bec1SHong Zhang   PetscFunctionReturn(0);
6945582bec1SHong Zhang }
6955582bec1SHong Zhang 
6965582bec1SHong Zhang #undef __FUNCT__
6975582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
6986ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
6995582bec1SHong Zhang {
7005582bec1SHong Zhang   PetscErrorCode  ierr;
7013751b4bdSBarry Smith   PetscInt        indx,PrintLevel;
7025582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
70301da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
70401da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
7055582bec1SHong Zhang 
7065582bec1SHong Zhang   PetscFunctionBegin;
7075582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
7085582bec1SHong Zhang   PrintLevel    = 0;
7095582bec1SHong Zhang   indx          = 0;
7105582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
7115582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
712573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
713573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
7143751b4bdSBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
7155582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
716573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
717573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
71840597110SHong Zhang   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
7195582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
7205582bec1SHong Zhang   PetscFunctionReturn(0);
7215582bec1SHong Zhang }
7225582bec1SHong Zhang 
7235582bec1SHong Zhang /* -------------------------------------------------------------------------- */
7245582bec1SHong Zhang /*
7255582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
7265582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
7275582bec1SHong Zhang    context, PC, that was created within PCCreate().
7285582bec1SHong Zhang 
7295582bec1SHong Zhang    Input Parameter:
7305582bec1SHong Zhang .  pc - the preconditioner context
7315582bec1SHong Zhang 
7325582bec1SHong Zhang    Application Interface Routine: PCCreate()
7335582bec1SHong Zhang */
7345582bec1SHong Zhang 
7355582bec1SHong Zhang /*MC
7361e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
7375582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
7386ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
7396ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
7405582bec1SHong Zhang 
7416ca4d86aSHong Zhang    Options Database Key:
7426ca4d86aSHong Zhang    Multigrid options(inherited)
7436ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
7446ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
7456ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
746f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
7476ca4d86aSHong Zhang 
74851f519a2SBarry Smith    ML options:
7496ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
7506ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
7516ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
752f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
7536ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
7546ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
7557ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
7565582bec1SHong Zhang 
7575582bec1SHong Zhang    Level: intermediate
7585582bec1SHong Zhang 
7595582bec1SHong Zhang   Concepts: multigrid
7605582bec1SHong Zhang 
7615582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
76297177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
76397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
76497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
76597177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
7665582bec1SHong Zhang M*/
7675582bec1SHong Zhang 
7685582bec1SHong Zhang EXTERN_C_BEGIN
7695582bec1SHong Zhang #undef __FUNCT__
7705582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
771dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
7725582bec1SHong Zhang {
7735582bec1SHong Zhang   PetscErrorCode  ierr;
7745582bec1SHong Zhang   PC_ML           *pc_ml;
77501da6913SBarry Smith   PC_MG           *mg;
7765582bec1SHong Zhang 
7775582bec1SHong Zhang   PetscFunctionBegin;
778573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
779c9e1c140SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
7805582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
7815582bec1SHong Zhang 
7825582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
78338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
78401da6913SBarry Smith   mg = (PC_MG*)pc->data;
78501da6913SBarry Smith   mg->innerctx = pc_ml;
7865582bec1SHong Zhang 
787573998d7SHong Zhang   pc_ml->ml_object     = 0;
788573998d7SHong Zhang   pc_ml->agg_object    = 0;
789573998d7SHong Zhang   pc_ml->gridctx       = 0;
790573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
791573998d7SHong Zhang   pc_ml->Nlevels       = -1;
792573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
793573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
7943751b4bdSBarry Smith   pc_ml->CoarsenScheme = 1;
795573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
796573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
797573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
798573998d7SHong Zhang   pc_ml->size          = 0;
799573998d7SHong Zhang 
8005582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
8015582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
8025582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
8035582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
8045582bec1SHong Zhang   PetscFunctionReturn(0);
8055582bec1SHong Zhang }
8065582bec1SHong Zhang EXTERN_C_END
807