15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4af0996ceSBarry Smith #include <petsc/private/matimpl.h> 5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 618c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/ 7f96513f1SMatthew G Knepley 8c9567895SMark #if defined(PETSC_HAVE_CUDA) 9c9567895SMark #include <cuda_runtime.h> 10c9567895SMark #endif 11c9567895SMark 12c9567895SMark #if defined(PETSC_HAVE_HIP) 13c9567895SMark #include <hip/hip_runtime.h> 14c9567895SMark #endif 15c9567895SMark 16849bee69SMark Adams PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET]; 174555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3]; 180cbbd2e1SMark F. Adams 19849bee69SMark Adams // #define GAMG_STAGES 204555aa8cSStefano Zampini #if defined(GAMG_STAGES) 2118c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS]; 22b4fbaa2aSMark F. Adams #endif 23f96513f1SMatthew G Knepley 240a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL; 253e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 269d5b6da9SMark F. Adams 27d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset_GAMG(PC pc) 28d71ae5a4SJacob Faibussowitsch { 29d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 30d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 31d3d6bff4SMark F. Adams 32d3d6bff4SMark F. Adams PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 341c1aac46SBarry Smith pc_gamg->data_sz = 0; 359566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->orig_data)); 365f80ce2aSJacob Faibussowitsch for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) { 3718c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 3818c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 3918c3aa7eSMark } 4018c3aa7eSMark pc_gamg->emin = 0; 4118c3aa7eSMark pc_gamg->emax = 0; 42a2f3521dSMark F. Adams PetscFunctionReturn(0); 43a2f3521dSMark F. Adams } 44a2f3521dSMark F. Adams 455b89ad90SMark F. Adams /* 46c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 47a147abb0SMark F. Adams of active processors. 485b89ad90SMark F. Adams 495b89ad90SMark F. Adams Input Parameter: 50a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 51a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 529d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 53c5bfad50SMark F. Adams . cr_bs - coarse block size 543530afc2SMark F. Adams In/Output Parameter: 55a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 56afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 5711e60469SMark F. Adams Output Parameter: 583530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 595b89ad90SMark F. Adams */ 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc, Mat Amat_fine, PetscInt cr_bs, Mat *a_P_inout, Mat *a_Amat_crs, PetscMPIInt *a_nactive_proc, IS *Pcolumnperm, PetscBool is_last) 61d71ae5a4SJacob Faibussowitsch { 629d5b6da9SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 63486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 64a2f3521dSMark F. Adams Mat Cmat, Pold = *a_P_inout; 653b4367a7SBarry Smith MPI_Comm comm; 66c5df96a5SBarry Smith PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc; 673ae0bb68SMark Adams PetscInt ncrs_eq, ncrs, f_bs; 685b89ad90SMark F. Adams 695b89ad90SMark F. Adams PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm)); 719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 739566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat_fine, &f_bs)); 74849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 759566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 769566063dSJacob Faibussowitsch PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat)); 779566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 78849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 79038e3b61SMark F. Adams 80ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 81ce7c7f2fSMark Adams 823ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 839566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL)); 843ae0bb68SMark Adams if (pc_gamg->data_cell_rows > 0) { 853ae0bb68SMark Adams ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows; 8673911c69SBarry Smith } else { 873ae0bb68SMark Adams PetscInt bs; 889566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Cmat, &bs)); 893ae0bb68SMark Adams ncrs = ncrs_eq / bs; 903ae0bb68SMark Adams } 91c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 92c9567895SMark if (pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0 && PetscDefined(HAVE_CUDA) && pc_gamg->current_level == 0) { /* 0 turns reducing to 1 process/device on; do for HIP, etc. */ 93c9567895SMark #if defined(PETSC_HAVE_CUDA) 94c9567895SMark PetscShmComm pshmcomm; 95c9567895SMark PetscMPIInt locrank; 96c9567895SMark MPI_Comm loccomm; 97c9567895SMark PetscInt s_nnodes, r_nnodes, new_new_size; 98c9567895SMark cudaError_t cerr; 99c9567895SMark int devCount; 1009566063dSJacob Faibussowitsch PetscCall(PetscShmCommGet(comm, &pshmcomm)); 1019566063dSJacob Faibussowitsch PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm)); 1029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(loccomm, &locrank)); 103c9567895SMark s_nnodes = !locrank; 1049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm)); 10563a3b9bcSJacob Faibussowitsch PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes); 106c9567895SMark devCount = 0; 107c9567895SMark cerr = cudaGetDeviceCount(&devCount); 108c9567895SMark cudaGetLastError(); /* Reset the last error */ 109c9567895SMark if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 110c9567895SMark new_new_size = r_nnodes * devCount; 111c9567895SMark new_size = new_new_size; 11263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Fine grid with Cuda. %" PetscInt_FMT " nodes. Change new active set size %d --> %d (devCount=%d #nodes=%" PetscInt_FMT ")\n", ((PetscObject)pc)->prefix, r_nnodes, nactive, new_size, devCount, r_nnodes)); 113c9567895SMark } else { 11463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.", ((PetscObject)pc)->prefix)); 115c9567895SMark goto HEURISTIC; 116c9567895SMark } 117c9567895SMark #else 118c9567895SMark SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here"); 119c9567895SMark #endif 120c9567895SMark } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 12163a3b9bcSJacob Faibussowitsch PetscCheck(nactive % pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of active process %d wrt reduction factor %" PetscInt_FMT, nactive, pc_gamg->level_reduction_factors[pc_gamg->current_level]); 122c9567895SMark new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level]; 12363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Manually setting reduction to %d active processes (%d/%" PetscInt_FMT ")\n", ((PetscObject)pc)->prefix, new_size, nactive, pc_gamg->level_reduction_factors[pc_gamg->current_level])); 124c9567895SMark } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 125c9567895SMark new_size = 1; 1269566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size)); 127c9567895SMark } else { 128472110cdSMark F. Adams PetscInt ncrs_eq_glob; 129c9567895SMark #if defined(PETSC_HAVE_CUDA) 130c9567895SMark HEURISTIC: 131c9567895SMark #endif 1329566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL)); 133a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1347f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 135c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 1369566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size)); 137a2f3521dSMark F. Adams } 1382e3501ffSMark Adams if (new_size == nactive) { 139ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 140ce7c7f2fSMark Adams if (new_size < size) { 141ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 1429566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n", ((PetscObject)pc)->prefix, nactive)); 143ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1449566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 1459566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 146ce7c7f2fSMark Adams } 147ce7c7f2fSMark Adams } 148ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1492e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 150192c0e8bSMark Adams PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1; 151885364a3SMark Adams IS is_eq_newproc, is_eq_num, is_eq_num_prim, new_eq_indices; 152849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 15371959b99SBarry Smith nloc_old = ncrs_eq / cr_bs; 15463a3b9bcSJacob Faibussowitsch PetscCheck(ncrs_eq % cr_bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ncrs_eq %" PetscInt_FMT " not divisible by cr_bs %" PetscInt_FMT, ncrs_eq, cr_bs); 155ce7c7f2fSMark Adams /* get new_size and rfactor */ 156ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 157ce7c7f2fSMark Adams /* find factor */ 158ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 159ce7c7f2fSMark Adams else { 160ce7c7f2fSMark Adams PetscReal best_fact = 0.; 161ce7c7f2fSMark Adams jj = -1; 162ce7c7f2fSMark Adams for (kk = 1; kk <= size; kk++) { 163ce7c7f2fSMark Adams if (!(size % kk)) { /* a candidate */ 164ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size; 165ce7c7f2fSMark Adams if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */ 166ce7c7f2fSMark Adams if (fact > best_fact) { 1679371c9d4SSatish Balay best_fact = fact; 1689371c9d4SSatish Balay jj = kk; 169ce7c7f2fSMark Adams } 170ce7c7f2fSMark Adams } 171ce7c7f2fSMark Adams } 172ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 173ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 174ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 175ce7c7f2fSMark Adams else expand_factor = rfactor; 176ce7c7f2fSMark Adams } 177ce7c7f2fSMark Adams new_size = size / rfactor; /* make new size one that is factor */ 1784cdfd227SMark if (new_size == nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1794cdfd227SMark *a_Amat_crs = Cmat; 18063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, new_size, ncrs_eq)); 181849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 182ce7c7f2fSMark Adams PetscFunctionReturn(0); 183ce7c7f2fSMark Adams } 184ce7c7f2fSMark Adams } 185a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 1869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &counts)); 187849bee69SMark Adams if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 1885a9b9e01SMark F. Adams Mat adj; 189849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 19063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Repartition: size (active): %d --> %d, %" PetscInt_FMT " local equations, using %s process layout\n", ((PetscObject)pc)->prefix, *a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread")); 191a2f3521dSMark F. Adams /* get 'adj' */ 192c5bfad50SMark F. Adams if (cr_bs == 1) { 1939566063dSJacob Faibussowitsch PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 194806fa848SBarry Smith } else { 195a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 196eb07cef2SMark F. Adams Mat tMat; 197a2f3521dSMark F. Adams PetscInt Istart_crs, Iend_crs, ncols, jj, Ii; 198b4fbaa2aSMark F. Adams const PetscScalar *vals; 199b4fbaa2aSMark F. Adams const PetscInt *idx; 200a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 20139d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 202d9558ea9SBarry Smith MatType mtype; 203b4fbaa2aSMark F. Adams 2049566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz)); 2059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs)); 2069566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &M, &N)); 207c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 2089566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat, Ii, &ncols, NULL, NULL)); 209c5bfad50SMark F. Adams d_nnz[jj] = ncols / cr_bs; 210c5bfad50SMark F. Adams o_nnz[jj] = ncols / cr_bs; 2119566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL)); 2123ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 2133ae0bb68SMark Adams if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs; 21458471d46SMark F. Adams } 2156876a03eSMark F. Adams 2169566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat_fine, &mtype)); 2179566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &tMat)); 2189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetType(tMat, mtype)); 2209566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(tMat, 0, d_nnz)); 2219566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz)); 2229566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz, o_nnz)); 223eb07cef2SMark F. Adams 224a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 225c5bfad50SMark F. Adams PetscInt dest_row = ii / cr_bs; 2269566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat, ii, &ncols, &idx, &vals)); 227eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 228c5bfad50SMark F. Adams PetscInt dest_col = idx[jj] / cr_bs; 229eb07cef2SMark F. Adams PetscScalar v = 1.0; 2309566063dSJacob Faibussowitsch PetscCall(MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES)); 231eb07cef2SMark F. Adams } 2329566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat, ii, &ncols, &idx, &vals)); 233eb07cef2SMark F. Adams } 2349566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY)); 2359566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY)); 236eb07cef2SMark F. Adams 237b4fbaa2aSMark F. Adams if (llev++ == -1) { 2389371c9d4SSatish Balay PetscViewer viewer; 2399371c9d4SSatish Balay char fname[32]; 24063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev)); 2413b4367a7SBarry Smith PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer); 2429566063dSJacob Faibussowitsch PetscCall(MatView(tMat, viewer)); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 244b4fbaa2aSMark F. Adams } 2459566063dSJacob Faibussowitsch PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 2469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tMat)); 247a2f3521dSMark F. Adams } /* create 'adj' */ 248f150b916SMark F. Adams 249a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2505a9b9e01SMark F. Adams char prefix[256]; 2515a9b9e01SMark F. Adams const char *pcpre; 252b4fbaa2aSMark F. Adams const PetscInt *is_idx; 253b4fbaa2aSMark F. Adams MatPartitioning mpart; 254a4b7d37bSMark F. Adams IS proc_is; 2552f03bc48SMark F. Adams 2569566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm, &mpart)); 2579566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 2589566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 2619566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 2629566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart, new_size)); 2639566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart, &proc_is)); 2649566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 2655a9b9e01SMark F. Adams 2665ef31b24SMark F. Adams /* collect IS info */ 2679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx)); 2689566063dSJacob Faibussowitsch PetscCall(ISGetIndices(proc_is, &is_idx)); 269a2f3521dSMark F. Adams for (kk = jj = 0; kk < nloc_old; kk++) { 2709371c9d4SSatish Balay for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ } 2715ef31b24SMark F. Adams } 2729566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(proc_is, &is_idx)); 2739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&proc_is)); 2745ef31b24SMark F. Adams } 2759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 2765a9b9e01SMark F. Adams 2779566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(newproc_idx)); 279849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 28031cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 281ce7c7f2fSMark Adams PetscInt targetPE; 28208401ef6SPierre Jolivet PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen"); 28363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq)); 284ce7c7f2fSMark Adams targetPE = (rank / rfactor) * expand_factor; 2859566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc)); 286a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 287e33ef3b1SMark F. Adams 28811e60469SMark F. Adams /* 289a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 29011e60469SMark F. Adams */ 2919566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 2927700e67bSMark Adams is_eq_num_prim = is_eq_num; 29311e60469SMark F. Adams /* 294a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 29511e60469SMark F. Adams */ 2969566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 297c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 2989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_newproc)); 299ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new / cr_bs; 300a2f3521dSMark F. Adams 3019566063dSJacob Faibussowitsch PetscCall(PetscFree(counts)); 3026aad120cSJose E. Roman /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxiliary data into some complex abstracted thing */ 303885364a3SMark Adams { 304885364a3SMark Adams Vec src_crd, dest_crd; 305885364a3SMark Adams const PetscInt *idx, ndata_rows = pc_gamg->data_cell_rows, ndata_cols = pc_gamg->data_cell_cols, node_data_sz = ndata_rows * ndata_cols; 306885364a3SMark Adams VecScatter vecscat; 307885364a3SMark Adams PetscScalar *array; 308885364a3SMark Adams IS isscat; 309a2f3521dSMark F. Adams /* move data (for primal equations only) */ 31022063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3119566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &dest_crd)); 3129566063dSJacob Faibussowitsch PetscCall(VecSetSizes(dest_crd, node_data_sz * ncrs_new, PETSC_DECIDE)); 3139566063dSJacob Faibussowitsch PetscCall(VecSetType(dest_crd, VECSTANDARD)); /* this is needed! */ 31411e60469SMark F. Adams /* 3159d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 316c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 31711e60469SMark F. Adams */ 3189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncrs * node_data_sz, &tidx)); 3199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_eq_num_prim, &idx)); 3203ae0bb68SMark Adams for (ii = 0, jj = 0; ii < ncrs; ii++) { 321c5bfad50SMark F. Adams PetscInt id = idx[ii * cr_bs] / cr_bs; /* get node back */ 322a2f3521dSMark F. Adams for (kk = 0; kk < node_data_sz; kk++, jj++) tidx[jj] = id * node_data_sz + kk; 32311e60469SMark F. Adams } 3249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_eq_num_prim, &idx)); 3259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, node_data_sz * ncrs, tidx, PETSC_COPY_VALUES, &isscat)); 3269566063dSJacob Faibussowitsch PetscCall(PetscFree(tidx)); 32711e60469SMark F. Adams /* 32811e60469SMark F. Adams Create a vector to contain the original vertex information for each element 32911e60469SMark F. Adams */ 3309566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz * ncrs, &src_crd)); 3319d5b6da9SMark F. Adams for (jj = 0; jj < ndata_cols; jj++) { 3323ae0bb68SMark Adams const PetscInt stride0 = ncrs * pc_gamg->data_cell_rows; 3333ae0bb68SMark Adams for (ii = 0; ii < ncrs; ii++) { 3349d5b6da9SMark F. Adams for (kk = 0; kk < ndata_rows; kk++) { 335a2f3521dSMark F. Adams PetscInt ix = ii * ndata_rows + kk + jj * stride0, jx = ii * node_data_sz + kk * ndata_cols + jj; 336c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 3379566063dSJacob Faibussowitsch PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES)); 338d3d6bff4SMark F. Adams } 339038e3b61SMark F. Adams } 340eb07cef2SMark F. Adams } 3419566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(src_crd)); 3429566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(src_crd)); 34311e60469SMark F. Adams /* 34411e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 34511e60469SMark F. Adams to the correct processor 34611e60469SMark F. Adams */ 3479566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat)); 3489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isscat)); 3499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 3509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 3519566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vecscat)); 3529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&src_crd)); 35311e60469SMark F. Adams /* 35411e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 35511e60469SMark F. Adams */ 3569566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(node_data_sz * ncrs_new, &pc_gamg->data)); 3582fa5cd67SKarl Rupp 3593ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz * ncrs_new; 3603ae0bb68SMark Adams strideNew = ncrs_new * ndata_rows; 3612fa5cd67SKarl Rupp 3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(dest_crd, &array)); 3639d5b6da9SMark F. Adams for (jj = 0; jj < ndata_cols; jj++) { 3643ae0bb68SMark Adams for (ii = 0; ii < ncrs_new; ii++) { 3659d5b6da9SMark F. Adams for (kk = 0; kk < ndata_rows; kk++) { 366a2f3521dSMark F. Adams PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj; 367c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 368d3d6bff4SMark F. Adams } 369038e3b61SMark F. Adams } 370038e3b61SMark F. Adams } 3719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(dest_crd, &array)); 3729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dest_crd)); 373885364a3SMark Adams } 374a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 375849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 37611e60469SMark F. Adams /* 3777dae84e0SHong Zhang Invert for MatCreateSubMatrix 37811e60469SMark F. Adams */ 3799566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices)); 3809566063dSJacob Faibussowitsch PetscCall(ISSort(new_eq_indices)); /* is this needed? */ 3819566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(new_eq_indices, cr_bs)); 3829371c9d4SSatish Balay if (is_eq_num != is_eq_num_prim) { PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ } 3833cb8563fSToby Isaac if (Pcolumnperm) { 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)new_eq_indices)); 3853cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3863cb8563fSToby Isaac } 3879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_num)); 388849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 389849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 390849bee69SMark Adams 391a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 392a2f3521dSMark F. Adams { 393a2f3521dSMark F. Adams Mat mat; 394b94d7dedSBarry Smith PetscBool isset, isspd, isher; 39590db8557SMark Adams #if !defined(PETSC_USE_COMPLEX) 396b94d7dedSBarry Smith PetscBool issym; 397b94d7dedSBarry Smith #endif 398b94d7dedSBarry Smith 399b94d7dedSBarry Smith PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 400b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ? 401b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd)); 402b94d7dedSBarry Smith else { 403b94d7dedSBarry Smith PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher)); 404b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher)); 405b94d7dedSBarry Smith else { 406b94d7dedSBarry Smith #if !defined(PETSC_USE_COMPLEX) 407b94d7dedSBarry Smith PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym)); 408b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym)); 40990db8557SMark Adams #endif 41090db8557SMark Adams } 41190db8557SMark Adams } 412a2f3521dSMark F. Adams *a_Amat_crs = mat; 413a2f3521dSMark F. Adams } 4149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cmat)); 415a2f3521dSMark F. Adams 416849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 41711e60469SMark F. Adams /* prolongator */ 41811e60469SMark F. Adams { 41911e60469SMark F. Adams IS findices; 420a2f3521dSMark F. Adams PetscInt Istart, Iend; 421a2f3521dSMark F. Adams Mat Pnew; 42262294041SBarry Smith 4239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 424849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4259566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices)); 4269566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(findices, f_bs)); 4279566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 4289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&findices)); 4299566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 430c5bfad50SMark F. Adams 431849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4329566063dSJacob Faibussowitsch PetscCall(MatDestroy(a_P_inout)); 433a2f3521dSMark F. Adams 434a2f3521dSMark F. Adams /* output - repartitioned */ 435a2f3521dSMark F. Adams *a_P_inout = Pnew; 436e33ef3b1SMark F. Adams } 4379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&new_eq_indices)); 4385b89ad90SMark F. Adams 439c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 440ce7c7f2fSMark Adams 441ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 442ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 443ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4448bca76a6SMark Adams static PetscInt llev = 2; 44563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Pinning level %" PetscInt_FMT " to the CPU\n", ((PetscObject)pc)->prefix, llev++)); 446ce7c7f2fSMark Adams #endif 4479566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 4489566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 449adf5291fSStefano Zampini if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 450ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 451ce7c7f2fSMark Adams PetscMPIInt size; 4529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 453ce7c7f2fSMark Adams if (size > 1) { 454ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *p = (Mat_MPIAIJ *)P->data; 4559566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(a->lvec, PETSC_TRUE)); 4569566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(p->lvec, PETSC_TRUE)); 457ce7c7f2fSMark Adams } 458ce7c7f2fSMark Adams } 459ce7c7f2fSMark Adams } 460849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 461849bee69SMark Adams } // processor reduce 4625b89ad90SMark F. Adams PetscFunctionReturn(0); 4635b89ad90SMark F. Adams } 4645b89ad90SMark F. Adams 465bae903cbSmarkadams4 // used in GEO 466d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2) 467d71ae5a4SJacob Faibussowitsch { 4684b1575e2SStefano Zampini const char *prefix; 4694b1575e2SStefano Zampini char addp[32]; 4704b1575e2SStefano Zampini PC_MG *mg = (PC_MG *)a_pc->data; 4714b1575e2SStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 4724b1575e2SStefano Zampini 4734b1575e2SStefano Zampini PetscFunctionBegin; 4749566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(a_pc, &prefix)); 47563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1)); 4769566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2)); 4779566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(*Gmat2, prefix)); 47863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level)); 4799566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(*Gmat2, addp)); 480b94d7dedSBarry Smith if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) { 4819566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB)); 482b4da3a1bSStefano Zampini } else { 4839566063dSJacob Faibussowitsch PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 4849566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB)); 485b4da3a1bSStefano Zampini } 4869566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(*Gmat2)); 4879566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 4889566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(*Gmat2)); 4899566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 4909566063dSJacob Faibussowitsch PetscCall(MatProductClear(*Gmat2)); 4914b1575e2SStefano Zampini /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 4924b1575e2SStefano Zampini (*Gmat2)->assembled = PETSC_TRUE; 4934b1575e2SStefano Zampini PetscFunctionReturn(0); 4944b1575e2SStefano Zampini } 4954b1575e2SStefano Zampini 4965b89ad90SMark F. Adams /* 4975b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4985b89ad90SMark F. Adams by setting data structures and options. 4995b89ad90SMark F. Adams 5005b89ad90SMark F. Adams Input Parameter: 5015b89ad90SMark F. Adams . pc - the preconditioner context 5025b89ad90SMark F. Adams 5035b89ad90SMark F. Adams */ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_GAMG(PC pc) 505d71ae5a4SJacob Faibussowitsch { 5069d5b6da9SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 5075b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5082adcac29SMark F. Adams Mat Pmat = pc->pmat; 50918c3aa7eSMark PetscInt fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS]; 5103b4367a7SBarry Smith MPI_Comm comm; 511c5df96a5SBarry Smith PetscMPIInt rank, size, nactivepe; 51218c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS]; 51318c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 514a2f3521dSMark F. Adams PetscLogDouble nnz0 = 0., nnztot = 0.; 515569f4572SMark Adams MatInfo info; 516171cca9aSMark Adams PetscBool is_last = PETSC_FALSE; 5175ef31b24SMark F. Adams 5185b89ad90SMark F. Adams PetscFunctionBegin; 5199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 5209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 522849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 5238abdc6daSStefano Zampini if (pc->setupcalled) { 5248abdc6daSStefano Zampini if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 525878e152fSMark F. Adams /* reset everything */ 5269566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 527878e152fSMark F. Adams pc->setupcalled = 0; 528806fa848SBarry Smith } else { 52984d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 53003a628feSMark F. Adams /* just do Galerkin grids */ 53158471d46SMark F. Adams Mat B, dA, dB; 5329d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 5334555aa8cSStefano Zampini PetscInt gl; 53458471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5359566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB)); 53658471d46SMark F. Adams /* (re)set to get dirty flag */ 5379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB)); 53858471d46SMark F. Adams 5394555aa8cSStefano Zampini for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) { 5408abdc6daSStefano Zampini MatReuse reuse = MAT_INITIAL_MATRIX; 541849bee69SMark Adams #if defined(GAMG_STAGES) 542849bee69SMark Adams PetscCall(PetscLogStagePush(gamg_stages[gl])); 543849bee69SMark Adams #endif 5448abdc6daSStefano Zampini /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 5459566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B)); 5468abdc6daSStefano Zampini if (B->product) { 547ad540459SPierre Jolivet if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX; 5488abdc6daSStefano Zampini } 5499566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 5508abdc6daSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 55163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 5528abdc6daSStefano Zampini } else { 55363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 5548abdc6daSStefano Zampini } 5559566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 5569566063dSJacob Faibussowitsch PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B)); 5579566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 55863b77682SMark Adams if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 5599566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B)); 56058471d46SMark F. Adams dB = B; 561849bee69SMark Adams #if defined(GAMG_STAGES) 562849bee69SMark Adams PetscCall(PetscLogStagePop()); 563849bee69SMark Adams #endif 56458471d46SMark F. Adams } 5655f8cf99dSMark F. Adams } 566d5280255SMark F. Adams 5679566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 568849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 56958471d46SMark F. Adams PetscFunctionReturn(0); 570eb07cef2SMark F. Adams } 571878e152fSMark F. Adams } 572f6536408SMark F. Adams 573878e152fSMark F. Adams if (!pc_gamg->data) { 574878e152fSMark F. Adams if (pc_gamg->orig_data) { 5759566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 5769566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 5772fa5cd67SKarl Rupp 578878e152fSMark F. Adams pc_gamg->data_sz = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols; 579878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 580878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5812fa5cd67SKarl Rupp 5829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 583878e152fSMark F. Adams for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 584806fa848SBarry Smith } else { 5855f80ce2aSJacob Faibussowitsch PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data"); 5869566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat)); 5879d5b6da9SMark F. Adams } 588878e152fSMark F. Adams } 589878e152fSMark F. Adams 590878e152fSMark F. Adams /* cache original data for reuse */ 5911c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 5929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 593878e152fSMark F. Adams for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 594878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 595878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 596878e152fSMark F. Adams } 597038e3b61SMark F. Adams 598302f38e8SMark F. Adams /* get basic dims */ 5999566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 6009566063dSJacob Faibussowitsch PetscCall(MatGetSize(Pmat, &M, &N)); 60184d3f75bSMark F. Adams 6029566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */ 603569f4572SMark Adams nnz0 = info.nz_used; 604569f4572SMark Adams nnztot = info.nz_used; 605bae903cbSmarkadams4 PetscCall(PetscInfo(pc, "%s: level %d) N=%" PetscInt_FMT ", n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", np=%d\n", ((PetscObject)pc)->prefix, 0, M, pc_gamg->data_cell_rows, pc_gamg->data_cell_cols, (PetscInt)(nnz0 / (PetscReal)M + 0.5), size)); 606569f4572SMark Adams 607a2f3521dSMark F. Adams /* Get A_i and R_i */ 60862294041SBarry Smith for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) { 6099ab59c8bSMark Adams pc_gamg->current_level = level; 61063a3b9bcSJacob Faibussowitsch PetscCheck(level < PETSC_MG_MAXLEVELS, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level); 6115b89ad90SMark F. Adams level1 = level + 1; 6124555aa8cSStefano Zampini #if defined(GAMG_STAGES) 613849bee69SMark Adams if (!gamg_stages[level]) { 614849bee69SMark Adams char str[32]; 615a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level)); 616849bee69SMark Adams PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 617849bee69SMark Adams } 6189566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(gamg_stages[level])); 619b4fbaa2aSMark F. Adams #endif 620c8b0795cSMark F. Adams { /* construct prolongator */ 621725b86d8SJed Brown Mat Gmat; 6220cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 6237700e67bSMark Adams Mat Prol11; 624c8b0795cSMark F. Adams 6252d776b49SBarry Smith PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat)); 6269566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 6279566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], Gmat, agg_lists, &Prol11)); 628c8b0795cSMark F. Adams 629a2f3521dSMark F. Adams /* could have failed to create new level */ 630a2f3521dSMark F. Adams if (Prol11) { 631f7df55f0SStefano Zampini const char *prefix; 632f7df55f0SStefano Zampini char addp[32]; 633f7df55f0SStefano Zampini 6349d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6359566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 636a2f3521dSMark F. Adams 637fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 638c8b0795cSMark F. Adams /* smooth */ 6399566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 640c8b0795cSMark F. Adams } 641c8b0795cSMark F. Adams 6420c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6431b18a24aSMark Adams PetscInt bs; 644849bee69SMark Adams PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 6459566063dSJacob Faibussowitsch PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 646ffc955d6SMark F. Adams } 647ffc955d6SMark F. Adams 6489566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 6499566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(Prol11, prefix)); 6509566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level)); 6519566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(Prol11, addp)); 65291f31d3dSStefano Zampini /* Always generate the transpose with CUDA 653f7df55f0SStefano Zampini Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 6549566063dSJacob Faibussowitsch PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 6559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Prol11)); 6564bde40a0SMark Adams Parr[level1] = Prol11; 6574bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 6584bde40a0SMark Adams 6599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat)); 6609566063dSJacob Faibussowitsch PetscCall(PetscCDDestroy(agg_lists)); 661a2f3521dSMark F. Adams } /* construct prolongator scope */ 6627f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 663171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 66463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 6654555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6669566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 667a90e85d9SMark Adams #endif 668c8b0795cSMark F. Adams break; 669c8b0795cSMark F. Adams } 6709566063dSJacob Faibussowitsch PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 6712472a847SBarry Smith PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?"); 672171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 6730e2909e1SMark Adams if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE; 674849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 6759566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 676849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 677a2f3521dSMark F. Adams 6789566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 6799566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 680569f4572SMark Adams nnztot += info.nz_used; 681bae903cbSmarkadams4 PetscCall(PetscInfo(pc, "%s: %d) N=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", %d active pes\n", ((PetscObject)pc)->prefix, (int)level1, M, pc_gamg->data_cell_cols, (PetscInt)(info.nz_used / (PetscReal)M), nactivepe)); 682569f4572SMark Adams 6834555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6849566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 685b4fbaa2aSMark F. Adams #endif 686a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6879ab59c8bSMark Adams if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) { 68863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: HARD stop of coarsening on level %" PetscInt_FMT ". Grid too small: %" PetscInt_FMT " block nodes\n", ((PetscObject)pc)->prefix, level, M / bs)); 689a90e85d9SMark Adams level++; 690a90e85d9SMark Adams break; 691a90e85d9SMark Adams } 692c8b0795cSMark F. Adams } /* levels */ 6939566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 694c8b0795cSMark F. Adams 695ba4ceb06Smarkadams4 PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0)); 6969d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6975b89ad90SMark F. Adams fine_level = level; 6989566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL)); 6995b89ad90SMark F. Adams 70062294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 7010ed2132dSStefano Zampini 702d5280255SMark F. Adams /* set default smoothers & set operators */ 70362294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) { 704ffc955d6SMark F. Adams KSP smoother; 705ffc955d6SMark F. Adams PC subpc; 706a2f3521dSMark F. Adams 7079566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7089566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 709ffc955d6SMark F. Adams 7109566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 711a2f3521dSMark F. Adams /* set ops */ 7129566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 7139566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1])); 714a2f3521dSMark F. Adams 715a2f3521dSMark F. Adams /* set defaults */ 7169566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 717a2f3521dSMark F. Adams 7180c3bc534SBarry Smith /* set blocks for ASM smoother that uses the 'aggregates' */ 7190c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 7202d3561bbSSatish Balay PetscInt sz; 7217a28f3e5SMark Adams IS *iss; 722a2f3521dSMark F. Adams 7232d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7247a28f3e5SMark Adams iss = ASMLocalIDsArr[level]; 7259566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCASM)); 7269566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(subpc, 0)); 7279566063dSJacob Faibussowitsch PetscCall(PCASMSetType(subpc, PC_ASM_BASIC)); 7287f66b68fSMark Adams if (!sz) { 729ffc955d6SMark F. Adams IS is; 7309566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 7319566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 7329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 733806fa848SBarry Smith } else { 734a94c3b12SMark F. Adams PetscInt kk; 735*2eab5cd7Smarkadams4 PetscCall(PCASMSetLocalSubdomains(subpc, sz, iss, NULL)); 73648a46eb9SPierre Jolivet for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk])); 7379566063dSJacob Faibussowitsch PetscCall(PetscFree(iss)); 738ffc955d6SMark F. Adams } 7390298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 740ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 741806fa848SBarry Smith } else { 7429566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCJACOBI)); 743ffc955d6SMark F. Adams } 744d5280255SMark F. Adams } 745d5280255SMark F. Adams { 746d5280255SMark F. Adams /* coarse grid */ 7479371c9d4SSatish Balay KSP smoother, *k2; 7489371c9d4SSatish Balay PC subpc, pc2; 7499371c9d4SSatish Balay PetscInt ii, first; 7509371c9d4SSatish Balay Mat Lmat = Aarr[(level = pc_gamg->Nlevels - 1)]; 7519371c9d4SSatish Balay lidx = 0; 7520ed2132dSStefano Zampini 7539566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7549566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 755cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 7569566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 7579566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 7589566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCBJACOBI)); 7599566063dSJacob Faibussowitsch PetscCall(PCSetUp(subpc)); 7609566063dSJacob Faibussowitsch PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2)); 76163a3b9bcSJacob Faibussowitsch PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii); 7629566063dSJacob Faibussowitsch PetscCall(KSPGetPC(k2[0], &pc2)); 7639566063dSJacob Faibussowitsch PetscCall(PCSetType(pc2, PCLU)); 7649566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS)); 7659566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 7669566063dSJacob Faibussowitsch PetscCall(KSPSetType(k2[0], KSPPREONLY)); 767d5280255SMark F. Adams } 768cf8ae1d3SMark Adams } 769d5280255SMark F. Adams 770d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 771d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)pc); 772dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); 773d0609cedSBarry Smith PetscOptionsEnd(); 7749566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 775d5280255SMark F. Adams 776f1580f4eSBarry Smith /* set cheby eigen estimates from SA to use in the solver */ 7777e6512fdSJed Brown if (pc_gamg->use_sa_esteig) { 77818c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) { 77918c3aa7eSMark KSP smoother; 78018c3aa7eSMark PetscBool ischeb; 7810ed2132dSStefano Zampini 7829566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb)); 78418c3aa7eSMark if (ischeb) { 78518c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data; 7860ed2132dSStefano Zampini 7872de708cbSJed Brown // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 788efe053fcSJed Brown if (mg->max_eigen_DinvA[level] > 0) { 7897e6512fdSJed Brown // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 7907e6512fdSJed Brown // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 79118c3aa7eSMark PetscReal emax, emin; 7920ed2132dSStefano Zampini 79318c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 79418c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 79563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %" PetscInt_FMT " (N=%" PetscInt_FMT ") with emax = %g emin = %g\n", ((PetscObject)pc)->prefix, level, Aarr[level]->rmap->N, (double)emax, (double)emin)); 7960a94a983SJed Brown cheb->emin_provided = emin; 7970a94a983SJed Brown cheb->emax_provided = emax; 79818c3aa7eSMark } 79918c3aa7eSMark } 80018c3aa7eSMark } 80118c3aa7eSMark } 8020ed2132dSStefano Zampini 8039566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8040ed2132dSStefano Zampini 805d5280255SMark F. Adams /* clean up */ 806d5280255SMark F. Adams for (level = 1; level < pc_gamg->Nlevels; level++) { 8079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Parr[level])); 8089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aarr[level])); 8095b89ad90SMark F. Adams } 810806fa848SBarry Smith } else { 8115f8cf99dSMark F. Adams KSP smoother; 8120ed2132dSStefano Zampini 81363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix)); 8149566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 8159566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 8169566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPPREONLY)); 8179566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8185f8cf99dSMark F. Adams } 819849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 8205b89ad90SMark F. Adams PetscFunctionReturn(0); 8215b89ad90SMark F. Adams } 8225b89ad90SMark F. Adams 8235b89ad90SMark F. Adams /* 8245b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 8255b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 8265b89ad90SMark F. Adams 8275b89ad90SMark F. Adams Input Parameter: 8285b89ad90SMark F. Adams . pc - the preconditioner context 8295b89ad90SMark F. Adams 8305b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 8315b89ad90SMark F. Adams */ 832d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_GAMG(PC pc) 833d71ae5a4SJacob Faibussowitsch { 8345b89ad90SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 8355b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8365b89ad90SMark F. Adams 8375b89ad90SMark F. Adams PetscFunctionBegin; 8389566063dSJacob Faibussowitsch PetscCall(PCReset_GAMG(pc)); 8391baa6e33SBarry Smith if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 8409566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->ops)); 8419566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 8429566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg)); 8432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL)); 8442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL)); 8452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL)); 8462e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL)); 8472e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL)); 8482e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL)); 8492e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL)); 8502e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", NULL)); 8512e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL)); 8522e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL)); 8532e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL)); 8542e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL)); 8552e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL)); 8562e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL)); 8572e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL)); 8582e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL)); 8599566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 8605b89ad90SMark F. Adams PetscFunctionReturn(0); 8615b89ad90SMark F. Adams } 8625b89ad90SMark F. Adams 863676e1743SMark F. Adams /*@ 864f1580f4eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG` 865676e1743SMark F. Adams 866c3339decSBarry Smith Logically Collective 867676e1743SMark F. Adams 868676e1743SMark F. Adams Input Parameters: 8691cc46a46SBarry Smith + pc - the preconditioner context 8701cc46a46SBarry Smith - n - the number of equations 871676e1743SMark F. Adams 872676e1743SMark F. Adams Options Database Key: 873147403d9SBarry Smith . -pc_gamg_process_eq_limit <limit> - set the limit 874676e1743SMark F. Adams 875f1580f4eSBarry Smith Note: 876f1580f4eSBarry Smith `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 877cab9ed1eSBarry Smith that has degrees of freedom 878cab9ed1eSBarry Smith 879676e1743SMark F. Adams Level: intermediate 880676e1743SMark F. Adams 881f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 882676e1743SMark F. Adams @*/ 883d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 884d71ae5a4SJacob Faibussowitsch { 885676e1743SMark F. Adams PetscFunctionBegin; 886676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 887cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n)); 888676e1743SMark F. Adams PetscFunctionReturn(0); 889676e1743SMark F. Adams } 890676e1743SMark F. Adams 891d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 892d71ae5a4SJacob Faibussowitsch { 893c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 894c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 895676e1743SMark F. Adams 896676e1743SMark F. Adams PetscFunctionBegin; 8979d5b6da9SMark F. Adams if (n > 0) pc_gamg->min_eq_proc = n; 898676e1743SMark F. Adams PetscFunctionReturn(0); 899676e1743SMark F. Adams } 900676e1743SMark F. Adams 901389730f3SMark F. Adams /*@ 902f1580f4eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG` 903389730f3SMark F. Adams 904c3339decSBarry Smith Collective 905389730f3SMark F. Adams 906389730f3SMark F. Adams Input Parameters: 9071cc46a46SBarry Smith + pc - the preconditioner context 9081cc46a46SBarry Smith - n - maximum number of equations to aim for 909389730f3SMark F. Adams 910389730f3SMark F. Adams Options Database Key: 911147403d9SBarry Smith . -pc_gamg_coarse_eq_limit <limit> - set the limit 912389730f3SMark F. Adams 913f1580f4eSBarry Smith Note: 914f1580f4eSBarry Smith For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 91574329af1SBarry Smith has less than 1000 unknowns. 91674329af1SBarry Smith 917389730f3SMark F. Adams Level: intermediate 918389730f3SMark F. Adams 919f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 920389730f3SMark F. Adams @*/ 921d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 922d71ae5a4SJacob Faibussowitsch { 923389730f3SMark F. Adams PetscFunctionBegin; 924389730f3SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 925cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n)); 926389730f3SMark F. Adams PetscFunctionReturn(0); 927389730f3SMark F. Adams } 928389730f3SMark F. Adams 929d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 930d71ae5a4SJacob Faibussowitsch { 931389730f3SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 932389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 933389730f3SMark F. Adams 934389730f3SMark F. Adams PetscFunctionBegin; 9359d5b6da9SMark F. Adams if (n > 0) pc_gamg->coarse_eq_limit = n; 936389730f3SMark F. Adams PetscFunctionReturn(0); 937389730f3SMark F. Adams } 938389730f3SMark F. Adams 939676e1743SMark F. Adams /*@ 940f1580f4eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use 941676e1743SMark F. Adams 942c3339decSBarry Smith Collective 943676e1743SMark F. Adams 944676e1743SMark F. Adams Input Parameters: 9451cc46a46SBarry Smith + pc - the preconditioner context 946f1580f4eSBarry Smith - n - `PETSC_TRUE` or `PETSC_FALSE` 947676e1743SMark F. Adams 948676e1743SMark F. Adams Options Database Key: 949147403d9SBarry Smith . -pc_gamg_repartition <true,false> - turn on the repartitioning 950676e1743SMark F. Adams 951f1580f4eSBarry Smith Note: 952f1580f4eSBarry Smith This will generally improve the loading balancing of the work on each level 953cab9ed1eSBarry Smith 954676e1743SMark F. Adams Level: intermediate 955676e1743SMark F. Adams 956f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 957676e1743SMark F. Adams @*/ 958d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 959d71ae5a4SJacob Faibussowitsch { 960676e1743SMark F. Adams PetscFunctionBegin; 961676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 962cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n)); 963676e1743SMark F. Adams PetscFunctionReturn(0); 964676e1743SMark F. Adams } 965676e1743SMark F. Adams 966d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 967d71ae5a4SJacob Faibussowitsch { 968c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 969c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 970676e1743SMark F. Adams 971676e1743SMark F. Adams PetscFunctionBegin; 9729d5b6da9SMark F. Adams pc_gamg->repart = n; 973676e1743SMark F. Adams PetscFunctionReturn(0); 974676e1743SMark F. Adams } 975676e1743SMark F. Adams 976dfd5c07aSMark F. Adams /*@ 977f1580f4eSBarry Smith PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process 97818c3aa7eSMark 979c3339decSBarry Smith Collective 98018c3aa7eSMark 98118c3aa7eSMark Input Parameters: 98218c3aa7eSMark + pc - the preconditioner context 98318c3aa7eSMark - n - number of its 98418c3aa7eSMark 98518c3aa7eSMark Options Database Key: 986147403d9SBarry Smith . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 98718c3aa7eSMark 98818c3aa7eSMark Notes: 9897e6512fdSJed Brown Smoothed aggregation constructs the smoothed prolongator $P = (I - \omega D^{-1} A) T$ where $T$ is the tentative prolongator and $D$ is the diagonal of $A$. 990f1580f4eSBarry Smith Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 991f1580f4eSBarry Smith If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused during the solution process 992f1580f4eSBarry Smith This option is only used when the smoother uses Jacobi, and should be turned off if a different `PCJacobiType` is used. 993efe053fcSJed Brown It became default in PETSc 3.17. 99418c3aa7eSMark 9957e6512fdSJed Brown Level: advanced 99618c3aa7eSMark 997f1580f4eSBarry Smith .seealso: `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 99818c3aa7eSMark @*/ 999d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 1000d71ae5a4SJacob Faibussowitsch { 100118c3aa7eSMark PetscFunctionBegin; 100218c3aa7eSMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1003cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, n)); 100418c3aa7eSMark PetscFunctionReturn(0); 100518c3aa7eSMark } 100618c3aa7eSMark 1007d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 1008d71ae5a4SJacob Faibussowitsch { 100918c3aa7eSMark PC_MG *mg = (PC_MG *)pc->data; 101018c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 101118c3aa7eSMark 101218c3aa7eSMark PetscFunctionBegin; 10137e6512fdSJed Brown pc_gamg->use_sa_esteig = n; 101418c3aa7eSMark PetscFunctionReturn(0); 101518c3aa7eSMark } 101618c3aa7eSMark 101718c3aa7eSMark /*@ 1018f1580f4eSBarry Smith PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY? 101918c3aa7eSMark 1020c3339decSBarry Smith Collective 102118c3aa7eSMark 102218c3aa7eSMark Input Parameters: 102318c3aa7eSMark + pc - the preconditioner context 102418c3aa7eSMark - emax - max eigenvalue 102518c3aa7eSMark . emin - min eigenvalue 102618c3aa7eSMark 102718c3aa7eSMark Options Database Key: 1028147403d9SBarry Smith . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 102918c3aa7eSMark 103018c3aa7eSMark Level: intermediate 103118c3aa7eSMark 1032f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetUseSAEstEig()` 103318c3aa7eSMark @*/ 1034d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin) 1035d71ae5a4SJacob Faibussowitsch { 103618c3aa7eSMark PetscFunctionBegin; 103718c3aa7eSMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1038cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin)); 103918c3aa7eSMark PetscFunctionReturn(0); 104018c3aa7eSMark } 104141ffd417SStefano Zampini 1042d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin) 1043d71ae5a4SJacob Faibussowitsch { 104418c3aa7eSMark PC_MG *mg = (PC_MG *)pc->data; 104518c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 104618c3aa7eSMark 104718c3aa7eSMark PetscFunctionBegin; 10485f80ce2aSJacob Faibussowitsch PetscCheck(emax > emin, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Maximum eigenvalue must be larger than minimum: max %g min %g", (double)emax, (double)emin); 10495f80ce2aSJacob Faibussowitsch PetscCheck(emax * emin > 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Both eigenvalues must be of the same sign: max %g min %g", (double)emax, (double)emin); 105018c3aa7eSMark pc_gamg->emax = emax; 105118c3aa7eSMark pc_gamg->emin = emin; 105218c3aa7eSMark PetscFunctionReturn(0); 105318c3aa7eSMark } 105418c3aa7eSMark 105518c3aa7eSMark /*@ 1056f1580f4eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner 1057dfd5c07aSMark F. Adams 1058c3339decSBarry Smith Collective 1059dfd5c07aSMark F. Adams 1060dfd5c07aSMark F. Adams Input Parameters: 10611cc46a46SBarry Smith + pc - the preconditioner context 1062f1580f4eSBarry Smith - n - `PETSC_TRUE` or `PETSC_FALSE` 1063dfd5c07aSMark F. Adams 1064dfd5c07aSMark F. Adams Options Database Key: 1065147403d9SBarry Smith . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1066dfd5c07aSMark F. Adams 1067dfd5c07aSMark F. Adams Level: intermediate 1068dfd5c07aSMark F. Adams 1069f1580f4eSBarry Smith Note: 1070147403d9SBarry Smith May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1071cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1072cab9ed1eSBarry Smith 1073f1580f4eSBarry Smith .seealso: `PCGAMG` 1074dfd5c07aSMark F. Adams @*/ 1075d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1076d71ae5a4SJacob Faibussowitsch { 1077dfd5c07aSMark F. Adams PetscFunctionBegin; 1078dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1079cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n)); 1080dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1081dfd5c07aSMark F. Adams } 1082dfd5c07aSMark F. Adams 1083d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1084d71ae5a4SJacob Faibussowitsch { 1085dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1086dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1087dfd5c07aSMark F. Adams 1088dfd5c07aSMark F. Adams PetscFunctionBegin; 1089dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1090dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1091dfd5c07aSMark F. Adams } 1092dfd5c07aSMark F. Adams 1093ffc955d6SMark F. Adams /*@ 1094f1580f4eSBarry Smith PCGAMGASMSetUseAggs - Have the `PCGAMG` smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner 1095f1580f4eSBarry Smith used as the smoother 1096ffc955d6SMark F. Adams 1097c3339decSBarry Smith Collective 1098ffc955d6SMark F. Adams 1099ffc955d6SMark F. Adams Input Parameters: 1100cab9ed1eSBarry Smith + pc - the preconditioner context 1101f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not 1102ffc955d6SMark F. Adams 1103ffc955d6SMark F. Adams Options Database Key: 1104147403d9SBarry Smith . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1105ffc955d6SMark F. Adams 1106ffc955d6SMark F. Adams Level: intermediate 1107ffc955d6SMark F. Adams 1108f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCASM`, `PCSetType` 1109ffc955d6SMark F. Adams @*/ 1110d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1111d71ae5a4SJacob Faibussowitsch { 1112ffc955d6SMark F. Adams PetscFunctionBegin; 1113ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1114cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg)); 1115ffc955d6SMark F. Adams PetscFunctionReturn(0); 1116ffc955d6SMark F. Adams } 1117ffc955d6SMark F. Adams 1118d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1119d71ae5a4SJacob Faibussowitsch { 1120ffc955d6SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1121ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1122ffc955d6SMark F. Adams 1123ffc955d6SMark F. Adams PetscFunctionBegin; 1124cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 1125ffc955d6SMark F. Adams PetscFunctionReturn(0); 1126ffc955d6SMark F. Adams } 1127ffc955d6SMark F. Adams 1128171cca9aSMark Adams /*@ 1129cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1130171cca9aSMark Adams 1131c3339decSBarry Smith Collective 1132171cca9aSMark Adams 1133171cca9aSMark Adams Input Parameters: 1134171cca9aSMark Adams + pc - the preconditioner context 1135f1580f4eSBarry Smith - flg - `PETSC_TRUE` to not force coarse grid onto one processor 1136171cca9aSMark Adams 1137171cca9aSMark Adams Options Database Key: 1138147403d9SBarry Smith . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1139171cca9aSMark Adams 1140171cca9aSMark Adams Level: intermediate 1141171cca9aSMark Adams 1142f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()` 1143171cca9aSMark Adams @*/ 1144d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1145d71ae5a4SJacob Faibussowitsch { 1146171cca9aSMark Adams PetscFunctionBegin; 1147171cca9aSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1148cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetUseParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg)); 1149171cca9aSMark Adams PetscFunctionReturn(0); 1150171cca9aSMark Adams } 1151171cca9aSMark Adams 1152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1153d71ae5a4SJacob Faibussowitsch { 1154171cca9aSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1155171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1156171cca9aSMark Adams 1157171cca9aSMark Adams PetscFunctionBegin; 1158171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 1159ffc955d6SMark F. Adams PetscFunctionReturn(0); 1160ffc955d6SMark F. Adams } 1161ffc955d6SMark F. Adams 11624ef23d27SMark F. Adams /*@ 1163f1580f4eSBarry Smith PCGAMGSetCpuPinCoarseGrids - pin the coarse grids created in `PCGAMG` to run only on the CPU since the problems may be too small to run efficiently on the GPUs 1164ce7c7f2fSMark Adams 1165c3339decSBarry Smith Collective 1166ce7c7f2fSMark Adams 1167ce7c7f2fSMark Adams Input Parameters: 1168ce7c7f2fSMark Adams + pc - the preconditioner context 1169f1580f4eSBarry Smith - flg - `PETSC_TRUE` to pin coarse grids to the CPU 1170ce7c7f2fSMark Adams 1171ce7c7f2fSMark Adams Options Database Key: 1172147403d9SBarry Smith . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1173ce7c7f2fSMark Adams 1174ce7c7f2fSMark Adams Level: intermediate 1175ce7c7f2fSMark Adams 1176f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()` 1177ce7c7f2fSMark Adams @*/ 1178d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1179d71ae5a4SJacob Faibussowitsch { 1180ce7c7f2fSMark Adams PetscFunctionBegin; 1181ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1182cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg)); 1183ce7c7f2fSMark Adams PetscFunctionReturn(0); 1184ce7c7f2fSMark Adams } 1185ce7c7f2fSMark Adams 1186d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1187d71ae5a4SJacob Faibussowitsch { 1188ce7c7f2fSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1189ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1190ce7c7f2fSMark Adams 1191ce7c7f2fSMark Adams PetscFunctionBegin; 1192ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 1193ce7c7f2fSMark Adams PetscFunctionReturn(0); 1194ce7c7f2fSMark Adams } 1195ce7c7f2fSMark Adams 1196ce7c7f2fSMark Adams /*@ 1197147403d9SBarry Smith PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1198ce7c7f2fSMark Adams 1199c3339decSBarry Smith Collective 1200ce7c7f2fSMark Adams 1201ce7c7f2fSMark Adams Input Parameters: 1202ce7c7f2fSMark Adams + pc - the preconditioner context 1203f1580f4eSBarry Smith - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD` 1204ce7c7f2fSMark Adams 1205ce7c7f2fSMark Adams Options Database Key: 1206147403d9SBarry Smith . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1207ce7c7f2fSMark Adams 1208ce7c7f2fSMark Adams Level: intermediate 1209ce7c7f2fSMark Adams 1210f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD` 1211ce7c7f2fSMark Adams @*/ 1212d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1213d71ae5a4SJacob Faibussowitsch { 1214ce7c7f2fSMark Adams PetscFunctionBegin; 1215ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1216cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg)); 1217ce7c7f2fSMark Adams PetscFunctionReturn(0); 1218ce7c7f2fSMark Adams } 1219ce7c7f2fSMark Adams 1220d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1221d71ae5a4SJacob Faibussowitsch { 1222ce7c7f2fSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1223ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1224ce7c7f2fSMark Adams 1225ce7c7f2fSMark Adams PetscFunctionBegin; 1226ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 1227ce7c7f2fSMark Adams PetscFunctionReturn(0); 1228ce7c7f2fSMark Adams } 1229ce7c7f2fSMark Adams 1230ce7c7f2fSMark Adams /*@ 1231f1580f4eSBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels `PCGAMG` will use 12324ef23d27SMark F. Adams 1233c3339decSBarry Smith Not collective 12344ef23d27SMark F. Adams 12354ef23d27SMark F. Adams Input Parameters: 12361cc46a46SBarry Smith + pc - the preconditioner 12371cc46a46SBarry Smith - n - the maximum number of levels to use 12384ef23d27SMark F. Adams 12394ef23d27SMark F. Adams Options Database Key: 1240147403d9SBarry Smith . -pc_mg_levels <n> - set the maximum number of levels to allow 12414ef23d27SMark F. Adams 12424ef23d27SMark F. Adams Level: intermediate 12434ef23d27SMark F. Adams 1244f1580f4eSBarry Smith Developer Note: 1245f1580f4eSBarry Smith Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG` 1246f1580f4eSBarry Smith 1247f1580f4eSBarry Smith .seealso: `PCGAMG` 12484ef23d27SMark F. Adams @*/ 1249d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1250d71ae5a4SJacob Faibussowitsch { 12514ef23d27SMark F. Adams PetscFunctionBegin; 12524ef23d27SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1253cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n)); 12544ef23d27SMark F. Adams PetscFunctionReturn(0); 12554ef23d27SMark F. Adams } 12564ef23d27SMark F. Adams 1257d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1258d71ae5a4SJacob Faibussowitsch { 12594ef23d27SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 12604ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 12614ef23d27SMark F. Adams 12624ef23d27SMark F. Adams PetscFunctionBegin; 12639d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12644ef23d27SMark F. Adams PetscFunctionReturn(0); 12654ef23d27SMark F. Adams } 12664ef23d27SMark F. Adams 12673542efc5SMark F. Adams /*@ 12683542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12693542efc5SMark F. Adams 1270c3339decSBarry Smith Not collective 12713542efc5SMark F. Adams 12723542efc5SMark F. Adams Input Parameters: 12731cc46a46SBarry Smith + pc - the preconditioner context 1274c9567895SMark . v - array of threshold values for finest n levels; 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph 1275055c8bd0SJed Brown - n - number of threshold values provided in array 12763542efc5SMark F. Adams 12773542efc5SMark F. Adams Options Database Key: 1278147403d9SBarry Smith . -pc_gamg_threshold <threshold> - the threshold to drop edges 12793542efc5SMark F. Adams 128095452b02SPatrick Sanan Notes: 1281af3c827dSMark Adams Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably. 1282f1580f4eSBarry Smith Before coarsening or aggregating the graph, `PCGAMG` removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points. 1283cab9ed1eSBarry Smith 1284f1580f4eSBarry Smith If n is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening. 1285f1580f4eSBarry Smith In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`. 1286055c8bd0SJed Brown If n is greater than the total number of levels, the excess entries in threshold will not be used. 1287055c8bd0SJed Brown 12883542efc5SMark F. Adams Level: intermediate 12893542efc5SMark F. Adams 1290*2eab5cd7Smarkadams4 .seealso: `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetThresholdScale()` 12913542efc5SMark F. Adams @*/ 1292d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 1293d71ae5a4SJacob Faibussowitsch { 12943542efc5SMark F. Adams PetscFunctionBegin; 12953542efc5SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1296055c8bd0SJed Brown if (n) PetscValidRealPointer(v, 2); 1297cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n)); 12983542efc5SMark F. Adams PetscFunctionReturn(0); 12993542efc5SMark F. Adams } 13003542efc5SMark F. Adams 1301d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 1302d71ae5a4SJacob Faibussowitsch { 1303c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1304c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1305c1eae691SMark Adams PetscInt i; 1306c1eae691SMark Adams PetscFunctionBegin; 1307055c8bd0SJed Brown for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1308055c8bd0SJed Brown for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1309c1eae691SMark Adams PetscFunctionReturn(0); 1310c1eae691SMark Adams } 1311c1eae691SMark Adams 1312c1eae691SMark Adams /*@ 1313f1580f4eSBarry Smith PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids 1314c9567895SMark 1315c3339decSBarry Smith Collective 1316c9567895SMark 1317c9567895SMark Input Parameters: 1318c9567895SMark + pc - the preconditioner context 1319f1580f4eSBarry Smith . v - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA 1320c9567895SMark - n - number of values provided in array 1321c9567895SMark 1322c9567895SMark Options Database Key: 1323147403d9SBarry Smith . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1324c9567895SMark 1325c9567895SMark Level: intermediate 1326c9567895SMark 1327f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()` 1328c9567895SMark @*/ 1329d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1330d71ae5a4SJacob Faibussowitsch { 1331c9567895SMark PetscFunctionBegin; 1332c9567895SMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1333c9567895SMark if (n) PetscValidIntPointer(v, 2); 1334cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n)); 1335c9567895SMark PetscFunctionReturn(0); 1336c9567895SMark } 1337c9567895SMark 1338d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1339d71ae5a4SJacob Faibussowitsch { 1340c9567895SMark PC_MG *mg = (PC_MG *)pc->data; 1341c9567895SMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1342c9567895SMark PetscInt i; 1343c9567895SMark PetscFunctionBegin; 1344c9567895SMark for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1345c9567895SMark for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1346c9567895SMark PetscFunctionReturn(0); 1347c9567895SMark } 1348c9567895SMark 1349c9567895SMark /*@ 1350c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1351c1eae691SMark Adams 1352c3339decSBarry Smith Not collective 1353c1eae691SMark Adams 1354c1eae691SMark Adams Input Parameters: 1355c1eae691SMark Adams + pc - the preconditioner context 1356147403d9SBarry Smith - scale - the threshold value reduction, usually < 1.0 1357c1eae691SMark Adams 1358c1eae691SMark Adams Options Database Key: 1359147403d9SBarry Smith . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1360c1eae691SMark Adams 1361f1580f4eSBarry Smith Note: 1362f1580f4eSBarry Smith The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`. 1363f1580f4eSBarry Smith This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`. 1364055c8bd0SJed Brown 1365c1eae691SMark Adams Level: advanced 1366c1eae691SMark Adams 1367f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetThreshold()` 1368c1eae691SMark Adams @*/ 1369d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1370d71ae5a4SJacob Faibussowitsch { 13713542efc5SMark F. Adams PetscFunctionBegin; 1372c1eae691SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1373cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v)); 1374c1eae691SMark Adams PetscFunctionReturn(0); 1375c1eae691SMark Adams } 1376c1eae691SMark Adams 1377d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1378d71ae5a4SJacob Faibussowitsch { 1379c1eae691SMark Adams PC_MG *mg = (PC_MG *)pc->data; 1380c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1381c1eae691SMark Adams PetscFunctionBegin; 1382c1eae691SMark Adams pc_gamg->threshold_scale = v; 13833542efc5SMark F. Adams PetscFunctionReturn(0); 13843542efc5SMark F. Adams } 13853542efc5SMark F. Adams 1386e20c40e8SBarry Smith /*@C 1387f1580f4eSBarry Smith PCGAMGSetType - Set the type of algorithm `PCGAMG` should use 1388676e1743SMark F. Adams 1389c3339decSBarry Smith Collective 1390676e1743SMark F. Adams 1391676e1743SMark F. Adams Input Parameters: 1392c60c7ad4SBarry Smith + pc - the preconditioner context 1393f1580f4eSBarry Smith - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL` 1394676e1743SMark F. Adams 1395676e1743SMark F. Adams Options Database Key: 1396cab9ed1eSBarry Smith . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1397676e1743SMark F. Adams 1398676e1743SMark F. Adams Level: intermediate 1399676e1743SMark F. Adams 1400db781477SPatrick Sanan .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1401676e1743SMark F. Adams @*/ 1402d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1403d71ae5a4SJacob Faibussowitsch { 1404676e1743SMark F. Adams PetscFunctionBegin; 1405676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1406cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type)); 1407676e1743SMark F. Adams PetscFunctionReturn(0); 1408676e1743SMark F. Adams } 1409676e1743SMark F. Adams 1410e20c40e8SBarry Smith /*@C 1411f1580f4eSBarry Smith PCGAMGGetType - Get the type of algorithm `PCGAMG` will use 1412c60c7ad4SBarry Smith 1413c3339decSBarry Smith Collective 1414c60c7ad4SBarry Smith 1415c60c7ad4SBarry Smith Input Parameter: 1416c60c7ad4SBarry Smith . pc - the preconditioner context 1417c60c7ad4SBarry Smith 1418c60c7ad4SBarry Smith Output Parameter: 1419c60c7ad4SBarry Smith . type - the type of algorithm used 1420c60c7ad4SBarry Smith 1421c60c7ad4SBarry Smith Level: intermediate 1422c60c7ad4SBarry Smith 1423f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType` 1424c60c7ad4SBarry Smith @*/ 1425d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1426d71ae5a4SJacob Faibussowitsch { 1427c60c7ad4SBarry Smith PetscFunctionBegin; 1428c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1429cac4c232SBarry Smith PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type)); 1430c60c7ad4SBarry Smith PetscFunctionReturn(0); 1431c60c7ad4SBarry Smith } 1432c60c7ad4SBarry Smith 1433d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1434d71ae5a4SJacob Faibussowitsch { 1435c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1436c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1437c60c7ad4SBarry Smith 1438c60c7ad4SBarry Smith PetscFunctionBegin; 1439c60c7ad4SBarry Smith *type = pc_gamg->type; 1440c60c7ad4SBarry Smith PetscFunctionReturn(0); 1441c60c7ad4SBarry Smith } 1442c60c7ad4SBarry Smith 1443d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1444d71ae5a4SJacob Faibussowitsch { 14451ab5ffc9SJed Brown PC_MG *mg = (PC_MG *)pc->data; 14461ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 14475f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(PC); 1448676e1743SMark F. Adams 1449676e1743SMark F. Adams PetscFunctionBegin; 1450c60c7ad4SBarry Smith pc_gamg->type = type; 14519566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(GAMGList, type, &r)); 14525f80ce2aSJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type); 14531ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 14549566063dSJacob Faibussowitsch PetscCall((*pc_gamg->ops->destroy)(pc)); 14559566063dSJacob Faibussowitsch PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps))); 1456e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 14573ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 14583ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 14593ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 14603ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 14613ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 14629566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 14633ae0bb68SMark Adams pc_gamg->data_sz = 0; 14641ab5ffc9SJed Brown } 14659566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 14669566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name)); 14679566063dSJacob Faibussowitsch PetscCall((*r)(pc)); 1468676e1743SMark F. Adams PetscFunctionReturn(0); 1469676e1743SMark F. Adams } 1470676e1743SMark F. Adams 1471d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer) 1472d71ae5a4SJacob Faibussowitsch { 14735adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14745adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1475e7d4b4cbSMark Adams PetscReal gc = 0, oc = 0; 147690db8557SMark Adams 14775adeb434SBarry Smith PetscFunctionBegin; 14789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GAMG specific options\n")); 14799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold for dropping small values in graph on each level =")); 14809566063dSJacob Faibussowitsch for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i])); 14819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 14829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale)); 148348a46eb9SPierre Jolivet if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, " Using aggregates from coarsening process to define subdomains for PCASM\n")); 148448a46eb9SPierre Jolivet if (pc_gamg->use_parallel_coarse_grid_solver) PetscCall(PetscViewerASCIIPrintf(viewer, " Using parallel coarse grid solver (all coarse grid equations not put on one process)\n")); 14851baa6e33SBarry Smith if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer)); 14869566063dSJacob Faibussowitsch PetscCall(PCMGGetGridComplexity(pc, &gc, &oc)); 148763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Complexity: grid = %g operator = %g\n", (double)gc, (double)oc)); 14885adeb434SBarry Smith PetscFunctionReturn(0); 14895adeb434SBarry Smith } 14905adeb434SBarry Smith 1491d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject) 1492d71ae5a4SJacob Faibussowitsch { 1493676e1743SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1494676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 14957e6512fdSJed Brown PetscBool flag; 14963b4367a7SBarry Smith MPI_Comm comm; 149718c3aa7eSMark char prefix[256], tname[32]; 1498c1eae691SMark Adams PetscInt i, n; 149914a9496bSBarry Smith const char *pcpre; 15000a545947SLisandro Dalcin static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL}; 15015b89ad90SMark F. Adams PetscFunctionBegin; 15029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 1503d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options"); 15049566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 15051baa6e33SBarry Smith if (flag) PetscCall(PCGAMGSetType(pc, tname)); 15069566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL)); 1507f1580f4eSBarry Smith PetscCall(PetscOptionsBool("-pc_gamg_use_sa_esteig", "Use eigen estimate from smoothed aggregation for smoother", "PCGAMGSetUseSAEstEig", pc_gamg->use_sa_esteig, &pc_gamg->use_sa_esteig, NULL)); 15089566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL)); 15099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_asm_use_agg", "Use aggregation aggregates for ASM smoother", "PCGAMGASMSetUseAggs", pc_gamg->use_aggs_in_asm, &pc_gamg->use_aggs_in_asm, NULL)); 15109566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver", "Use parallel coarse grid solver (otherwise put last grid on one process)", "PCGAMGSetUseParallelCoarseGridSolve", pc_gamg->use_parallel_coarse_grid_solver, &pc_gamg->use_parallel_coarse_grid_solver, NULL)); 15119566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids", "Pin coarse grids to the CPU", "PCGAMGSetCpuPinCoarseGrids", pc_gamg->cpu_pin_coarse_grids, &pc_gamg->cpu_pin_coarse_grids, NULL)); 15129371c9d4SSatish Balay PetscCall(PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type", "compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth", "PCGAMGSetCoarseGridLayoutType", LayoutTypes, 15139371c9d4SSatish Balay (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL)); 15149566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_process_eq_limit", "Limit (goal) on number of equations per process on coarse grids", "PCGAMGSetProcEqLim", pc_gamg->min_eq_proc, &pc_gamg->min_eq_proc, NULL)); 15159566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_coarse_eq_limit", "Limit on number of equations for the coarse grid", "PCGAMGSetCoarseEqLim", pc_gamg->coarse_eq_limit, &pc_gamg->coarse_eq_limit, NULL)); 15169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale", "Scaling of threshold for each level not specified", "PCGAMGSetThresholdScale", pc_gamg->threshold_scale, &pc_gamg->threshold_scale, NULL)); 151718c3aa7eSMark n = PETSC_MG_MAXLEVELS; 15189566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag)); 151918c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1520efd3c5ceSMark Adams if (!flag) n = 1; 1521c1eae691SMark Adams i = n; 1522d71ae5a4SJacob Faibussowitsch do { 1523d71ae5a4SJacob Faibussowitsch pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1524d71ae5a4SJacob Faibussowitsch } while (++i < PETSC_MG_MAXLEVELS); 1525c1eae691SMark Adams } 1526c9567895SMark n = PETSC_MG_MAXLEVELS; 15279566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-pc_gamg_rank_reduction_factors", "Manual schedule of coarse grid reduction factors that overrides internal heuristics (0 for first reduction puts one process/device)", "PCGAMGSetRankReductionFactors", pc_gamg->level_reduction_factors, &n, &flag)); 1528c9567895SMark if (!flag) i = 0; 1529c9567895SMark else i = n; 1530d71ae5a4SJacob Faibussowitsch do { 1531d71ae5a4SJacob Faibussowitsch pc_gamg->level_reduction_factors[i] = -1; 1532d71ae5a4SJacob Faibussowitsch } while (++i < PETSC_MG_MAXLEVELS); 15339566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL)); 153418c3aa7eSMark { 153518c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 153618c3aa7eSMark n = 2; 15379566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag)); 153818c3aa7eSMark if (flag) { 153908401ef6SPierre Jolivet PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 15409566063dSJacob Faibussowitsch PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 154118c3aa7eSMark } 154218c3aa7eSMark } 1543b7cbab4eSMark Adams /* set options for subtype */ 1544dbbe0bcdSBarry Smith PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject)); 154518c3aa7eSMark 15469566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 15479566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 1548d0609cedSBarry Smith PetscOptionsHeadEnd(); 15495b89ad90SMark F. Adams PetscFunctionReturn(0); 15505b89ad90SMark F. Adams } 15515b89ad90SMark F. Adams 15525b89ad90SMark F. Adams /*MC 15531cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 15545b89ad90SMark F. Adams 1555280d9858SJed Brown Options Database Keys: 155621d928e4Smarkadams4 + -pc_gamg_type <type,default=agg> - one of agg, geo, or classical 155721d928e4Smarkadams4 . -pc_gamg_repartition <bool,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 155821d928e4Smarkadams4 . -pc_gamg_asm_use_agg <bool,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother 1559f1580f4eSBarry Smith . -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> 1560cab9ed1eSBarry Smith equations on each process that has degrees of freedom 15612d776b49SBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 156221d928e4Smarkadams4 . -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true) 156321d928e4Smarkadams4 . -pc_gamg_threshold[] <thresh,default=[-1,...]> - Before aggregating the graph `PCGAMG` will remove small values from the graph on each level (< 0 does no filtering) 15642d776b49SBarry Smith - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1565cab9ed1eSBarry Smith 1566f1580f4eSBarry Smith Options Database Keys for Aggregation: 1567cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 15682d776b49SBarry Smith . -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated) 1569bae903cbSmarkadams4 - -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1570cab9ed1eSBarry Smith 1571f1580f4eSBarry Smith Options Database Keys for Multigrid: 1572f1580f4eSBarry Smith + -pc_mg_cycles <v> - v or w, see `PCMGSetCycleType()` 1573db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1574db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 157521d928e4Smarkadams4 - -pc_mg_levels <levels> - Number of levels of multigrid to use. GAMG has a heuristic so pc_mg_levels is not usually used with GAMG 15765b89ad90SMark F. Adams 157795452b02SPatrick Sanan Notes: 1578f1580f4eSBarry Smith To obtain good performance for `PCGAMG` for vector valued problems you must 1579f1580f4eSBarry Smith call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point 1580f1580f4eSBarry Smith call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1581f1580f4eSBarry Smith 1582f1580f4eSBarry Smith See [the Users Manual section on PCGAMG](sec_amg) for more details. 15831cc46a46SBarry Smith 15845b89ad90SMark F. Adams Level: intermediate 1585280d9858SJed Brown 1586db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1587db781477SPatrick Sanan `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()` 15885b89ad90SMark F. Adams M*/ 1589b2573a8aSBarry Smith 1590d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 1591d71ae5a4SJacob Faibussowitsch { 15925b89ad90SMark F. Adams PC_GAMG *pc_gamg; 15935b89ad90SMark F. Adams PC_MG *mg; 15945b89ad90SMark F. Adams 15955b89ad90SMark F. Adams PetscFunctionBegin; 15961c1aac46SBarry Smith /* register AMG type */ 15979566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 15981c1aac46SBarry Smith 15995b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 16009566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); 16019566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 16025b89ad90SMark F. Adams 16035b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 16044dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg)); 16059566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 16065b89ad90SMark F. Adams mg = (PC_MG *)pc->data; 16075b89ad90SMark F. Adams mg->innerctx = pc_gamg; 16085b89ad90SMark F. Adams 16094dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg->ops)); 16101ab5ffc9SJed Brown 16119d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 16129d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 16130a545947SLisandro Dalcin pc_gamg->data = NULL; 16145b89ad90SMark F. Adams 16159d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 16165b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 16175b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 16185b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 16195b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 16205adeb434SBarry Smith mg->view = PCView_GAMG; 16215b89ad90SMark F. Adams 16229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG)); 16239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG)); 16249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG)); 16259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG)); 16269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG)); 16279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG)); 16289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG)); 16299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 16309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG)); 16319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG)); 16329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG)); 16339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG)); 16349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG)); 16359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG)); 16369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG)); 16379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG)); 16389d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 163921d928e4Smarkadams4 pc_gamg->reuse_prol = PETSC_TRUE; 16400c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1641171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1642a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1643a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1644038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 164525a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 164653134ebeSMark Adams for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1; 1647c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 164818c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 16499ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 16507e6512fdSJed Brown pc_gamg->use_sa_esteig = PETSC_TRUE; 165118c3aa7eSMark pc_gamg->emin = 0; 165218c3aa7eSMark pc_gamg->emax = 0; 165318c3aa7eSMark 1654c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 16559d5b6da9SMark F. Adams 1656bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 16579566063dSJacob Faibussowitsch PetscCall(PCGAMGSetType(pc, PCGAMGAGG)); 16585b89ad90SMark F. Adams PetscFunctionReturn(0); 16595b89ad90SMark F. Adams } 16603e3471ccSMark Adams 16613e3471ccSMark Adams /*@C 1662f1580f4eSBarry Smith PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called 1663f1580f4eSBarry Smith from `PCInitializePackage()`. 16643e3471ccSMark Adams 16653e3471ccSMark Adams Level: developer 16663e3471ccSMark Adams 1667db781477SPatrick Sanan .seealso: `PetscInitialize()` 16683e3471ccSMark Adams @*/ 1669d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGInitializePackage(void) 1670d71ae5a4SJacob Faibussowitsch { 16714555aa8cSStefano Zampini PetscInt l; 16723e3471ccSMark Adams 16733e3471ccSMark Adams PetscFunctionBegin; 16743e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 16753e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 16769566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO)); 16779566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG)); 16789566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical)); 16799566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1680c1c463dbSMark Adams 1681c1c463dbSMark Adams /* general events */ 1682849bee69SMark Adams PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1683849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1684849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1685849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1686849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1687849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1688849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1689849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1690849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1691849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1692849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1693849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1694849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1695849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1696849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1697849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 16984555aa8cSStefano Zampini for (l = 0; l < PETSC_MG_MAXLEVELS; l++) { 16994555aa8cSStefano Zampini char ename[32]; 17005b89ad90SMark F. Adams 170163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l)); 17029566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 170363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l)); 17049566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 170563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l)); 17069566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 17074555aa8cSStefano Zampini } 17084555aa8cSStefano Zampini #if defined(GAMG_STAGES) 1709849bee69SMark Adams { /* create timer stages */ 17105b89ad90SMark F. Adams char str[32]; 1711a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0)); 17129566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 17135b89ad90SMark F. Adams } 17145b89ad90SMark F. Adams #endif 17153e3471ccSMark Adams PetscFunctionReturn(0); 17163e3471ccSMark Adams } 17173e3471ccSMark Adams 17183e3471ccSMark Adams /*@C 1719f1580f4eSBarry Smith PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is 1720f1580f4eSBarry Smith called from `PetscFinalize()` automatically. 17213e3471ccSMark Adams 17223e3471ccSMark Adams Level: developer 17233e3471ccSMark Adams 1724db781477SPatrick Sanan .seealso: `PetscFinalize()` 17253e3471ccSMark Adams @*/ 1726d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGFinalizePackage(void) 1727d71ae5a4SJacob Faibussowitsch { 17283e3471ccSMark Adams PetscFunctionBegin; 17293e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 17309566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&GAMGList)); 17313e3471ccSMark Adams PetscFunctionReturn(0); 17323e3471ccSMark Adams } 1733a36cf38bSToby Isaac 1734a36cf38bSToby Isaac /*@C 1735f1580f4eSBarry Smith PCGAMGRegister - Register a `PCGAMG` implementation. 1736a36cf38bSToby Isaac 1737a36cf38bSToby Isaac Input Parameters: 1738f1580f4eSBarry Smith + type - string that will be used as the name of the `PCGAMG` type. 1739a36cf38bSToby Isaac - create - function for creating the gamg context. 1740a36cf38bSToby Isaac 1741f1580f4eSBarry Smith Level: developer 1742a36cf38bSToby Isaac 1743db781477SPatrick Sanan .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 1744a36cf38bSToby Isaac @*/ 1745d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1746d71ae5a4SJacob Faibussowitsch { 1747a36cf38bSToby Isaac PetscFunctionBegin; 17489566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 17499566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, type, create)); 1750a36cf38bSToby Isaac PetscFunctionReturn(0); 1751a36cf38bSToby Isaac } 17522d776b49SBarry Smith 17532d776b49SBarry Smith /*@C 17542d776b49SBarry Smith PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process 17552d776b49SBarry Smith 17562d776b49SBarry Smith Input Parameters: 17572d776b49SBarry Smith + pc - the `PCGAMG` 17582d776b49SBarry Smith - A - the matrix, for any level 17592d776b49SBarry Smith 17602d776b49SBarry Smith Output Parameter: 17612d776b49SBarry Smith . G - the graph 17622d776b49SBarry Smith 17632d776b49SBarry Smith Level: advanced 17642d776b49SBarry Smith 17652d776b49SBarry Smith .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 17662d776b49SBarry Smith @*/ 17672d776b49SBarry Smith PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G) 17682d776b49SBarry Smith { 17692d776b49SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 17702d776b49SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 17712d776b49SBarry Smith 17722d776b49SBarry Smith PetscFunctionBegin; 17732d776b49SBarry Smith PetscCall(pc_gamg->ops->creategraph(pc, A, G)); 17742d776b49SBarry Smith PetscFunctionReturn(0); 17752d776b49SBarry Smith } 1776