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 229371c9d4SSatish Balay typedef enum { 239371c9d4SSatish Balay PCML_NULLSPACE_AUTO, 249371c9d4SSatish Balay PCML_NULLSPACE_USER, 259371c9d4SSatish Balay PCML_NULLSPACE_BLOCK, 269371c9d4SSatish Balay PCML_NULLSPACE_SCALAR 279371c9d4SSatish Balay } PCMLNullSpaceType; 28fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO", "USER", "BLOCK", "SCALAR", "PCMLNullSpaceType", "PCML_NULLSPACE_", 0}; 29fb6a8e6dSJed Brown 305582bec1SHong Zhang /* The context (data structure) at each grid level */ 315582bec1SHong Zhang typedef struct { 325582bec1SHong Zhang Vec x, b, r; /* global vectors */ 335582bec1SHong Zhang Mat A, P, R; 345582bec1SHong Zhang KSP ksp; 3539381ba2SJed Brown Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */ 365582bec1SHong Zhang } GridCtx; 375582bec1SHong Zhang 385582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 395582bec1SHong Zhang typedef struct { 40573998d7SHong Zhang Mat A; /* Petsc matrix in aij format */ 41573998d7SHong Zhang Mat Aloc; /* local portion of A to be used by ML */ 4224a42b14SHong Zhang Vec x, y; 435582bec1SHong Zhang ML_Operator *mlmat; 445582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 455582bec1SHong Zhang } FineGridCtx; 465582bec1SHong Zhang 475582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 485582bec1SHong Zhang typedef struct { 495582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 505582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 515582bec1SHong Zhang } Mat_MLShell; 525582bec1SHong Zhang 535582bec1SHong Zhang /* Private context for the ML preconditioner */ 545582bec1SHong Zhang typedef struct { 555582bec1SHong Zhang ML *ml_object; 565582bec1SHong Zhang ML_Aggregate *agg_object; 575582bec1SHong Zhang GridCtx *gridctx; 585582bec1SHong Zhang FineGridCtx *PetscMLdata; 5939381ba2SJed Brown PetscInt Nlevels, MaxNlevels, MaxCoarseSize, CoarsenScheme, EnergyMinimization, MinPerProc, PutOnSingleProc, RepartitionType, ZoltanScheme; 6039381ba2SJed Brown PetscReal Threshold, DampingFactor, EnergyMinimizationDropTol, MaxMinRatio, AuxThreshold; 6139381ba2SJed Brown PetscBool SpectralNormScheme_Anorm, BlockScaling, EnergyMinimizationCheap, Symmetrize, OldHierarchy, KeepAggInfo, Reusable, Repartition, Aux; 6248268eb4SJed Brown PetscBool reuse_interpolation; 63fb6a8e6dSJed Brown PCMLNullSpaceType nulltype; 64573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 6539381ba2SJed Brown PetscInt dim; /* data from PCSetCoordinates(_ML) */ 6639381ba2SJed Brown PetscInt nloc; 6739381ba2SJed Brown PetscReal *coords; /* ML has a grid object for each level: the finest grid will point into coords */ 685582bec1SHong Zhang } PC_ML; 6941ca0015SHong Zhang 709371c9d4SSatish Balay 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[]) { 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 765f80ce2aSJacob Faibussowitsch if (MatGetSize(ml->Aloc, &m, NULL)) 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 939371c9d4SSatish Balay static PetscErrorCode PetscML_comm(double p[], void *ML_data) { 946562c4e1SBarry Smith FineGridCtx *ml = (FineGridCtx *)ML_data; 956562c4e1SBarry Smith Mat A = ml->A; 966562c4e1SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data; 976562c4e1SBarry Smith PetscMPIInt size; 986562c4e1SBarry Smith PetscInt i, in_length = A->rmap->n, out_length = ml->Aloc->cmap->n; 99d9ca1df4SBarry Smith const PetscScalar *array; 1006562c4e1SBarry Smith 1016562c4e1SBarry Smith PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 103708418deSStefano Zampini if (size == 1) PetscFunctionReturn(0); 1046562c4e1SBarry Smith 1059566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ml->y, p)); 1069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD)); 1079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD)); 1089566063dSJacob Faibussowitsch PetscCall(VecResetArray(ml->y)); 1099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(a->lvec, &array)); 1102fa5cd67SKarl Rupp for (i = in_length; i < out_length; i++) p[i] = array[i - in_length]; 1119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(a->lvec, &array)); 1126562c4e1SBarry Smith PetscFunctionReturn(0); 1136562c4e1SBarry Smith } 1146562c4e1SBarry Smith 1159371c9d4SSatish Balay static int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length, double ap[]) { 1166562c4e1SBarry Smith FineGridCtx *ml = (FineGridCtx *)ML_Get_MyMatvecData(ML_data); 1176562c4e1SBarry Smith Mat A = ml->A, Aloc = ml->Aloc; 1186562c4e1SBarry Smith PetscMPIInt size; 1196562c4e1SBarry Smith PetscScalar *pwork = ml->pwork; 1206562c4e1SBarry Smith PetscInt i; 1216562c4e1SBarry Smith 1226562c4e1SBarry Smith PetscFunctionBegin; 1239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 1246562c4e1SBarry Smith if (size == 1) { 1259566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ml->x, p)); 1266562c4e1SBarry Smith } else { 1276562c4e1SBarry Smith for (i = 0; i < in_length; i++) pwork[i] = p[i]; 1289566063dSJacob Faibussowitsch PetscCall(PetscML_comm(pwork, ml)); 1299566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ml->x, pwork)); 1306562c4e1SBarry Smith } 1319566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ml->y, ap)); 1329566063dSJacob Faibussowitsch PetscCall(MatMult(Aloc, ml->x, ml->y)); 1339566063dSJacob Faibussowitsch PetscCall(VecResetArray(ml->x)); 1349566063dSJacob Faibussowitsch PetscCall(VecResetArray(ml->y)); 1356562c4e1SBarry Smith PetscFunctionReturn(0); 1366562c4e1SBarry Smith } 1376562c4e1SBarry Smith 1389371c9d4SSatish Balay static PetscErrorCode MatMult_ML(Mat A, Vec x, Vec y) { 1396562c4e1SBarry Smith Mat_MLShell *shell; 140d9ca1df4SBarry Smith PetscScalar *yarray; 141d9ca1df4SBarry Smith const PetscScalar *xarray; 1426562c4e1SBarry Smith PetscInt x_length, y_length; 1436562c4e1SBarry Smith 1446562c4e1SBarry Smith PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &shell)); 1469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xarray)); 1479566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yarray)); 1486562c4e1SBarry Smith x_length = shell->mlmat->invec_leng; 1496562c4e1SBarry Smith y_length = shell->mlmat->outvec_leng; 150e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Operator_Apply", ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar *)xarray, y_length, yarray)); 1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xarray)); 1529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yarray)); 1536562c4e1SBarry Smith PetscFunctionReturn(0); 1546562c4e1SBarry Smith } 1556562c4e1SBarry Smith 15679d04de1SBarry Smith /* newtype is ignored since only handles one case */ 1579371c9d4SSatish Balay static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A, MatType newtype, MatReuse scall, Mat *Aloc) { 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; 1665f80ce2aSJacob Faibussowitsch PetscCheck(am == an, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A must have a square diagonal portion, am: %d != an: %d", am, an); 1679566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(mpimat->A, (const PetscScalar **)&aa)); 1689566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(mpimat->B, (const PetscScalar **)&ba)); 1696562c4e1SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1709566063dSJacob Faibussowitsch PetscCall(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]); 1739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + ci[am], &cj)); 1749566063dSJacob Faibussowitsch PetscCall(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++) { 1879371c9d4SSatish Balay cj[k] = an + (*bj); 1889371c9d4SSatish Balay bj++; 1896562c4e1SBarry Smith ca[k++] = *ba++; 1906562c4e1SBarry Smith } 1916562c4e1SBarry Smith } 1925f80ce2aSJacob Faibussowitsch PetscCheck(k == ci[am], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "k: %d != ci[am]: %d", k, ci[am]); 1936562c4e1SBarry Smith 1946562c4e1SBarry Smith /* put together the new matrix */ 1956562c4e1SBarry Smith an = mpimat->A->cmap->n + mpimat->B->cmap->n; 1969566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, am, an, ci, cj, ca, Aloc)); 1976562c4e1SBarry Smith 1986562c4e1SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 1996562c4e1SBarry Smith /* Since these are PETSc arrays, change flags to free them as necessary. */ 2006562c4e1SBarry Smith mat = (Mat_SeqAIJ *)(*Aloc)->data; 2016562c4e1SBarry Smith mat->free_a = PETSC_TRUE; 2026562c4e1SBarry Smith mat->free_ij = PETSC_TRUE; 2036562c4e1SBarry Smith 2046562c4e1SBarry Smith mat->nonew = 0; 2056562c4e1SBarry Smith } else if (scall == MAT_REUSE_MATRIX) { 2066562c4e1SBarry Smith mat = (Mat_SeqAIJ *)(*Aloc)->data; 2079371c9d4SSatish Balay ci = mat->i; 2089371c9d4SSatish Balay cj = mat->j; 2099371c9d4SSatish Balay ca = mat->a; 2106562c4e1SBarry Smith for (i = 0; i < am; i++) { 2116562c4e1SBarry Smith /* diagonal portion of A */ 2126562c4e1SBarry Smith ncols = ai[i + 1] - ai[i]; 2136562c4e1SBarry Smith for (j = 0; j < ncols; j++) *ca++ = *aa++; 2146562c4e1SBarry Smith /* off-diagonal portion of A */ 2156562c4e1SBarry Smith ncols = bi[i + 1] - bi[i]; 2166562c4e1SBarry Smith for (j = 0; j < ncols; j++) *ca++ = *ba++; 2176562c4e1SBarry Smith } 21898921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid MatReuse %d", (int)scall); 2199566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArrayRead(mpimat->A, (const PetscScalar **)&aa)); 2209566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArrayRead(mpimat->B, (const PetscScalar **)&ba)); 2216562c4e1SBarry Smith PetscFunctionReturn(0); 2226562c4e1SBarry Smith } 2236562c4e1SBarry Smith 2249371c9d4SSatish Balay static PetscErrorCode MatDestroy_ML(Mat A) { 2256562c4e1SBarry Smith Mat_MLShell *shell; 2266562c4e1SBarry Smith 2276562c4e1SBarry Smith PetscFunctionBegin; 2289566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &shell)); 2299566063dSJacob Faibussowitsch PetscCall(PetscFree(shell)); 2306562c4e1SBarry Smith PetscFunctionReturn(0); 2316562c4e1SBarry Smith } 2326562c4e1SBarry Smith 2339371c9d4SSatish Balay static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat) { 2346562c4e1SBarry Smith struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 2350298fd71SBarry Smith PetscInt m = mlmat->outvec_leng, n = mlmat->invec_leng, *nnz = NULL, nz_max; 23639381ba2SJed Brown PetscInt *ml_cols = matdata->columns, *ml_rowptr = matdata->rowptr, *aj, i; 2376562c4e1SBarry Smith PetscScalar *ml_vals = matdata->values, *aa; 2386562c4e1SBarry Smith 2396562c4e1SBarry Smith PetscFunctionBegin; 2405f80ce2aSJacob Faibussowitsch PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL"); 2416562c4e1SBarry Smith if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 2426562c4e1SBarry Smith if (reuse) { 2436562c4e1SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)(*newmat)->data; 2446562c4e1SBarry Smith aij->i = ml_rowptr; 2456562c4e1SBarry Smith aij->j = ml_cols; 2466562c4e1SBarry Smith aij->a = ml_vals; 2476562c4e1SBarry Smith } else { 2486562c4e1SBarry Smith /* sort ml_cols and ml_vals */ 2499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m + 1, &nnz)); 2502fa5cd67SKarl Rupp for (i = 0; i < m; i++) nnz[i] = ml_rowptr[i + 1] - ml_rowptr[i]; 2519371c9d4SSatish Balay aj = ml_cols; 2529371c9d4SSatish Balay aa = ml_vals; 2536562c4e1SBarry Smith for (i = 0; i < m; i++) { 2549566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithScalarArray(nnz[i], aj, aa)); 2559371c9d4SSatish Balay aj += nnz[i]; 2569371c9d4SSatish Balay aa += nnz[i]; 2576562c4e1SBarry Smith } 2589566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, m, n, ml_rowptr, ml_cols, ml_vals, newmat)); 2599566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 2606562c4e1SBarry Smith } 2619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY)); 2629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY)); 2636562c4e1SBarry Smith PetscFunctionReturn(0); 2646562c4e1SBarry Smith } 2656562c4e1SBarry Smith 26639381ba2SJed Brown nz_max = PetscMax(1, mlmat->max_nz_per_row); 2679566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj)); 26839381ba2SJed Brown if (!reuse) { 2699566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, newmat)); 2709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*newmat, m, n, PETSC_DECIDE, PETSC_DECIDE)); 2719566063dSJacob Faibussowitsch PetscCall(MatSetType(*newmat, MATSEQAIJ)); 27239381ba2SJed Brown /* keep track of block size for A matrices */ 2739566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(*newmat, mlmat->num_PDEs)); 2746562c4e1SBarry Smith 2759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &nnz)); 276ad540459SPierre Jolivet for (i = 0; i < m; i++) PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i])); 2779566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*newmat, 0, nnz)); 278ae7fe62dSJed Brown } 2796562c4e1SBarry Smith for (i = 0; i < m; i++) { 280ae7fe62dSJed Brown PetscInt ncols; 28139381ba2SJed Brown 282e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols)); 2839566063dSJacob Faibussowitsch PetscCall(MatSetValues(*newmat, 1, &i, ncols, aj, aa, INSERT_VALUES)); 2846562c4e1SBarry Smith } 2859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY)); 2869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY)); 2876562c4e1SBarry Smith 2889566063dSJacob Faibussowitsch PetscCall(PetscFree2(aa, aj)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 2906562c4e1SBarry Smith PetscFunctionReturn(0); 2916562c4e1SBarry Smith } 2926562c4e1SBarry Smith 2939371c9d4SSatish Balay static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat, MatReuse reuse, Mat *newmat) { 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) { 3039566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(*newmat, &shellctx)); 3046562c4e1SBarry Smith shellctx->mlmat = mlmat; 3056562c4e1SBarry Smith PetscFunctionReturn(0); 3066562c4e1SBarry Smith } 3076562c4e1SBarry Smith 3086562c4e1SBarry Smith MLcomm = mlmat->comm; 3092fa5cd67SKarl Rupp 3109566063dSJacob Faibussowitsch PetscCall(PetscNew(&shellctx)); 3119566063dSJacob Faibussowitsch PetscCall(MatCreateShell(MLcomm->USR_comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, shellctx, newmat)); 3129566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*newmat, MATOP_MULT, (void (*)(void))MatMult_ML)); 3139566063dSJacob Faibussowitsch PetscCall(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 3209371c9d4SSatish Balay static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat) { 32139381ba2SJed Brown PetscInt *aj; 32239381ba2SJed Brown PetscScalar *aa; 32339381ba2SJed Brown PetscInt i, j, *gordering; 324ae7fe62dSJed Brown PetscInt m = mlmat->outvec_leng, n, nz_max, row; 3256562c4e1SBarry Smith Mat A; 3266562c4e1SBarry Smith 3276562c4e1SBarry Smith PetscFunctionBegin; 3285f80ce2aSJacob Faibussowitsch PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL"); 3296562c4e1SBarry Smith n = mlmat->invec_leng; 3305f80ce2aSJacob Faibussowitsch PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %d must equal to n %d", m, n); 3316562c4e1SBarry Smith 3327be6b909SBarry Smith /* create global row numbering for a ML_Operator */ 333e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_build_global_numbering", ML_build_global_numbering(mlmat, &gordering, "rows")); 3347be6b909SBarry Smith 3351d94bf15SBarry Smith nz_max = PetscMax(1, mlmat->max_nz_per_row) + 1; 3369566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj)); 3377be6b909SBarry Smith if (reuse) { 3387be6b909SBarry Smith A = *newmat; 3397be6b909SBarry Smith } else { 340ae7fe62dSJed Brown PetscInt *nnzA, *nnzB, *nnz; 3417be6b909SBarry Smith PetscInt rstart; 3429566063dSJacob Faibussowitsch PetscCall(MatCreate(mlmat->comm->USR_comm, &A)); 3439566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, PETSC_DECIDE, PETSC_DECIDE)); 3449566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ)); 34539381ba2SJed Brown /* keep track of block size for A matrices */ 3469566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, mlmat->num_PDEs)); 3479566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(m, &nnzA, m, &nnzB, m, &nnz)); 3489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&m, &rstart, 1, MPIU_INT, MPI_SUM, mlmat->comm->USR_comm)); 3497be6b909SBarry Smith rstart -= m; 3506562c4e1SBarry Smith 3516562c4e1SBarry Smith for (i = 0; i < m; i++) { 3527be6b909SBarry Smith row = gordering[i] - rstart; 353e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i])); 3547be6b909SBarry Smith nnzA[row] = 0; 35539381ba2SJed Brown for (j = 0; j < nnz[i]; j++) { 3567be6b909SBarry Smith if (aj[j] < m) nnzA[row]++; 3576562c4e1SBarry Smith } 3587be6b909SBarry Smith nnzB[row] = nnz[i] - nnzA[row]; 3596562c4e1SBarry Smith } 3609566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 0, nnzA, 0, nnzB)); 3619566063dSJacob Faibussowitsch PetscCall(PetscFree3(nnzA, nnzB, nnz)); 362ae7fe62dSJed Brown } 3636562c4e1SBarry Smith for (i = 0; i < m; i++) { 364ae7fe62dSJed Brown PetscInt ncols; 3656562c4e1SBarry Smith row = gordering[i]; 36639381ba2SJed Brown 367e77caa6dSBarry Smith PetscStackCallExternalVoid(",ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols)); 3682fa5cd67SKarl Rupp for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]]; 3699566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, ncols, aj, aa, INSERT_VALUES)); 3706562c4e1SBarry Smith } 371e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_free", ML_free(gordering)); 3729566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3739566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 3746562c4e1SBarry Smith *newmat = A; 3756562c4e1SBarry Smith 3769566063dSJacob Faibussowitsch PetscCall(PetscFree2(aa, aj)); 3776562c4e1SBarry Smith PetscFunctionReturn(0); 3786562c4e1SBarry Smith } 3796562c4e1SBarry Smith 38039381ba2SJed Brown /* -------------------------------------------------------------------------- */ 38139381ba2SJed Brown /* 38239381ba2SJed Brown PCSetCoordinates_ML 38339381ba2SJed Brown 38439381ba2SJed Brown Input Parameter: 38539381ba2SJed Brown . pc - the preconditioner context 38639381ba2SJed Brown */ 3879371c9d4SSatish Balay static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) { 38839381ba2SJed Brown PC_MG *mg = (PC_MG *)pc->data; 38939381ba2SJed Brown PC_ML *pc_ml = (PC_ML *)mg->innerctx; 39090fbc344SStefano Zampini PetscInt arrsz, oldarrsz, bs, my0, kk, ii, nloc, Iend, aloc; 39139381ba2SJed Brown Mat Amat = pc->pmat; 39239381ba2SJed Brown 39339381ba2SJed Brown /* this function copied and modified from PCSetCoordinates_GEO -TGI */ 39439381ba2SJed Brown PetscFunctionBegin; 39539381ba2SJed Brown PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1); 3969566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 39739381ba2SJed Brown 3989566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend)); 39990fbc344SStefano Zampini aloc = (Iend - my0); 40039381ba2SJed Brown nloc = (Iend - my0) / bs; 40139381ba2SJed Brown 40263a3b9bcSJacob Faibussowitsch PetscCheck((nloc == a_nloc) || (aloc == a_nloc), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of local blocks %" PetscInt_FMT " must be %" PetscInt_FMT " or %" PetscInt_FMT ".", a_nloc, nloc, aloc); 40339381ba2SJed Brown 40439381ba2SJed Brown oldarrsz = pc_ml->dim * pc_ml->nloc; 40539381ba2SJed Brown pc_ml->dim = ndm; 40690fbc344SStefano Zampini pc_ml->nloc = nloc; 40790fbc344SStefano Zampini arrsz = ndm * nloc; 40839381ba2SJed Brown 40939381ba2SJed Brown /* create data - syntactic sugar that should be refactored at some point */ 41039381ba2SJed Brown if (pc_ml->coords == 0 || (oldarrsz != arrsz)) { 4119566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml->coords)); 4129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz, &pc_ml->coords)); 41339381ba2SJed Brown } 41439381ba2SJed Brown for (kk = 0; kk < arrsz; kk++) pc_ml->coords[kk] = -999.; 41539381ba2SJed Brown /* copy data in - column oriented */ 41690fbc344SStefano Zampini if (nloc == a_nloc) { 41739381ba2SJed Brown for (kk = 0; kk < nloc; kk++) { 418ad540459SPierre Jolivet for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[kk * ndm + ii]; 41939381ba2SJed Brown } 42090fbc344SStefano Zampini } else { /* assumes the coordinates are blocked */ 42190fbc344SStefano Zampini for (kk = 0; kk < nloc; kk++) { 422ad540459SPierre Jolivet for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[bs * kk * ndm + ii]; 42390fbc344SStefano Zampini } 42490fbc344SStefano Zampini } 42539381ba2SJed Brown PetscFunctionReturn(0); 42639381ba2SJed Brown } 42739381ba2SJed Brown 4286562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/ 429e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC); 4309371c9d4SSatish Balay PetscErrorCode PCReset_ML(PC pc) { 431e0262f48SMatthew G Knepley PC_MG *mg = (PC_MG *)pc->data; 432e0262f48SMatthew G Knepley PC_ML *pc_ml = (PC_ML *)mg->innerctx; 43339381ba2SJed Brown PetscInt level, fine_level = pc_ml->Nlevels - 1, dim = pc_ml->dim; 43401da6913SBarry Smith 43501da6913SBarry Smith PetscFunctionBegin; 43639381ba2SJed Brown if (dim) { 43748a46eb9SPierre Jolivet for (level = 0; level <= fine_level; level++) PetscCall(VecDestroy(&pc_ml->gridctx[level].coords)); 438448f31a9SStefano Zampini if (pc_ml->ml_object && pc_ml->ml_object->Grid) { 439448f31a9SStefano Zampini ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)pc_ml->ml_object->Grid[0].Grid; 44039381ba2SJed Brown grid_info->x = 0; /* do this so ML doesn't try to free coordinates */ 44139381ba2SJed Brown grid_info->y = 0; 44239381ba2SJed Brown grid_info->z = 0; 443e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object)); 44439381ba2SJed Brown } 445448f31a9SStefano Zampini } 446e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Aggregate_Destroy(&pc_ml->agg_object)); 447e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Destroy(&pc_ml->ml_object)); 44801da6913SBarry Smith 44901da6913SBarry Smith if (pc_ml->PetscMLdata) { 4509566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml->PetscMLdata->pwork)); 4519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc_ml->PetscMLdata->Aloc)); 4529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ml->PetscMLdata->x)); 4539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ml->PetscMLdata->y)); 45401da6913SBarry Smith } 4559566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml->PetscMLdata)); 45601da6913SBarry Smith 457f5a5dd59SJed Brown if (pc_ml->gridctx) { 45801da6913SBarry Smith for (level = 0; level < fine_level; level++) { 4599566063dSJacob Faibussowitsch if (pc_ml->gridctx[level].A) PetscCall(MatDestroy(&pc_ml->gridctx[level].A)); 4609566063dSJacob Faibussowitsch if (pc_ml->gridctx[level].P) PetscCall(MatDestroy(&pc_ml->gridctx[level].P)); 4619566063dSJacob Faibussowitsch if (pc_ml->gridctx[level].R) PetscCall(MatDestroy(&pc_ml->gridctx[level].R)); 4629566063dSJacob Faibussowitsch if (pc_ml->gridctx[level].x) PetscCall(VecDestroy(&pc_ml->gridctx[level].x)); 4639566063dSJacob Faibussowitsch if (pc_ml->gridctx[level].b) PetscCall(VecDestroy(&pc_ml->gridctx[level].b)); 4649566063dSJacob Faibussowitsch if (pc_ml->gridctx[level + 1].r) PetscCall(VecDestroy(&pc_ml->gridctx[level + 1].r)); 46501da6913SBarry Smith } 466f5a5dd59SJed Brown } 4679566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml->gridctx)); 4689566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml->coords)); 4692fa5cd67SKarl Rupp 47039381ba2SJed Brown pc_ml->dim = 0; 47139381ba2SJed Brown pc_ml->nloc = 0; 4729566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 47301da6913SBarry Smith PetscFunctionReturn(0); 47401da6913SBarry Smith } 4755582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4765582bec1SHong Zhang /* 4775582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 4785582bec1SHong Zhang by setting data structures and options. 4795582bec1SHong Zhang 4805582bec1SHong Zhang Input Parameter: 4815582bec1SHong Zhang . pc - the preconditioner context 4825582bec1SHong Zhang 4835582bec1SHong Zhang Application Interface Routine: PCSetUp() 4845582bec1SHong Zhang 485*f1580f4eSBarry Smith Note: 4865582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 4875582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 4885582bec1SHong Zhang */ 489dbbe0bcdSBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PC, PetscOptionItems *PetscOptionsObject); 490a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC); 491c07bf074SBarry Smith 4929371c9d4SSatish Balay PetscErrorCode PCSetUp_ML(PC pc) { 493eef31507SHong Zhang PetscMPIInt size; 4945582bec1SHong Zhang FineGridCtx *PetscMLdata; 4955582bec1SHong Zhang ML *ml_object; 4965582bec1SHong Zhang ML_Aggregate *agg_object; 4975582bec1SHong Zhang ML_Operator *mlmat; 4984f8eab3cSJed Brown PetscInt nlocal_allcols, Nlevels, mllevel, level, level1, m, fine_level, bs; 4995582bec1SHong Zhang Mat A, Aloc; 5005582bec1SHong Zhang GridCtx *gridctx; 50101da6913SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 50201da6913SBarry Smith PC_ML *pc_ml = (PC_ML *)mg->innerctx; 503ace3abfcSBarry Smith PetscBool isSeq, isMPI; 504c07bf074SBarry Smith KSP smoother; 505c07bf074SBarry Smith PC subpc; 50648268eb4SJed Brown PetscInt mesh_level, old_mesh_level; 5078a62b701SToby Isaac MatInfo info; 5081f817a21SBarry Smith static PetscBool cite = PETSC_FALSE; 50948268eb4SJed Brown 5105582bec1SHong Zhang PetscFunctionBegin; 5119371c9d4SSatish Balay PetscCall(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 = " 5129371c9d4SSatish Balay "{SAND2004-4821},\n year = 2004\n}\n", 5139371c9d4SSatish Balay &cite)); 51448268eb4SJed Brown A = pc->pmat; 5159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 51648268eb4SJed Brown 517573998d7SHong Zhang if (pc->setupcalled) { 51848268eb4SJed Brown if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) { 51948268eb4SJed Brown /* 52048268eb4SJed Brown Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for 52148268eb4SJed Brown multiple solves in which the matrix is not changing too quickly. 52248268eb4SJed Brown */ 52348268eb4SJed Brown ml_object = pc_ml->ml_object; 52448268eb4SJed Brown gridctx = pc_ml->gridctx; 52548268eb4SJed Brown Nlevels = pc_ml->Nlevels; 52648268eb4SJed Brown fine_level = Nlevels - 1; 52748268eb4SJed Brown gridctx[fine_level].A = A; 52848268eb4SJed Brown 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq)); 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI)); 53148268eb4SJed Brown if (isMPI) { 5329566063dSJacob Faibussowitsch PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc)); 53348268eb4SJed Brown } else if (isSeq) { 53448268eb4SJed Brown Aloc = A; 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Aloc)); 53698921bdaSJacob 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); 53748268eb4SJed Brown 5389566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols)); 53948268eb4SJed Brown PetscMLdata = pc_ml->PetscMLdata; 5409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PetscMLdata->Aloc)); 54148268eb4SJed Brown PetscMLdata->A = A; 54248268eb4SJed Brown PetscMLdata->Aloc = Aloc; 543e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata)); 544e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec)); 54548268eb4SJed Brown 54648268eb4SJed Brown mesh_level = ml_object->ML_finest_level; 54748268eb4SJed Brown while (ml_object->SingleLevel[mesh_level].Rmat->to) { 54848268eb4SJed Brown old_mesh_level = mesh_level; 54948268eb4SJed Brown mesh_level = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum; 55048268eb4SJed Brown 55148268eb4SJed Brown /* clean and regenerate A */ 55248268eb4SJed Brown mlmat = &(ml_object->Amat[mesh_level]); 553e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Operator_Clean", ML_Operator_Clean(mlmat)); 554e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Operator_Init", ML_Operator_Init(mlmat, ml_object->comm)); 555e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Gen_AmatrixRAP", ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level)); 55648268eb4SJed Brown } 55748268eb4SJed Brown 55848268eb4SJed Brown level = fine_level - 1; 55948268eb4SJed Brown if (size == 1) { /* convert ML P, R and A into seqaij format */ 56048268eb4SJed Brown for (mllevel = 1; mllevel < Nlevels; mllevel++) { 56148268eb4SJed Brown mlmat = &(ml_object->Amat[mllevel]); 5629566063dSJacob Faibussowitsch PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A)); 56348268eb4SJed Brown level--; 56448268eb4SJed Brown } 56548268eb4SJed Brown } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 56648268eb4SJed Brown for (mllevel = 1; mllevel < Nlevels; mllevel++) { 56748268eb4SJed Brown mlmat = &(ml_object->Amat[mllevel]); 5689566063dSJacob Faibussowitsch PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A)); 56948268eb4SJed Brown level--; 57048268eb4SJed Brown } 57148268eb4SJed Brown } 57248268eb4SJed Brown 57348268eb4SJed Brown for (level = 0; level < fine_level; level++) { 57448a46eb9SPierre Jolivet if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A)); 5759566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A)); 57648268eb4SJed Brown } 5779566063dSJacob Faibussowitsch PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A)); 5789566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A)); 57948268eb4SJed Brown 5809566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 58148268eb4SJed Brown PetscFunctionReturn(0); 58248268eb4SJed Brown } else { 583c07bf074SBarry Smith /* since ML can change the size of vectors/matrices at any level we must destroy everything */ 5849566063dSJacob Faibussowitsch PetscCall(PCReset_ML(pc)); 585573998d7SHong Zhang } 58648268eb4SJed Brown } 587573998d7SHong Zhang 5885582bec1SHong Zhang /* setup special features of PCML */ 5895582bec1SHong Zhang /*--------------------------------*/ 5905582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 5915582bec1SHong Zhang pc_ml->size = size; 5929566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq)); 5939566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI)); 594864b637dSMatthew Knepley if (isMPI) { 5959566063dSJacob Faibussowitsch PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc)); 596864b637dSMatthew Knepley } else if (isSeq) { 5975582bec1SHong Zhang Aloc = A; 5989566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Aloc)); 59998921bdaSJacob 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); 6005582bec1SHong Zhang 6015582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 6029566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &PetscMLdata)); 6035582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 6049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Aloc->cmap->n + 1, &PetscMLdata->pwork)); 6055582bec1SHong Zhang 6069566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Aloc, &PetscMLdata->x, &PetscMLdata->y)); 60724a42b14SHong Zhang 608573998d7SHong Zhang PetscMLdata->A = A; 609573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 61039381ba2SJed Brown if (pc_ml->dim) { /* create vecs around the coordinate data given */ 61139381ba2SJed Brown PetscInt i, j, dim = pc_ml->dim; 61239381ba2SJed Brown PetscInt nloc = pc_ml->nloc, nlocghost; 61339381ba2SJed Brown PetscReal *ghostedcoords; 61439381ba2SJed Brown 6159566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 61639381ba2SJed Brown nlocghost = Aloc->cmap->n / bs; 6179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nlocghost, &ghostedcoords)); 61839381ba2SJed Brown for (i = 0; i < dim; i++) { 61939381ba2SJed Brown /* copy coordinate values into first component of pwork */ 620ad540459SPierre Jolivet for (j = 0; j < nloc; j++) PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j]; 62139381ba2SJed Brown /* get the ghost values */ 6229566063dSJacob Faibussowitsch PetscCall(PetscML_comm(PetscMLdata->pwork, PetscMLdata)); 62339381ba2SJed Brown /* write into the vector */ 624ad540459SPierre Jolivet for (j = 0; j < nlocghost; j++) ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j]; 62539381ba2SJed Brown } 62639381ba2SJed Brown /* replace the original coords with the ghosted coords, because these are 62739381ba2SJed Brown * what ML needs */ 6289566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml->coords)); 62939381ba2SJed Brown pc_ml->coords = ghostedcoords; 63039381ba2SJed Brown } 63124a42b14SHong Zhang 6325582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 63345cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 6349566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols)); 6359566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 636e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Create", ML_Create(&ml_object, pc_ml->MaxNlevels)); 637e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Comm_Set_UsrComm", ML_Comm_Set_UsrComm(ml_object->comm, PetscObjectComm((PetscObject)A))); 638573998d7SHong Zhang pc_ml->ml_object = ml_object; 639e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Init_Amatrix", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata)); 640e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Set_Amatrix_Getrow", ML_Set_Amatrix_Getrow(ml_object, 0, PetscML_getrow, PetscML_comm, nlocal_allcols)); 641e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec)); 6425582bec1SHong Zhang 643e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Set_Symmetrize", ML_Set_Symmetrize(ml_object, pc_ml->Symmetrize ? ML_YES : ML_NO)); 644b5c8bdf8SJed Brown 6455582bec1SHong Zhang /* aggregation */ 646e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Create", ML_Aggregate_Create(&agg_object)); 647573998d7SHong Zhang pc_ml->agg_object = agg_object; 648573998d7SHong Zhang 649fb6a8e6dSJed Brown { 650fb6a8e6dSJed Brown MatNullSpace mnull; 6519566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(A, &mnull)); 652fb6a8e6dSJed Brown if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) { 653fb6a8e6dSJed Brown if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER; 654fb6a8e6dSJed Brown else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK; 655fb6a8e6dSJed Brown else pc_ml->nulltype = PCML_NULLSPACE_SCALAR; 656fb6a8e6dSJed Brown } 657fb6a8e6dSJed Brown switch (pc_ml->nulltype) { 658fb6a8e6dSJed Brown case PCML_NULLSPACE_USER: { 659fb6a8e6dSJed Brown PetscScalar *nullvec; 660fb6a8e6dSJed Brown const PetscScalar *v; 661fb6a8e6dSJed Brown PetscBool has_const; 6621c547e14SJed Brown PetscInt i, j, mlocal, nvec, M; 663fb6a8e6dSJed Brown const Vec *vecs; 6642fa5cd67SKarl Rupp 6655f80ce2aSJacob Faibussowitsch PetscCheck(mnull, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space"); 6669566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL)); 6679566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Aloc, &mlocal, NULL)); 6689566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 6699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 6709371c9d4SSatish Balay if (has_const) 6719371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0 / M; 672fb6a8e6dSJed Brown for (i = 0; i < nvec; i++) { 6739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 674fb6a8e6dSJed Brown for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = v[j]; 6759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 676fb6a8e6dSJed Brown } 677e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Create", PetscCall(ML_Aggregate_Set_NullSpace(agg_object, bs, nvec + !!has_const, nullvec, mlocal))); 6789566063dSJacob Faibussowitsch PetscCall(PetscFree(nullvec)); 679fb6a8e6dSJed Brown } break; 6809371c9d4SSatish Balay case PCML_NULLSPACE_BLOCK: PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", PetscCall(ML_Aggregate_Set_NullSpace(agg_object, bs, bs, 0, 0))); break; 6819371c9d4SSatish Balay case PCML_NULLSPACE_SCALAR: break; 682ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unknown null space type"); 683fb6a8e6dSJed Brown } 684fb6a8e6dSJed Brown } 685e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Set_MaxCoarseSize", ML_Aggregate_Set_MaxCoarseSize(agg_object, pc_ml->MaxCoarseSize)); 6865582bec1SHong Zhang /* set options */ 6875582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 6889371c9d4SSatish Balay case 1: PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_Coupled", ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object)); break; 6899371c9d4SSatish Balay case 2: PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_MIS", ML_Aggregate_Set_CoarsenScheme_MIS(agg_object)); break; 6909371c9d4SSatish Balay case 3: PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_METIS", ML_Aggregate_Set_CoarsenScheme_METIS(agg_object)); break; 6915582bec1SHong Zhang } 692e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Set_Threshold", ML_Aggregate_Set_Threshold(agg_object, pc_ml->Threshold)); 693e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Set_DampingFactor", ML_Aggregate_Set_DampingFactor(agg_object, pc_ml->DampingFactor)); 694ad540459SPierre Jolivet if (pc_ml->SpectralNormScheme_Anorm) PetscStackCallExternalVoid("ML_Set_SpectralNormScheme_Anorm", ML_Set_SpectralNormScheme_Anorm(ml_object)); 695b5c8bdf8SJed Brown agg_object->keep_agg_information = (int)pc_ml->KeepAggInfo; 696b5c8bdf8SJed Brown agg_object->keep_P_tentative = (int)pc_ml->Reusable; 697b5c8bdf8SJed Brown agg_object->block_scaled_SA = (int)pc_ml->BlockScaling; 698b5c8bdf8SJed Brown agg_object->minimizing_energy = (int)pc_ml->EnergyMinimization; 699b5c8bdf8SJed Brown agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol; 700b5c8bdf8SJed Brown agg_object->cheap_minimizing_energy = (int)pc_ml->EnergyMinimizationCheap; 7015582bec1SHong Zhang 70239381ba2SJed Brown if (pc_ml->Aux) { 7035f80ce2aSJacob Faibussowitsch PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Auxiliary matrix requires coordinates"); 70439381ba2SJed Brown ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold; 70539381ba2SJed Brown ml_object->Amat[0].aux_data->enable = 1; 70639381ba2SJed Brown ml_object->Amat[0].aux_data->max_level = 10; 70739381ba2SJed Brown ml_object->Amat[0].num_PDEs = bs; 70839381ba2SJed Brown } 70939381ba2SJed Brown 7109566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A, MAT_LOCAL, &info)); 7118a62b701SToby Isaac ml_object->Amat[0].N_nonzeros = (int)info.nz_used; 7128a62b701SToby Isaac 71339381ba2SJed Brown if (pc_ml->dim) { 71439381ba2SJed Brown PetscInt i, dim = pc_ml->dim; 71539381ba2SJed Brown ML_Aggregate_Viz_Stats *grid_info; 71639381ba2SJed Brown PetscInt nlocghost; 71739381ba2SJed Brown 7189566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 71939381ba2SJed Brown nlocghost = Aloc->cmap->n / bs; 72039381ba2SJed Brown 721e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_VizAndStats_Setup(", ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */ 72239381ba2SJed Brown grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[0].Grid; 72339381ba2SJed Brown for (i = 0; i < dim; i++) { 72439381ba2SJed Brown /* set the finest level coordinates to point to the column-order array 72539381ba2SJed Brown * in pc_ml */ 72639381ba2SJed Brown /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */ 72739381ba2SJed Brown switch (i) { 72839381ba2SJed Brown case 0: grid_info->x = pc_ml->coords + nlocghost * i; break; 72939381ba2SJed Brown case 1: grid_info->y = pc_ml->coords + nlocghost * i; break; 73039381ba2SJed Brown case 2: grid_info->z = pc_ml->coords + nlocghost * i; break; 731ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3"); 73239381ba2SJed Brown } 73339381ba2SJed Brown } 73439381ba2SJed Brown grid_info->Ndim = dim; 73539381ba2SJed Brown } 73639381ba2SJed Brown 73739381ba2SJed Brown /* repartitioning */ 73839381ba2SJed Brown if (pc_ml->Repartition) { 739e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Activate", ML_Repartition_Activate(ml_object)); 740e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Set_LargestMinMaxRatio", ML_Repartition_Set_LargestMinMaxRatio(ml_object, pc_ml->MaxMinRatio)); 741e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Set_MinPerProc", ML_Repartition_Set_MinPerProc(ml_object, pc_ml->MinPerProc)); 742e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Set_PutOnSingleProc", ML_Repartition_Set_PutOnSingleProc(ml_object, pc_ml->PutOnSingleProc)); 74339381ba2SJed Brown #if 0 /* Function not yet defined in ml-6.2 */ 74439381ba2SJed Brown /* I'm not sure what compatibility issues might crop up if we partitioned 74539381ba2SJed Brown * on the finest level, so to be safe repartition starts on the next 74639381ba2SJed Brown * finest level (reflection default behavior in 74739381ba2SJed Brown * ml_MultiLevelPreconditioner) */ 748e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1)); 74939381ba2SJed Brown #endif 75039381ba2SJed Brown 75139381ba2SJed Brown if (!pc_ml->RepartitionType) { 75239381ba2SJed Brown PetscInt i; 75339381ba2SJed Brown 7545f80ce2aSJacob Faibussowitsch PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "ML Zoltan repartitioning requires coordinates"); 755e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEZOLTAN)); 756e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Aggregate_Set_Dimensions", ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim)); 75739381ba2SJed Brown 75839381ba2SJed Brown for (i = 0; i < ml_object->ML_num_levels; i++) { 75939381ba2SJed Brown ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[i].Grid; 76039381ba2SJed Brown grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */ 76139381ba2SJed Brown /* defaults from ml_agg_info.c */ 76239381ba2SJed Brown grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */ 76339381ba2SJed Brown grid_info->zoltan_timers = 0; 76439381ba2SJed Brown grid_info->smoothing_steps = 4; /* only relevant to hypergraph / fast hypergraph */ 76539381ba2SJed Brown } 7662fa5cd67SKarl Rupp } else { 767e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEPARMETIS)); 76839381ba2SJed Brown } 76939381ba2SJed Brown } 77039381ba2SJed Brown 771b5c8bdf8SJed Brown if (pc_ml->OldHierarchy) { 772e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Gen_MGHierarchy_UsingAggregation", Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object)); 773b5c8bdf8SJed Brown } else { 774e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Gen_MultiLevelHierarchy_UsingAggregation", Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object)); 775b5c8bdf8SJed Brown } 7765f80ce2aSJacob Faibussowitsch PetscCheck(Nlevels > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Nlevels %d must > 0", Nlevels); 777573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 778aa85bbbfSHong Zhang fine_level = Nlevels - 1; 779c07bf074SBarry Smith 7809566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, Nlevels, NULL)); 781aa85bbbfSHong Zhang /* set default smoothers */ 782aa85bbbfSHong Zhang for (level = 1; level <= fine_level; level++) { 7839566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, level, &smoother)); 7849566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPRICHARDSON)); 7859566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 7869566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCSOR)); 787aa85bbbfSHong Zhang } 788d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)pc); 789dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 790d0609cedSBarry Smith PetscOptionsEnd(); 7915582bec1SHong Zhang 7929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nlevels, &gridctx)); 7932fa5cd67SKarl Rupp 7945582bec1SHong Zhang pc_ml->gridctx = gridctx; 7955582bec1SHong Zhang 7965582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 7975582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 798e14861a4SHong Zhang gridctx[fine_level].A = A; 799573998d7SHong Zhang 800e14861a4SHong Zhang level = fine_level - 1; 80143ef1857SStefano Zampini /* TODO: support for GPUs */ 802ab718edeSHong Zhang if (size == 1) { /* convert ML P, R and A into seqaij format */ 8035582bec1SHong Zhang for (mllevel = 1; mllevel < Nlevels; mllevel++) { 804e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 8059566063dSJacob Faibussowitsch PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P)); 806e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel - 1]); 8079566063dSJacob Faibussowitsch PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R)); 808573998d7SHong Zhang 809573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 8109566063dSJacob Faibussowitsch PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A)); 8115582bec1SHong Zhang level--; 8125582bec1SHong Zhang } 813ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 8145582bec1SHong Zhang for (mllevel = 1; mllevel < Nlevels; mllevel++) { 8155582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 8169566063dSJacob Faibussowitsch PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P)); 817ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel - 1]); 8189566063dSJacob Faibussowitsch PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R)); 819573998d7SHong Zhang 8205582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 8219566063dSJacob Faibussowitsch PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A)); 8225582bec1SHong Zhang level--; 8235582bec1SHong Zhang } 8245582bec1SHong Zhang } 8255582bec1SHong Zhang 826573998d7SHong Zhang /* create vectors and ksp at all levels */ 827ac346b81SHong Zhang for (level = 0; level < fine_level; level++) { 828573998d7SHong Zhang level1 = level + 1; 829708418deSStefano Zampini 8309566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(gridctx[level].A, &gridctx[level].x, &gridctx[level].b)); 8319566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(gridctx[level1].A, NULL, &gridctx[level1].r)); 8329566063dSJacob Faibussowitsch PetscCall(PCMGSetX(pc, level, gridctx[level].x)); 8339566063dSJacob Faibussowitsch PetscCall(PCMGSetRhs(pc, level, gridctx[level].b)); 8349566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, level1, gridctx[level1].r)); 835ac346b81SHong Zhang 8365582bec1SHong Zhang if (level == 0) { 8379566063dSJacob Faibussowitsch PetscCall(PCMGGetCoarseSolve(pc, &gridctx[level].ksp)); 8385582bec1SHong Zhang } else { 8399566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, level, &gridctx[level].ksp)); 840573998d7SHong Zhang } 841573998d7SHong Zhang } 8429566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, fine_level, &gridctx[fine_level].ksp)); 843573998d7SHong Zhang 844573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 845573998d7SHong Zhang for (level = 0; level < fine_level; level++) { 846573998d7SHong Zhang level1 = level + 1; 847708418deSStefano Zampini 8489566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, level1, gridctx[level].P)); 8499566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc, level1, gridctx[level].R)); 85048a46eb9SPierre Jolivet if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A)); 8519566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A)); 8525582bec1SHong Zhang } 8539566063dSJacob Faibussowitsch PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A)); 8549566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A)); 8555582bec1SHong Zhang 85639381ba2SJed Brown /* put coordinate info in levels */ 85739381ba2SJed Brown if (pc_ml->dim) { 85839381ba2SJed Brown PetscInt i, j, dim = pc_ml->dim; 85939381ba2SJed Brown PetscInt bs, nloc; 86039381ba2SJed Brown PC subpc; 86139381ba2SJed Brown PetscReal *array; 86239381ba2SJed Brown 86339381ba2SJed Brown level = fine_level; 86439381ba2SJed Brown for (mllevel = 0; mllevel < Nlevels; mllevel++) { 865ebbbbe33SJed Brown ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Amat[mllevel].to->Grid->Grid; 86639381ba2SJed Brown MPI_Comm comm = ((PetscObject)gridctx[level].A)->comm; 86739381ba2SJed Brown 8689566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(gridctx[level].A, &bs)); 8699566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(gridctx[level].A, NULL, &nloc)); 87039381ba2SJed Brown nloc /= bs; /* number of local nodes */ 87139381ba2SJed Brown 8729566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &gridctx[level].coords)); 8739566063dSJacob Faibussowitsch PetscCall(VecSetSizes(gridctx[level].coords, dim * nloc, PETSC_DECIDE)); 8749566063dSJacob Faibussowitsch PetscCall(VecSetType(gridctx[level].coords, VECMPI)); 8759566063dSJacob Faibussowitsch PetscCall(VecGetArray(gridctx[level].coords, &array)); 87639381ba2SJed Brown for (j = 0; j < nloc; j++) { 87739381ba2SJed Brown for (i = 0; i < dim; i++) { 87839381ba2SJed Brown switch (i) { 87939381ba2SJed Brown case 0: array[dim * j + i] = grid_info->x[j]; break; 88039381ba2SJed Brown case 1: array[dim * j + i] = grid_info->y[j]; break; 88139381ba2SJed Brown case 2: array[dim * j + i] = grid_info->z[j]; break; 882ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3"); 88339381ba2SJed Brown } 88439381ba2SJed Brown } 88539381ba2SJed Brown } 88639381ba2SJed Brown 88739381ba2SJed Brown /* passing coordinates to smoothers/coarse solver, should they need them */ 8889566063dSJacob Faibussowitsch PetscCall(KSPGetPC(gridctx[level].ksp, &subpc)); 8899566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates(subpc, dim, nloc, array)); 8909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gridctx[level].coords, &array)); 89139381ba2SJed Brown level--; 89239381ba2SJed Brown } 89339381ba2SJed Brown } 89439381ba2SJed Brown 895c07bf074SBarry Smith /* setupcalled is set to 0 so that MG is setup from scratch */ 896c07bf074SBarry Smith pc->setupcalled = 0; 8979566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8985582bec1SHong Zhang PetscFunctionReturn(0); 8995582bec1SHong Zhang } 9005582bec1SHong Zhang 9015582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 9025582bec1SHong Zhang /* 9035582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 9045582bec1SHong Zhang that was created with PCCreate_ML(). 9055582bec1SHong Zhang 9065582bec1SHong Zhang Input Parameter: 9075582bec1SHong Zhang . pc - the preconditioner context 9085582bec1SHong Zhang 9095582bec1SHong Zhang Application Interface Routine: PCDestroy() 9105582bec1SHong Zhang */ 9119371c9d4SSatish Balay PetscErrorCode PCDestroy_ML(PC pc) { 91201da6913SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 91301da6913SBarry Smith PC_ML *pc_ml = (PC_ML *)mg->innerctx; 9145582bec1SHong Zhang 9155582bec1SHong Zhang PetscFunctionBegin; 9169566063dSJacob Faibussowitsch PetscCall(PCReset_ML(pc)); 9179566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ml)); 9189566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 9199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 9205582bec1SHong Zhang PetscFunctionReturn(0); 9215582bec1SHong Zhang } 9225582bec1SHong Zhang 9239371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_ML(PC pc, PetscOptionItems *PetscOptionsObject) { 92439381ba2SJed Brown PetscInt indx, PrintLevel, partindx; 9255582bec1SHong Zhang const char *scheme[] = {"Uncoupled", "Coupled", "MIS", "METIS"}; 92639381ba2SJed Brown const char *part[] = {"Zoltan", "ParMETIS"}; 92739381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN) 92839381ba2SJed Brown const char *zscheme[] = {"RCB", "hypergraph", "fast_hypergraph"}; 92939381ba2SJed Brown #endif 93001da6913SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 93101da6913SBarry Smith PC_ML *pc_ml = (PC_ML *)mg->innerctx; 932b5c8bdf8SJed Brown PetscMPIInt size; 933ce94432eSBarry Smith MPI_Comm comm; 9345582bec1SHong Zhang 9355582bec1SHong Zhang PetscFunctionBegin; 9369566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 9379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 938d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ML options"); 9392fa5cd67SKarl Rupp 9405582bec1SHong Zhang PrintLevel = 0; 9415582bec1SHong Zhang indx = 0; 94239381ba2SJed Brown partindx = 0; 9432fa5cd67SKarl Rupp 9449566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_ml_PrintLevel", "Print level", "ML_Set_PrintLevel", PrintLevel, &PrintLevel, NULL)); 945e77caa6dSBarry Smith PetscStackCallExternalVoid("ML_Set_PrintLevel", ML_Set_PrintLevel(PrintLevel)); 9469566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_ml_maxNlevels", "Maximum number of levels", "None", pc_ml->MaxNlevels, &pc_ml->MaxNlevels, NULL)); 9479566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_ml_maxCoarseSize", "Maximum coarsest mesh size", "ML_Aggregate_Set_MaxCoarseSize", pc_ml->MaxCoarseSize, &pc_ml->MaxCoarseSize, NULL)); 9489566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-pc_ml_CoarsenScheme", "Aggregate Coarsen Scheme", "ML_Aggregate_Set_CoarsenScheme_*", scheme, 4, scheme[0], &indx, NULL)); 9492fa5cd67SKarl Rupp 9505582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 9512fa5cd67SKarl Rupp 9529566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_ml_DampingFactor", "P damping factor", "ML_Aggregate_Set_DampingFactor", pc_ml->DampingFactor, &pc_ml->DampingFactor, NULL)); 9539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_ml_Threshold", "Smoother drop tol", "ML_Aggregate_Set_Threshold", pc_ml->Threshold, &pc_ml->Threshold, NULL)); 9549566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm", "Method used for estimating spectral radius", "ML_Set_SpectralNormScheme_Anorm", pc_ml->SpectralNormScheme_Anorm, &pc_ml->SpectralNormScheme_Anorm, NULL)); 9559566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_Symmetrize", "Symmetrize aggregation", "ML_Set_Symmetrize", pc_ml->Symmetrize, &pc_ml->Symmetrize, NULL)); 9569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_BlockScaling", "Scale all dofs at each node together", "None", pc_ml->BlockScaling, &pc_ml->BlockScaling, NULL)); 9579566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_ml_nullspace", "Which type of null space information to use", "None", PCMLNullSpaceTypes, (PetscEnum)pc_ml->nulltype, (PetscEnum *)&pc_ml->nulltype, NULL)); 9589566063dSJacob Faibussowitsch PetscCall(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)); 9599566063dSJacob Faibussowitsch PetscCall(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)); 960b5c8bdf8SJed Brown /* 961b5c8bdf8SJed Brown The following checks a number of conditions. If we let this stuff slip by, then ML's error handling will take over. 962b5c8bdf8SJed Brown This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions. 963b5c8bdf8SJed Brown 964b5c8bdf8SJed Brown We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working 965b5c8bdf8SJed Brown combination of options and ML's exit(1) explanations don't help matters. 966b5c8bdf8SJed Brown */ 9672472a847SBarry Smith PetscCheck(pc_ml->EnergyMinimization >= -1 && pc_ml->EnergyMinimization <= 4, comm, PETSC_ERR_ARG_OUTOFRANGE, "EnergyMinimization must be in range -1..4"); 9682472a847SBarry Smith PetscCheck(pc_ml->EnergyMinimization != 4 || size == 1, comm, PETSC_ERR_SUP, "Energy minimization type 4 does not work in parallel"); 9699566063dSJacob Faibussowitsch if (pc_ml->EnergyMinimization == 4) PetscCall(PetscInfo(pc, "Mandel's energy minimization scheme is experimental and broken in ML-6.2\n")); 97048a46eb9SPierre Jolivet if (pc_ml->EnergyMinimization) PetscCall(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol", "Energy minimization drop tolerance", "None", pc_ml->EnergyMinimizationDropTol, &pc_ml->EnergyMinimizationDropTol, NULL)); 971b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == 2) { 972b5c8bdf8SJed Brown /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */ 9739566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap", "Use cheaper variant of norm type 2", "None", pc_ml->EnergyMinimizationCheap, &pc_ml->EnergyMinimizationCheap, NULL)); 974b5c8bdf8SJed Brown } 975b5c8bdf8SJed Brown /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */ 976b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE; 9779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_KeepAggInfo", "Allows the preconditioner to be reused, or auxilliary matrices to be generated", "None", pc_ml->KeepAggInfo, &pc_ml->KeepAggInfo, NULL)); 978b5c8bdf8SJed Brown /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */ 979b5c8bdf8SJed Brown if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE; 9809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_Reusable", "Store intermedaiate data structures so that the multilevel hierarchy is reusable", "None", pc_ml->Reusable, &pc_ml->Reusable, NULL)); 981b5c8bdf8SJed Brown /* 982b5c8bdf8SJed Brown ML's C API is severely underdocumented and lacks significant functionality. The C++ API calls 983b5c8bdf8SJed Brown ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function 984b5c8bdf8SJed Brown ML_Gen_MGHierarchy_UsingAggregation(). This modification, however, does not provide a strict superset of the 985b5c8bdf8SJed Brown functionality in the old function, so some users may still want to use it. Note that many options are ignored in 986b5c8bdf8SJed Brown this context, but ML doesn't provide a way to find out which ones. 987b5c8bdf8SJed Brown */ 9889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_OldHierarchy", "Use old routine to generate hierarchy", "None", pc_ml->OldHierarchy, &pc_ml->OldHierarchy, NULL)); 9899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy", "ML_Repartition_Activate", pc_ml->Repartition, &pc_ml->Repartition, NULL)); 99039381ba2SJed Brown if (pc_ml->Repartition) { 9919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes", "ML_Repartition_Set_LargestMinMaxRatio", pc_ml->MaxMinRatio, &pc_ml->MaxMinRatio, NULL)); 9929566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size", "ML_Repartition_Set_MinPerProc", pc_ml->MinPerProc, &pc_ml->MinPerProc, NULL)); 9939566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor", "ML_Repartition_Set_PutOnSingleProc", pc_ml->PutOnSingleProc, &pc_ml->PutOnSingleProc, NULL)); 99439381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN) 99539381ba2SJed Brown partindx = 0; 9969566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[0], &partindx, NULL)); 9972fa5cd67SKarl Rupp 99839381ba2SJed Brown pc_ml->RepartitionType = partindx; 99939381ba2SJed Brown if (!partindx) { 10005572b5bbSJed Brown PetscInt zindx = 0; 10012fa5cd67SKarl Rupp 10029566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use", "None", zscheme, 3, zscheme[0], &zindx, NULL)); 10032fa5cd67SKarl Rupp 100439381ba2SJed Brown pc_ml->ZoltanScheme = zindx; 100539381ba2SJed Brown } 100639381ba2SJed Brown #else 100739381ba2SJed Brown partindx = 1; 10089566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[1], &partindx, NULL)); 1009e6b1cc6bSSatish Balay pc_ml->RepartitionType = partindx; 10105f80ce2aSJacob Faibussowitsch PetscCheck(partindx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP_SYS, "ML not compiled with Zoltan"); 101139381ba2SJed Brown #endif 10129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_ml_Aux", "Aggregate using auxiliary coordinate-based laplacian", "None", pc_ml->Aux, &pc_ml->Aux, NULL)); 10139566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_ml_AuxThreshold", "Auxiliary smoother drop tol", "None", pc_ml->AuxThreshold, &pc_ml->AuxThreshold, NULL)); 101439381ba2SJed Brown } 1015d0609cedSBarry Smith PetscOptionsHeadEnd(); 10165582bec1SHong Zhang PetscFunctionReturn(0); 10175582bec1SHong Zhang } 10185582bec1SHong Zhang 10195582bec1SHong Zhang /* 10205582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 10215582bec1SHong Zhang and sets this as the private data within the generic preconditioning 10225582bec1SHong Zhang context, PC, that was created within PCCreate(). 10235582bec1SHong Zhang 10245582bec1SHong Zhang Input Parameter: 10255582bec1SHong Zhang . pc - the preconditioner context 10265582bec1SHong Zhang 10275582bec1SHong Zhang Application Interface Routine: PCCreate() 10285582bec1SHong Zhang */ 10295582bec1SHong Zhang 10305582bec1SHong Zhang /*MC 1031*f1580f4eSBarry Smith PCML - Use the SNL ML algebraic multigrid preconditioner. 10325582bec1SHong Zhang 1033*f1580f4eSBarry Smith Options Database Keys: 10342612397fSMatthew G. Knepley Multigrid options(inherited): 1035*f1580f4eSBarry Smith + -pc_mg_cycles <1> - 1 for V cycle, 2 for W-cycle (`MGSetCycles()`) 1036*f1580f4eSBarry Smith . -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (`PCMGSetDistinctSmoothUp()`) 1037a2b725a8SWilliam Gropp - -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1038a2b725a8SWilliam Gropp 1039*f1580f4eSBarry Smith ML Options Database Key: 1040*f1580f4eSBarry Smith + -pc_ml_PrintLevel <0> - Print level (`ML_Set_PrintLevel()`) 1041a2b725a8SWilliam Gropp . -pc_ml_maxNlevels <10> - Maximum number of levels (None) 1042*f1580f4eSBarry Smith . -pc_ml_maxCoarseSize <1> - Maximum coarsest mesh size (`ML_Aggregate_Set_MaxCoarseSize()`) 1043a2b725a8SWilliam Gropp . -pc_ml_CoarsenScheme <Uncoupled> - (one of) Uncoupled Coupled MIS METIS 1044*f1580f4eSBarry Smith . -pc_ml_DampingFactor <1.33333> - P damping factor (`ML_Aggregate_Set_DampingFactor()`) 1045*f1580f4eSBarry Smith . -pc_ml_Threshold <0> - Smoother drop tol (`ML_Aggregate_Set_Threshold()`) 1046*f1580f4eSBarry Smith . -pc_ml_SpectralNormScheme_Anorm <false> - Method used for estimating spectral radius (`ML_Set_SpectralNormScheme_Anorm()`) 1047*f1580f4eSBarry Smith . -pc_ml_repartition <false> - Allow ML to repartition levels of the hierarchy (`ML_Repartition_Activate()`) 1048*f1580f4eSBarry Smith . -pc_ml_repartitionMaxMinRatio <1.3> - Acceptable ratio of repartitioned sizes (`ML_Repartition_Set_LargestMinMaxRatio()`) 1049*f1580f4eSBarry Smith . -pc_ml_repartitionMinPerProc <512> - Smallest repartitioned size (`ML_Repartition_Set_MinPerProc()`) 1050*f1580f4eSBarry Smith . -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (`ML_Repartition_Set_PutOnSingleProc()`) 1051*f1580f4eSBarry Smith . -pc_ml_repartitionType <Zoltan> - Repartitioning library to use (`ML_Repartition_Set_Partitioner()`) 1052a2b725a8SWilliam Gropp . -pc_ml_repartitionZoltanScheme <RCB> - Repartitioning scheme to use (None) 1053147403d9SBarry Smith . -pc_ml_Aux <false> - Aggregate using auxiliary coordinate-based Laplacian (None) 1054a2b725a8SWilliam Gropp - -pc_ml_AuxThreshold <0.0> - Auxiliary smoother drop tol (None) 10555582bec1SHong Zhang 10565582bec1SHong Zhang Level: intermediate 10575582bec1SHong Zhang 1058*f1580f4eSBarry Smith Developer Note: 1059*f1580f4eSBarry Smith The coarser grid matrices and restriction/interpolation 1060*f1580f4eSBarry Smith operators are computed by ML, with the matrices coverted to PETSc matrices in `MATAIJ` format 1061*f1580f4eSBarry Smith and the restriction/interpolation operators wrapped as PETSc shell matrices. 1062*f1580f4eSBarry Smith 1063*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCMG`, `PCHYPRE`, `PCGAMG`, 1064db781477SPatrick Sanan `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `MPSetCycles()`, `PCMGSetDistinctSmoothUp()`, 1065db781477SPatrick Sanan `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`, 1066db781477SPatrick Sanan `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, 1067db781477SPatrick Sanan `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()` 10685582bec1SHong Zhang M*/ 10695582bec1SHong Zhang 10709371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc) { 10715582bec1SHong Zhang PC_ML *pc_ml; 107201da6913SBarry Smith PC_MG *mg; 10735582bec1SHong Zhang 10745582bec1SHong Zhang PetscFunctionBegin; 1075573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 10769566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */ 10779566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCML)); 1078e0f5d30fSBarry Smith /* Since PCMG tries to use DM assocated with PC must delete it */ 10799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&pc->dm)); 10809566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 1081e0f5d30fSBarry Smith mg = (PC_MG *)pc->data; 10825582bec1SHong Zhang 10835582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 10849566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &pc_ml)); 108501da6913SBarry Smith mg->innerctx = pc_ml; 10865582bec1SHong Zhang 1087573998d7SHong Zhang pc_ml->ml_object = 0; 1088573998d7SHong Zhang pc_ml->agg_object = 0; 1089573998d7SHong Zhang pc_ml->gridctx = 0; 1090573998d7SHong Zhang pc_ml->PetscMLdata = 0; 1091573998d7SHong Zhang pc_ml->Nlevels = -1; 1092573998d7SHong Zhang pc_ml->MaxNlevels = 10; 1093573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 10943751b4bdSBarry Smith pc_ml->CoarsenScheme = 1; 1095573998d7SHong Zhang pc_ml->Threshold = 0.0; 1096573998d7SHong Zhang pc_ml->DampingFactor = 4.0 / 3.0; 1097573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 1098573998d7SHong Zhang pc_ml->size = 0; 109939381ba2SJed Brown pc_ml->dim = 0; 110039381ba2SJed Brown pc_ml->nloc = 0; 110139381ba2SJed Brown pc_ml->coords = 0; 110239381ba2SJed Brown pc_ml->Repartition = PETSC_FALSE; 110339381ba2SJed Brown pc_ml->MaxMinRatio = 1.3; 110439381ba2SJed Brown pc_ml->MinPerProc = 512; 110539381ba2SJed Brown pc_ml->PutOnSingleProc = 5000; 110639381ba2SJed Brown pc_ml->RepartitionType = 0; 110739381ba2SJed Brown pc_ml->ZoltanScheme = 0; 110839381ba2SJed Brown pc_ml->Aux = PETSC_FALSE; 110939381ba2SJed Brown pc_ml->AuxThreshold = 0.0; 111039381ba2SJed Brown 111139381ba2SJed Brown /* allow for coordinates to be passed */ 11129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_ML)); 1113573998d7SHong Zhang 11145582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 11155582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 11165582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 1117a06653b4SBarry Smith pc->ops->reset = PCReset_ML; 11185582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 11195582bec1SHong Zhang PetscFunctionReturn(0); 11205582bec1SHong Zhang } 1121