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 */ 7af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 8af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 9c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 10c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 111e25c274SJed Brown #include <petscdm.h> /* for DMDestroy(&pc->mg) hack */ 12cb5d8e9eSHong Zhang 132cf39c26SSatish Balay EXTERN_C_BEGIN 1468210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */ 1568210224SSatish Balay #if !defined(HAVE_CONFIG_H) 1668210224SSatish Balay #define HAVE_CONFIG_H 1768210224SSatish Balay #endif 18c6db04a5SJed Brown #include <ml_include.h> 1939381ba2SJed Brown #include <ml_viz_stats.h> 205582bec1SHong Zhang EXTERN_C_END 215582bec1SHong Zhang 22fb6a8e6dSJed Brown typedef enum {PCML_NULLSPACE_AUTO,PCML_NULLSPACE_USER,PCML_NULLSPACE_BLOCK,PCML_NULLSPACE_SCALAR} PCMLNullSpaceType; 23fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO","USER","BLOCK","SCALAR","PCMLNullSpaceType","PCML_NULLSPACE_",0}; 24fb6a8e6dSJed Brown 255582bec1SHong Zhang /* The context (data structure) at each grid level */ 265582bec1SHong Zhang typedef struct { 275582bec1SHong Zhang Vec x,b,r; /* global vectors */ 285582bec1SHong Zhang Mat A,P,R; 295582bec1SHong Zhang KSP ksp; 3039381ba2SJed Brown Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */ 315582bec1SHong Zhang } GridCtx; 325582bec1SHong Zhang 335582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 345582bec1SHong Zhang typedef struct { 35573998d7SHong Zhang Mat A; /* Petsc matrix in aij format */ 36573998d7SHong Zhang Mat Aloc; /* local portion of A to be used by ML */ 3724a42b14SHong Zhang Vec x,y; 385582bec1SHong Zhang ML_Operator *mlmat; 395582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 405582bec1SHong Zhang } FineGridCtx; 415582bec1SHong Zhang 425582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 435582bec1SHong Zhang typedef struct { 445582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 455582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 465582bec1SHong Zhang } Mat_MLShell; 475582bec1SHong Zhang 485582bec1SHong Zhang /* Private context for the ML preconditioner */ 495582bec1SHong Zhang typedef struct { 505582bec1SHong Zhang ML *ml_object; 515582bec1SHong Zhang ML_Aggregate *agg_object; 525582bec1SHong Zhang GridCtx *gridctx; 535582bec1SHong Zhang FineGridCtx *PetscMLdata; 5439381ba2SJed Brown PetscInt Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization,MinPerProc,PutOnSingleProc,RepartitionType,ZoltanScheme; 5539381ba2SJed Brown PetscReal Threshold,DampingFactor,EnergyMinimizationDropTol,MaxMinRatio,AuxThreshold; 5639381ba2SJed Brown PetscBool SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable,Repartition,Aux; 5748268eb4SJed Brown PetscBool reuse_interpolation; 58fb6a8e6dSJed Brown PCMLNullSpaceType nulltype; 59573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 6039381ba2SJed Brown PetscInt dim; /* data from PCSetCoordinates(_ML) */ 6139381ba2SJed Brown PetscInt nloc; 6239381ba2SJed Brown PetscReal *coords; /* ML has a grid object for each level: the finest grid will point into coords */ 635582bec1SHong Zhang } PC_ML; 6441ca0015SHong Zhang 656562c4e1SBarry 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[]) 666562c4e1SBarry Smith { 676562c4e1SBarry Smith PetscInt m,i,j,k=0,row,*aj; 686562c4e1SBarry Smith PetscScalar *aa; 696562c4e1SBarry Smith FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 706562c4e1SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)ml->Aloc->data; 715582bec1SHong Zhang 72*5f80ce2aSJacob Faibussowitsch if (MatGetSize(ml->Aloc,&m,NULL)) return(0); 736562c4e1SBarry Smith for (i = 0; i<N_requested_rows; i++) { 746562c4e1SBarry Smith row = requested_rows[i]; 756562c4e1SBarry Smith row_lengths[i] = a->ilen[row]; 766562c4e1SBarry Smith if (allocated_space < k+row_lengths[i]) return(0); 776562c4e1SBarry Smith if ((row >= 0) || (row <= (m-1))) { 786562c4e1SBarry Smith aj = a->j + a->i[row]; 796562c4e1SBarry Smith aa = a->a + a->i[row]; 806562c4e1SBarry Smith for (j=0; j<row_lengths[i]; j++) { 816562c4e1SBarry Smith columns[k] = aj[j]; 826562c4e1SBarry Smith values[k++] = aa[j]; 836562c4e1SBarry Smith } 846562c4e1SBarry Smith } 856562c4e1SBarry Smith } 866562c4e1SBarry Smith return(1); 876562c4e1SBarry Smith } 886562c4e1SBarry Smith 896562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data) 906562c4e1SBarry Smith { 916562c4e1SBarry Smith FineGridCtx *ml = (FineGridCtx*)ML_data; 926562c4e1SBarry Smith Mat A = ml->A; 936562c4e1SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 946562c4e1SBarry Smith PetscMPIInt size; 956562c4e1SBarry Smith PetscInt i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n; 96d9ca1df4SBarry Smith const PetscScalar *array; 976562c4e1SBarry Smith 986562c4e1SBarry Smith PetscFunctionBegin; 99*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 100708418deSStefano Zampini if (size == 1) PetscFunctionReturn(0); 1016562c4e1SBarry Smith 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(ml->y,p)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(ml->y)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(a->lvec,&array)); 1072fa5cd67SKarl Rupp for (i=in_length; i<out_length; i++) p[i] = array[i-in_length]; 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(a->lvec,&array)); 1096562c4e1SBarry Smith PetscFunctionReturn(0); 1106562c4e1SBarry Smith } 1116562c4e1SBarry Smith 1126562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 1136562c4e1SBarry Smith { 1146562c4e1SBarry Smith FineGridCtx *ml = (FineGridCtx*)ML_Get_MyMatvecData(ML_data); 1156562c4e1SBarry Smith Mat A = ml->A, Aloc=ml->Aloc; 1166562c4e1SBarry Smith PetscMPIInt size; 1176562c4e1SBarry Smith PetscScalar *pwork=ml->pwork; 1186562c4e1SBarry Smith PetscInt i; 1196562c4e1SBarry Smith 1206562c4e1SBarry Smith PetscFunctionBegin; 121*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 1226562c4e1SBarry Smith if (size == 1) { 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(ml->x,p)); 1246562c4e1SBarry Smith } else { 1256562c4e1SBarry Smith for (i=0; i<in_length; i++) pwork[i] = p[i]; 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscML_comm(pwork,ml)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(ml->x,pwork)); 1286562c4e1SBarry Smith } 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(ml->y,ap)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Aloc,ml->x,ml->y)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(ml->x)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(ml->y)); 1336562c4e1SBarry Smith PetscFunctionReturn(0); 1346562c4e1SBarry Smith } 1356562c4e1SBarry Smith 1366562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 1376562c4e1SBarry Smith { 1386562c4e1SBarry Smith Mat_MLShell *shell; 139d9ca1df4SBarry Smith PetscScalar *yarray; 140d9ca1df4SBarry Smith const PetscScalar *xarray; 1416562c4e1SBarry Smith PetscInt x_length,y_length; 1426562c4e1SBarry Smith 1436562c4e1SBarry Smith PetscFunctionBegin; 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&shell)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xarray)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(y,&yarray)); 1476562c4e1SBarry Smith x_length = shell->mlmat->invec_leng; 1486562c4e1SBarry Smith y_length = shell->mlmat->outvec_leng; 149d9ca1df4SBarry Smith PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat,x_length,(PetscScalar*)xarray,y_length,yarray)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xarray)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(y,&yarray)); 1526562c4e1SBarry Smith PetscFunctionReturn(0); 1536562c4e1SBarry Smith } 1546562c4e1SBarry Smith 15579d04de1SBarry Smith /* newtype is ignored since only handles one case */ 1566562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 1576562c4e1SBarry Smith { 1586562c4e1SBarry Smith Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 1596562c4e1SBarry Smith Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 1606562c4e1SBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 161708418deSStefano Zampini PetscScalar *aa,*ba,*ca; 1626562c4e1SBarry Smith PetscInt am =A->rmap->n,an=A->cmap->n,i,j,k; 1636562c4e1SBarry Smith PetscInt *ci,*cj,ncols; 1646562c4e1SBarry Smith 1656562c4e1SBarry Smith PetscFunctionBegin; 166*5f80ce2aSJacob Faibussowitsch PetscCheck(am == an,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArrayRead(mpimat->A,(const PetscScalar**)&aa)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArrayRead(mpimat->B,(const PetscScalar**)&ba)); 1696562c4e1SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(1+am,&ci)); 1716562c4e1SBarry Smith ci[0] = 0; 1722fa5cd67SKarl Rupp for (i=0; i<am; i++) ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(1+ci[am],&cj)); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(1+ci[am],&ca)); 1756562c4e1SBarry Smith 1766562c4e1SBarry Smith k = 0; 1776562c4e1SBarry Smith for (i=0; i<am; i++) { 1786562c4e1SBarry Smith /* diagonal portion of A */ 1796562c4e1SBarry Smith ncols = ai[i+1] - ai[i]; 1806562c4e1SBarry Smith for (j=0; j<ncols; j++) { 1816562c4e1SBarry Smith cj[k] = *aj++; 1826562c4e1SBarry Smith ca[k++] = *aa++; 1836562c4e1SBarry Smith } 1846562c4e1SBarry Smith /* off-diagonal portion of A */ 1856562c4e1SBarry Smith ncols = bi[i+1] - bi[i]; 1866562c4e1SBarry Smith for (j=0; j<ncols; j++) { 1876562c4e1SBarry Smith cj[k] = an + (*bj); bj++; 1886562c4e1SBarry Smith ca[k++] = *ba++; 1896562c4e1SBarry Smith } 1906562c4e1SBarry Smith } 191*5f80ce2aSJacob Faibussowitsch PetscCheck(k == ci[am],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 1926562c4e1SBarry Smith 1936562c4e1SBarry Smith /* put together the new matrix */ 1946562c4e1SBarry Smith an = mpimat->A->cmap->n+mpimat->B->cmap->n; 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc)); 1966562c4e1SBarry Smith 1976562c4e1SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 1986562c4e1SBarry Smith /* Since these are PETSc arrays, change flags to free them as necessary. */ 1996562c4e1SBarry Smith mat = (Mat_SeqAIJ*)(*Aloc)->data; 2006562c4e1SBarry Smith mat->free_a = PETSC_TRUE; 2016562c4e1SBarry Smith mat->free_ij = PETSC_TRUE; 2026562c4e1SBarry Smith 2036562c4e1SBarry Smith mat->nonew = 0; 2046562c4e1SBarry Smith } else if (scall == MAT_REUSE_MATRIX) { 2056562c4e1SBarry Smith mat=(Mat_SeqAIJ*)(*Aloc)->data; 2066562c4e1SBarry Smith ci = mat->i; cj = mat->j; ca = mat->a; 2076562c4e1SBarry Smith for (i=0; i<am; i++) { 2086562c4e1SBarry Smith /* diagonal portion of A */ 2096562c4e1SBarry Smith ncols = ai[i+1] - ai[i]; 2106562c4e1SBarry Smith for (j=0; j<ncols; j++) *ca++ = *aa++; 2116562c4e1SBarry Smith /* off-diagonal portion of A */ 2126562c4e1SBarry Smith ncols = bi[i+1] - bi[i]; 2136562c4e1SBarry Smith for (j=0; j<ncols; j++) *ca++ = *ba++; 2146562c4e1SBarry Smith } 21598921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArrayRead(mpimat->A,(const PetscScalar**)&aa)); 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArrayRead(mpimat->B,(const PetscScalar**)&ba)); 2186562c4e1SBarry Smith PetscFunctionReturn(0); 2196562c4e1SBarry Smith } 2206562c4e1SBarry Smith 2216562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A) 2226562c4e1SBarry Smith { 2236562c4e1SBarry Smith Mat_MLShell *shell; 2246562c4e1SBarry Smith 2256562c4e1SBarry Smith PetscFunctionBegin; 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&shell)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(shell)); 2286562c4e1SBarry Smith PetscFunctionReturn(0); 2296562c4e1SBarry Smith } 2306562c4e1SBarry Smith 2316562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 2326562c4e1SBarry Smith { 2336562c4e1SBarry Smith struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata*)mlmat->data; 2340298fd71SBarry Smith PetscInt m =mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = NULL,nz_max; 23539381ba2SJed Brown PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i; 2366562c4e1SBarry Smith PetscScalar *ml_vals=matdata->values,*aa; 2376562c4e1SBarry Smith 2386562c4e1SBarry Smith PetscFunctionBegin; 239*5f80ce2aSJacob Faibussowitsch PetscCheck(mlmat->getrow,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 2406562c4e1SBarry Smith if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 2416562c4e1SBarry Smith if (reuse) { 2426562c4e1SBarry Smith Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 2436562c4e1SBarry Smith aij->i = ml_rowptr; 2446562c4e1SBarry Smith aij->j = ml_cols; 2456562c4e1SBarry Smith aij->a = ml_vals; 2466562c4e1SBarry Smith } else { 2476562c4e1SBarry Smith /* sort ml_cols and ml_vals */ 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m+1,&nnz)); 2492fa5cd67SKarl Rupp for (i=0; i<m; i++) nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 2506562c4e1SBarry Smith aj = ml_cols; aa = ml_vals; 2516562c4e1SBarry Smith for (i=0; i<m; i++) { 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortIntWithScalarArray(nnz[i],aj,aa)); 2536562c4e1SBarry Smith aj += nnz[i]; aa += nnz[i]; 2546562c4e1SBarry Smith } 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat)); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(nnz)); 2576562c4e1SBarry Smith } 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY)); 259*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY)); 2606562c4e1SBarry Smith PetscFunctionReturn(0); 2616562c4e1SBarry Smith } 2626562c4e1SBarry Smith 26339381ba2SJed Brown nz_max = PetscMax(1,mlmat->max_nz_per_row); 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nz_max,&aa,nz_max,&aj)); 26539381ba2SJed Brown if (!reuse) { 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,newmat)); 267*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE)); 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(*newmat,MATSEQAIJ)); 26939381ba2SJed Brown /* keep track of block size for A matrices */ 270*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize (*newmat, mlmat->num_PDEs)); 2716562c4e1SBarry Smith 272*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m,&nnz)); 2736562c4e1SBarry Smith for (i=0; i<m; i++) { 274815d23e5SBarry Smith PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i])); 2756562c4e1SBarry Smith } 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(*newmat,0,nnz)); 277ae7fe62dSJed Brown } 2786562c4e1SBarry Smith for (i=0; i<m; i++) { 279ae7fe62dSJed Brown PetscInt ncols; 28039381ba2SJed Brown 281815d23e5SBarry Smith PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols)); 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES)); 2836562c4e1SBarry Smith } 284*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY)); 285*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY)); 2866562c4e1SBarry Smith 287*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(aa,aj)); 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(nnz)); 2896562c4e1SBarry Smith PetscFunctionReturn(0); 2906562c4e1SBarry Smith } 2916562c4e1SBarry Smith 2926562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 2936562c4e1SBarry Smith { 2946562c4e1SBarry Smith PetscInt m,n; 2956562c4e1SBarry Smith ML_Comm *MLcomm; 2966562c4e1SBarry Smith Mat_MLShell *shellctx; 2976562c4e1SBarry Smith 2986562c4e1SBarry Smith PetscFunctionBegin; 2996562c4e1SBarry Smith m = mlmat->outvec_leng; 3006562c4e1SBarry Smith n = mlmat->invec_leng; 3016562c4e1SBarry Smith 3026562c4e1SBarry Smith if (reuse) { 303*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(*newmat,&shellctx)); 3046562c4e1SBarry Smith shellctx->mlmat = mlmat; 3056562c4e1SBarry Smith PetscFunctionReturn(0); 3066562c4e1SBarry Smith } 3076562c4e1SBarry Smith 3086562c4e1SBarry Smith MLcomm = mlmat->comm; 3092fa5cd67SKarl Rupp 310*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&shellctx)); 311*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat)); 312*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML)); 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*newmat,MATOP_DESTROY,(void(*)(void))MatDestroy_ML)); 3142fa5cd67SKarl Rupp 3156562c4e1SBarry Smith shellctx->A = *newmat; 3166562c4e1SBarry Smith shellctx->mlmat = mlmat; 3176562c4e1SBarry Smith PetscFunctionReturn(0); 3186562c4e1SBarry Smith } 3196562c4e1SBarry Smith 320ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 3216562c4e1SBarry Smith { 32239381ba2SJed Brown PetscInt *aj; 32339381ba2SJed Brown PetscScalar *aa; 32439381ba2SJed Brown PetscInt i,j,*gordering; 325ae7fe62dSJed Brown PetscInt m=mlmat->outvec_leng,n,nz_max,row; 3266562c4e1SBarry Smith Mat A; 3276562c4e1SBarry Smith 3286562c4e1SBarry Smith PetscFunctionBegin; 329*5f80ce2aSJacob Faibussowitsch PetscCheck(mlmat->getrow,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 3306562c4e1SBarry Smith n = mlmat->invec_leng; 331*5f80ce2aSJacob Faibussowitsch PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 3326562c4e1SBarry Smith 3337be6b909SBarry Smith /* create global row numbering for a ML_Operator */ 3347be6b909SBarry Smith PetscStackCall("ML_build_global_numbering",ML_build_global_numbering(mlmat,&gordering,"rows")); 3357be6b909SBarry Smith 3361d94bf15SBarry Smith nz_max = PetscMax(1,mlmat->max_nz_per_row) + 1; 337*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nz_max,&aa,nz_max,&aj)); 3387be6b909SBarry Smith if (reuse) { 3397be6b909SBarry Smith A = *newmat; 3407be6b909SBarry Smith } else { 341ae7fe62dSJed Brown PetscInt *nnzA,*nnzB,*nnz; 3427be6b909SBarry Smith PetscInt rstart; 343*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(mlmat->comm->USR_comm,&A)); 344*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); 345*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATMPIAIJ)); 34639381ba2SJed Brown /* keep track of block size for A matrices */ 347*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize (A,mlmat->num_PDEs)); 348*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(m,&nnzA,m,&nnzB,m,&nnz)); 349*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Scan(&m,&rstart,1,MPIU_INT,MPI_SUM,mlmat->comm->USR_comm)); 3507be6b909SBarry Smith rstart -= m; 3516562c4e1SBarry Smith 3526562c4e1SBarry Smith for (i=0; i<m; i++) { 3537be6b909SBarry Smith row = gordering[i] - rstart; 354815d23e5SBarry Smith PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i])); 3557be6b909SBarry Smith nnzA[row] = 0; 35639381ba2SJed Brown for (j=0; j<nnz[i]; j++) { 3577be6b909SBarry Smith if (aj[j] < m) nnzA[row]++; 3586562c4e1SBarry Smith } 3597be6b909SBarry Smith nnzB[row] = nnz[i] - nnzA[row]; 3606562c4e1SBarry Smith } 361*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB)); 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(nnzA,nnzB,nnz)); 363ae7fe62dSJed Brown } 3646562c4e1SBarry Smith for (i=0; i<m; i++) { 365ae7fe62dSJed Brown PetscInt ncols; 3666562c4e1SBarry Smith row = gordering[i]; 36739381ba2SJed Brown 368815d23e5SBarry Smith PetscStackCall(",ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols)); 3692fa5cd67SKarl Rupp for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]]; 370*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES)); 3716562c4e1SBarry Smith } 3727be6b909SBarry Smith PetscStackCall("ML_free",ML_free(gordering)); 373*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 3756562c4e1SBarry Smith *newmat = A; 3766562c4e1SBarry Smith 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(aa,aj)); 3786562c4e1SBarry Smith PetscFunctionReturn(0); 3796562c4e1SBarry Smith } 3806562c4e1SBarry Smith 38139381ba2SJed Brown /* -------------------------------------------------------------------------- */ 38239381ba2SJed Brown /* 38339381ba2SJed Brown PCSetCoordinates_ML 38439381ba2SJed Brown 38539381ba2SJed Brown Input Parameter: 38639381ba2SJed Brown . pc - the preconditioner context 38739381ba2SJed Brown */ 388f7a08781SBarry Smith static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 38939381ba2SJed Brown { 39039381ba2SJed Brown PC_MG *mg = (PC_MG*)pc->data; 39139381ba2SJed Brown PC_ML *pc_ml = (PC_ML*)mg->innerctx; 39290fbc344SStefano Zampini PetscInt arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend,aloc; 39339381ba2SJed Brown Mat Amat = pc->pmat; 39439381ba2SJed Brown 39539381ba2SJed Brown /* this function copied and modified from PCSetCoordinates_GEO -TGI */ 39639381ba2SJed Brown PetscFunctionBegin; 39739381ba2SJed Brown PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1); 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(Amat, &bs)); 39939381ba2SJed Brown 400*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(Amat, &my0, &Iend)); 40190fbc344SStefano Zampini aloc = (Iend-my0); 40239381ba2SJed Brown nloc = (Iend-my0)/bs; 40339381ba2SJed Brown 404*5f80ce2aSJacob Faibussowitsch PetscCheck((nloc == a_nloc) || (aloc == a_nloc),PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of local blocks %D must be %D or %D.",a_nloc,nloc,aloc); 40539381ba2SJed Brown 40639381ba2SJed Brown oldarrsz = pc_ml->dim * pc_ml->nloc; 40739381ba2SJed Brown pc_ml->dim = ndm; 40890fbc344SStefano Zampini pc_ml->nloc = nloc; 40990fbc344SStefano Zampini arrsz = ndm * nloc; 41039381ba2SJed Brown 41139381ba2SJed Brown /* create data - syntactic sugar that should be refactored at some point */ 41239381ba2SJed Brown if (pc_ml->coords==0 || (oldarrsz != arrsz)) { 413*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml->coords)); 414*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(arrsz, &pc_ml->coords)); 41539381ba2SJed Brown } 41639381ba2SJed Brown for (kk=0; kk<arrsz; kk++) pc_ml->coords[kk] = -999.; 41739381ba2SJed Brown /* copy data in - column oriented */ 41890fbc344SStefano Zampini if (nloc == a_nloc) { 41939381ba2SJed Brown for (kk = 0; kk < nloc; kk++) { 42039381ba2SJed Brown for (ii = 0; ii < ndm; ii++) { 42139381ba2SJed Brown pc_ml->coords[ii*nloc + kk] = coords[kk*ndm + ii]; 42239381ba2SJed Brown } 42339381ba2SJed Brown } 42490fbc344SStefano Zampini } else { /* assumes the coordinates are blocked */ 42590fbc344SStefano Zampini for (kk = 0; kk < nloc; kk++) { 42690fbc344SStefano Zampini for (ii = 0; ii < ndm; ii++) { 42790fbc344SStefano Zampini pc_ml->coords[ii*nloc + kk] = coords[bs*kk*ndm + ii]; 42890fbc344SStefano Zampini } 42990fbc344SStefano Zampini } 43090fbc344SStefano Zampini } 43139381ba2SJed Brown PetscFunctionReturn(0); 43239381ba2SJed Brown } 43339381ba2SJed Brown 4346562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/ 435e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC); 43616336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc) 43701da6913SBarry Smith { 438e0262f48SMatthew G Knepley PC_MG *mg = (PC_MG*)pc->data; 439e0262f48SMatthew G Knepley PC_ML *pc_ml = (PC_ML*)mg->innerctx; 44039381ba2SJed Brown PetscInt level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim; 44101da6913SBarry Smith 44201da6913SBarry Smith PetscFunctionBegin; 44339381ba2SJed Brown if (dim) { 44439381ba2SJed Brown for (level=0; level<=fine_level; level++) { 445*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&pc_ml->gridctx[level].coords)); 44639381ba2SJed Brown } 447448f31a9SStefano Zampini if (pc_ml->ml_object && pc_ml->ml_object->Grid) { 448448f31a9SStefano Zampini ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats*) pc_ml->ml_object->Grid[0].Grid; 44939381ba2SJed Brown grid_info->x = 0; /* do this so ML doesn't try to free coordinates */ 45039381ba2SJed Brown grid_info->y = 0; 45139381ba2SJed Brown grid_info->z = 0; 452815d23e5SBarry Smith PetscStackCall("ML_Operator_Getrow",ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object)); 45339381ba2SJed Brown } 454448f31a9SStefano Zampini } 455815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Destroy",ML_Aggregate_Destroy(&pc_ml->agg_object)); 456815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Destroy",ML_Destroy(&pc_ml->ml_object)); 45701da6913SBarry Smith 45801da6913SBarry Smith if (pc_ml->PetscMLdata) { 459*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml->PetscMLdata->pwork)); 460*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&pc_ml->PetscMLdata->Aloc)); 461*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&pc_ml->PetscMLdata->x)); 462*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&pc_ml->PetscMLdata->y)); 46301da6913SBarry Smith } 464*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml->PetscMLdata)); 46501da6913SBarry Smith 466f5a5dd59SJed Brown if (pc_ml->gridctx) { 46701da6913SBarry Smith for (level=0; level<fine_level; level++) { 468*5f80ce2aSJacob Faibussowitsch if (pc_ml->gridctx[level].A) CHKERRQ(MatDestroy(&pc_ml->gridctx[level].A)); 469*5f80ce2aSJacob Faibussowitsch if (pc_ml->gridctx[level].P) CHKERRQ(MatDestroy(&pc_ml->gridctx[level].P)); 470*5f80ce2aSJacob Faibussowitsch if (pc_ml->gridctx[level].R) CHKERRQ(MatDestroy(&pc_ml->gridctx[level].R)); 471*5f80ce2aSJacob Faibussowitsch if (pc_ml->gridctx[level].x) CHKERRQ(VecDestroy(&pc_ml->gridctx[level].x)); 472*5f80ce2aSJacob Faibussowitsch if (pc_ml->gridctx[level].b) CHKERRQ(VecDestroy(&pc_ml->gridctx[level].b)); 473*5f80ce2aSJacob Faibussowitsch if (pc_ml->gridctx[level+1].r) CHKERRQ(VecDestroy(&pc_ml->gridctx[level+1].r)); 47401da6913SBarry Smith } 475f5a5dd59SJed Brown } 476*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml->gridctx)); 477*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml->coords)); 4782fa5cd67SKarl Rupp 47939381ba2SJed Brown pc_ml->dim = 0; 48039381ba2SJed Brown pc_ml->nloc = 0; 481*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_MG(pc)); 48201da6913SBarry Smith PetscFunctionReturn(0); 48301da6913SBarry Smith } 4845582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4855582bec1SHong Zhang /* 4865582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 4875582bec1SHong Zhang by setting data structures and options. 4885582bec1SHong Zhang 4895582bec1SHong Zhang Input Parameter: 4905582bec1SHong Zhang . pc - the preconditioner context 4915582bec1SHong Zhang 4925582bec1SHong Zhang Application Interface Routine: PCSetUp() 4935582bec1SHong Zhang 4945582bec1SHong Zhang Notes: 4955582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 4965582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 4975582bec1SHong Zhang */ 4984416b707SBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC); 499a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC); 500c07bf074SBarry Smith 5016ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 5025582bec1SHong Zhang { 5035582bec1SHong Zhang PetscErrorCode ierr; 504eef31507SHong Zhang PetscMPIInt size; 5055582bec1SHong Zhang FineGridCtx *PetscMLdata; 5065582bec1SHong Zhang ML *ml_object; 5075582bec1SHong Zhang ML_Aggregate *agg_object; 5085582bec1SHong Zhang ML_Operator *mlmat; 5094f8eab3cSJed Brown PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs; 5105582bec1SHong Zhang Mat A,Aloc; 5115582bec1SHong Zhang GridCtx *gridctx; 51201da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 51301da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 514ace3abfcSBarry Smith PetscBool isSeq, isMPI; 515c07bf074SBarry Smith KSP smoother; 516c07bf074SBarry Smith PC subpc; 51748268eb4SJed Brown PetscInt mesh_level, old_mesh_level; 5188a62b701SToby Isaac MatInfo info; 5191f817a21SBarry Smith static PetscBool cite = PETSC_FALSE; 52048268eb4SJed Brown 5215582bec1SHong Zhang PetscFunctionBegin; 522*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCitationsRegister("@TechReport{ml_users_guide,\n author = {M. Sala and J.J. Hu and R.S. Tuminaro},\n title = {{ML}3.1 {S}moothed {A}ggregation {U}ser's {G}uide},\n institution = {Sandia National Laboratories},\n number = {SAND2004-4821},\n year = 2004\n}\n",&cite)); 52348268eb4SJed Brown A = pc->pmat; 524*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 52548268eb4SJed Brown 526573998d7SHong Zhang if (pc->setupcalled) { 52748268eb4SJed Brown if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) { 52848268eb4SJed Brown /* 52948268eb4SJed Brown Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for 53048268eb4SJed Brown multiple solves in which the matrix is not changing too quickly. 53148268eb4SJed Brown */ 53248268eb4SJed Brown ml_object = pc_ml->ml_object; 53348268eb4SJed Brown gridctx = pc_ml->gridctx; 53448268eb4SJed Brown Nlevels = pc_ml->Nlevels; 53548268eb4SJed Brown fine_level = Nlevels - 1; 53648268eb4SJed Brown gridctx[fine_level].A = A; 53748268eb4SJed Brown 538*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq)); 539*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI)); 54048268eb4SJed Brown if (isMPI) { 541*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc)); 54248268eb4SJed Brown } else if (isSeq) { 54348268eb4SJed Brown Aloc = A; 544*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)Aloc)); 54598921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name); 54648268eb4SJed Brown 547*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(Aloc,&m,&nlocal_allcols)); 54848268eb4SJed Brown PetscMLdata = pc_ml->PetscMLdata; 549*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&PetscMLdata->Aloc)); 55048268eb4SJed Brown PetscMLdata->A = A; 55148268eb4SJed Brown PetscMLdata->Aloc = Aloc; 552815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Destroy",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata)); 553815d23e5SBarry Smith PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec)); 55448268eb4SJed Brown 55548268eb4SJed Brown mesh_level = ml_object->ML_finest_level; 55648268eb4SJed Brown while (ml_object->SingleLevel[mesh_level].Rmat->to) { 55748268eb4SJed Brown old_mesh_level = mesh_level; 55848268eb4SJed Brown mesh_level = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum; 55948268eb4SJed Brown 56048268eb4SJed Brown /* clean and regenerate A */ 56148268eb4SJed Brown mlmat = &(ml_object->Amat[mesh_level]); 562815d23e5SBarry Smith PetscStackCall("ML_Operator_Clean",ML_Operator_Clean(mlmat)); 563815d23e5SBarry Smith PetscStackCall("ML_Operator_Init",ML_Operator_Init(mlmat,ml_object->comm)); 564815d23e5SBarry Smith PetscStackCall("ML_Gen_AmatrixRAP",ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level)); 56548268eb4SJed Brown } 56648268eb4SJed Brown 56748268eb4SJed Brown level = fine_level - 1; 56848268eb4SJed Brown if (size == 1) { /* convert ML P, R and A into seqaij format */ 56948268eb4SJed Brown for (mllevel=1; mllevel<Nlevels; mllevel++) { 57048268eb4SJed Brown mlmat = &(ml_object->Amat[mllevel]); 571*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A)); 57248268eb4SJed Brown level--; 57348268eb4SJed Brown } 57448268eb4SJed Brown } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 57548268eb4SJed Brown for (mllevel=1; mllevel<Nlevels; mllevel++) { 57648268eb4SJed Brown mlmat = &(ml_object->Amat[mllevel]); 577*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A)); 57848268eb4SJed Brown level--; 57948268eb4SJed Brown } 58048268eb4SJed Brown } 58148268eb4SJed Brown 58248268eb4SJed Brown for (level=0; level<fine_level; level++) { 58348268eb4SJed Brown if (level > 0) { 584*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A)); 58548268eb4SJed Brown } 586*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A)); 58748268eb4SJed Brown } 588*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A)); 589*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A)); 59048268eb4SJed Brown 591*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetUp_MG(pc)); 59248268eb4SJed Brown PetscFunctionReturn(0); 59348268eb4SJed Brown } else { 594c07bf074SBarry Smith /* since ML can change the size of vectors/matrices at any level we must destroy everything */ 595*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_ML(pc)); 596573998d7SHong Zhang } 59748268eb4SJed Brown } 598573998d7SHong Zhang 5995582bec1SHong Zhang /* setup special features of PCML */ 6005582bec1SHong Zhang /*--------------------------------*/ 6015582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 6025582bec1SHong Zhang pc_ml->size = size; 603*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq)); 604*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI)); 605864b637dSMatthew Knepley if (isMPI) { 606*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc)); 607864b637dSMatthew Knepley } else if (isSeq) { 6085582bec1SHong Zhang Aloc = A; 609*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)Aloc)); 61098921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name); 6115582bec1SHong Zhang 6125582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 613*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&PetscMLdata)); 6145582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 615*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Aloc->cmap->n+1,&PetscMLdata->pwork)); 6165582bec1SHong Zhang 617*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(Aloc,&PetscMLdata->x,&PetscMLdata->y)); 61824a42b14SHong Zhang 619573998d7SHong Zhang PetscMLdata->A = A; 620573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 62139381ba2SJed Brown if (pc_ml->dim) { /* create vecs around the coordinate data given */ 62239381ba2SJed Brown PetscInt i,j,dim=pc_ml->dim; 62339381ba2SJed Brown PetscInt nloc = pc_ml->nloc,nlocghost; 62439381ba2SJed Brown PetscReal *ghostedcoords; 62539381ba2SJed Brown 626*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 62739381ba2SJed Brown nlocghost = Aloc->cmap->n / bs; 628*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(dim*nlocghost,&ghostedcoords)); 62939381ba2SJed Brown for (i = 0; i < dim; i++) { 63039381ba2SJed Brown /* copy coordinate values into first component of pwork */ 63139381ba2SJed Brown for (j = 0; j < nloc; j++) { 63239381ba2SJed Brown PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j]; 63339381ba2SJed Brown } 63439381ba2SJed Brown /* get the ghost values */ 635*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscML_comm(PetscMLdata->pwork,PetscMLdata)); 63639381ba2SJed Brown /* write into the vector */ 63739381ba2SJed Brown for (j = 0; j < nlocghost; j++) { 63839381ba2SJed Brown ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j]; 63939381ba2SJed Brown } 64039381ba2SJed Brown } 64139381ba2SJed Brown /* replace the original coords with the ghosted coords, because these are 64239381ba2SJed Brown * what ML needs */ 643*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml->coords)); 64439381ba2SJed Brown pc_ml->coords = ghostedcoords; 64539381ba2SJed Brown } 64624a42b14SHong Zhang 6475582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 64845cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 649*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(Aloc,&m,&nlocal_allcols)); 650*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 651815d23e5SBarry Smith PetscStackCall("ML_Create",ML_Create(&ml_object,pc_ml->MaxNlevels)); 652ce94432eSBarry Smith PetscStackCall("ML_Comm_Set_UsrComm",ML_Comm_Set_UsrComm(ml_object->comm,PetscObjectComm((PetscObject)A))); 653573998d7SHong Zhang pc_ml->ml_object = ml_object; 654815d23e5SBarry Smith PetscStackCall("ML_Init_Amatrix",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata)); 655815d23e5SBarry Smith PetscStackCall("ML_Set_Amatrix_Getrow",ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols)); 656815d23e5SBarry Smith PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec)); 6575582bec1SHong Zhang 658815d23e5SBarry Smith PetscStackCall("ML_Set_Symmetrize",ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO)); 659b5c8bdf8SJed Brown 6605582bec1SHong Zhang /* aggregation */ 661815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Create",ML_Aggregate_Create(&agg_object)); 662573998d7SHong Zhang pc_ml->agg_object = agg_object; 663573998d7SHong Zhang 664fb6a8e6dSJed Brown { 665fb6a8e6dSJed Brown MatNullSpace mnull; 666*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetNearNullSpace(A,&mnull)); 667fb6a8e6dSJed Brown if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) { 668fb6a8e6dSJed Brown if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER; 669fb6a8e6dSJed Brown else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK; 670fb6a8e6dSJed Brown else pc_ml->nulltype = PCML_NULLSPACE_SCALAR; 671fb6a8e6dSJed Brown } 672fb6a8e6dSJed Brown switch (pc_ml->nulltype) { 673fb6a8e6dSJed Brown case PCML_NULLSPACE_USER: { 674fb6a8e6dSJed Brown PetscScalar *nullvec; 675fb6a8e6dSJed Brown const PetscScalar *v; 676fb6a8e6dSJed Brown PetscBool has_const; 6771c547e14SJed Brown PetscInt i,j,mlocal,nvec,M; 678fb6a8e6dSJed Brown const Vec *vecs; 6792fa5cd67SKarl Rupp 680*5f80ce2aSJacob Faibussowitsch PetscCheck(mnull,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space"); 681*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,NULL)); 682*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(Aloc,&mlocal,NULL)); 683*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs)); 684*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec)); 6851c547e14SJed Brown if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M; 686fb6a8e6dSJed Brown for (i=0; i<nvec; i++) { 687*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(vecs[i],&v)); 688fb6a8e6dSJed Brown for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j]; 689*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(vecs[i],&v)); 690fb6a8e6dSJed Brown } 691*5f80ce2aSJacob Faibussowitsch PetscStackCall("ML_Aggregate_Create",CHKERRQ(ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal))); 692*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(nullvec)); 693fb6a8e6dSJed Brown } break; 694fb6a8e6dSJed Brown case PCML_NULLSPACE_BLOCK: 695*5f80ce2aSJacob Faibussowitsch PetscStackCall("ML_Aggregate_Set_NullSpace",CHKERRQ(ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0))); 696fb6a8e6dSJed Brown break; 697fb6a8e6dSJed Brown case PCML_NULLSPACE_SCALAR: 698fb6a8e6dSJed Brown break; 699ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unknown null space type"); 700fb6a8e6dSJed Brown } 701fb6a8e6dSJed Brown } 702815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_MaxCoarseSize",ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize)); 7035582bec1SHong Zhang /* set options */ 7045582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 7055582bec1SHong Zhang case 1: 706815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_CoarsenScheme_Coupled",ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));break; 7075582bec1SHong Zhang case 2: 708815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_CoarsenScheme_MIS",ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));break; 7095582bec1SHong Zhang case 3: 710815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_CoarsenScheme_METIS",ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));break; 7115582bec1SHong Zhang } 712815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_Threshold",ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold)); 713815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_DampingFactor",ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor)); 7145582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm) { 715815d23e5SBarry Smith PetscStackCall("ML_Set_SpectralNormScheme_Anorm",ML_Set_SpectralNormScheme_Anorm(ml_object)); 7165582bec1SHong Zhang } 717b5c8bdf8SJed Brown agg_object->keep_agg_information = (int)pc_ml->KeepAggInfo; 718b5c8bdf8SJed Brown agg_object->keep_P_tentative = (int)pc_ml->Reusable; 719b5c8bdf8SJed Brown agg_object->block_scaled_SA = (int)pc_ml->BlockScaling; 720b5c8bdf8SJed Brown agg_object->minimizing_energy = (int)pc_ml->EnergyMinimization; 721b5c8bdf8SJed Brown agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol; 722b5c8bdf8SJed Brown agg_object->cheap_minimizing_energy = (int)pc_ml->EnergyMinimizationCheap; 7235582bec1SHong Zhang 72439381ba2SJed Brown if (pc_ml->Aux) { 725*5f80ce2aSJacob Faibussowitsch PetscCheck(pc_ml->dim,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Auxiliary matrix requires coordinates"); 72639381ba2SJed Brown ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold; 72739381ba2SJed Brown ml_object->Amat[0].aux_data->enable = 1; 72839381ba2SJed Brown ml_object->Amat[0].aux_data->max_level = 10; 72939381ba2SJed Brown ml_object->Amat[0].num_PDEs = bs; 73039381ba2SJed Brown } 73139381ba2SJed Brown 732*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(A,MAT_LOCAL,&info)); 7338a62b701SToby Isaac ml_object->Amat[0].N_nonzeros = (int) info.nz_used; 7348a62b701SToby Isaac 73539381ba2SJed Brown if (pc_ml->dim) { 73639381ba2SJed Brown PetscInt i,dim = pc_ml->dim; 73739381ba2SJed Brown ML_Aggregate_Viz_Stats *grid_info; 73839381ba2SJed Brown PetscInt nlocghost; 73939381ba2SJed Brown 740*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 74139381ba2SJed Brown nlocghost = Aloc->cmap->n / bs; 74239381ba2SJed Brown 743815d23e5SBarry Smith PetscStackCall("ML_Aggregate_VizAndStats_Setup(",ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */ 74439381ba2SJed Brown grid_info = (ML_Aggregate_Viz_Stats*) ml_object->Grid[0].Grid; 74539381ba2SJed Brown for (i = 0; i < dim; i++) { 74639381ba2SJed Brown /* set the finest level coordinates to point to the column-order array 74739381ba2SJed Brown * in pc_ml */ 74839381ba2SJed Brown /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */ 74939381ba2SJed Brown switch (i) { 75039381ba2SJed Brown case 0: grid_info->x = pc_ml->coords + nlocghost * i; break; 75139381ba2SJed Brown case 1: grid_info->y = pc_ml->coords + nlocghost * i; break; 75239381ba2SJed Brown case 2: grid_info->z = pc_ml->coords + nlocghost * i; break; 753ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3"); 75439381ba2SJed Brown } 75539381ba2SJed Brown } 75639381ba2SJed Brown grid_info->Ndim = dim; 75739381ba2SJed Brown } 75839381ba2SJed Brown 75939381ba2SJed Brown /* repartitioning */ 76039381ba2SJed Brown if (pc_ml->Repartition) { 761815d23e5SBarry Smith PetscStackCall("ML_Repartition_Activate",ML_Repartition_Activate(ml_object)); 762815d23e5SBarry Smith PetscStackCall("ML_Repartition_Set_LargestMinMaxRatio",ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio)); 763815d23e5SBarry Smith PetscStackCall("ML_Repartition_Set_MinPerProc",ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc)); 764815d23e5SBarry Smith PetscStackCall("ML_Repartition_Set_PutOnSingleProc",ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc)); 76539381ba2SJed Brown #if 0 /* Function not yet defined in ml-6.2 */ 76639381ba2SJed Brown /* I'm not sure what compatibility issues might crop up if we partitioned 76739381ba2SJed Brown * on the finest level, so to be safe repartition starts on the next 76839381ba2SJed Brown * finest level (reflection default behavior in 76939381ba2SJed Brown * ml_MultiLevelPreconditioner) */ 770815d23e5SBarry Smith PetscStackCall("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1)); 77139381ba2SJed Brown #endif 77239381ba2SJed Brown 77339381ba2SJed Brown if (!pc_ml->RepartitionType) { 77439381ba2SJed Brown PetscInt i; 77539381ba2SJed Brown 776*5f80ce2aSJacob Faibussowitsch PetscCheck(pc_ml->dim,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates"); 777815d23e5SBarry Smith PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN)); 778815d23e5SBarry Smith PetscStackCall("ML_Aggregate_Set_Dimensions",ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim)); 77939381ba2SJed Brown 78039381ba2SJed Brown for (i = 0; i < ml_object->ML_num_levels; i++) { 78139381ba2SJed Brown ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Grid[i].Grid; 78239381ba2SJed Brown grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */ 78339381ba2SJed Brown /* defaults from ml_agg_info.c */ 78439381ba2SJed Brown grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */ 78539381ba2SJed Brown grid_info->zoltan_timers = 0; 78639381ba2SJed Brown grid_info->smoothing_steps = 4; /* only relevant to hypergraph / fast hypergraph */ 78739381ba2SJed Brown } 7882fa5cd67SKarl Rupp } else { 789815d23e5SBarry Smith PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS)); 79039381ba2SJed Brown } 79139381ba2SJed Brown } 79239381ba2SJed Brown 793b5c8bdf8SJed Brown if (pc_ml->OldHierarchy) { 794815d23e5SBarry Smith PetscStackCall("ML_Gen_MGHierarchy_UsingAggregation",Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object)); 795b5c8bdf8SJed Brown } else { 796815d23e5SBarry Smith PetscStackCall("ML_Gen_MultiLevelHierarchy_UsingAggregation",Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object)); 797b5c8bdf8SJed Brown } 798*5f80ce2aSJacob Faibussowitsch PetscCheck(Nlevels>0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 799573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 800aa85bbbfSHong Zhang fine_level = Nlevels - 1; 801c07bf074SBarry Smith 802*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetLevels(pc,Nlevels,NULL)); 803aa85bbbfSHong Zhang /* set default smoothers */ 804aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++) { 805*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetSmoother(pc,level,&smoother)); 806*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(smoother,KSPRICHARDSON)); 807*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(smoother,&subpc)); 808*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(subpc,PCSOR)); 809aa85bbbfSHong Zhang } 810f2e59741SMatthew G Knepley ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 811*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetFromOptions_MG(PetscOptionsObject,pc)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 812f2e59741SMatthew G Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 8135582bec1SHong Zhang 814*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nlevels,&gridctx)); 8152fa5cd67SKarl Rupp 8165582bec1SHong Zhang pc_ml->gridctx = gridctx; 8175582bec1SHong Zhang 8185582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 8195582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 820e14861a4SHong Zhang gridctx[fine_level].A = A; 821573998d7SHong Zhang 822e14861a4SHong Zhang level = fine_level - 1; 82343ef1857SStefano Zampini /* TODO: support for GPUs */ 824ab718edeSHong Zhang if (size == 1) { /* convert ML P, R and A into seqaij format */ 8255582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++) { 826e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 827*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P)); 828e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 829*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R)); 830573998d7SHong Zhang 831573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 832*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A)); 8335582bec1SHong Zhang level--; 8345582bec1SHong Zhang } 835ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 8365582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++) { 8375582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 838*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P)); 839ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 840*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R)); 841573998d7SHong Zhang 8425582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 843*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A)); 8445582bec1SHong Zhang level--; 8455582bec1SHong Zhang } 8465582bec1SHong Zhang } 8475582bec1SHong Zhang 848573998d7SHong Zhang /* create vectors and ksp at all levels */ 849ac346b81SHong Zhang for (level=0; level<fine_level; level++) { 850573998d7SHong Zhang level1 = level + 1; 851708418deSStefano Zampini 852*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(gridctx[level].A,&gridctx[level].x,&gridctx[level].b)); 853*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(gridctx[level1].A,NULL,&gridctx[level1].r)); 854*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetX(pc,level,gridctx[level].x)); 855*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetRhs(pc,level,gridctx[level].b)); 856*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetR(pc,level1,gridctx[level1].r)); 857ac346b81SHong Zhang 8585582bec1SHong Zhang if (level == 0) { 859*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetCoarseSolve(pc,&gridctx[level].ksp)); 8605582bec1SHong Zhang } else { 861*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetSmoother(pc,level,&gridctx[level].ksp)); 862573998d7SHong Zhang } 863573998d7SHong Zhang } 864*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp)); 865573998d7SHong Zhang 866573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 867573998d7SHong Zhang for (level=0; level<fine_level; level++) { 868573998d7SHong Zhang level1 = level + 1; 869708418deSStefano Zampini 870*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetInterpolation(pc,level1,gridctx[level].P)); 871*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetRestriction(pc,level1,gridctx[level].R)); 872573998d7SHong Zhang if (level > 0) { 873*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A)); 8745582bec1SHong Zhang } 875*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A)); 8765582bec1SHong Zhang } 877*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A)); 878*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A)); 8795582bec1SHong Zhang 88039381ba2SJed Brown /* put coordinate info in levels */ 88139381ba2SJed Brown if (pc_ml->dim) { 88239381ba2SJed Brown PetscInt i,j,dim = pc_ml->dim; 88339381ba2SJed Brown PetscInt bs, nloc; 88439381ba2SJed Brown PC subpc; 88539381ba2SJed Brown PetscReal *array; 88639381ba2SJed Brown 88739381ba2SJed Brown level = fine_level; 88839381ba2SJed Brown for (mllevel = 0; mllevel < Nlevels; mllevel++) { 889ebbbbe33SJed Brown ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Amat[mllevel].to->Grid->Grid; 89039381ba2SJed Brown MPI_Comm comm = ((PetscObject)gridctx[level].A)->comm; 89139381ba2SJed Brown 892*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize (gridctx[level].A, &bs)); 893*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize (gridctx[level].A, NULL, &nloc)); 89439381ba2SJed Brown nloc /= bs; /* number of local nodes */ 89539381ba2SJed Brown 896*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(comm,&gridctx[level].coords)); 897*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE)); 898*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetType(gridctx[level].coords,VECMPI)); 899*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(gridctx[level].coords,&array)); 90039381ba2SJed Brown for (j = 0; j < nloc; j++) { 90139381ba2SJed Brown for (i = 0; i < dim; i++) { 90239381ba2SJed Brown switch (i) { 90339381ba2SJed Brown case 0: array[dim * j + i] = grid_info->x[j]; break; 90439381ba2SJed Brown case 1: array[dim * j + i] = grid_info->y[j]; break; 90539381ba2SJed Brown case 2: array[dim * j + i] = grid_info->z[j]; break; 906ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3"); 90739381ba2SJed Brown } 90839381ba2SJed Brown } 90939381ba2SJed Brown } 91039381ba2SJed Brown 91139381ba2SJed Brown /* passing coordinates to smoothers/coarse solver, should they need them */ 912*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(gridctx[level].ksp,&subpc)); 913*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetCoordinates(subpc,dim,nloc,array)); 914*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(gridctx[level].coords,&array)); 91539381ba2SJed Brown level--; 91639381ba2SJed Brown } 91739381ba2SJed Brown } 91839381ba2SJed Brown 919c07bf074SBarry Smith /* setupcalled is set to 0 so that MG is setup from scratch */ 920c07bf074SBarry Smith pc->setupcalled = 0; 921*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetUp_MG(pc)); 9225582bec1SHong Zhang PetscFunctionReturn(0); 9235582bec1SHong Zhang } 9245582bec1SHong Zhang 9255582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 9265582bec1SHong Zhang /* 9275582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 9285582bec1SHong Zhang that was created with PCCreate_ML(). 9295582bec1SHong Zhang 9305582bec1SHong Zhang Input Parameter: 9315582bec1SHong Zhang . pc - the preconditioner context 9325582bec1SHong Zhang 9335582bec1SHong Zhang Application Interface Routine: PCDestroy() 9345582bec1SHong Zhang */ 9356ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 9365582bec1SHong Zhang { 93701da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 93801da6913SBarry Smith PC_ML *pc_ml= (PC_ML*)mg->innerctx; 9395582bec1SHong Zhang 9405582bec1SHong Zhang PetscFunctionBegin; 941*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_ML(pc)); 942*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc_ml)); 943*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCDestroy_MG(pc)); 944*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL)); 9455582bec1SHong Zhang PetscFunctionReturn(0); 9465582bec1SHong Zhang } 9475582bec1SHong Zhang 9484416b707SBarry Smith PetscErrorCode PCSetFromOptions_ML(PetscOptionItems *PetscOptionsObject,PC pc) 9495582bec1SHong Zhang { 95039381ba2SJed Brown PetscInt indx,PrintLevel,partindx; 9515582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 95239381ba2SJed Brown const char *part[] = {"Zoltan","ParMETIS"}; 95339381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN) 95439381ba2SJed Brown const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"}; 95539381ba2SJed Brown #endif 95601da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 95701da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 958b5c8bdf8SJed Brown PetscMPIInt size; 959ce94432eSBarry Smith MPI_Comm comm; 9605582bec1SHong Zhang 9615582bec1SHong Zhang PetscFunctionBegin; 962*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)pc,&comm)); 963*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 964*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"ML options")); 9652fa5cd67SKarl Rupp 9665582bec1SHong Zhang PrintLevel = 0; 9675582bec1SHong Zhang indx = 0; 96839381ba2SJed Brown partindx = 0; 9692fa5cd67SKarl Rupp 970*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,NULL)); 971448f31a9SStefano Zampini PetscStackCall("ML_Set_PrintLevel",ML_Set_PrintLevel(PrintLevel)); 972*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,NULL)); 973*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,NULL)); 974*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,NULL)); 9752fa5cd67SKarl Rupp 9765582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 9772fa5cd67SKarl Rupp 978*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,NULL)); 979*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,NULL)); 980*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,NULL)); 981*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,NULL)); 982*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,NULL)); 983*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-pc_ml_nullspace","Which type of null space information to use","None",PCMLNullSpaceTypes,(PetscEnum)pc_ml->nulltype,(PetscEnum*)&pc_ml->nulltype,NULL)); 984*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_ml_EnergyMinimization","Energy minimization norm type (0=no minimization; see ML manual for 1,2,3; -1 and 4 undocumented)","None",pc_ml->EnergyMinimization,&pc_ml->EnergyMinimization,NULL)); 985*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_reuse_interpolation","Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)","None",pc_ml->reuse_interpolation,&pc_ml->reuse_interpolation,NULL)); 986b5c8bdf8SJed Brown /* 987b5c8bdf8SJed Brown The following checks a number of conditions. If we let this stuff slip by, then ML's error handling will take over. 988b5c8bdf8SJed Brown This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions. 989b5c8bdf8SJed Brown 990b5c8bdf8SJed Brown We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working 991b5c8bdf8SJed Brown combination of options and ML's exit(1) explanations don't help matters. 992b5c8bdf8SJed Brown */ 9932c71b3e2SJacob Faibussowitsch PetscCheckFalse(pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4,comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4"); 9942c71b3e2SJacob Faibussowitsch PetscCheckFalse(pc_ml->EnergyMinimization == 4 && size > 1,comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel"); 995*5f80ce2aSJacob Faibussowitsch if (pc_ml->EnergyMinimization == 4) CHKERRQ(PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2\n")); 996b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) { 997*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,NULL)); 998b5c8bdf8SJed Brown } 999b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == 2) { 1000b5c8bdf8SJed Brown /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */ 1001*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,NULL)); 1002b5c8bdf8SJed Brown } 1003b5c8bdf8SJed Brown /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */ 1004b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE; 1005*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_KeepAggInfo","Allows the preconditioner to be reused, or auxilliary matrices to be generated","None",pc_ml->KeepAggInfo,&pc_ml->KeepAggInfo,NULL)); 1006b5c8bdf8SJed Brown /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */ 1007b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE; 1008*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_Reusable","Store intermedaiate data structures so that the multilevel hierarchy is reusable","None",pc_ml->Reusable,&pc_ml->Reusable,NULL)); 1009b5c8bdf8SJed Brown /* 1010b5c8bdf8SJed Brown ML's C API is severely underdocumented and lacks significant functionality. The C++ API calls 1011b5c8bdf8SJed Brown ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function 1012b5c8bdf8SJed Brown ML_Gen_MGHierarchy_UsingAggregation(). This modification, however, does not provide a strict superset of the 1013b5c8bdf8SJed Brown functionality in the old function, so some users may still want to use it. Note that many options are ignored in 1014b5c8bdf8SJed Brown this context, but ML doesn't provide a way to find out which ones. 1015b5c8bdf8SJed Brown */ 1016*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,NULL)); 1017*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy","ML_Repartition_Activate",pc_ml->Repartition,&pc_ml->Repartition,NULL)); 101839381ba2SJed Brown if (pc_ml->Repartition) { 1019*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,NULL)); 1020*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,NULL)); 1021*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor","ML_Repartition_Set_PutOnSingleProc",pc_ml->PutOnSingleProc,&pc_ml->PutOnSingleProc,NULL)); 102239381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN) 102339381ba2SJed Brown partindx = 0; 1024*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,NULL)); 10252fa5cd67SKarl Rupp 102639381ba2SJed Brown pc_ml->RepartitionType = partindx; 102739381ba2SJed Brown if (!partindx) { 10285572b5bbSJed Brown PetscInt zindx = 0; 10292fa5cd67SKarl Rupp 1030*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,NULL)); 10312fa5cd67SKarl Rupp 103239381ba2SJed Brown pc_ml->ZoltanScheme = zindx; 103339381ba2SJed Brown } 103439381ba2SJed Brown #else 103539381ba2SJed Brown partindx = 1; 1036*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,NULL)); 1037e6b1cc6bSSatish Balay pc_ml->RepartitionType = partindx; 1038*5f80ce2aSJacob Faibussowitsch PetscCheck(partindx,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan"); 103939381ba2SJed Brown #endif 1040*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,NULL)); 1041*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,NULL)); 104239381ba2SJed Brown } 1043*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 10445582bec1SHong Zhang PetscFunctionReturn(0); 10455582bec1SHong Zhang } 10465582bec1SHong Zhang 10475582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 10485582bec1SHong Zhang /* 10495582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 10505582bec1SHong Zhang and sets this as the private data within the generic preconditioning 10515582bec1SHong Zhang context, PC, that was created within PCCreate(). 10525582bec1SHong Zhang 10535582bec1SHong Zhang Input Parameter: 10545582bec1SHong Zhang . pc - the preconditioner context 10555582bec1SHong Zhang 10565582bec1SHong Zhang Application Interface Routine: PCCreate() 10575582bec1SHong Zhang */ 10585582bec1SHong Zhang 10595582bec1SHong Zhang /*MC 10601e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 10615582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 10626ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 10636ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 10645582bec1SHong Zhang 10656ca4d86aSHong Zhang Options Database Key: 10662612397fSMatthew G. Knepley Multigrid options(inherited): 1067a2b725a8SWilliam Gropp + -pc_mg_cycles <1> - 1 for V cycle, 2 for W-cycle (MGSetCycles) 1068a2b725a8SWilliam Gropp . -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (PCMGSetDistinctSmoothUp) 1069a2b725a8SWilliam Gropp - -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1070a2b725a8SWilliam Gropp 107151f519a2SBarry Smith ML options: 1072a2b725a8SWilliam Gropp + -pc_ml_PrintLevel <0> - Print level (ML_Set_PrintLevel) 1073a2b725a8SWilliam Gropp . -pc_ml_maxNlevels <10> - Maximum number of levels (None) 1074a2b725a8SWilliam Gropp . -pc_ml_maxCoarseSize <1> - Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 1075a2b725a8SWilliam Gropp . -pc_ml_CoarsenScheme <Uncoupled> - (one of) Uncoupled Coupled MIS METIS 1076a2b725a8SWilliam Gropp . -pc_ml_DampingFactor <1.33333> - P damping factor (ML_Aggregate_Set_DampingFactor) 1077a2b725a8SWilliam Gropp . -pc_ml_Threshold <0> - Smoother drop tol (ML_Aggregate_Set_Threshold) 1078a2b725a8SWilliam Gropp . -pc_ml_SpectralNormScheme_Anorm <false> - Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm) 1079a5b23f4aSJose E. Roman . -pc_ml_repartition <false> - Allow ML to repartition levels of the hierarchy (ML_Repartition_Activate) 1080a2b725a8SWilliam Gropp . -pc_ml_repartitionMaxMinRatio <1.3> - Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio) 1081147403d9SBarry Smith . -pc_ml_repartitionMinPerProc <512> - Smallest repartitioned size (ML_Repartition_Set_MinPerProc) 1082a2b725a8SWilliam Gropp . -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc) 1083a2b725a8SWilliam Gropp . -pc_ml_repartitionType <Zoltan> - Repartitioning library to use (ML_Repartition_Set_Partitioner) 1084a2b725a8SWilliam Gropp . -pc_ml_repartitionZoltanScheme <RCB> - Repartitioning scheme to use (None) 1085147403d9SBarry Smith . -pc_ml_Aux <false> - Aggregate using auxiliary coordinate-based Laplacian (None) 1086a2b725a8SWilliam Gropp - -pc_ml_AuxThreshold <0.0> - Auxiliary smoother drop tol (None) 10875582bec1SHong Zhang 10885582bec1SHong Zhang Level: intermediate 10895582bec1SHong Zhang 10905582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1091710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetDistinctSmoothUp(), 1092710315b6SLawrence Mitchell PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 109397177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 109410167fecSBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 10955582bec1SHong Zhang M*/ 10965582bec1SHong Zhang 10978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc) 10985582bec1SHong Zhang { 10995582bec1SHong Zhang PC_ML *pc_ml; 110001da6913SBarry Smith PC_MG *mg; 11015582bec1SHong Zhang 11025582bec1SHong Zhang PetscFunctionBegin; 1103573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 1104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */ 1105*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)pc,PCML)); 1106e0f5d30fSBarry Smith /* Since PCMG tries to use DM assocated with PC must delete it */ 1107*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&pc->dm)); 1108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 1109e0f5d30fSBarry Smith mg = (PC_MG*)pc->data; 11105582bec1SHong Zhang 11115582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 1112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&pc_ml)); 111301da6913SBarry Smith mg->innerctx = pc_ml; 11145582bec1SHong Zhang 1115573998d7SHong Zhang pc_ml->ml_object = 0; 1116573998d7SHong Zhang pc_ml->agg_object = 0; 1117573998d7SHong Zhang pc_ml->gridctx = 0; 1118573998d7SHong Zhang pc_ml->PetscMLdata = 0; 1119573998d7SHong Zhang pc_ml->Nlevels = -1; 1120573998d7SHong Zhang pc_ml->MaxNlevels = 10; 1121573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 11223751b4bdSBarry Smith pc_ml->CoarsenScheme = 1; 1123573998d7SHong Zhang pc_ml->Threshold = 0.0; 1124573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 1125573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 1126573998d7SHong Zhang pc_ml->size = 0; 112739381ba2SJed Brown pc_ml->dim = 0; 112839381ba2SJed Brown pc_ml->nloc = 0; 112939381ba2SJed Brown pc_ml->coords = 0; 113039381ba2SJed Brown pc_ml->Repartition = PETSC_FALSE; 113139381ba2SJed Brown pc_ml->MaxMinRatio = 1.3; 113239381ba2SJed Brown pc_ml->MinPerProc = 512; 113339381ba2SJed Brown pc_ml->PutOnSingleProc = 5000; 113439381ba2SJed Brown pc_ml->RepartitionType = 0; 113539381ba2SJed Brown pc_ml->ZoltanScheme = 0; 113639381ba2SJed Brown pc_ml->Aux = PETSC_FALSE; 113739381ba2SJed Brown pc_ml->AuxThreshold = 0.0; 113839381ba2SJed Brown 113939381ba2SJed Brown /* allow for coordinates to be passed */ 1140*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_ML)); 1141573998d7SHong Zhang 11425582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 11435582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 11445582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 1145a06653b4SBarry Smith pc->ops->reset = PCReset_ML; 11465582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 11475582bec1SHong Zhang PetscFunctionReturn(0); 11485582bec1SHong Zhang } 1149