1ab718edeSHong Zhang 25582bec1SHong Zhang /* 32dccc152SHong Zhang Provides an interface to the ML smoothed Aggregation 47ffd031bSHong Zhang Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs 57ffd031bSHong Zhang Jed Brown, see [PETSC #18321, #18449]. 65582bec1SHong Zhang */ 7b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 8c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 9c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 10c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 11cb5d8e9eSHong Zhang 125582bec1SHong Zhang #include <math.h> 132cf39c26SSatish Balay EXTERN_C_BEGIN 1468210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */ 1568210224SSatish Balay #if !defined(HAVE_CONFIG_H) 1668210224SSatish Balay #define HAVE_CONFIG_H 1768210224SSatish Balay #endif 18c6db04a5SJed Brown #include <ml_include.h> 1939381ba2SJed Brown #include <ml_viz_stats.h> 205582bec1SHong Zhang EXTERN_C_END 215582bec1SHong Zhang 22fb6a8e6dSJed Brown typedef enum {PCML_NULLSPACE_AUTO,PCML_NULLSPACE_USER,PCML_NULLSPACE_BLOCK,PCML_NULLSPACE_SCALAR} PCMLNullSpaceType; 23fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO","USER","BLOCK","SCALAR","PCMLNullSpaceType","PCML_NULLSPACE_",0}; 24fb6a8e6dSJed Brown 255582bec1SHong Zhang /* The context (data structure) at each grid level */ 265582bec1SHong Zhang typedef struct { 275582bec1SHong Zhang Vec x,b,r; /* global vectors */ 285582bec1SHong Zhang Mat A,P,R; 295582bec1SHong Zhang KSP ksp; 3039381ba2SJed Brown Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */ 315582bec1SHong Zhang } GridCtx; 325582bec1SHong Zhang 335582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 345582bec1SHong Zhang typedef struct { 35573998d7SHong Zhang Mat A; /* Petsc matrix in aij format */ 36573998d7SHong Zhang Mat Aloc; /* local portion of A to be used by ML */ 3724a42b14SHong Zhang Vec x,y; 385582bec1SHong Zhang ML_Operator *mlmat; 395582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 405582bec1SHong Zhang } FineGridCtx; 415582bec1SHong Zhang 425582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 435582bec1SHong Zhang typedef struct { 445582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 455582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 4667d6f150SMatthew G Knepley Vec y, work; 475582bec1SHong Zhang } Mat_MLShell; 485582bec1SHong Zhang 495582bec1SHong Zhang /* Private context for the ML preconditioner */ 505582bec1SHong Zhang typedef struct { 515582bec1SHong Zhang ML *ml_object; 525582bec1SHong Zhang ML_Aggregate *agg_object; 535582bec1SHong Zhang GridCtx *gridctx; 545582bec1SHong Zhang FineGridCtx *PetscMLdata; 5539381ba2SJed Brown PetscInt Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization,MinPerProc,PutOnSingleProc,RepartitionType,ZoltanScheme; 5639381ba2SJed Brown PetscReal Threshold,DampingFactor,EnergyMinimizationDropTol,MaxMinRatio,AuxThreshold; 5739381ba2SJed Brown PetscBool SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable,Repartition,Aux; 5848268eb4SJed Brown PetscBool reuse_interpolation; 59fb6a8e6dSJed Brown PCMLNullSpaceType nulltype; 60573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 6139381ba2SJed Brown PetscInt dim; /* data from PCSetCoordinates(_ML) */ 6239381ba2SJed Brown PetscInt nloc; 6339381ba2SJed Brown PetscReal *coords; /* ML has a grid object for each level: the finest grid will point into coords */ 645582bec1SHong Zhang } PC_ML; 6541ca0015SHong Zhang 666562c4e1SBarry Smith #undef __FUNCT__ 676562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow" 686562c4e1SBarry Smith static int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],int allocated_space, int columns[], double values[], int row_lengths[]) 696562c4e1SBarry Smith { 706562c4e1SBarry Smith PetscErrorCode ierr; 716562c4e1SBarry Smith PetscInt m,i,j,k=0,row,*aj; 726562c4e1SBarry Smith PetscScalar *aa; 736562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 746562c4e1SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)ml->Aloc->data; 755582bec1SHong Zhang 766562c4e1SBarry Smith ierr = MatGetSize(ml->Aloc,&m,PETSC_NULL); if (ierr) return(0); 776562c4e1SBarry Smith for (i = 0; i<N_requested_rows; i++) { 786562c4e1SBarry Smith row = requested_rows[i]; 796562c4e1SBarry Smith row_lengths[i] = a->ilen[row]; 806562c4e1SBarry Smith if (allocated_space < k+row_lengths[i]) return(0); 816562c4e1SBarry Smith if ((row >= 0) || (row <= (m-1))) { 826562c4e1SBarry Smith aj = a->j + a->i[row]; 836562c4e1SBarry Smith aa = a->a + a->i[row]; 846562c4e1SBarry Smith for (j=0; j<row_lengths[i]; j++) { 856562c4e1SBarry Smith columns[k] = aj[j]; 866562c4e1SBarry Smith values[k++] = aa[j]; 876562c4e1SBarry Smith } 886562c4e1SBarry Smith } 896562c4e1SBarry Smith } 906562c4e1SBarry Smith return(1); 916562c4e1SBarry Smith } 926562c4e1SBarry Smith 936562c4e1SBarry Smith #undef __FUNCT__ 946562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm" 956562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data) 966562c4e1SBarry Smith { 976562c4e1SBarry Smith PetscErrorCode ierr; 986562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_data; 996562c4e1SBarry Smith Mat A=ml->A; 1006562c4e1SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 1016562c4e1SBarry Smith PetscMPIInt size; 1026562c4e1SBarry Smith PetscInt i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n; 1036562c4e1SBarry Smith PetscScalar *array; 1046562c4e1SBarry Smith 1056562c4e1SBarry Smith PetscFunctionBegin; 1066562c4e1SBarry Smith ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1076562c4e1SBarry Smith if (size == 1) return 0; 1086562c4e1SBarry Smith 1096562c4e1SBarry Smith ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 1106562c4e1SBarry Smith ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1116562c4e1SBarry Smith ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1126562c4e1SBarry Smith ierr = VecResetArray(ml->y);CHKERRQ(ierr); 1136562c4e1SBarry Smith ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 1146562c4e1SBarry Smith for (i=in_length; i<out_length; i++) { 1156562c4e1SBarry Smith p[i] = array[i-in_length]; 1166562c4e1SBarry Smith } 1176562c4e1SBarry Smith ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 1186562c4e1SBarry Smith PetscFunctionReturn(0); 1196562c4e1SBarry Smith } 1206562c4e1SBarry Smith 1216562c4e1SBarry Smith #undef __FUNCT__ 1226562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec" 1236562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 1246562c4e1SBarry Smith { 1256562c4e1SBarry Smith PetscErrorCode ierr; 1266562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 1276562c4e1SBarry Smith Mat A=ml->A, Aloc=ml->Aloc; 1286562c4e1SBarry Smith PetscMPIInt size; 1296562c4e1SBarry Smith PetscScalar *pwork=ml->pwork; 1306562c4e1SBarry Smith PetscInt i; 1316562c4e1SBarry Smith 1326562c4e1SBarry Smith PetscFunctionBegin; 1336562c4e1SBarry Smith ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1346562c4e1SBarry Smith if (size == 1) { 1356562c4e1SBarry Smith ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 1366562c4e1SBarry Smith } else { 1376562c4e1SBarry Smith for (i=0; i<in_length; i++) pwork[i] = p[i]; 1386562c4e1SBarry Smith PetscML_comm(pwork,ml); 1396562c4e1SBarry Smith ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 1406562c4e1SBarry Smith } 1416562c4e1SBarry Smith ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 1426562c4e1SBarry Smith ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 1436562c4e1SBarry Smith ierr = VecResetArray(ml->x);CHKERRQ(ierr); 1446562c4e1SBarry Smith ierr = VecResetArray(ml->y);CHKERRQ(ierr); 1456562c4e1SBarry Smith PetscFunctionReturn(0); 1466562c4e1SBarry Smith } 1476562c4e1SBarry Smith 1486562c4e1SBarry Smith #undef __FUNCT__ 1496562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML" 1506562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 1516562c4e1SBarry Smith { 1526562c4e1SBarry Smith PetscErrorCode ierr; 1536562c4e1SBarry Smith Mat_MLShell *shell; 1546562c4e1SBarry Smith PetscScalar *xarray,*yarray; 1556562c4e1SBarry Smith PetscInt x_length,y_length; 1566562c4e1SBarry Smith 1576562c4e1SBarry Smith PetscFunctionBegin; 1586562c4e1SBarry Smith ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 1596562c4e1SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 1606562c4e1SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 1616562c4e1SBarry Smith x_length = shell->mlmat->invec_leng; 1626562c4e1SBarry Smith y_length = shell->mlmat->outvec_leng; 1636562c4e1SBarry Smith ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 1646562c4e1SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 1656562c4e1SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 1666562c4e1SBarry Smith PetscFunctionReturn(0); 1676562c4e1SBarry Smith } 1686562c4e1SBarry Smith 1696562c4e1SBarry Smith #undef __FUNCT__ 1706562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML" 17167d6f150SMatthew G Knepley /* Computes y = w + A * x 17267d6f150SMatthew G Knepley It is possible that w == y, but not x == y 17367d6f150SMatthew G Knepley */ 1746562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 1756562c4e1SBarry Smith { 1766562c4e1SBarry Smith Mat_MLShell *shell; 1776562c4e1SBarry Smith PetscScalar *xarray,*yarray; 1786562c4e1SBarry Smith PetscInt x_length,y_length; 17967d6f150SMatthew G Knepley PetscErrorCode ierr; 1806562c4e1SBarry Smith 1816562c4e1SBarry Smith PetscFunctionBegin; 1826562c4e1SBarry Smith ierr = MatShellGetContext(A, (void **) &shell);CHKERRQ(ierr); 18367d6f150SMatthew G Knepley if (y == w) { 18467d6f150SMatthew G Knepley if (!shell->work) { 18567d6f150SMatthew G Knepley ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr); 18667d6f150SMatthew G Knepley } 18767d6f150SMatthew G Knepley ierr = VecGetArray(x, &xarray);CHKERRQ(ierr); 18867d6f150SMatthew G Knepley ierr = VecGetArray(shell->work, &yarray);CHKERRQ(ierr); 18967d6f150SMatthew G Knepley x_length = shell->mlmat->invec_leng; 19067d6f150SMatthew G Knepley y_length = shell->mlmat->outvec_leng; 19167d6f150SMatthew G Knepley ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray); 19267d6f150SMatthew G Knepley ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr); 19367d6f150SMatthew G Knepley ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr); 1943ba3408dSMatthew G Knepley ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr); 19567d6f150SMatthew G Knepley } else { 1966562c4e1SBarry Smith ierr = VecGetArray(x, &xarray);CHKERRQ(ierr); 1976562c4e1SBarry Smith ierr = VecGetArray(y, &yarray);CHKERRQ(ierr); 1986562c4e1SBarry Smith x_length = shell->mlmat->invec_leng; 1996562c4e1SBarry Smith y_length = shell->mlmat->outvec_leng; 2006562c4e1SBarry Smith ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray); 2016562c4e1SBarry Smith ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr); 2026562c4e1SBarry Smith ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr); 2036562c4e1SBarry Smith ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr); 20467d6f150SMatthew G Knepley } 2056562c4e1SBarry Smith PetscFunctionReturn(0); 2066562c4e1SBarry Smith } 2076562c4e1SBarry Smith 20879d04de1SBarry Smith /* newtype is ignored since only handles one case */ 2096562c4e1SBarry Smith #undef __FUNCT__ 2106562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML" 2116562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 2126562c4e1SBarry Smith { 2136562c4e1SBarry Smith PetscErrorCode ierr; 2146562c4e1SBarry Smith Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 2156562c4e1SBarry Smith Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 2166562c4e1SBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 2176562c4e1SBarry Smith PetscScalar *aa=a->a,*ba=b->a,*ca; 2186562c4e1SBarry Smith PetscInt am=A->rmap->n,an=A->cmap->n,i,j,k; 2196562c4e1SBarry Smith PetscInt *ci,*cj,ncols; 2206562c4e1SBarry Smith 2216562c4e1SBarry Smith PetscFunctionBegin; 222e32f2f54SBarry 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); 2236562c4e1SBarry Smith 2246562c4e1SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 2256562c4e1SBarry Smith ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2266562c4e1SBarry Smith ci[0] = 0; 2276562c4e1SBarry Smith for (i=0; i<am; i++) { 2286562c4e1SBarry Smith ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 2296562c4e1SBarry Smith } 2306562c4e1SBarry Smith ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 2316562c4e1SBarry Smith ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 2326562c4e1SBarry Smith 2336562c4e1SBarry Smith k = 0; 2346562c4e1SBarry Smith for (i=0; i<am; i++) { 2356562c4e1SBarry Smith /* diagonal portion of A */ 2366562c4e1SBarry Smith ncols = ai[i+1] - ai[i]; 2376562c4e1SBarry Smith for (j=0; j<ncols; j++) { 2386562c4e1SBarry Smith cj[k] = *aj++; 2396562c4e1SBarry Smith ca[k++] = *aa++; 2406562c4e1SBarry Smith } 2416562c4e1SBarry Smith /* off-diagonal portion of A */ 2426562c4e1SBarry Smith ncols = bi[i+1] - bi[i]; 2436562c4e1SBarry Smith for (j=0; j<ncols; j++) { 2446562c4e1SBarry Smith cj[k] = an + (*bj); bj++; 2456562c4e1SBarry Smith ca[k++] = *ba++; 2466562c4e1SBarry Smith } 2476562c4e1SBarry Smith } 248e32f2f54SBarry Smith if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 2496562c4e1SBarry Smith 2506562c4e1SBarry Smith /* put together the new matrix */ 2516562c4e1SBarry Smith an = mpimat->A->cmap->n+mpimat->B->cmap->n; 2526562c4e1SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 2536562c4e1SBarry Smith 2546562c4e1SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 2556562c4e1SBarry Smith /* Since these are PETSc arrays, change flags to free them as necessary. */ 2566562c4e1SBarry Smith mat = (Mat_SeqAIJ*)(*Aloc)->data; 2576562c4e1SBarry Smith mat->free_a = PETSC_TRUE; 2586562c4e1SBarry Smith mat->free_ij = PETSC_TRUE; 2596562c4e1SBarry Smith 2606562c4e1SBarry Smith mat->nonew = 0; 2616562c4e1SBarry Smith } else if (scall == MAT_REUSE_MATRIX) { 2626562c4e1SBarry Smith mat=(Mat_SeqAIJ*)(*Aloc)->data; 2636562c4e1SBarry Smith ci = mat->i; cj = mat->j; ca = mat->a; 2646562c4e1SBarry Smith for (i=0; i<am; i++) { 2656562c4e1SBarry Smith /* diagonal portion of A */ 2666562c4e1SBarry Smith ncols = ai[i+1] - ai[i]; 2676562c4e1SBarry Smith for (j=0; j<ncols; j++) *ca++ = *aa++; 2686562c4e1SBarry Smith /* off-diagonal portion of A */ 2696562c4e1SBarry Smith ncols = bi[i+1] - bi[i]; 2706562c4e1SBarry Smith for (j=0; j<ncols; j++) *ca++ = *ba++; 2716562c4e1SBarry Smith } 27279d04de1SBarry Smith } else SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 2736562c4e1SBarry Smith PetscFunctionReturn(0); 2746562c4e1SBarry Smith } 2756562c4e1SBarry Smith 2766562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat); 2776562c4e1SBarry Smith #undef __FUNCT__ 2786562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML" 2796562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A) 2806562c4e1SBarry Smith { 2816562c4e1SBarry Smith PetscErrorCode ierr; 2826562c4e1SBarry Smith Mat_MLShell *shell; 2836562c4e1SBarry Smith 2846562c4e1SBarry Smith PetscFunctionBegin; 2856562c4e1SBarry Smith ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 286601cad40SBrad Aagaard ierr = VecDestroy(&shell->y);CHKERRQ(ierr); 287601cad40SBrad Aagaard if (shell->work) {ierr = VecDestroy(&shell->work);CHKERRQ(ierr);} 2886562c4e1SBarry Smith ierr = PetscFree(shell);CHKERRQ(ierr); 2896562c4e1SBarry Smith ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 2906562c4e1SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 2916562c4e1SBarry Smith PetscFunctionReturn(0); 2926562c4e1SBarry Smith } 2936562c4e1SBarry Smith 2946562c4e1SBarry Smith #undef __FUNCT__ 2956562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ" 2966562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 2976562c4e1SBarry Smith { 2986562c4e1SBarry Smith struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 2996562c4e1SBarry Smith PetscErrorCode ierr; 300ae7fe62dSJed Brown PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = PETSC_NULL,nz_max; 30139381ba2SJed Brown PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i; 3026562c4e1SBarry Smith PetscScalar *ml_vals=matdata->values,*aa; 3036562c4e1SBarry Smith 3046562c4e1SBarry Smith PetscFunctionBegin; 305e7e72b3dSBarry Smith if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 3066562c4e1SBarry Smith if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 3076562c4e1SBarry Smith if (reuse) { 3086562c4e1SBarry Smith Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 3096562c4e1SBarry Smith aij->i = ml_rowptr; 3106562c4e1SBarry Smith aij->j = ml_cols; 3116562c4e1SBarry Smith aij->a = ml_vals; 3126562c4e1SBarry Smith } else { 3136562c4e1SBarry Smith /* sort ml_cols and ml_vals */ 3146562c4e1SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 3156562c4e1SBarry Smith for (i=0; i<m; i++) { 3166562c4e1SBarry Smith nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 3176562c4e1SBarry Smith } 3186562c4e1SBarry Smith aj = ml_cols; aa = ml_vals; 3196562c4e1SBarry Smith for (i=0; i<m; i++) { 3206562c4e1SBarry Smith ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 3216562c4e1SBarry Smith aj += nnz[i]; aa += nnz[i]; 3226562c4e1SBarry Smith } 3236562c4e1SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 3246562c4e1SBarry Smith ierr = PetscFree(nnz);CHKERRQ(ierr); 3256562c4e1SBarry Smith } 3266562c4e1SBarry Smith PetscFunctionReturn(0); 3276562c4e1SBarry Smith } 3286562c4e1SBarry Smith 32939381ba2SJed Brown nz_max = PetscMax(1,mlmat->max_nz_per_row); 33039381ba2SJed Brown ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 33139381ba2SJed Brown if (!reuse) { 3326562c4e1SBarry Smith ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 3336562c4e1SBarry Smith ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 3346562c4e1SBarry Smith ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 33539381ba2SJed Brown /* keep track of block size for A matrices */ 33639381ba2SJed Brown ierr = MatSetBlockSize (*newmat, mlmat->num_PDEs);CHKERRQ(ierr); 3376562c4e1SBarry Smith 33839381ba2SJed Brown ierr = PetscMalloc(m*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 3396562c4e1SBarry Smith for (i=0; i<m; i++) { 34039381ba2SJed Brown ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]); 3416562c4e1SBarry Smith } 3426562c4e1SBarry Smith ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 343ae7fe62dSJed Brown } 3446562c4e1SBarry Smith for (i=0; i<m; i++) { 345ae7fe62dSJed Brown PetscInt ncols; 34639381ba2SJed Brown 34739381ba2SJed Brown ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols); 348ae7fe62dSJed Brown ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr); 3496562c4e1SBarry Smith } 3506562c4e1SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3516562c4e1SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3526562c4e1SBarry Smith 3536562c4e1SBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 3546562c4e1SBarry Smith ierr = PetscFree(nnz);CHKERRQ(ierr); 3556562c4e1SBarry Smith PetscFunctionReturn(0); 3566562c4e1SBarry Smith } 3576562c4e1SBarry Smith 3586562c4e1SBarry Smith #undef __FUNCT__ 3596562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL" 3606562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 3616562c4e1SBarry Smith { 3626562c4e1SBarry Smith PetscErrorCode ierr; 3636562c4e1SBarry Smith PetscInt m,n; 3646562c4e1SBarry Smith ML_Comm *MLcomm; 3656562c4e1SBarry Smith Mat_MLShell *shellctx; 3666562c4e1SBarry Smith 3676562c4e1SBarry Smith PetscFunctionBegin; 3686562c4e1SBarry Smith m = mlmat->outvec_leng; 3696562c4e1SBarry Smith n = mlmat->invec_leng; 3706562c4e1SBarry Smith 3716562c4e1SBarry Smith if (reuse) { 3726562c4e1SBarry Smith ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 3736562c4e1SBarry Smith shellctx->mlmat = mlmat; 3746562c4e1SBarry Smith PetscFunctionReturn(0); 3756562c4e1SBarry Smith } 3766562c4e1SBarry Smith 3776562c4e1SBarry Smith MLcomm = mlmat->comm; 3786562c4e1SBarry Smith ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 3796562c4e1SBarry Smith ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 3806562c4e1SBarry Smith ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 3816562c4e1SBarry Smith ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 3826562c4e1SBarry Smith shellctx->A = *newmat; 3836562c4e1SBarry Smith shellctx->mlmat = mlmat; 38467d6f150SMatthew G Knepley shellctx->work = PETSC_NULL; 3859bb5392cSJed Brown ierr = VecCreate(MLcomm->USR_comm,&shellctx->y);CHKERRQ(ierr); 3866562c4e1SBarry Smith ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 3876562c4e1SBarry Smith ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 3886562c4e1SBarry Smith (*newmat)->ops->destroy = MatDestroy_ML; 3896562c4e1SBarry Smith PetscFunctionReturn(0); 3906562c4e1SBarry Smith } 3916562c4e1SBarry Smith 3926562c4e1SBarry Smith #undef __FUNCT__ 3936562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ" 394ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 3956562c4e1SBarry Smith { 39639381ba2SJed Brown PetscInt *aj; 39739381ba2SJed Brown PetscScalar *aa; 3986562c4e1SBarry Smith PetscErrorCode ierr; 39939381ba2SJed Brown PetscInt i,j,*gordering; 400ae7fe62dSJed Brown PetscInt m=mlmat->outvec_leng,n,nz_max,row; 4016562c4e1SBarry Smith Mat A; 4026562c4e1SBarry Smith 4036562c4e1SBarry Smith PetscFunctionBegin; 404e7e72b3dSBarry Smith if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 4056562c4e1SBarry Smith n = mlmat->invec_leng; 406e32f2f54SBarry Smith if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 4076562c4e1SBarry Smith 40839381ba2SJed Brown nz_max = PetscMax(1,mlmat->max_nz_per_row); 40939381ba2SJed Brown ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 410ae7fe62dSJed Brown if (reuse) { 411ae7fe62dSJed Brown A = *newmat; 412ae7fe62dSJed Brown } else { 413ae7fe62dSJed Brown PetscInt *nnzA,*nnzB,*nnz; 4146562c4e1SBarry Smith ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 4156562c4e1SBarry Smith ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 4166562c4e1SBarry Smith ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 41739381ba2SJed Brown /* keep track of block size for A matrices */ 41839381ba2SJed Brown ierr = MatSetBlockSize (A,mlmat->num_PDEs);CHKERRQ(ierr); 4196562c4e1SBarry Smith ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 4206562c4e1SBarry Smith 4216562c4e1SBarry Smith for (i=0; i<m; i++) { 42239381ba2SJed Brown ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]); 42339381ba2SJed Brown nnzA[i] = 0; 42439381ba2SJed Brown for (j=0; j<nnz[i]; j++) { 42539381ba2SJed Brown if (aj[j] < m) nnzA[i]++; 4266562c4e1SBarry Smith } 4276562c4e1SBarry Smith nnzB[i] = nnz[i] - nnzA[i]; 4286562c4e1SBarry Smith } 4296562c4e1SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 430ae7fe62dSJed Brown ierr = PetscFree3(nnzA,nnzB,nnz); 431ae7fe62dSJed Brown } 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++) { 435ae7fe62dSJed Brown PetscInt ncols; 4366562c4e1SBarry Smith row = gordering[i]; 43739381ba2SJed Brown 43839381ba2SJed Brown ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols); 43939381ba2SJed Brown for (j = 0; j < ncols; j++) { 44039381ba2SJed Brown aj[j] = gordering[aj[j]]; 4416562c4e1SBarry Smith } 442ae7fe62dSJed Brown ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr); 4436562c4e1SBarry Smith } 444eb4736aaSBarry Smith ML_free(gordering); 4456562c4e1SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4466562c4e1SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4476562c4e1SBarry Smith *newmat = A; 4486562c4e1SBarry Smith 4496562c4e1SBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 4506562c4e1SBarry Smith PetscFunctionReturn(0); 4516562c4e1SBarry Smith } 4526562c4e1SBarry Smith 45339381ba2SJed Brown /* -------------------------------------------------------------------------- */ 45439381ba2SJed Brown /* 45539381ba2SJed Brown PCSetCoordinates_ML 45639381ba2SJed Brown 45739381ba2SJed Brown Input Parameter: 45839381ba2SJed Brown . pc - the preconditioner context 45939381ba2SJed Brown */ 46039381ba2SJed Brown #undef __FUNCT__ 46139381ba2SJed Brown #define __FUNCT__ "PCSetCoordinates_ML" 46239381ba2SJed Brown PETSC_EXTERN_C PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 46339381ba2SJed Brown { 46439381ba2SJed Brown PC_MG *mg = (PC_MG*)pc->data; 46539381ba2SJed Brown PC_ML *pc_ml = (PC_ML*)mg->innerctx; 46639381ba2SJed Brown PetscErrorCode ierr; 46739381ba2SJed Brown PetscInt arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend; 46839381ba2SJed Brown Mat Amat = pc->pmat; 46939381ba2SJed Brown 47039381ba2SJed Brown /* this function copied and modified from PCSetCoordinates_GEO -TGI */ 47139381ba2SJed Brown PetscFunctionBegin; 47239381ba2SJed Brown PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1); 47339381ba2SJed Brown ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 47439381ba2SJed Brown 47539381ba2SJed Brown ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr); 47639381ba2SJed Brown nloc = (Iend-my0)/bs; 47739381ba2SJed Brown 47839381ba2SJed Brown if (nloc!=a_nloc)SETERRQ2(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Number of local blocks must locations = %d %d.",a_nloc,nloc); 47939381ba2SJed Brown if ((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 48039381ba2SJed Brown 48139381ba2SJed Brown oldarrsz = pc_ml->dim * pc_ml->nloc; 48239381ba2SJed Brown pc_ml->dim = ndm; 48339381ba2SJed Brown pc_ml->nloc = a_nloc; 48439381ba2SJed Brown arrsz = ndm * a_nloc; 48539381ba2SJed Brown 48639381ba2SJed Brown /* create data - syntactic sugar that should be refactored at some point */ 48739381ba2SJed Brown if (pc_ml->coords==0 || (oldarrsz != arrsz)) { 48839381ba2SJed Brown ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr); 48939381ba2SJed Brown ierr = PetscMalloc((arrsz)*sizeof(PetscReal), &pc_ml->coords);CHKERRQ(ierr); 49039381ba2SJed Brown } 49139381ba2SJed Brown for (kk=0;kk<arrsz;kk++)pc_ml->coords[kk] = -999.; 49239381ba2SJed Brown /* copy data in - column oriented */ 49339381ba2SJed Brown for (kk = 0; kk < nloc; kk++) { 49439381ba2SJed Brown for (ii = 0; ii < ndm; ii++) { 49539381ba2SJed Brown pc_ml->coords[ii*nloc + kk] = coords[kk*ndm + ii]; 49639381ba2SJed Brown } 49739381ba2SJed Brown } 49839381ba2SJed Brown 49939381ba2SJed Brown PetscFunctionReturn(0); 50039381ba2SJed Brown } 50139381ba2SJed Brown 5026562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/ 50301da6913SBarry Smith #undef __FUNCT__ 504a06653b4SBarry Smith #define __FUNCT__ "PCReset_ML" 50516336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc) 50601da6913SBarry Smith { 50701da6913SBarry Smith PetscErrorCode ierr; 508e0262f48SMatthew G Knepley PC_MG *mg = (PC_MG*)pc->data; 509e0262f48SMatthew G Knepley PC_ML *pc_ml = (PC_ML*)mg->innerctx; 51039381ba2SJed Brown PetscInt level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim; 51101da6913SBarry Smith 51201da6913SBarry Smith PetscFunctionBegin; 51339381ba2SJed Brown if (dim) { 51439381ba2SJed Brown ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats *) pc_ml->ml_object->Grid[0].Grid; 51539381ba2SJed Brown 51639381ba2SJed Brown for (level=0; level<=fine_level; level++) { 51739381ba2SJed Brown ierr = VecDestroy(&pc_ml->gridctx[level].coords);CHKERRQ(ierr); 51839381ba2SJed Brown } 51939381ba2SJed Brown 52039381ba2SJed Brown grid_info->x = 0; /* do this so ML doesn't try to free coordinates */ 52139381ba2SJed Brown grid_info->y = 0; 52239381ba2SJed Brown grid_info->z = 0; 52339381ba2SJed Brown 52439381ba2SJed Brown ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object); 52539381ba2SJed Brown } 52601da6913SBarry Smith ML_Aggregate_Destroy(&pc_ml->agg_object); 52701da6913SBarry Smith ML_Destroy(&pc_ml->ml_object); 52801da6913SBarry Smith 52901da6913SBarry Smith if (pc_ml->PetscMLdata) { 53001da6913SBarry Smith ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 531ae7fe62dSJed Brown ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr); 532ae7fe62dSJed Brown ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr); 533ae7fe62dSJed Brown ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr); 53401da6913SBarry Smith } 53501da6913SBarry Smith ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 53601da6913SBarry Smith 537f5a5dd59SJed Brown if (pc_ml->gridctx) { 53801da6913SBarry Smith for (level=0; level<fine_level; level++) { 539601cad40SBrad Aagaard if (pc_ml->gridctx[level].A) {ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);} 540601cad40SBrad Aagaard if (pc_ml->gridctx[level].P) {ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);} 541601cad40SBrad Aagaard if (pc_ml->gridctx[level].R) {ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);} 542601cad40SBrad Aagaard if (pc_ml->gridctx[level].x) {ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);} 543601cad40SBrad Aagaard if (pc_ml->gridctx[level].b) {ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);} 544601cad40SBrad Aagaard if (pc_ml->gridctx[level+1].r) {ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 54501da6913SBarry Smith } 546f5a5dd59SJed Brown } 54701da6913SBarry Smith ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 54839381ba2SJed Brown ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr); 54939381ba2SJed Brown pc_ml->dim = 0; 55039381ba2SJed Brown pc_ml->nloc = 0; 55101da6913SBarry Smith PetscFunctionReturn(0); 55201da6913SBarry Smith } 5535582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 5545582bec1SHong Zhang /* 5555582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 5565582bec1SHong Zhang by setting data structures and options. 5575582bec1SHong Zhang 5585582bec1SHong Zhang Input Parameter: 5595582bec1SHong Zhang . pc - the preconditioner context 5605582bec1SHong Zhang 5615582bec1SHong Zhang Application Interface Routine: PCSetUp() 5625582bec1SHong Zhang 5635582bec1SHong Zhang Notes: 5645582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 5655582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 5665582bec1SHong Zhang */ 5676ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC); 568a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC); 569c07bf074SBarry Smith 5705582bec1SHong Zhang #undef __FUNCT__ 5715582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 5726ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 5735582bec1SHong Zhang { 5745582bec1SHong Zhang PetscErrorCode ierr; 575eef31507SHong Zhang PetscMPIInt size; 5765582bec1SHong Zhang FineGridCtx *PetscMLdata; 5775582bec1SHong Zhang ML *ml_object; 5785582bec1SHong Zhang ML_Aggregate *agg_object; 5795582bec1SHong Zhang ML_Operator *mlmat; 5804f8eab3cSJed Brown PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs; 5815582bec1SHong Zhang Mat A,Aloc; 5825582bec1SHong Zhang GridCtx *gridctx; 58301da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 58401da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 585ace3abfcSBarry Smith PetscBool isSeq, isMPI; 586c07bf074SBarry Smith KSP smoother; 587c07bf074SBarry Smith PC subpc; 58848268eb4SJed Brown PetscInt mesh_level, old_mesh_level; 58948268eb4SJed Brown 5905582bec1SHong Zhang PetscFunctionBegin; 59148268eb4SJed Brown A = pc->pmat; 59248268eb4SJed Brown ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 59348268eb4SJed Brown 594573998d7SHong Zhang if (pc->setupcalled) { 59548268eb4SJed Brown if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) { 59648268eb4SJed Brown /* 59748268eb4SJed Brown Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for 59848268eb4SJed Brown multiple solves in which the matrix is not changing too quickly. 59948268eb4SJed Brown */ 60048268eb4SJed Brown ml_object = pc_ml->ml_object; 60148268eb4SJed Brown gridctx = pc_ml->gridctx; 60248268eb4SJed Brown Nlevels = pc_ml->Nlevels; 60348268eb4SJed Brown fine_level = Nlevels - 1; 60448268eb4SJed Brown gridctx[fine_level].A = A; 60548268eb4SJed Brown 606251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 607251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 60848268eb4SJed Brown if (isMPI) { 60948268eb4SJed Brown ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr); 61048268eb4SJed Brown } else if (isSeq) { 61148268eb4SJed Brown Aloc = A; 612ae7fe62dSJed Brown ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr); 61348268eb4SJed Brown } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name); 61448268eb4SJed Brown 61548268eb4SJed Brown ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 61648268eb4SJed Brown PetscMLdata = pc_ml->PetscMLdata; 617ae7fe62dSJed Brown ierr = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr); 61848268eb4SJed Brown PetscMLdata->A = A; 61948268eb4SJed Brown PetscMLdata->Aloc = Aloc; 62048268eb4SJed Brown ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 62148268eb4SJed Brown ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 62248268eb4SJed Brown 62348268eb4SJed Brown mesh_level = ml_object->ML_finest_level; 62448268eb4SJed Brown while (ml_object->SingleLevel[mesh_level].Rmat->to) { 62548268eb4SJed Brown old_mesh_level = mesh_level; 62648268eb4SJed Brown mesh_level = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum; 62748268eb4SJed Brown 62848268eb4SJed Brown /* clean and regenerate A */ 62948268eb4SJed Brown mlmat = &(ml_object->Amat[mesh_level]); 63048268eb4SJed Brown ML_Operator_Clean(mlmat); 63148268eb4SJed Brown ML_Operator_Init(mlmat,ml_object->comm); 63248268eb4SJed Brown ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level); 63348268eb4SJed Brown } 63448268eb4SJed Brown 63548268eb4SJed Brown level = fine_level - 1; 63648268eb4SJed Brown if (size == 1) { /* convert ML P, R and A into seqaij format */ 63748268eb4SJed Brown for (mllevel=1; mllevel<Nlevels; mllevel++) { 63848268eb4SJed Brown mlmat = &(ml_object->Amat[mllevel]); 639ae7fe62dSJed Brown ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 64048268eb4SJed Brown level--; 64148268eb4SJed Brown } 64248268eb4SJed Brown } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 64348268eb4SJed Brown for (mllevel=1; mllevel<Nlevels; mllevel++) { 64448268eb4SJed Brown mlmat = &(ml_object->Amat[mllevel]); 645ae7fe62dSJed Brown ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 64648268eb4SJed Brown level--; 64748268eb4SJed Brown } 64848268eb4SJed Brown } 64948268eb4SJed Brown 65048268eb4SJed Brown for (level=0; level<fine_level; level++) { 65148268eb4SJed Brown if (level > 0) { 65248268eb4SJed Brown ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 65348268eb4SJed Brown } 6547721a10bSJed Brown ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 65548268eb4SJed Brown } 65648268eb4SJed Brown ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 6577721a10bSJed Brown ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 65848268eb4SJed Brown 65948268eb4SJed Brown ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 66048268eb4SJed Brown PetscFunctionReturn(0); 66148268eb4SJed Brown } else { 662c07bf074SBarry Smith /* since ML can change the size of vectors/matrices at any level we must destroy everything */ 66316336fedSMatthew G Knepley ierr = PCReset_ML(pc);CHKERRQ(ierr); 664a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 665573998d7SHong Zhang } 66648268eb4SJed Brown } 667573998d7SHong Zhang 6685582bec1SHong Zhang /* setup special features of PCML */ 6695582bec1SHong Zhang /*--------------------------------*/ 6705582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 6715582bec1SHong Zhang pc_ml->size = size; 672251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 673251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 674864b637dSMatthew Knepley if (isMPI) { 675db571536SBarry Smith ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr); 676864b637dSMatthew Knepley } else if (isSeq) { 6775582bec1SHong Zhang Aloc = A; 678ae7fe62dSJed Brown ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr); 679c5b2ea42SJed Brown } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name); 6805582bec1SHong Zhang 6815582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 68238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 6835582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 684d0f46423SBarry Smith ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 6855582bec1SHong Zhang 68624a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 687d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr); 68824a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 68924a42b14SHong Zhang 69024a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 691d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 69224a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 693573998d7SHong Zhang PetscMLdata->A = A; 694573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 69539381ba2SJed Brown if (pc_ml->dim) { /* create vecs around the coordinate data given */ 69639381ba2SJed Brown PetscInt i,j,dim=pc_ml->dim; 69739381ba2SJed Brown PetscInt nloc = pc_ml->nloc,nlocghost; 69839381ba2SJed Brown PetscReal *ghostedcoords; 69939381ba2SJed Brown 70039381ba2SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 70139381ba2SJed Brown nlocghost = Aloc->cmap->n / bs; 70239381ba2SJed Brown ierr = PetscMalloc(dim*nlocghost*sizeof(PetscReal),&ghostedcoords);CHKERRQ(ierr); 70339381ba2SJed Brown for (i = 0; i < dim; i++) { 70439381ba2SJed Brown /* copy coordinate values into first component of pwork */ 70539381ba2SJed Brown for (j = 0; j < nloc; j++) { 70639381ba2SJed Brown PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j]; 70739381ba2SJed Brown } 70839381ba2SJed Brown /* get the ghost values */ 70939381ba2SJed Brown ierr = PetscML_comm(PetscMLdata->pwork,PetscMLdata);CHKERRQ(ierr); 71039381ba2SJed Brown /* write into the vector */ 71139381ba2SJed Brown for (j = 0; j < nlocghost; j++) { 71239381ba2SJed Brown ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j]; 71339381ba2SJed Brown } 71439381ba2SJed Brown } 71539381ba2SJed Brown /* replace the original coords with the ghosted coords, because these are 71639381ba2SJed Brown * what ML needs */ 71739381ba2SJed Brown ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr); 71839381ba2SJed Brown pc_ml->coords = ghostedcoords; 71939381ba2SJed Brown } 72024a42b14SHong Zhang 7215582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 72245cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 7235582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 7244f8eab3cSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 7255582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 726ead7dcbeSHong Zhang ML_Comm_Set_UsrComm(ml_object->comm,((PetscObject)A)->comm); 727573998d7SHong Zhang pc_ml->ml_object = ml_object; 7285582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 7295582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 7305582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 7315582bec1SHong Zhang 732b5c8bdf8SJed Brown ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO); 733b5c8bdf8SJed Brown 7345582bec1SHong Zhang /* aggregation */ 7355582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 736573998d7SHong Zhang pc_ml->agg_object = agg_object; 737573998d7SHong Zhang 738fb6a8e6dSJed Brown { 739fb6a8e6dSJed Brown MatNullSpace mnull; 740fb6a8e6dSJed Brown ierr = MatGetNearNullSpace(A,&mnull);CHKERRQ(ierr); 741fb6a8e6dSJed Brown if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) { 742fb6a8e6dSJed Brown if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER; 743fb6a8e6dSJed Brown else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK; 744fb6a8e6dSJed Brown else pc_ml->nulltype = PCML_NULLSPACE_SCALAR; 745fb6a8e6dSJed Brown } 746fb6a8e6dSJed Brown switch (pc_ml->nulltype) { 747fb6a8e6dSJed Brown case PCML_NULLSPACE_USER: { 748fb6a8e6dSJed Brown PetscScalar *nullvec; 749fb6a8e6dSJed Brown const PetscScalar *v; 750fb6a8e6dSJed Brown PetscBool has_const; 7511c547e14SJed Brown PetscInt i,j,mlocal,nvec,M; 752fb6a8e6dSJed Brown const Vec *vecs; 753fb6a8e6dSJed Brown if (!mnull) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space"); 7541c547e14SJed Brown ierr = MatGetSize(A,&M,PETSC_NULL);CHKERRQ(ierr); 755fb6a8e6dSJed Brown ierr = MatGetLocalSize(Aloc,&mlocal,PETSC_NULL);CHKERRQ(ierr); 756fb6a8e6dSJed Brown ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr); 7578caf3d72SBarry Smith ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr); 7581c547e14SJed Brown if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M; 759fb6a8e6dSJed Brown for (i=0; i<nvec; i++) { 760fb6a8e6dSJed Brown ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 761fb6a8e6dSJed Brown for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j]; 762fb6a8e6dSJed Brown ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 763fb6a8e6dSJed Brown } 764fb6a8e6dSJed Brown ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal);CHKERRQ(ierr); 765fb6a8e6dSJed Brown ierr = PetscFree(nullvec);CHKERRQ(ierr); 766fb6a8e6dSJed Brown } break; 767fb6a8e6dSJed Brown case PCML_NULLSPACE_BLOCK: 768fb6a8e6dSJed Brown ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr); 769fb6a8e6dSJed Brown break; 770fb6a8e6dSJed Brown case PCML_NULLSPACE_SCALAR: 771fb6a8e6dSJed Brown break; 772fb6a8e6dSJed Brown default: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unknown null space type"); 773fb6a8e6dSJed Brown } 774fb6a8e6dSJed Brown } 7755582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 7765582bec1SHong Zhang /* set options */ 7775582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 7785582bec1SHong Zhang case 1: 7795582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 7805582bec1SHong Zhang case 2: 7815582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 7825582bec1SHong Zhang case 3: 7835582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 7845582bec1SHong Zhang } 7855582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 7865582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 7875582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm) { 7887ffd031bSHong Zhang ML_Set_SpectralNormScheme_Anorm(ml_object); 7895582bec1SHong Zhang } 790b5c8bdf8SJed Brown agg_object->keep_agg_information = (int)pc_ml->KeepAggInfo; 791b5c8bdf8SJed Brown agg_object->keep_P_tentative = (int)pc_ml->Reusable; 792b5c8bdf8SJed Brown agg_object->block_scaled_SA = (int)pc_ml->BlockScaling; 793b5c8bdf8SJed Brown agg_object->minimizing_energy = (int)pc_ml->EnergyMinimization; 794b5c8bdf8SJed Brown agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol; 795b5c8bdf8SJed Brown agg_object->cheap_minimizing_energy = (int)pc_ml->EnergyMinimizationCheap; 7965582bec1SHong Zhang 79739381ba2SJed Brown if (pc_ml->Aux) { 798*f23aa3ddSBarry Smith if (!pc_ml->dim) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_USER,"Auxiliary matrix requires coordinates"); 79939381ba2SJed Brown ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold; 80039381ba2SJed Brown ml_object->Amat[0].aux_data->enable = 1; 80139381ba2SJed Brown ml_object->Amat[0].aux_data->max_level = 10; 80239381ba2SJed Brown ml_object->Amat[0].num_PDEs = bs; 80339381ba2SJed Brown } 80439381ba2SJed Brown 80539381ba2SJed Brown if (pc_ml->dim) { 80639381ba2SJed Brown PetscInt i,dim = pc_ml->dim; 80739381ba2SJed Brown ML_Aggregate_Viz_Stats *grid_info; 80839381ba2SJed Brown PetscInt nlocghost; 80939381ba2SJed Brown 81039381ba2SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 81139381ba2SJed Brown nlocghost = Aloc->cmap->n / bs; 81239381ba2SJed Brown 81339381ba2SJed Brown ML_Aggregate_VizAndStats_Setup(ml_object); /* create ml info for coords */ 81439381ba2SJed Brown grid_info = (ML_Aggregate_Viz_Stats *) ml_object->Grid[0].Grid; 81539381ba2SJed Brown for (i = 0; i < dim; i++) { 81639381ba2SJed Brown /* set the finest level coordinates to point to the column-order array 81739381ba2SJed Brown * in pc_ml */ 81839381ba2SJed Brown /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */ 81939381ba2SJed Brown switch (i) { 82039381ba2SJed Brown case 0: grid_info->x = pc_ml->coords + nlocghost * i; break; 82139381ba2SJed Brown case 1: grid_info->y = pc_ml->coords + nlocghost * i; break; 82239381ba2SJed Brown case 2: grid_info->z = pc_ml->coords + nlocghost * i; break; 82339381ba2SJed Brown default: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3"); 82439381ba2SJed Brown } 82539381ba2SJed Brown } 82639381ba2SJed Brown grid_info->Ndim = dim; 82739381ba2SJed Brown } 82839381ba2SJed Brown 82939381ba2SJed Brown /* repartitioning */ 83039381ba2SJed Brown if (pc_ml->Repartition) { 83139381ba2SJed Brown ML_Repartition_Activate (ml_object); 83239381ba2SJed Brown ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio); 83339381ba2SJed Brown ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc); 83439381ba2SJed Brown ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc); 83539381ba2SJed Brown #if 0 /* Function not yet defined in ml-6.2 */ 83639381ba2SJed Brown /* I'm not sure what compatibility issues might crop up if we partitioned 83739381ba2SJed Brown * on the finest level, so to be safe repartition starts on the next 83839381ba2SJed Brown * finest level (reflection default behavior in 83939381ba2SJed Brown * ml_MultiLevelPreconditioner) */ 84039381ba2SJed Brown ML_Repartition_Set_StartLevel(ml_object,1); 84139381ba2SJed Brown #endif 84239381ba2SJed Brown 84339381ba2SJed Brown if (!pc_ml->RepartitionType) { 84439381ba2SJed Brown PetscInt i; 84539381ba2SJed Brown 846*f23aa3ddSBarry Smith if (!pc_ml->dim) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates"); 84739381ba2SJed Brown ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN); 84839381ba2SJed Brown ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim); 84939381ba2SJed Brown 85039381ba2SJed Brown for (i = 0;i < ml_object->ML_num_levels;i++) { 85139381ba2SJed Brown ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[i].Grid; 85239381ba2SJed Brown grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */ 85339381ba2SJed Brown /* defaults from ml_agg_info.c */ 85439381ba2SJed Brown grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */ 85539381ba2SJed Brown grid_info->zoltan_timers = 0; 85639381ba2SJed Brown grid_info->smoothing_steps = 4; /* only relevant to hypergraph / fast hypergraph */ 85739381ba2SJed Brown } 85839381ba2SJed Brown } 85939381ba2SJed Brown else { 86039381ba2SJed Brown ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS); 86139381ba2SJed Brown } 86239381ba2SJed Brown } 86339381ba2SJed Brown 864b5c8bdf8SJed Brown if (pc_ml->OldHierarchy) { 8655582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 866b5c8bdf8SJed Brown } else { 867b5c8bdf8SJed Brown Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 868b5c8bdf8SJed Brown } 86965e19b50SBarry Smith if (Nlevels<=0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 870573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 871aa85bbbfSHong Zhang fine_level = Nlevels - 1; 872c07bf074SBarry Smith 87397177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 874aa85bbbfSHong Zhang /* set default smoothers */ 875aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++) { 876aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 877aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 878aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 879aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 880aa85bbbfSHong Zhang } 881f2e59741SMatthew G Knepley ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 88297177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 883f2e59741SMatthew G Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 8845582bec1SHong Zhang 8855582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 8865582bec1SHong Zhang pc_ml->gridctx = gridctx; 8875582bec1SHong Zhang 8885582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 8895582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 890e14861a4SHong Zhang gridctx[fine_level].A = A; 891573998d7SHong Zhang 892e14861a4SHong Zhang level = fine_level - 1; 893ab718edeSHong Zhang if (size == 1) { /* convert ML P, R and A into seqaij format */ 8945582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++) { 895e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 896db571536SBarry Smith ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr); 897e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 898db571536SBarry Smith ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr); 899573998d7SHong Zhang 900573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 901573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 9025582bec1SHong Zhang level--; 9035582bec1SHong Zhang } 904ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 9055582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++) { 9065582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 907db571536SBarry Smith ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr); 908ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 909db571536SBarry Smith ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr); 910573998d7SHong Zhang 9115582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 912ae7fe62dSJed Brown ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 9135582bec1SHong Zhang level--; 9145582bec1SHong Zhang } 9155582bec1SHong Zhang } 9165582bec1SHong Zhang 917573998d7SHong Zhang /* create vectors and ksp at all levels */ 918ac346b81SHong Zhang for (level=0; level<fine_level; level++) { 919573998d7SHong Zhang level1 = level + 1; 920e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr); 921d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr); 9225582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 92397177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 9245582bec1SHong Zhang 925e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr); 926d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 9275582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 92897177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 929ac346b81SHong Zhang 930e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr); 931d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 932ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 93397177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 934ac346b81SHong Zhang 9355582bec1SHong Zhang if (level == 0) { 93697177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 9375582bec1SHong Zhang } else { 93897177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 939573998d7SHong Zhang } 940573998d7SHong Zhang } 941573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 942573998d7SHong Zhang 943573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 944573998d7SHong Zhang for (level=0; level<fine_level; level++) { 945573998d7SHong Zhang level1 = level + 1; 946aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr); 947573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 948573998d7SHong Zhang if (level > 0) { 94997177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 9505582bec1SHong Zhang } 9515582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 9525582bec1SHong Zhang } 95397177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 954ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 9555582bec1SHong Zhang 95639381ba2SJed Brown /* put coordinate info in levels */ 95739381ba2SJed Brown if (pc_ml->dim) { 95839381ba2SJed Brown PetscInt i,j,dim = pc_ml->dim; 95939381ba2SJed Brown PetscInt bs, nloc; 96039381ba2SJed Brown PC subpc; 96139381ba2SJed Brown PetscReal *array; 96239381ba2SJed Brown 96339381ba2SJed Brown level = fine_level; 96439381ba2SJed Brown for (mllevel = 0; mllevel < Nlevels; mllevel++) { 965ebbbbe33SJed Brown ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Amat[mllevel].to->Grid->Grid; 96639381ba2SJed Brown MPI_Comm comm = ((PetscObject)gridctx[level].A)->comm; 96739381ba2SJed Brown 96839381ba2SJed Brown ierr = MatGetBlockSize (gridctx[level].A, &bs);CHKERRQ(ierr); 96939381ba2SJed Brown ierr = MatGetLocalSize (gridctx[level].A, PETSC_NULL, &nloc);CHKERRQ(ierr); 97039381ba2SJed Brown nloc /= bs; /* number of local nodes */ 97139381ba2SJed Brown 97239381ba2SJed Brown ierr = VecCreate(comm,&gridctx[level].coords);CHKERRQ(ierr); 97339381ba2SJed Brown ierr = VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE);CHKERRQ(ierr); 97439381ba2SJed Brown ierr = VecSetType(gridctx[level].coords,VECMPI);CHKERRQ(ierr); 97539381ba2SJed Brown ierr = VecGetArray(gridctx[level].coords,&array);CHKERRQ(ierr); 97639381ba2SJed Brown for (j = 0; j < nloc; j++) { 97739381ba2SJed Brown for (i = 0; i < dim; i++) { 97839381ba2SJed Brown switch (i) { 97939381ba2SJed Brown case 0: array[dim * j + i] = grid_info->x[j]; break; 98039381ba2SJed Brown case 1: array[dim * j + i] = grid_info->y[j]; break; 98139381ba2SJed Brown case 2: array[dim * j + i] = grid_info->z[j]; break; 98239381ba2SJed Brown default: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3"); 98339381ba2SJed Brown } 98439381ba2SJed Brown } 98539381ba2SJed Brown } 98639381ba2SJed Brown 98739381ba2SJed Brown /* passing coordinates to smoothers/coarse solver, should they need them */ 98839381ba2SJed Brown ierr = KSPGetPC(gridctx[level].ksp,&subpc);CHKERRQ(ierr); 98939381ba2SJed Brown ierr = PCSetCoordinates(subpc,dim,nloc,array);CHKERRQ(ierr); 99039381ba2SJed Brown ierr = VecRestoreArray(gridctx[level].coords,&array);CHKERRQ(ierr); 99139381ba2SJed Brown level--; 99239381ba2SJed Brown } 99339381ba2SJed Brown } 99439381ba2SJed Brown 995c07bf074SBarry Smith /* setupcalled is set to 0 so that MG is setup from scratch */ 996c07bf074SBarry Smith pc->setupcalled = 0; 9973751b4bdSBarry Smith ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 9985582bec1SHong Zhang PetscFunctionReturn(0); 9995582bec1SHong Zhang } 10005582bec1SHong Zhang 10015582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 10025582bec1SHong Zhang /* 10035582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 10045582bec1SHong Zhang that was created with PCCreate_ML(). 10055582bec1SHong Zhang 10065582bec1SHong Zhang Input Parameter: 10075582bec1SHong Zhang . pc - the preconditioner context 10085582bec1SHong Zhang 10095582bec1SHong Zhang Application Interface Routine: PCDestroy() 10105582bec1SHong Zhang */ 10115582bec1SHong Zhang #undef __FUNCT__ 10125582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 10136ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 10145582bec1SHong Zhang { 10155582bec1SHong Zhang PetscErrorCode ierr; 101601da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 101701da6913SBarry Smith PC_ML *pc_ml= (PC_ML*)mg->innerctx; 10185582bec1SHong Zhang 10195582bec1SHong Zhang PetscFunctionBegin; 102016336fedSMatthew G Knepley ierr = PCReset_ML(pc);CHKERRQ(ierr); 102101da6913SBarry Smith ierr = PetscFree(pc_ml);CHKERRQ(ierr); 102201da6913SBarry Smith ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 102339381ba2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","",PETSC_NULL);CHKERRQ(ierr); 10245582bec1SHong Zhang PetscFunctionReturn(0); 10255582bec1SHong Zhang } 10265582bec1SHong Zhang 10275582bec1SHong Zhang #undef __FUNCT__ 10285582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 10296ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 10305582bec1SHong Zhang { 10315582bec1SHong Zhang PetscErrorCode ierr; 103239381ba2SJed Brown PetscInt indx,PrintLevel,partindx; 10335582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 103439381ba2SJed Brown const char *part[] = {"Zoltan","ParMETIS"}; 103539381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN) 103639381ba2SJed Brown PetscInt zidx; 103739381ba2SJed Brown const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"}; 103839381ba2SJed Brown #endif 103901da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 104001da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 1041b5c8bdf8SJed Brown PetscMPIInt size; 104288ff4cc7SJed Brown MPI_Comm comm = ((PetscObject)pc)->comm; 10435582bec1SHong Zhang 10445582bec1SHong Zhang PetscFunctionBegin; 104588ff4cc7SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 10465582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 10475582bec1SHong Zhang PrintLevel = 0; 10485582bec1SHong Zhang indx = 0; 104939381ba2SJed Brown partindx = 0; 10505582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 10515582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 1052573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 1053573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 10543751b4bdSBarry Smith ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 10555582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 1056573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 1057573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 1058acfcf0e5SJed 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); 1059acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,PETSC_NULL);CHKERRQ(ierr); 1060acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,PETSC_NULL);CHKERRQ(ierr); 1061fb6a8e6dSJed Brown ierr = PetscOptionsEnum("-pc_ml_nullspace","Which type of null space information to use","None",PCMLNullSpaceTypes,(PetscEnum)pc_ml->nulltype,(PetscEnum*)&pc_ml->nulltype,PETSC_NULL);CHKERRQ(ierr); 1062b5c8bdf8SJed 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); 106348268eb4SJed Brown ierr = PetscOptionsBool("-pc_ml_reuse_interpolation","Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)","None",pc_ml->reuse_interpolation,&pc_ml->reuse_interpolation,PETSC_NULL);CHKERRQ(ierr); 1064b5c8bdf8SJed Brown /* 1065b5c8bdf8SJed Brown The following checks a number of conditions. If we let this stuff slip by, then ML's error handling will take over. 1066b5c8bdf8SJed Brown This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions. 1067b5c8bdf8SJed Brown 1068b5c8bdf8SJed Brown We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working 1069b5c8bdf8SJed Brown combination of options and ML's exit(1) explanations don't help matters. 1070b5c8bdf8SJed Brown */ 107188ff4cc7SJed Brown if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4"); 107288ff4cc7SJed Brown if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel"); 1073b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2");CHKERRQ(ierr);} 1074b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) { 1075b5c8bdf8SJed Brown ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,PETSC_NULL);CHKERRQ(ierr); 1076b5c8bdf8SJed Brown } 1077b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == 2) { 1078b5c8bdf8SJed Brown /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */ 1079acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,PETSC_NULL);CHKERRQ(ierr); 1080b5c8bdf8SJed Brown } 1081b5c8bdf8SJed Brown /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */ 1082b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE; 1083acfcf0e5SJed 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); 1084b5c8bdf8SJed Brown /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */ 1085b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE; 1086acfcf0e5SJed 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); 1087b5c8bdf8SJed Brown /* 1088b5c8bdf8SJed Brown ML's C API is severely underdocumented and lacks significant functionality. The C++ API calls 1089b5c8bdf8SJed Brown ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function 1090b5c8bdf8SJed Brown ML_Gen_MGHierarchy_UsingAggregation(). This modification, however, does not provide a strict superset of the 1091b5c8bdf8SJed Brown functionality in the old function, so some users may still want to use it. Note that many options are ignored in 1092b5c8bdf8SJed Brown this context, but ML doesn't provide a way to find out which ones. 1093b5c8bdf8SJed Brown */ 1094acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,PETSC_NULL);CHKERRQ(ierr); 109539381ba2SJed Brown ierr = PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the heirarchy","ML_Repartition_Activate",pc_ml->Repartition,&pc_ml->Repartition,PETSC_NULL);CHKERRQ(ierr); 109639381ba2SJed Brown if (pc_ml->Repartition) { 109739381ba2SJed Brown ierr = PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,PETSC_NULL);CHKERRQ(ierr); 109839381ba2SJed Brown ierr = PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,PETSC_NULL);CHKERRQ(ierr); 109939381ba2SJed Brown ierr = PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor","ML_Repartition_Set_PutOnSingleProc",pc_ml->PutOnSingleProc,&pc_ml->PutOnSingleProc,PETSC_NULL);CHKERRQ(ierr); 110039381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN) 110139381ba2SJed Brown partindx = 0; 110239381ba2SJed Brown ierr = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,PETSC_NULL);CHKERRQ(ierr); 110339381ba2SJed Brown pc_ml->RepartitionType = partindx; 110439381ba2SJed Brown if (!partindx) { 11055572b5bbSJed Brown PetscInt zindx = 0; 110639381ba2SJed Brown ierr = PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,PETSC_NULL);CHKERRQ(ierr); 110739381ba2SJed Brown pc_ml->ZoltanScheme = zindx; 110839381ba2SJed Brown } 110939381ba2SJed Brown #else 111039381ba2SJed Brown partindx = 1; 111139381ba2SJed Brown ierr = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,PETSC_NULL);CHKERRQ(ierr); 1112*f23aa3ddSBarry Smith if (!partindx) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan"); 111339381ba2SJed Brown #endif 111439381ba2SJed Brown ierr = PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,PETSC_NULL);CHKERRQ(ierr); 111539381ba2SJed Brown ierr = PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,PETSC_NULL);CHKERRQ(ierr); 111639381ba2SJed Brown } 11175582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 11185582bec1SHong Zhang PetscFunctionReturn(0); 11195582bec1SHong Zhang } 11205582bec1SHong Zhang 11215582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 11225582bec1SHong Zhang /* 11235582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 11245582bec1SHong Zhang and sets this as the private data within the generic preconditioning 11255582bec1SHong Zhang context, PC, that was created within PCCreate(). 11265582bec1SHong Zhang 11275582bec1SHong Zhang Input Parameter: 11285582bec1SHong Zhang . pc - the preconditioner context 11295582bec1SHong Zhang 11305582bec1SHong Zhang Application Interface Routine: PCCreate() 11315582bec1SHong Zhang */ 11325582bec1SHong Zhang 11335582bec1SHong Zhang /*MC 11341e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 11355582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 11366ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 11376ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 11385582bec1SHong Zhang 11396ca4d86aSHong Zhang Options Database Key: 11406ca4d86aSHong Zhang Multigrid options(inherited) 11416ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 11426ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 11436ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 11448c1c2452SJed Brown -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 114551f519a2SBarry Smith ML options: 1146ccd75124SBarry Smith . -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 11476ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 11486ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 1149f41ab451SVictor Eijkhout . -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS 11506ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 11516ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 115239381ba2SJed Brown . -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm) 115339381ba2SJed Brown . -pc_ml_repartition <false>: Allow ML to repartition levels of the heirarchy (ML_Repartition_Activate) 115439381ba2SJed Brown . -pc_ml_repartitionMaxMinRatio <1.3>: Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio) 115539381ba2SJed Brown . -pc_ml_repartitionMinPerProc <512>: Smallest repartitioned size (ML_Repartition_Set_MinPerProc) 115639381ba2SJed Brown . -pc_ml_repartitionPutOnSingleProc <5000>: Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc) 115739381ba2SJed Brown . -pc_ml_repartitionType <Zoltan>: Repartitioning library to use (ML_Repartition_Set_Partitioner) 115839381ba2SJed Brown . -pc_ml_repartitionZoltanScheme <RCB>: Repartitioning scheme to use (None) 115939381ba2SJed Brown . -pc_ml_Aux <false>: Aggregate using auxiliary coordinate-based laplacian (None) 116039381ba2SJed Brown - -pc_ml_AuxThreshold <0.0>: Auxiliary smoother drop tol (None) 11615582bec1SHong Zhang 11625582bec1SHong Zhang Level: intermediate 11635582bec1SHong Zhang 11645582bec1SHong Zhang Concepts: multigrid 11655582bec1SHong Zhang 11665582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 116797177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 116897177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 116997177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 117097177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 11715582bec1SHong Zhang M*/ 11725582bec1SHong Zhang 11735582bec1SHong Zhang EXTERN_C_BEGIN 11745582bec1SHong Zhang #undef __FUNCT__ 11755582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 11767087cfbeSBarry Smith PetscErrorCode PCCreate_ML(PC pc) 11775582bec1SHong Zhang { 11785582bec1SHong Zhang PetscErrorCode ierr; 11795582bec1SHong Zhang PC_ML *pc_ml; 118001da6913SBarry Smith PC_MG *mg; 11815582bec1SHong Zhang 11825582bec1SHong Zhang PetscFunctionBegin; 1183573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 11845582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 118503bfa161SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr); 1186e0f5d30fSBarry Smith /* Since PCMG tries to use DM assocated with PC must delete it */ 1187e0f5d30fSBarry Smith ierr = DMDestroy(&pc->dm);CHKERRQ(ierr); 1188e0f5d30fSBarry Smith mg = (PC_MG*)pc->data; 1189c91913d3SJed Brown mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 11905582bec1SHong Zhang 11915582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 119238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr); 119301da6913SBarry Smith mg->innerctx = pc_ml; 11945582bec1SHong Zhang 1195573998d7SHong Zhang pc_ml->ml_object = 0; 1196573998d7SHong Zhang pc_ml->agg_object = 0; 1197573998d7SHong Zhang pc_ml->gridctx = 0; 1198573998d7SHong Zhang pc_ml->PetscMLdata = 0; 1199573998d7SHong Zhang pc_ml->Nlevels = -1; 1200573998d7SHong Zhang pc_ml->MaxNlevels = 10; 1201573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 12023751b4bdSBarry Smith pc_ml->CoarsenScheme = 1; 1203573998d7SHong Zhang pc_ml->Threshold = 0.0; 1204573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 1205573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 1206573998d7SHong Zhang pc_ml->size = 0; 120739381ba2SJed Brown pc_ml->dim = 0; 120839381ba2SJed Brown pc_ml->nloc = 0; 120939381ba2SJed Brown pc_ml->coords = 0; 121039381ba2SJed Brown pc_ml->Repartition = PETSC_FALSE; 121139381ba2SJed Brown pc_ml->MaxMinRatio = 1.3; 121239381ba2SJed Brown pc_ml->MinPerProc = 512; 121339381ba2SJed Brown pc_ml->PutOnSingleProc = 5000; 121439381ba2SJed Brown pc_ml->RepartitionType = 0; 121539381ba2SJed Brown pc_ml->ZoltanScheme = 0; 121639381ba2SJed Brown pc_ml->Aux = PETSC_FALSE; 121739381ba2SJed Brown pc_ml->AuxThreshold = 0.0; 121839381ba2SJed Brown 121939381ba2SJed Brown /* allow for coordinates to be passed */ 122039381ba2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","PCSetCoordinates_ML",PCSetCoordinates_ML);CHKERRQ(ierr); 1221573998d7SHong Zhang 12225582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 12235582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 12245582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 1225a06653b4SBarry Smith pc->ops->reset = PCReset_ML; 12265582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 12275582bec1SHong Zhang PetscFunctionReturn(0); 12285582bec1SHong Zhang } 12295582bec1SHong Zhang EXTERN_C_END 1230