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 */ 4267d6f150SMatthew G Knepley Vec y, work; 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; 51b5c8bdf8SJed Brown PetscInt Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization; 52b5c8bdf8SJed Brown PetscReal Threshold,DampingFactor,EnergyMinimizationDropTol; 53ace3abfcSBarry Smith PetscBool SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable; 54573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 555582bec1SHong Zhang } PC_ML; 5641ca0015SHong Zhang 576562c4e1SBarry Smith #undef __FUNCT__ 586562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow" 596562c4e1SBarry 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[]) 606562c4e1SBarry Smith { 616562c4e1SBarry Smith PetscErrorCode ierr; 626562c4e1SBarry Smith PetscInt m,i,j,k=0,row,*aj; 636562c4e1SBarry Smith PetscScalar *aa; 646562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 656562c4e1SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)ml->Aloc->data; 665582bec1SHong Zhang 676562c4e1SBarry Smith 686562c4e1SBarry Smith ierr = MatGetSize(ml->Aloc,&m,PETSC_NULL); if (ierr) return(0); 696562c4e1SBarry Smith for (i = 0; i<N_requested_rows; i++) { 706562c4e1SBarry Smith row = requested_rows[i]; 716562c4e1SBarry Smith row_lengths[i] = a->ilen[row]; 726562c4e1SBarry Smith if (allocated_space < k+row_lengths[i]) return(0); 736562c4e1SBarry Smith if ( (row >= 0) || (row <= (m-1)) ) { 746562c4e1SBarry Smith aj = a->j + a->i[row]; 756562c4e1SBarry Smith aa = a->a + a->i[row]; 766562c4e1SBarry Smith for (j=0; j<row_lengths[i]; j++){ 776562c4e1SBarry Smith columns[k] = aj[j]; 786562c4e1SBarry Smith values[k++] = aa[j]; 796562c4e1SBarry Smith } 806562c4e1SBarry Smith } 816562c4e1SBarry Smith } 826562c4e1SBarry Smith return(1); 836562c4e1SBarry Smith } 846562c4e1SBarry Smith 856562c4e1SBarry Smith #undef __FUNCT__ 866562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm" 876562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data) 886562c4e1SBarry Smith { 896562c4e1SBarry Smith PetscErrorCode ierr; 906562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_data; 916562c4e1SBarry Smith Mat A=ml->A; 926562c4e1SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 936562c4e1SBarry Smith PetscMPIInt size; 946562c4e1SBarry Smith PetscInt i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n; 956562c4e1SBarry Smith PetscScalar *array; 966562c4e1SBarry Smith 976562c4e1SBarry Smith PetscFunctionBegin; 986562c4e1SBarry Smith ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 996562c4e1SBarry Smith if (size == 1) return 0; 1006562c4e1SBarry Smith 1016562c4e1SBarry Smith ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 1026562c4e1SBarry Smith ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1036562c4e1SBarry Smith ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1046562c4e1SBarry Smith ierr = VecResetArray(ml->y);CHKERRQ(ierr); 1056562c4e1SBarry Smith ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 1066562c4e1SBarry Smith for (i=in_length; i<out_length; i++){ 1076562c4e1SBarry Smith p[i] = array[i-in_length]; 1086562c4e1SBarry Smith } 1096562c4e1SBarry Smith ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 1106562c4e1SBarry Smith PetscFunctionReturn(0); 1116562c4e1SBarry Smith } 1126562c4e1SBarry Smith 1136562c4e1SBarry Smith #undef __FUNCT__ 1146562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec" 1156562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 1166562c4e1SBarry Smith { 1176562c4e1SBarry Smith PetscErrorCode ierr; 1186562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 1196562c4e1SBarry Smith Mat A=ml->A, Aloc=ml->Aloc; 1206562c4e1SBarry Smith PetscMPIInt size; 1216562c4e1SBarry Smith PetscScalar *pwork=ml->pwork; 1226562c4e1SBarry Smith PetscInt i; 1236562c4e1SBarry Smith 1246562c4e1SBarry Smith PetscFunctionBegin; 1256562c4e1SBarry Smith ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1266562c4e1SBarry Smith if (size == 1){ 1276562c4e1SBarry Smith ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 1286562c4e1SBarry Smith } else { 1296562c4e1SBarry Smith for (i=0; i<in_length; i++) pwork[i] = p[i]; 1306562c4e1SBarry Smith PetscML_comm(pwork,ml); 1316562c4e1SBarry Smith ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 1326562c4e1SBarry Smith } 1336562c4e1SBarry Smith ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 1346562c4e1SBarry Smith ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 1356562c4e1SBarry Smith ierr = VecResetArray(ml->x);CHKERRQ(ierr); 1366562c4e1SBarry Smith ierr = VecResetArray(ml->y);CHKERRQ(ierr); 1376562c4e1SBarry Smith PetscFunctionReturn(0); 1386562c4e1SBarry Smith } 1396562c4e1SBarry Smith 1406562c4e1SBarry Smith #undef __FUNCT__ 1416562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML" 1426562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 1436562c4e1SBarry Smith { 1446562c4e1SBarry Smith PetscErrorCode ierr; 1456562c4e1SBarry Smith Mat_MLShell *shell; 1466562c4e1SBarry Smith PetscScalar *xarray,*yarray; 1476562c4e1SBarry Smith PetscInt x_length,y_length; 1486562c4e1SBarry Smith 1496562c4e1SBarry Smith PetscFunctionBegin; 1506562c4e1SBarry Smith ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 1516562c4e1SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 1526562c4e1SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 1536562c4e1SBarry Smith x_length = shell->mlmat->invec_leng; 1546562c4e1SBarry Smith y_length = shell->mlmat->outvec_leng; 1556562c4e1SBarry Smith ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 1566562c4e1SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 1576562c4e1SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 1586562c4e1SBarry Smith PetscFunctionReturn(0); 1596562c4e1SBarry Smith } 1606562c4e1SBarry Smith 1616562c4e1SBarry Smith #undef __FUNCT__ 1626562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML" 16367d6f150SMatthew G Knepley /* Computes y = w + A * x 16467d6f150SMatthew G Knepley It is possible that w == y, but not x == y 16567d6f150SMatthew G Knepley */ 1666562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 1676562c4e1SBarry Smith { 1686562c4e1SBarry Smith Mat_MLShell *shell; 1696562c4e1SBarry Smith PetscScalar *xarray,*yarray; 1706562c4e1SBarry Smith PetscInt x_length,y_length; 17167d6f150SMatthew G Knepley PetscErrorCode ierr; 1726562c4e1SBarry Smith 1736562c4e1SBarry Smith PetscFunctionBegin; 1746562c4e1SBarry Smith ierr = MatShellGetContext(A, (void **) &shell);CHKERRQ(ierr); 17567d6f150SMatthew G Knepley if (y == w) { 17667d6f150SMatthew G Knepley if (!shell->work) { 17767d6f150SMatthew G Knepley ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr); 17867d6f150SMatthew G Knepley } 17967d6f150SMatthew G Knepley ierr = VecGetArray(x, &xarray);CHKERRQ(ierr); 18067d6f150SMatthew G Knepley ierr = VecGetArray(shell->work, &yarray);CHKERRQ(ierr); 18167d6f150SMatthew G Knepley x_length = shell->mlmat->invec_leng; 18267d6f150SMatthew G Knepley y_length = shell->mlmat->outvec_leng; 18367d6f150SMatthew G Knepley ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray); 18467d6f150SMatthew G Knepley ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr); 18567d6f150SMatthew G Knepley ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr); 1863ba3408dSMatthew G Knepley ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr); 18767d6f150SMatthew G Knepley } else { 1886562c4e1SBarry Smith ierr = VecGetArray(x, &xarray);CHKERRQ(ierr); 1896562c4e1SBarry Smith ierr = VecGetArray(y, &yarray);CHKERRQ(ierr); 1906562c4e1SBarry Smith x_length = shell->mlmat->invec_leng; 1916562c4e1SBarry Smith y_length = shell->mlmat->outvec_leng; 1926562c4e1SBarry Smith ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray); 1936562c4e1SBarry Smith ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr); 1946562c4e1SBarry Smith ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr); 1956562c4e1SBarry Smith ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr); 19667d6f150SMatthew G Knepley } 1976562c4e1SBarry Smith PetscFunctionReturn(0); 1986562c4e1SBarry Smith } 1996562c4e1SBarry Smith 2006562c4e1SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */ 2016562c4e1SBarry Smith #undef __FUNCT__ 2026562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML" 2036562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 2046562c4e1SBarry Smith { 2056562c4e1SBarry Smith PetscErrorCode ierr; 2066562c4e1SBarry Smith Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 2076562c4e1SBarry Smith Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 2086562c4e1SBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 2096562c4e1SBarry Smith PetscScalar *aa=a->a,*ba=b->a,*ca; 2106562c4e1SBarry Smith PetscInt am=A->rmap->n,an=A->cmap->n,i,j,k; 2116562c4e1SBarry Smith PetscInt *ci,*cj,ncols; 2126562c4e1SBarry Smith 2136562c4e1SBarry Smith PetscFunctionBegin; 214e32f2f54SBarry 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); 2156562c4e1SBarry Smith 2166562c4e1SBarry Smith if (scall == MAT_INITIAL_MATRIX){ 2176562c4e1SBarry Smith ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2186562c4e1SBarry Smith ci[0] = 0; 2196562c4e1SBarry Smith for (i=0; i<am; i++){ 2206562c4e1SBarry Smith ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 2216562c4e1SBarry Smith } 2226562c4e1SBarry Smith ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 2236562c4e1SBarry Smith ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 2246562c4e1SBarry Smith 2256562c4e1SBarry Smith k = 0; 2266562c4e1SBarry Smith for (i=0; i<am; i++){ 2276562c4e1SBarry Smith /* diagonal portion of A */ 2286562c4e1SBarry Smith ncols = ai[i+1] - ai[i]; 2296562c4e1SBarry Smith for (j=0; j<ncols; j++) { 2306562c4e1SBarry Smith cj[k] = *aj++; 2316562c4e1SBarry Smith ca[k++] = *aa++; 2326562c4e1SBarry Smith } 2336562c4e1SBarry Smith /* off-diagonal portion of A */ 2346562c4e1SBarry Smith ncols = bi[i+1] - bi[i]; 2356562c4e1SBarry Smith for (j=0; j<ncols; j++) { 2366562c4e1SBarry Smith cj[k] = an + (*bj); bj++; 2376562c4e1SBarry Smith ca[k++] = *ba++; 2386562c4e1SBarry Smith } 2396562c4e1SBarry Smith } 240e32f2f54SBarry Smith if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 2416562c4e1SBarry Smith 2426562c4e1SBarry Smith /* put together the new matrix */ 2436562c4e1SBarry Smith an = mpimat->A->cmap->n+mpimat->B->cmap->n; 2446562c4e1SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 2456562c4e1SBarry Smith 2466562c4e1SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 2476562c4e1SBarry Smith /* Since these are PETSc arrays, change flags to free them as necessary. */ 2486562c4e1SBarry Smith mat = (Mat_SeqAIJ*)(*Aloc)->data; 2496562c4e1SBarry Smith mat->free_a = PETSC_TRUE; 2506562c4e1SBarry Smith mat->free_ij = PETSC_TRUE; 2516562c4e1SBarry Smith 2526562c4e1SBarry Smith mat->nonew = 0; 2536562c4e1SBarry Smith } else if (scall == MAT_REUSE_MATRIX){ 2546562c4e1SBarry Smith mat=(Mat_SeqAIJ*)(*Aloc)->data; 2556562c4e1SBarry Smith ci = mat->i; cj = mat->j; ca = mat->a; 2566562c4e1SBarry Smith for (i=0; i<am; i++) { 2576562c4e1SBarry Smith /* diagonal portion of A */ 2586562c4e1SBarry Smith ncols = ai[i+1] - ai[i]; 2596562c4e1SBarry Smith for (j=0; j<ncols; j++) *ca++ = *aa++; 2606562c4e1SBarry Smith /* off-diagonal portion of A */ 2616562c4e1SBarry Smith ncols = bi[i+1] - bi[i]; 2626562c4e1SBarry Smith for (j=0; j<ncols; j++) *ca++ = *ba++; 2636562c4e1SBarry Smith } 2646562c4e1SBarry Smith } else { 2650005702aSJed Brown SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 2666562c4e1SBarry Smith } 2676562c4e1SBarry Smith PetscFunctionReturn(0); 2686562c4e1SBarry Smith } 2696562c4e1SBarry Smith 2706562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat); 2716562c4e1SBarry Smith #undef __FUNCT__ 2726562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML" 2736562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A) 2746562c4e1SBarry Smith { 2756562c4e1SBarry Smith PetscErrorCode ierr; 2766562c4e1SBarry Smith Mat_MLShell *shell; 2776562c4e1SBarry Smith 2786562c4e1SBarry Smith PetscFunctionBegin; 2796562c4e1SBarry Smith ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 2806562c4e1SBarry Smith ierr = VecDestroy(shell->y);CHKERRQ(ierr); 28167d6f150SMatthew G Knepley if (shell->work) {ierr = VecDestroy(shell->work);CHKERRQ(ierr);} 2826562c4e1SBarry Smith ierr = PetscFree(shell);CHKERRQ(ierr); 2836562c4e1SBarry Smith ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 2846562c4e1SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 2856562c4e1SBarry Smith PetscFunctionReturn(0); 2866562c4e1SBarry Smith } 2876562c4e1SBarry Smith 2886562c4e1SBarry Smith #undef __FUNCT__ 2896562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ" 2906562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 2916562c4e1SBarry Smith { 2926562c4e1SBarry Smith struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 2936562c4e1SBarry Smith PetscErrorCode ierr; 2946562c4e1SBarry Smith PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 2956562c4e1SBarry Smith PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k; 2966562c4e1SBarry Smith PetscScalar *ml_vals=matdata->values,*aa; 2976562c4e1SBarry Smith 2986562c4e1SBarry Smith PetscFunctionBegin; 299e7e72b3dSBarry Smith if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 3006562c4e1SBarry Smith if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 3016562c4e1SBarry Smith if (reuse){ 3026562c4e1SBarry Smith Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 3036562c4e1SBarry Smith aij->i = ml_rowptr; 3046562c4e1SBarry Smith aij->j = ml_cols; 3056562c4e1SBarry Smith aij->a = ml_vals; 3066562c4e1SBarry Smith } else { 3076562c4e1SBarry Smith /* sort ml_cols and ml_vals */ 3086562c4e1SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 3096562c4e1SBarry Smith for (i=0; i<m; i++) { 3106562c4e1SBarry Smith nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 3116562c4e1SBarry Smith } 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 3236562c4e1SBarry Smith /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 3246562c4e1SBarry Smith ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 3256562c4e1SBarry Smith ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 3266562c4e1SBarry Smith ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 3276562c4e1SBarry Smith 3286562c4e1SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 3296562c4e1SBarry Smith nz_max = 1; 3306562c4e1SBarry Smith for (i=0; i<m; i++) { 3316562c4e1SBarry Smith nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 3326562c4e1SBarry Smith if (nnz[i] > nz_max) nz_max += nnz[i]; 3336562c4e1SBarry Smith } 3346562c4e1SBarry Smith 3356562c4e1SBarry Smith ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 3366562c4e1SBarry Smith ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 3376562c4e1SBarry Smith for (i=0; i<m; i++){ 3386562c4e1SBarry Smith k = 0; 3396562c4e1SBarry Smith /* diagonal entry */ 3406562c4e1SBarry Smith aj[k] = i; aa[k++] = ml_vals[i]; 3416562c4e1SBarry Smith /* off diagonal entries */ 3426562c4e1SBarry Smith for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 3436562c4e1SBarry Smith aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 3446562c4e1SBarry Smith } 3456562c4e1SBarry Smith /* sort aj and aa */ 3466562c4e1SBarry Smith ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 3476562c4e1SBarry Smith ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 3486562c4e1SBarry Smith } 3496562c4e1SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3506562c4e1SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3516562c4e1SBarry Smith 3526562c4e1SBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 3536562c4e1SBarry Smith ierr = PetscFree(nnz);CHKERRQ(ierr); 3546562c4e1SBarry Smith PetscFunctionReturn(0); 3556562c4e1SBarry Smith } 3566562c4e1SBarry Smith 3576562c4e1SBarry Smith #undef __FUNCT__ 3586562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL" 3596562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 3606562c4e1SBarry Smith { 3616562c4e1SBarry Smith PetscErrorCode ierr; 3626562c4e1SBarry Smith PetscInt m,n; 3636562c4e1SBarry Smith ML_Comm *MLcomm; 3646562c4e1SBarry Smith Mat_MLShell *shellctx; 3656562c4e1SBarry Smith 3666562c4e1SBarry Smith PetscFunctionBegin; 3676562c4e1SBarry Smith m = mlmat->outvec_leng; 3686562c4e1SBarry Smith n = mlmat->invec_leng; 3696562c4e1SBarry Smith if (!m || !n){ 3706562c4e1SBarry Smith newmat = PETSC_NULL; 3716562c4e1SBarry Smith PetscFunctionReturn(0); 3726562c4e1SBarry Smith } 3736562c4e1SBarry Smith 3746562c4e1SBarry Smith if (reuse){ 3756562c4e1SBarry Smith ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 3766562c4e1SBarry Smith shellctx->mlmat = mlmat; 3776562c4e1SBarry Smith PetscFunctionReturn(0); 3786562c4e1SBarry Smith } 3796562c4e1SBarry Smith 3806562c4e1SBarry Smith MLcomm = mlmat->comm; 3816562c4e1SBarry Smith ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 3826562c4e1SBarry Smith ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 3836562c4e1SBarry Smith ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 3846562c4e1SBarry Smith ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 3856562c4e1SBarry Smith shellctx->A = *newmat; 3866562c4e1SBarry Smith shellctx->mlmat = mlmat; 38767d6f150SMatthew G Knepley shellctx->work = PETSC_NULL; 3886562c4e1SBarry Smith ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 3896562c4e1SBarry Smith ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 3906562c4e1SBarry Smith ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 3916562c4e1SBarry Smith (*newmat)->ops->destroy = MatDestroy_ML; 3926562c4e1SBarry Smith PetscFunctionReturn(0); 3936562c4e1SBarry Smith } 3946562c4e1SBarry Smith 3956562c4e1SBarry Smith #undef __FUNCT__ 3966562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ" 3976562c4e1SBarry Smith static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 3986562c4e1SBarry Smith { 3996562c4e1SBarry Smith struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 4006562c4e1SBarry Smith PetscInt *ml_cols=matdata->columns,*aj; 4016562c4e1SBarry Smith PetscScalar *ml_vals=matdata->values,*aa; 4026562c4e1SBarry Smith PetscErrorCode ierr; 4036562c4e1SBarry Smith PetscInt i,j,k,*gordering; 4046562c4e1SBarry Smith PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 4056562c4e1SBarry Smith Mat A; 4066562c4e1SBarry Smith 4076562c4e1SBarry Smith PetscFunctionBegin; 408e7e72b3dSBarry Smith if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 4096562c4e1SBarry Smith n = mlmat->invec_leng; 410e32f2f54SBarry Smith if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 4116562c4e1SBarry Smith 4126562c4e1SBarry Smith ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 4136562c4e1SBarry Smith ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 4146562c4e1SBarry Smith ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 4156562c4e1SBarry Smith ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 4166562c4e1SBarry Smith 4176562c4e1SBarry Smith nz_max = 0; 4186562c4e1SBarry Smith for (i=0; i<m; i++){ 4196562c4e1SBarry Smith nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 4206562c4e1SBarry Smith if (nz_max < nnz[i]) nz_max = nnz[i]; 4216562c4e1SBarry Smith nnzA[i] = 1; /* diag */ 4226562c4e1SBarry Smith for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 4236562c4e1SBarry Smith if (ml_cols[j] < m) nnzA[i]++; 4246562c4e1SBarry Smith } 4256562c4e1SBarry Smith nnzB[i] = nnz[i] - nnzA[i]; 4266562c4e1SBarry Smith } 4276562c4e1SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 4286562c4e1SBarry Smith 4296562c4e1SBarry Smith /* insert mat values -- remap row and column indices */ 4306562c4e1SBarry Smith nz_max++; 4316562c4e1SBarry Smith ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 4326562c4e1SBarry Smith /* create global row numbering for a ML_Operator */ 4336562c4e1SBarry Smith ML_build_global_numbering(mlmat,&gordering,"rows"); 4346562c4e1SBarry Smith for (i=0; i<m; i++){ 4356562c4e1SBarry Smith row = gordering[i]; 4366562c4e1SBarry Smith k = 0; 4376562c4e1SBarry Smith /* diagonal entry */ 4386562c4e1SBarry Smith aj[k] = row; aa[k++] = ml_vals[i]; 4396562c4e1SBarry Smith /* off diagonal entries */ 4406562c4e1SBarry Smith for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 4416562c4e1SBarry Smith aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 4426562c4e1SBarry Smith } 4436562c4e1SBarry Smith ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 4446562c4e1SBarry Smith } 445*eb4736aaSBarry Smith ML_free(gordering); 4466562c4e1SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4476562c4e1SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4486562c4e1SBarry Smith *newmat = A; 4496562c4e1SBarry Smith 4506562c4e1SBarry Smith ierr = PetscFree3(nnzA,nnzB,nnz); 4516562c4e1SBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 4526562c4e1SBarry Smith PetscFunctionReturn(0); 4536562c4e1SBarry Smith } 4546562c4e1SBarry Smith 4556562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/ 45601da6913SBarry Smith #undef __FUNCT__ 4573751b4bdSBarry Smith #define __FUNCT__ "PCDestroy_ML_Private" 4583751b4bdSBarry Smith PetscErrorCode PCDestroy_ML_Private(void *ptr) 45901da6913SBarry Smith { 46001da6913SBarry Smith PetscErrorCode ierr; 46101da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)ptr; 46201da6913SBarry Smith PetscInt level,fine_level=pc_ml->Nlevels-1; 46301da6913SBarry Smith 46401da6913SBarry Smith PetscFunctionBegin; 46501da6913SBarry Smith ML_Aggregate_Destroy(&pc_ml->agg_object); 46601da6913SBarry Smith ML_Destroy(&pc_ml->ml_object); 46701da6913SBarry Smith 46801da6913SBarry Smith if (pc_ml->PetscMLdata) { 46901da6913SBarry Smith ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 47001da6913SBarry Smith if (pc_ml->size > 1) {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 47101da6913SBarry Smith if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);} 47201da6913SBarry Smith if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);} 47301da6913SBarry Smith } 47401da6913SBarry Smith ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 47501da6913SBarry Smith 47601da6913SBarry Smith for (level=0; level<fine_level; level++){ 47701da6913SBarry Smith if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);} 47801da6913SBarry Smith if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);} 47901da6913SBarry Smith if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);} 48001da6913SBarry Smith if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);} 48101da6913SBarry Smith if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);} 48201da6913SBarry Smith if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 48301da6913SBarry Smith } 48401da6913SBarry Smith ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 48501da6913SBarry Smith PetscFunctionReturn(0); 48601da6913SBarry Smith } 4875582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4885582bec1SHong Zhang /* 4895582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 4905582bec1SHong Zhang by setting data structures and options. 4915582bec1SHong Zhang 4925582bec1SHong Zhang Input Parameter: 4935582bec1SHong Zhang . pc - the preconditioner context 4945582bec1SHong Zhang 4955582bec1SHong Zhang Application Interface Routine: PCSetUp() 4965582bec1SHong Zhang 4975582bec1SHong Zhang Notes: 4985582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 4995582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 5005582bec1SHong Zhang */ 5016ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC); 502c07bf074SBarry Smith extern PetscErrorCode PCDestroy_MG_Private(PC); 503c07bf074SBarry Smith 5045582bec1SHong Zhang #undef __FUNCT__ 5055582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 5066ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 5075582bec1SHong Zhang { 5085582bec1SHong Zhang PetscErrorCode ierr; 509eef31507SHong Zhang PetscMPIInt size; 5105582bec1SHong Zhang FineGridCtx *PetscMLdata; 5115582bec1SHong Zhang ML *ml_object; 5125582bec1SHong Zhang ML_Aggregate *agg_object; 5135582bec1SHong Zhang ML_Operator *mlmat; 5144f8eab3cSJed Brown PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs; 5155582bec1SHong Zhang Mat A,Aloc; 5165582bec1SHong Zhang GridCtx *gridctx; 51701da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 51801da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 519ace3abfcSBarry Smith PetscBool isSeq, isMPI; 520c07bf074SBarry Smith KSP smoother; 521c07bf074SBarry Smith PC subpc; 5225582bec1SHong Zhang 5235582bec1SHong Zhang PetscFunctionBegin; 524573998d7SHong Zhang if (pc->setupcalled){ 525c07bf074SBarry Smith /* since ML can change the size of vectors/matrices at any level we must destroy everything */ 5263751b4bdSBarry Smith ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr); 527c07bf074SBarry Smith ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr); 528573998d7SHong Zhang } 529573998d7SHong Zhang 5305582bec1SHong Zhang /* setup special features of PCML */ 5315582bec1SHong Zhang /*--------------------------------*/ 5325582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 5335582bec1SHong Zhang A = pc->pmat; 5347adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5355582bec1SHong Zhang pc_ml->size = size; 536864b637dSMatthew Knepley ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 537864b637dSMatthew Knepley ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 538864b637dSMatthew Knepley if (isMPI){ 539db571536SBarry Smith ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr); 540864b637dSMatthew Knepley } else if (isSeq) { 5415582bec1SHong Zhang Aloc = A; 542e7e72b3dSBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices."); 5435582bec1SHong Zhang 5445582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 54538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 5465582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 547d0f46423SBarry Smith ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 5485582bec1SHong Zhang 54924a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 550d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr); 55124a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 55224a42b14SHong Zhang 55324a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 554d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 55524a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 556573998d7SHong Zhang PetscMLdata->A = A; 557573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 55824a42b14SHong Zhang 5595582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 56045cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 5615582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 5624f8eab3cSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 5635582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 564ead7dcbeSHong Zhang ML_Comm_Set_UsrComm(ml_object->comm,((PetscObject)A)->comm); 565573998d7SHong Zhang pc_ml->ml_object = ml_object; 5665582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 5675582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 5685582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 5695582bec1SHong Zhang 570b5c8bdf8SJed Brown ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO); 571b5c8bdf8SJed Brown 5725582bec1SHong Zhang /* aggregation */ 5735582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 574573998d7SHong Zhang pc_ml->agg_object = agg_object; 575573998d7SHong Zhang 5764f8eab3cSJed Brown ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr); 5775582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 5785582bec1SHong Zhang /* set options */ 5795582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 5805582bec1SHong Zhang case 1: 5815582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 5825582bec1SHong Zhang case 2: 5835582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 5845582bec1SHong Zhang case 3: 5855582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 5865582bec1SHong Zhang } 5875582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 5885582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 5895582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 5907ffd031bSHong Zhang ML_Set_SpectralNormScheme_Anorm(ml_object); 5915582bec1SHong Zhang } 592b5c8bdf8SJed Brown agg_object->keep_agg_information = (int)pc_ml->KeepAggInfo; 593b5c8bdf8SJed Brown agg_object->keep_P_tentative = (int)pc_ml->Reusable; 594b5c8bdf8SJed Brown agg_object->block_scaled_SA = (int)pc_ml->BlockScaling; 595b5c8bdf8SJed Brown agg_object->minimizing_energy = (int)pc_ml->EnergyMinimization; 596b5c8bdf8SJed Brown agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol; 597b5c8bdf8SJed Brown agg_object->cheap_minimizing_energy = (int)pc_ml->EnergyMinimizationCheap; 5985582bec1SHong Zhang 599b5c8bdf8SJed Brown if (pc_ml->OldHierarchy) { 6005582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 601b5c8bdf8SJed Brown } else { 602b5c8bdf8SJed Brown Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 603b5c8bdf8SJed Brown } 60465e19b50SBarry Smith if (Nlevels<=0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 605573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 606aa85bbbfSHong Zhang fine_level = Nlevels - 1; 607c07bf074SBarry Smith 60897177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 609aa85bbbfSHong Zhang /* set default smoothers */ 610aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++){ 611aa85bbbfSHong Zhang if (size == 1){ 612aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 613aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 614aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 615aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 616aa85bbbfSHong Zhang } else { 617aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 618aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 619aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 620aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 621aa85bbbfSHong Zhang } 622aa85bbbfSHong Zhang } 62397177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 6245582bec1SHong Zhang 6255582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 6265582bec1SHong Zhang pc_ml->gridctx = gridctx; 6275582bec1SHong Zhang 6285582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 6295582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 630e14861a4SHong Zhang gridctx[fine_level].A = A; 631573998d7SHong Zhang 632e14861a4SHong Zhang level = fine_level - 1; 633ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 6345582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 635e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 636db571536SBarry Smith ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr); 637e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 638db571536SBarry Smith ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr); 639573998d7SHong Zhang 640573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 641573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 6425582bec1SHong Zhang level--; 6435582bec1SHong Zhang } 644ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 6455582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 6465582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 647db571536SBarry Smith ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr); 648ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 649db571536SBarry Smith ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr); 650573998d7SHong Zhang 6515582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 652eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 6535582bec1SHong Zhang level--; 6545582bec1SHong Zhang } 6555582bec1SHong Zhang } 6565582bec1SHong Zhang 657573998d7SHong Zhang /* create vectors and ksp at all levels */ 658ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 659573998d7SHong Zhang level1 = level + 1; 660e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr); 661d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr); 6625582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 66397177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 6645582bec1SHong Zhang 665e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr); 666d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 6675582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 66897177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 669ac346b81SHong Zhang 670e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr); 671d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 672ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 67397177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 674ac346b81SHong Zhang 6755582bec1SHong Zhang if (level == 0){ 67697177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 6775582bec1SHong Zhang } else { 67897177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 679573998d7SHong Zhang } 680573998d7SHong Zhang } 681573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 682573998d7SHong Zhang 683573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 684573998d7SHong Zhang for (level=0; level<fine_level; level++){ 685573998d7SHong Zhang level1 = level + 1; 686aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr); 687573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 688573998d7SHong Zhang if (level > 0){ 68997177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 6905582bec1SHong Zhang } 6915582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 6925582bec1SHong Zhang } 69397177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 694ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 6955582bec1SHong Zhang 696c07bf074SBarry Smith /* setupcalled is set to 0 so that MG is setup from scratch */ 697c07bf074SBarry Smith pc->setupcalled = 0; 6983751b4bdSBarry Smith ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 6995582bec1SHong Zhang PetscFunctionReturn(0); 7005582bec1SHong Zhang } 7015582bec1SHong Zhang 7025582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 7035582bec1SHong Zhang /* 7045582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 7055582bec1SHong Zhang that was created with PCCreate_ML(). 7065582bec1SHong Zhang 7075582bec1SHong Zhang Input Parameter: 7085582bec1SHong Zhang . pc - the preconditioner context 7095582bec1SHong Zhang 7105582bec1SHong Zhang Application Interface Routine: PCDestroy() 7115582bec1SHong Zhang */ 7125582bec1SHong Zhang #undef __FUNCT__ 7135582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 7146ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 7155582bec1SHong Zhang { 7165582bec1SHong Zhang PetscErrorCode ierr; 71701da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 71801da6913SBarry Smith PC_ML *pc_ml= (PC_ML*)mg->innerctx; 7195582bec1SHong Zhang 7205582bec1SHong Zhang PetscFunctionBegin; 7213751b4bdSBarry Smith ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr); 72201da6913SBarry Smith ierr = PetscFree(pc_ml);CHKERRQ(ierr); 72301da6913SBarry Smith ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 7245582bec1SHong Zhang PetscFunctionReturn(0); 7255582bec1SHong Zhang } 7265582bec1SHong Zhang 7275582bec1SHong Zhang #undef __FUNCT__ 7285582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 7296ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 7305582bec1SHong Zhang { 7315582bec1SHong Zhang PetscErrorCode ierr; 7323751b4bdSBarry Smith PetscInt indx,PrintLevel; 7335582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 73401da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 73501da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 736b5c8bdf8SJed Brown PetscMPIInt size; 73788ff4cc7SJed Brown MPI_Comm comm = ((PetscObject)pc)->comm; 7385582bec1SHong Zhang 7395582bec1SHong Zhang PetscFunctionBegin; 74088ff4cc7SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 7415582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 7425582bec1SHong Zhang PrintLevel = 0; 7435582bec1SHong Zhang indx = 0; 7445582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 7455582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 746573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 747573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 7483751b4bdSBarry Smith ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 7495582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 750573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 751573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 752acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);CHKERRQ(ierr); 753acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,PETSC_NULL);CHKERRQ(ierr); 754acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,PETSC_NULL);CHKERRQ(ierr); 755b5c8bdf8SJed Brown ierr = PetscOptionsInt("-pc_ml_EnergyMinimization","Energy minimization norm type (0=no minimization; see ML manual for 1,2,3; -1 and 4 undocumented)","None",pc_ml->EnergyMinimization,&pc_ml->EnergyMinimization,PETSC_NULL);CHKERRQ(ierr); 756b5c8bdf8SJed Brown /* 757b5c8bdf8SJed Brown The following checks a number of conditions. If we let this stuff slip by, then ML's error handling will take over. 758b5c8bdf8SJed Brown This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions. 759b5c8bdf8SJed Brown 760b5c8bdf8SJed Brown We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working 761b5c8bdf8SJed Brown combination of options and ML's exit(1) explanations don't help matters. 762b5c8bdf8SJed Brown */ 76388ff4cc7SJed Brown if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4"); 76488ff4cc7SJed Brown if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel"); 765b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2");CHKERRQ(ierr);} 766b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) { 767b5c8bdf8SJed Brown ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,PETSC_NULL);CHKERRQ(ierr); 768b5c8bdf8SJed Brown } 769b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == 2) { 770b5c8bdf8SJed Brown /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */ 771acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,PETSC_NULL);CHKERRQ(ierr); 772b5c8bdf8SJed Brown } 773b5c8bdf8SJed Brown /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */ 774b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE; 775acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_KeepAggInfo","Allows the preconditioner to be reused, or auxilliary matrices to be generated","None",pc_ml->KeepAggInfo,&pc_ml->KeepAggInfo,PETSC_NULL);CHKERRQ(ierr); 776b5c8bdf8SJed Brown /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */ 777b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE; 778acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_Reusable","Store intermedaiate data structures so that the multilevel hierarchy is reusable","None",pc_ml->Reusable,&pc_ml->Reusable,PETSC_NULL);CHKERRQ(ierr); 779b5c8bdf8SJed Brown /* 780b5c8bdf8SJed Brown ML's C API is severely underdocumented and lacks significant functionality. The C++ API calls 781b5c8bdf8SJed Brown ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function 782b5c8bdf8SJed Brown ML_Gen_MGHierarchy_UsingAggregation(). This modification, however, does not provide a strict superset of the 783b5c8bdf8SJed Brown functionality in the old function, so some users may still want to use it. Note that many options are ignored in 784b5c8bdf8SJed Brown this context, but ML doesn't provide a way to find out which ones. 785b5c8bdf8SJed Brown */ 786acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,PETSC_NULL);CHKERRQ(ierr); 7875582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 7885582bec1SHong Zhang PetscFunctionReturn(0); 7895582bec1SHong Zhang } 7905582bec1SHong Zhang 7915582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 7925582bec1SHong Zhang /* 7935582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 7945582bec1SHong Zhang and sets this as the private data within the generic preconditioning 7955582bec1SHong Zhang context, PC, that was created within PCCreate(). 7965582bec1SHong Zhang 7975582bec1SHong Zhang Input Parameter: 7985582bec1SHong Zhang . pc - the preconditioner context 7995582bec1SHong Zhang 8005582bec1SHong Zhang Application Interface Routine: PCCreate() 8015582bec1SHong Zhang */ 8025582bec1SHong Zhang 8035582bec1SHong Zhang /*MC 8041e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 8055582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 8066ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 8076ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 8085582bec1SHong Zhang 8096ca4d86aSHong Zhang Options Database Key: 8106ca4d86aSHong Zhang Multigrid options(inherited) 8116ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 8126ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 8136ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 814f41ab451SVictor Eijkhout - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 8156ca4d86aSHong Zhang 81651f519a2SBarry Smith ML options: 8176ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 8186ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 8196ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 820f41ab451SVictor Eijkhout . -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS 8216ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 8226ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 8237ffd031bSHong Zhang - -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm) 8245582bec1SHong Zhang 8255582bec1SHong Zhang Level: intermediate 8265582bec1SHong Zhang 8275582bec1SHong Zhang Concepts: multigrid 8285582bec1SHong Zhang 8295582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 83097177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 83197177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 83297177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 83397177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 8345582bec1SHong Zhang M*/ 8355582bec1SHong Zhang 8365582bec1SHong Zhang EXTERN_C_BEGIN 8375582bec1SHong Zhang #undef __FUNCT__ 8385582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 8397087cfbeSBarry Smith PetscErrorCode PCCreate_ML(PC pc) 8405582bec1SHong Zhang { 8415582bec1SHong Zhang PetscErrorCode ierr; 8425582bec1SHong Zhang PC_ML *pc_ml; 84301da6913SBarry Smith PC_MG *mg; 8445582bec1SHong Zhang 8455582bec1SHong Zhang PetscFunctionBegin; 846573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 847c9e1c140SHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr); 8485582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 8495582bec1SHong Zhang 8505582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 85138f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr); 85201da6913SBarry Smith mg = (PC_MG*)pc->data; 85301da6913SBarry Smith mg->innerctx = pc_ml; 8545582bec1SHong Zhang 855573998d7SHong Zhang pc_ml->ml_object = 0; 856573998d7SHong Zhang pc_ml->agg_object = 0; 857573998d7SHong Zhang pc_ml->gridctx = 0; 858573998d7SHong Zhang pc_ml->PetscMLdata = 0; 859573998d7SHong Zhang pc_ml->Nlevels = -1; 860573998d7SHong Zhang pc_ml->MaxNlevels = 10; 861573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 8623751b4bdSBarry Smith pc_ml->CoarsenScheme = 1; 863573998d7SHong Zhang pc_ml->Threshold = 0.0; 864573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 865573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 866573998d7SHong Zhang pc_ml->size = 0; 867573998d7SHong Zhang 8685582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 8695582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 8705582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 8715582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 8725582bec1SHong Zhang PetscFunctionReturn(0); 8735582bec1SHong Zhang } 8745582bec1SHong Zhang EXTERN_C_END 875