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 27d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 28d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 29d3d6bff4SMark F. Adams { 30d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 31d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 32d3d6bff4SMark F. Adams 33d3d6bff4SMark F. Adams PetscFunctionBegin; 349566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 351c1aac46SBarry Smith pc_gamg->data_sz = 0; 369566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->orig_data)); 375f80ce2aSJacob Faibussowitsch for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS ; level++) { 3818c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 3918c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 4018c3aa7eSMark } 4118c3aa7eSMark pc_gamg->emin = 0; 4218c3aa7eSMark pc_gamg->emax = 0; 43a2f3521dSMark F. Adams PetscFunctionReturn(0); 44a2f3521dSMark F. Adams } 45a2f3521dSMark F. Adams 465b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 475b89ad90SMark F. Adams /* 48c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 49a147abb0SMark F. Adams of active processors. 505b89ad90SMark F. Adams 515b89ad90SMark F. Adams Input Parameter: 52a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 53a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 549d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 55c5bfad50SMark F. Adams . cr_bs - coarse block size 563530afc2SMark F. Adams In/Output Parameter: 57a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 58afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 5911e60469SMark F. Adams Output Parameter: 603530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 615b89ad90SMark F. Adams */ 62171cca9aSMark Adams 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) 635b89ad90SMark F. Adams { 649d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 65486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 66a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 673b4367a7SBarry Smith MPI_Comm comm; 68c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 693ae0bb68SMark Adams PetscInt ncrs_eq,ncrs,f_bs; 705b89ad90SMark F. Adams 715b89ad90SMark F. Adams PetscFunctionBegin; 729566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat_fine,&comm)); 739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 759566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat_fine, &f_bs)); 76849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP],0,0,0,0)); 779566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0)); 789566063dSJacob Faibussowitsch PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat)); 799566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0)); 80849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP],0,0,0,0)); 81038e3b61SMark F. Adams 82ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 83ce7c7f2fSMark Adams 843ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 859566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL)); 863ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 873ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 8873911c69SBarry Smith } else { 893ae0bb68SMark Adams PetscInt bs; 909566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Cmat, &bs)); 913ae0bb68SMark Adams ncrs = ncrs_eq/bs; 923ae0bb68SMark Adams } 93c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 94c9567895SMark 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. */ 95c9567895SMark #if defined(PETSC_HAVE_CUDA) 96c9567895SMark PetscShmComm pshmcomm; 97c9567895SMark PetscMPIInt locrank; 98c9567895SMark MPI_Comm loccomm; 99c9567895SMark PetscInt s_nnodes,r_nnodes, new_new_size; 100c9567895SMark cudaError_t cerr; 101c9567895SMark int devCount; 1029566063dSJacob Faibussowitsch PetscCall(PetscShmCommGet(comm,&pshmcomm)); 1039566063dSJacob Faibussowitsch PetscCall(PetscShmCommGetMpiShmComm(pshmcomm,&loccomm)); 1049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(loccomm, &locrank)); 105c9567895SMark s_nnodes = !locrank; 1069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm)); 1077827d75bSBarry Smith PetscCheck((size%r_nnodes) == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes); 108c9567895SMark devCount = 0; 109c9567895SMark cerr = cudaGetDeviceCount(&devCount); 110c9567895SMark cudaGetLastError(); /* Reset the last error */ 111c9567895SMark if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 112c9567895SMark new_new_size = r_nnodes * devCount; 113c9567895SMark new_size = new_new_size; 1149566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Fine grid with Cuda. %D nodes. Change new active set size %d --> %d (devCount=%d #nodes=%D)\n",((PetscObject)pc)->prefix,r_nnodes,nactive,new_size,devCount,r_nnodes)); 115c9567895SMark } else { 1169566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: With Cuda but no device. Use heuristics.")); 117c9567895SMark goto HEURISTIC; 118c9567895SMark } 119c9567895SMark #else 120c9567895SMark SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here"); 121c9567895SMark #endif 122c9567895SMark } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 123*08401ef6SPierre Jolivet 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 %D",nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]); 124c9567895SMark new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level]; 1259566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Manually setting reduction to %d active processes (%d/%D)\n",((PetscObject)pc)->prefix,new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level])); 126c9567895SMark } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 127c9567895SMark new_size = 1; 1289566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size)); 129c9567895SMark } else { 130472110cdSMark F. Adams PetscInt ncrs_eq_glob; 131c9567895SMark #if defined(PETSC_HAVE_CUDA) 132c9567895SMark HEURISTIC: 133c9567895SMark #endif 1349566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL)); 135a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1367f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 137c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 1389566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size)); 139a2f3521dSMark F. Adams } 1402e3501ffSMark Adams if (new_size==nactive) { 141ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 142ce7c7f2fSMark Adams if (new_size < size) { 143ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 1449566063dSJacob 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)); 145ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1469566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE)); 1479566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE)); 148ce7c7f2fSMark Adams } 149ce7c7f2fSMark Adams } 150ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1512e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 152192c0e8bSMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1; 153885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 154849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 15571959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 156*08401ef6SPierre Jolivet PetscCheck(ncrs_eq % cr_bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 157ce7c7f2fSMark Adams /* get new_size and rfactor */ 158ce7c7f2fSMark Adams if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 159ce7c7f2fSMark Adams /* find factor */ 160ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 161ce7c7f2fSMark Adams else { 162ce7c7f2fSMark Adams PetscReal best_fact = 0.; 163ce7c7f2fSMark Adams jj = -1; 164ce7c7f2fSMark Adams for (kk = 1 ; kk <= size ; kk++) { 165ce7c7f2fSMark Adams if (!(size%kk)) { /* a candidate */ 166ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 167ce7c7f2fSMark Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 168ce7c7f2fSMark Adams if (fact > best_fact) { 169ce7c7f2fSMark Adams best_fact = fact; jj = kk; 170ce7c7f2fSMark Adams } 171ce7c7f2fSMark Adams } 172ce7c7f2fSMark Adams } 173ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 174ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 175ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 176ce7c7f2fSMark Adams else expand_factor = rfactor; 177ce7c7f2fSMark Adams } 178ce7c7f2fSMark Adams new_size = size/rfactor; /* make new size one that is factor */ 1794cdfd227SMark if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1804cdfd227SMark *a_Amat_crs = Cmat; 1819566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",((PetscObject)pc)->prefix,new_size,ncrs_eq)); 182849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 183ce7c7f2fSMark Adams PetscFunctionReturn(0); 184ce7c7f2fSMark Adams } 185ce7c7f2fSMark Adams } 186a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 1879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &counts)); 188849bee69SMark Adams if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 1895a9b9e01SMark F. Adams Mat adj; 190849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART],0,0,0,0)); 1919566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Repartition: size (active): %d --> %d, %D 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")); 192a2f3521dSMark F. Adams /* get 'adj' */ 193c5bfad50SMark F. Adams if (cr_bs == 1) { 1949566063dSJacob Faibussowitsch PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 195806fa848SBarry Smith } else { 196a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 197eb07cef2SMark F. Adams Mat tMat; 198a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 199b4fbaa2aSMark F. Adams const PetscScalar *vals; 200b4fbaa2aSMark F. Adams const PetscInt *idx; 201a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 20239d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 203d9558ea9SBarry Smith MatType mtype; 204b4fbaa2aSMark F. Adams 2059566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz)); 2069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs)); 2079566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &M, &N)); 208c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 2099566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat,Ii,&ncols,NULL,NULL)); 210c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 211c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 2129566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL)); 2133ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 2143ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 21558471d46SMark F. Adams } 2166876a03eSMark F. Adams 2179566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat_fine,&mtype)); 2189566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &tMat)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE)); 2209566063dSJacob Faibussowitsch PetscCall(MatSetType(tMat,mtype)); 2219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(tMat,0,d_nnz)); 2229566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz)); 2239566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz,o_nnz)); 224eb07cef2SMark F. Adams 225a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 226c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 2279566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat,ii,&ncols,&idx,&vals)); 228eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 229c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 230eb07cef2SMark F. Adams PetscScalar v = 1.0; 2319566063dSJacob Faibussowitsch PetscCall(MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES)); 232eb07cef2SMark F. Adams } 2339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat,ii,&ncols,&idx,&vals)); 234eb07cef2SMark F. Adams } 2359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY)); 2369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY)); 237eb07cef2SMark F. Adams 238b4fbaa2aSMark F. Adams if (llev++ == -1) { 239b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2409566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.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++) { 270c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 271ce7c7f2fSMark Adams newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ 272eb07cef2SMark F. Adams } 2735ef31b24SMark F. Adams } 2749566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(proc_is, &is_idx)); 2759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&proc_is)); 2765ef31b24SMark F. Adams } 2779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 2785a9b9e01SMark F. Adams 2799566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc)); 2809566063dSJacob Faibussowitsch PetscCall(PetscFree(newproc_idx)); 281849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART],0,0,0,0)); 28231cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 283ce7c7f2fSMark Adams PetscInt targetPE; 284*08401ef6SPierre Jolivet PetscCheck(new_size!=nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen"); 2859566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Number of equations (loc) %D with simple aggregation\n",((PetscObject)pc)->prefix,ncrs_eq)); 286ce7c7f2fSMark Adams targetPE = (rank/rfactor)*expand_factor; 2879566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc)); 288a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 289e33ef3b1SMark F. Adams 29011e60469SMark F. Adams /* 291a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 29211e60469SMark F. Adams */ 2939566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 2947700e67bSMark Adams is_eq_num_prim = is_eq_num; 29511e60469SMark F. Adams /* 296a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 29711e60469SMark F. Adams */ 2989566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 299c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 3009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_newproc)); 301ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new/cr_bs; 302a2f3521dSMark F. Adams 3039566063dSJacob Faibussowitsch PetscCall(PetscFree(counts)); 304885364a3SMark Adams /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */ 305885364a3SMark Adams { 306885364a3SMark Adams Vec src_crd, dest_crd; 307885364a3SMark 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; 308885364a3SMark Adams VecScatter vecscat; 309885364a3SMark Adams PetscScalar *array; 310885364a3SMark Adams IS isscat; 311a2f3521dSMark F. Adams /* move data (for primal equations only) */ 31222063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3139566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &dest_crd)); 3149566063dSJacob Faibussowitsch PetscCall(VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE)); 3159566063dSJacob Faibussowitsch PetscCall(VecSetType(dest_crd,VECSTANDARD)); /* this is needed! */ 31611e60469SMark F. Adams /* 3179d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 318c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 31911e60469SMark F. Adams */ 3209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncrs*node_data_sz, &tidx)); 3219566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_eq_num_prim, &idx)); 3223ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 323c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 324a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 32511e60469SMark F. Adams } 3269566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_eq_num_prim, &idx)); 3279566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat)); 3289566063dSJacob Faibussowitsch PetscCall(PetscFree(tidx)); 32911e60469SMark F. Adams /* 33011e60469SMark F. Adams Create a vector to contain the original vertex information for each element 33111e60469SMark F. Adams */ 3329566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd)); 3339d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3343ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3353ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3369d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 337a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 338c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 3399566063dSJacob Faibussowitsch PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES)); 340d3d6bff4SMark F. Adams } 341038e3b61SMark F. Adams } 342eb07cef2SMark F. Adams } 3439566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(src_crd)); 3449566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(src_crd)); 34511e60469SMark F. Adams /* 34611e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 34711e60469SMark F. Adams to the correct processor 34811e60469SMark F. Adams */ 3499566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat)); 3509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isscat)); 3519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD)); 3529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD)); 3539566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vecscat)); 3549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&src_crd)); 35511e60469SMark F. Adams /* 35611e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 35711e60469SMark F. Adams */ 3589566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data)); 3602fa5cd67SKarl Rupp 3613ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3623ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3632fa5cd67SKarl Rupp 3649566063dSJacob Faibussowitsch PetscCall(VecGetArray(dest_crd, &array)); 3659d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3663ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3679d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 368a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 369c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 370d3d6bff4SMark F. Adams } 371038e3b61SMark F. Adams } 372038e3b61SMark F. Adams } 3739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(dest_crd, &array)); 3749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dest_crd)); 375885364a3SMark Adams } 376a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 377849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 37811e60469SMark F. Adams /* 3797dae84e0SHong Zhang Invert for MatCreateSubMatrix 38011e60469SMark F. Adams */ 3819566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices)); 3829566063dSJacob Faibussowitsch PetscCall(ISSort(new_eq_indices)); /* is this needed? */ 3839566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(new_eq_indices, cr_bs)); 384a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 3859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ 386a2f3521dSMark F. Adams } 3873cb8563fSToby Isaac if (Pcolumnperm) { 3889566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)new_eq_indices)); 3893cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3903cb8563fSToby Isaac } 3919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_num)); 392849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 393849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 394849bee69SMark Adams 395a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 396a2f3521dSMark F. Adams { 397a2f3521dSMark F. Adams Mat mat; 39890db8557SMark Adams PetscBool flg; 3999566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 4009566063dSJacob Faibussowitsch PetscCall(MatGetOption(Cmat, MAT_SPD, &flg)); 40190db8557SMark Adams if (flg) { 4029566063dSJacob Faibussowitsch PetscCall(MatSetOption(mat, MAT_SPD,PETSC_TRUE)); 40390db8557SMark Adams } else { 4049566063dSJacob Faibussowitsch PetscCall(MatGetOption(Cmat, MAT_HERMITIAN, &flg)); 40590db8557SMark Adams if (flg) { 4069566063dSJacob Faibussowitsch PetscCall(MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE)); 40790db8557SMark Adams } else { 40890db8557SMark Adams #if !defined(PETSC_USE_COMPLEX) 4099566063dSJacob Faibussowitsch PetscCall(MatGetOption(Cmat, MAT_SYMMETRIC, &flg)); 41090db8557SMark Adams if (flg) { 4119566063dSJacob Faibussowitsch PetscCall(MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE)); 41290db8557SMark Adams } 41390db8557SMark Adams #endif 41490db8557SMark Adams } 41590db8557SMark Adams } 416a2f3521dSMark F. Adams *a_Amat_crs = mat; 417a2f3521dSMark F. Adams } 4189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cmat)); 419a2f3521dSMark F. Adams 420849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 42111e60469SMark F. Adams /* prolongator */ 42211e60469SMark F. Adams { 42311e60469SMark F. Adams IS findices; 424a2f3521dSMark F. Adams PetscInt Istart,Iend; 425a2f3521dSMark F. Adams Mat Pnew; 42662294041SBarry Smith 4279566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 428849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4299566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm,Iend-Istart,Istart,1,&findices)); 4309566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(findices,f_bs)); 4319566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 4329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&findices)); 4339566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 434c5bfad50SMark F. Adams 435849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4369566063dSJacob Faibussowitsch PetscCall(MatDestroy(a_P_inout)); 437a2f3521dSMark F. Adams 438a2f3521dSMark F. Adams /* output - repartitioned */ 439a2f3521dSMark F. Adams *a_P_inout = Pnew; 440e33ef3b1SMark F. Adams } 4419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&new_eq_indices)); 4425b89ad90SMark F. Adams 443c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 444ce7c7f2fSMark Adams 445ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 446ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 447ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4488bca76a6SMark Adams static PetscInt llev = 2; 4499566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Pinning level %D to the CPU\n",((PetscObject)pc)->prefix,llev++)); 450ce7c7f2fSMark Adams #endif 4519566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE)); 4529566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE)); 453adf5291fSStefano Zampini if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 454ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 455ce7c7f2fSMark Adams PetscMPIInt size; 4569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 457ce7c7f2fSMark Adams if (size > 1) { 458ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data; 4599566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(a->lvec,PETSC_TRUE)); 4609566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(p->lvec,PETSC_TRUE)); 461ce7c7f2fSMark Adams } 462ce7c7f2fSMark Adams } 463ce7c7f2fSMark Adams } 464849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 465849bee69SMark Adams } // processor reduce 4665b89ad90SMark F. Adams PetscFunctionReturn(0); 4675b89ad90SMark F. Adams } 4685b89ad90SMark F. Adams 4694b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2) 4704b1575e2SStefano Zampini { 4714b1575e2SStefano Zampini const char *prefix; 4724b1575e2SStefano Zampini char addp[32]; 4734b1575e2SStefano Zampini PC_MG *mg = (PC_MG*)a_pc->data; 4744b1575e2SStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4754b1575e2SStefano Zampini 4764b1575e2SStefano Zampini PetscFunctionBegin; 4779566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(a_pc,&prefix)); 4789566063dSJacob Faibussowitsch PetscCall(PetscInfo(a_pc,"%s: Square Graph on level %D\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1)); 4799566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Gmat1,Gmat1,NULL,Gmat2)); 4809566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(*Gmat2,prefix)); 4819566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level)); 4829566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(*Gmat2,addp)); 483b4da3a1bSStefano Zampini if ((*Gmat2)->structurally_symmetric) { 4849566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AB)); 485b4da3a1bSStefano Zampini } else { 4869566063dSJacob Faibussowitsch PetscCall(MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 4879566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AtB)); 488b4da3a1bSStefano Zampini } 4899566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(*Gmat2)); 4909566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0)); 4919566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(*Gmat2)); 4929566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0)); 4939566063dSJacob Faibussowitsch PetscCall(MatProductClear(*Gmat2)); 4944b1575e2SStefano Zampini /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 4954b1575e2SStefano Zampini (*Gmat2)->assembled = PETSC_TRUE; 4964b1575e2SStefano Zampini PetscFunctionReturn(0); 4974b1575e2SStefano Zampini } 4984b1575e2SStefano Zampini 4995b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 5005b89ad90SMark F. Adams /* 5015b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 5025b89ad90SMark F. Adams by setting data structures and options. 5035b89ad90SMark F. Adams 5045b89ad90SMark F. Adams Input Parameter: 5055b89ad90SMark F. Adams . pc - the preconditioner context 5065b89ad90SMark F. Adams 5075b89ad90SMark F. Adams */ 5089d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 5095b89ad90SMark F. Adams { 5109d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 5115b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 5122adcac29SMark F. Adams Mat Pmat = pc->pmat; 51318c3aa7eSMark PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS]; 5143b4367a7SBarry Smith MPI_Comm comm; 515c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 51618c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS]; 51718c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 518a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 519569f4572SMark Adams MatInfo info; 520171cca9aSMark Adams PetscBool is_last = PETSC_FALSE; 5215ef31b24SMark F. Adams 5225b89ad90SMark F. Adams PetscFunctionBegin; 5239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 5249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 5259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 526849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 5278abdc6daSStefano Zampini if (pc->setupcalled) { 5288abdc6daSStefano Zampini if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 529878e152fSMark F. Adams /* reset everything */ 5309566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 531878e152fSMark F. Adams pc->setupcalled = 0; 532806fa848SBarry Smith } else { 53384d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 53403a628feSMark F. Adams /* just do Galerkin grids */ 53558471d46SMark F. Adams Mat B,dA,dB; 5369d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 5374555aa8cSStefano Zampini PetscInt gl; 53858471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5399566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB)); 54058471d46SMark F. Adams /* (re)set to get dirty flag */ 5419566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB)); 54258471d46SMark F. Adams 5434555aa8cSStefano Zampini for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) { 5448abdc6daSStefano Zampini MatReuse reuse = MAT_INITIAL_MATRIX ; 545849bee69SMark Adams #if defined(GAMG_STAGES) 546849bee69SMark Adams PetscCall(PetscLogStagePush(gamg_stages[gl])); 547849bee69SMark Adams #endif 5488abdc6daSStefano Zampini /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 5499566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[level]->smoothd,NULL,&B)); 5508abdc6daSStefano Zampini if (B->product) { 5518abdc6daSStefano Zampini if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) { 5528abdc6daSStefano Zampini reuse = MAT_REUSE_MATRIX; 55303a628feSMark F. Adams } 5548abdc6daSStefano Zampini } 5559566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 5568abdc6daSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 5579566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %D\n",((PetscObject)pc)->prefix,level)); 5588abdc6daSStefano Zampini } else { 5599566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: RAP after first solve, new matrix level %D\n",((PetscObject)pc)->prefix,level)); 5608abdc6daSStefano Zampini } 5619566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0)); 5629566063dSJacob Faibussowitsch PetscCall(MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B)); 5639566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0)); 56463b77682SMark Adams if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 5659566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[level]->smoothd,B,B)); 56658471d46SMark F. Adams dB = B; 567849bee69SMark Adams #if defined(GAMG_STAGES) 568849bee69SMark Adams PetscCall(PetscLogStagePop()); 569849bee69SMark Adams #endif 57058471d46SMark F. Adams } 5715f8cf99dSMark F. Adams } 572d5280255SMark F. Adams 5739566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 574849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 57558471d46SMark F. Adams PetscFunctionReturn(0); 576eb07cef2SMark F. Adams } 577878e152fSMark F. Adams } 578f6536408SMark F. Adams 579878e152fSMark F. Adams if (!pc_gamg->data) { 580878e152fSMark F. Adams if (pc_gamg->orig_data) { 5819566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 5829566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 5832fa5cd67SKarl Rupp 584878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 585878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 586878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5872fa5cd67SKarl Rupp 5889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 589878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 590806fa848SBarry Smith } else { 5915f80ce2aSJacob Faibussowitsch PetscCheck(pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 5929566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createdefaultdata(pc,Pmat)); 5939d5b6da9SMark F. Adams } 594878e152fSMark F. Adams } 595878e152fSMark F. Adams 596878e152fSMark F. Adams /* cache original data for reuse */ 5971c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 5989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 599878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 600878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 601878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 602878e152fSMark F. Adams } 603038e3b61SMark F. Adams 604302f38e8SMark F. Adams /* get basic dims */ 6059566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 6069566063dSJacob Faibussowitsch PetscCall(MatGetSize(Pmat, &M, &N)); 60784d3f75bSMark F. Adams 6089566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info)); /* global reduction */ 609569f4572SMark Adams nnz0 = info.nz_used; 610569f4572SMark Adams nnztot = info.nz_used; 6119566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: level %D) N=%D, n data rows=%D, n data cols=%D, nnz/row (ave)=%d, np=%D\n",((PetscObject)pc)->prefix,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size)); 612569f4572SMark Adams 613a2f3521dSMark F. Adams /* Get A_i and R_i */ 61462294041SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 6159ab59c8bSMark Adams pc_gamg->current_level = level; 6165f80ce2aSJacob Faibussowitsch PetscCheck(level < PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level); 6175b89ad90SMark F. Adams level1 = level + 1; 6184555aa8cSStefano Zampini #if defined(GAMG_STAGES) 619849bee69SMark Adams if (!gamg_stages[level]) { 620849bee69SMark Adams char str[32]; 621849bee69SMark Adams sprintf(str,"GAMG Level %d",(int)level); 622849bee69SMark Adams PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 623849bee69SMark Adams } 6249566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(gamg_stages[level])); 625b4fbaa2aSMark F. Adams #endif 626c8b0795cSMark F. Adams { /* construct prolongator */ 627725b86d8SJed Brown Mat Gmat; 6280cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 6297700e67bSMark Adams Mat Prol11; 630c8b0795cSMark F. Adams 6319566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->graph(pc,Aarr[level], &Gmat)); 6329566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 6339566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11)); 634c8b0795cSMark F. Adams 635a2f3521dSMark F. Adams /* could have failed to create new level */ 636a2f3521dSMark F. Adams if (Prol11) { 637f7df55f0SStefano Zampini const char *prefix; 638f7df55f0SStefano Zampini char addp[32]; 639f7df55f0SStefano Zampini 6409d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6419566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 642a2f3521dSMark F. Adams 643fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 644c8b0795cSMark F. Adams /* smooth */ 6459566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 646c8b0795cSMark F. Adams } 647c8b0795cSMark F. Adams 6480c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6491b18a24aSMark Adams PetscInt bs; 650849bee69SMark Adams PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 6519566063dSJacob Faibussowitsch PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 652ffc955d6SMark F. Adams } 653ffc955d6SMark F. Adams 6549566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 6559566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(Prol11,prefix)); 6569566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level)); 6579566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(Prol11,addp)); 65891f31d3dSStefano Zampini /* Always generate the transpose with CUDA 659f7df55f0SStefano Zampini Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 6609566063dSJacob Faibussowitsch PetscCall(MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 6619566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Prol11)); 6624bde40a0SMark Adams Parr[level1] = Prol11; 6634bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 6644bde40a0SMark Adams 6659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat)); 6669566063dSJacob Faibussowitsch PetscCall(PetscCDDestroy(agg_lists)); 667a2f3521dSMark F. Adams } /* construct prolongator scope */ 6687f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 669171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 6709566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: Stop gridding, level %D\n",((PetscObject)pc)->prefix,level)); 6714555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6729566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 673a90e85d9SMark Adams #endif 674c8b0795cSMark F. Adams break; 675c8b0795cSMark F. Adams } 6769566063dSJacob Faibussowitsch PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 677849bee69SMark Adams PetscCheckFalse(is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?"); 678171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 6790e2909e1SMark Adams if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE; 680849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0)); 6819566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 682849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0)); 683a2f3521dSMark F. Adams 6849566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 6859566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 686569f4572SMark Adams nnztot += info.nz_used; 6879566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: %D) N=%D, n data cols=%D, nnz/row (ave)=%d, %D active pes\n",((PetscObject)pc)->prefix,level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe)); 688569f4572SMark Adams 6894555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6909566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 691b4fbaa2aSMark F. Adams #endif 692a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6939ab59c8bSMark Adams if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) { 6949566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",((PetscObject)pc)->prefix,level,M/bs)); 695a90e85d9SMark Adams level++; 696a90e85d9SMark Adams break; 697a90e85d9SMark Adams } 698c8b0795cSMark F. Adams } /* levels */ 6999566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 700c8b0795cSMark F. Adams 7019566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: %D levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0)); 7029d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 7035b89ad90SMark F. Adams fine_level = level; 7049566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc,pc_gamg->Nlevels,NULL)); 7055b89ad90SMark F. Adams 70662294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 7070ed2132dSStefano Zampini 708d5280255SMark F. Adams /* set default smoothers & set operators */ 70962294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 710ffc955d6SMark F. Adams KSP smoother; 711ffc955d6SMark F. Adams PC subpc; 712a2f3521dSMark F. Adams 7139566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7149566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 715ffc955d6SMark F. Adams 7169566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 717a2f3521dSMark F. Adams /* set ops */ 7189566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 7199566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level+1])); 720a2f3521dSMark F. Adams 721a2f3521dSMark F. Adams /* set defaults */ 7229566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 723a2f3521dSMark F. Adams 7240c3bc534SBarry Smith /* set blocks for ASM smoother that uses the 'aggregates' */ 7250c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 7262d3561bbSSatish Balay PetscInt sz; 7277a28f3e5SMark Adams IS *iss; 728a2f3521dSMark F. Adams 7292d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7307a28f3e5SMark Adams iss = ASMLocalIDsArr[level]; 7319566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCASM)); 7329566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(subpc, 0)); 7339566063dSJacob Faibussowitsch PetscCall(PCASMSetType(subpc,PC_ASM_BASIC)); 7347f66b68fSMark Adams if (!sz) { 735ffc955d6SMark F. Adams IS is; 7369566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 7379566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 7389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 739806fa848SBarry Smith } else { 740a94c3b12SMark F. Adams PetscInt kk; 7419566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss)); 742a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 7439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iss[kk])); 744a94c3b12SMark F. Adams } 7459566063dSJacob Faibussowitsch PetscCall(PetscFree(iss)); 746ffc955d6SMark F. Adams } 7470298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 748ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 749806fa848SBarry Smith } else { 7509566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCJACOBI)); 751ffc955d6SMark F. Adams } 752d5280255SMark F. Adams } 753d5280255SMark F. Adams { 754d5280255SMark F. Adams /* coarse grid */ 755d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 756d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 7570ed2132dSStefano Zampini 7589566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7599566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 760cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 7619566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 7629566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 7639566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCBJACOBI)); 7649566063dSJacob Faibussowitsch PetscCall(PCSetUp(subpc)); 7659566063dSJacob Faibussowitsch PetscCall(PCBJacobiGetSubKSP(subpc,&ii,&first,&k2)); 7665f80ce2aSJacob Faibussowitsch PetscCheck(ii == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 7679566063dSJacob Faibussowitsch PetscCall(KSPGetPC(k2[0],&pc2)); 7689566063dSJacob Faibussowitsch PetscCall(PCSetType(pc2, PCLU)); 7699566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS)); 7709566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1)); 7719566063dSJacob Faibussowitsch PetscCall(KSPSetType(k2[0], KSPPREONLY)); 772d5280255SMark F. Adams } 773cf8ae1d3SMark Adams } 774d5280255SMark F. Adams 775d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 776849bee69SMark Adams PetscErrorCode ierr; 7779566063dSJacob Faibussowitsch ierr = PetscObjectOptionsBegin((PetscObject)pc);PetscCall(ierr); 7789566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc)); 7799566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 7809566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 781d5280255SMark F. Adams 78218c3aa7eSMark /* setup cheby eigen estimates from SA */ 7837e6512fdSJed Brown if (pc_gamg->use_sa_esteig) { 78418c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 78518c3aa7eSMark KSP smoother; 78618c3aa7eSMark PetscBool ischeb; 7870ed2132dSStefano Zampini 7889566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb)); 79018c3aa7eSMark if (ischeb) { 79118c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data; 7920ed2132dSStefano Zampini 7932de708cbSJed Brown // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 794efe053fcSJed Brown if (mg->max_eigen_DinvA[level] > 0) { 7957e6512fdSJed Brown // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 7967e6512fdSJed Brown // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 79718c3aa7eSMark PetscReal emax,emin; 7980ed2132dSStefano Zampini 79918c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 80018c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 8019566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",((PetscObject)pc)->prefix,level,Aarr[level]->rmap->N,(double)emax,(double)emin)); 8020a94a983SJed Brown cheb->emin_provided = emin; 8030a94a983SJed Brown cheb->emax_provided = emax; 80418c3aa7eSMark } 80518c3aa7eSMark } 80618c3aa7eSMark } 80718c3aa7eSMark } 8080ed2132dSStefano Zampini 8099566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8100ed2132dSStefano Zampini 811d5280255SMark F. Adams /* clean up */ 812d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 8139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Parr[level])); 8149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aarr[level])); 8155b89ad90SMark F. Adams } 816806fa848SBarry Smith } else { 8175f8cf99dSMark F. Adams KSP smoother; 8180ed2132dSStefano Zampini 8199566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n")); 8209566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 8219566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 8229566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPPREONLY)); 8239566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8245f8cf99dSMark F. Adams } 825849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 8265b89ad90SMark F. Adams PetscFunctionReturn(0); 8275b89ad90SMark F. Adams } 8285b89ad90SMark F. Adams 829eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 8305b89ad90SMark F. Adams /* 8315b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 8325b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 8335b89ad90SMark F. Adams 8345b89ad90SMark F. Adams Input Parameter: 8355b89ad90SMark F. Adams . pc - the preconditioner context 8365b89ad90SMark F. Adams 8375b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 8385b89ad90SMark F. Adams */ 8395b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 8405b89ad90SMark F. Adams { 8415b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 8425b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 8435b89ad90SMark F. Adams 8445b89ad90SMark F. Adams PetscFunctionBegin; 8459566063dSJacob Faibussowitsch PetscCall(PCReset_GAMG(pc)); 8469b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 8479566063dSJacob Faibussowitsch PetscCall((*pc_gamg->ops->destroy)(pc)); 8489b8ffb57SJed Brown } 8499566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->ops)); 8509566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 8519566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg)); 8529566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 8535b89ad90SMark F. Adams PetscFunctionReturn(0); 8545b89ad90SMark F. Adams } 8555b89ad90SMark F. Adams 856676e1743SMark F. Adams /*@ 857cab9ed1eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 858676e1743SMark F. Adams 8591cc46a46SBarry Smith Logically Collective on PC 860676e1743SMark F. Adams 861676e1743SMark F. Adams Input Parameters: 8621cc46a46SBarry Smith + pc - the preconditioner context 8631cc46a46SBarry Smith - n - the number of equations 864676e1743SMark F. Adams 865676e1743SMark F. Adams Options Database Key: 866147403d9SBarry Smith . -pc_gamg_process_eq_limit <limit> - set the limit 867676e1743SMark F. Adams 86895452b02SPatrick Sanan Notes: 86995452b02SPatrick Sanan GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 870cab9ed1eSBarry Smith that has degrees of freedom 871cab9ed1eSBarry Smith 872676e1743SMark F. Adams Level: intermediate 873676e1743SMark F. Adams 874c9567895SMark .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors() 875676e1743SMark F. Adams @*/ 876676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 877676e1743SMark F. Adams { 878676e1743SMark F. Adams PetscFunctionBegin; 879676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 880cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n)); 881676e1743SMark F. Adams PetscFunctionReturn(0); 882676e1743SMark F. Adams } 883676e1743SMark F. Adams 8841e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 885676e1743SMark F. Adams { 886c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 887c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 888676e1743SMark F. Adams 889676e1743SMark F. Adams PetscFunctionBegin; 8909d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 891676e1743SMark F. Adams PetscFunctionReturn(0); 892676e1743SMark F. Adams } 893676e1743SMark F. Adams 894389730f3SMark F. Adams /*@ 895cab9ed1eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 896389730f3SMark F. Adams 897389730f3SMark F. Adams Collective on PC 898389730f3SMark F. Adams 899389730f3SMark F. Adams Input Parameters: 9001cc46a46SBarry Smith + pc - the preconditioner context 9011cc46a46SBarry Smith - n - maximum number of equations to aim for 902389730f3SMark F. Adams 903389730f3SMark F. Adams Options Database Key: 904147403d9SBarry Smith . -pc_gamg_coarse_eq_limit <limit> - set the limit 905389730f3SMark F. Adams 90674329af1SBarry Smith Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 90774329af1SBarry Smith has less than 1000 unknowns. 90874329af1SBarry Smith 909389730f3SMark F. Adams Level: intermediate 910389730f3SMark F. Adams 911c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors() 912389730f3SMark F. Adams @*/ 913389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 914389730f3SMark F. Adams { 915389730f3SMark F. Adams PetscFunctionBegin; 916389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 917cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n)); 918389730f3SMark F. Adams PetscFunctionReturn(0); 919389730f3SMark F. Adams } 920389730f3SMark F. Adams 9211e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 922389730f3SMark F. Adams { 923389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 924389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 925389730f3SMark F. Adams 926389730f3SMark F. Adams PetscFunctionBegin; 9279d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 928389730f3SMark F. Adams PetscFunctionReturn(0); 929389730f3SMark F. Adams } 930389730f3SMark F. Adams 931676e1743SMark F. Adams /*@ 932cab9ed1eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 933676e1743SMark F. Adams 934676e1743SMark F. Adams Collective on PC 935676e1743SMark F. Adams 936676e1743SMark F. Adams Input Parameters: 9371cc46a46SBarry Smith + pc - the preconditioner context 9381cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 939676e1743SMark F. Adams 940676e1743SMark F. Adams Options Database Key: 941147403d9SBarry Smith . -pc_gamg_repartition <true,false> - turn on the repartitioning 942676e1743SMark F. Adams 94395452b02SPatrick Sanan Notes: 94495452b02SPatrick Sanan this will generally improve the loading balancing of the work on each level 945cab9ed1eSBarry Smith 946676e1743SMark F. Adams Level: intermediate 947676e1743SMark F. Adams 948676e1743SMark F. Adams @*/ 949cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 950676e1743SMark F. Adams { 951676e1743SMark F. Adams PetscFunctionBegin; 952676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 953cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n)); 954676e1743SMark F. Adams PetscFunctionReturn(0); 955676e1743SMark F. Adams } 956676e1743SMark F. Adams 957cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 958676e1743SMark F. Adams { 959c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 960c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 961676e1743SMark F. Adams 962676e1743SMark F. Adams PetscFunctionBegin; 9639d5b6da9SMark F. Adams pc_gamg->repart = n; 964676e1743SMark F. Adams PetscFunctionReturn(0); 965676e1743SMark F. Adams } 966676e1743SMark F. Adams 967dfd5c07aSMark F. Adams /*@ 9687e6512fdSJed Brown PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother 96918c3aa7eSMark 97018c3aa7eSMark Collective on PC 97118c3aa7eSMark 97218c3aa7eSMark Input Parameters: 97318c3aa7eSMark + pc - the preconditioner context 97418c3aa7eSMark - n - number of its 97518c3aa7eSMark 97618c3aa7eSMark Options Database Key: 977147403d9SBarry Smith . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 97818c3aa7eSMark 97918c3aa7eSMark Notes: 9807e6512fdSJed 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$. 9817e6512fdSJed Brown Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 9827e6512fdSJed Brown If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused. 983efe053fcSJed Brown This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used. 984efe053fcSJed Brown It became default in PETSc 3.17. 98518c3aa7eSMark 9867e6512fdSJed Brown Level: advanced 98718c3aa7eSMark 98873f7197eSJed Brown .seealso: KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet() 98918c3aa7eSMark @*/ 99018c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 99118c3aa7eSMark { 99218c3aa7eSMark PetscFunctionBegin; 99318c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 994cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n)); 99518c3aa7eSMark PetscFunctionReturn(0); 99618c3aa7eSMark } 99718c3aa7eSMark 9980ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 99918c3aa7eSMark { 100018c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 100118c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 100218c3aa7eSMark 100318c3aa7eSMark PetscFunctionBegin; 10047e6512fdSJed Brown pc_gamg->use_sa_esteig = n; 100518c3aa7eSMark PetscFunctionReturn(0); 100618c3aa7eSMark } 100718c3aa7eSMark 100818c3aa7eSMark /*@ 100918c3aa7eSMark PCGAMGSetEigenvalues - Set eigenvalues 101018c3aa7eSMark 101118c3aa7eSMark Collective on PC 101218c3aa7eSMark 101318c3aa7eSMark Input Parameters: 101418c3aa7eSMark + pc - the preconditioner context 101518c3aa7eSMark - emax - max eigenvalue 101618c3aa7eSMark . emin - min eigenvalue 101718c3aa7eSMark 101818c3aa7eSMark Options Database Key: 1019147403d9SBarry Smith . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 102018c3aa7eSMark 102118c3aa7eSMark Level: intermediate 102218c3aa7eSMark 102373f7197eSJed Brown .seealso: PCGAMGSetUseSAEstEig() 102418c3aa7eSMark @*/ 102518c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin) 102618c3aa7eSMark { 102718c3aa7eSMark PetscFunctionBegin; 102818c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1029cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin)); 103018c3aa7eSMark PetscFunctionReturn(0); 103118c3aa7eSMark } 103241ffd417SStefano Zampini 103318c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin) 103418c3aa7eSMark { 103518c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 103618c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 103718c3aa7eSMark 103818c3aa7eSMark PetscFunctionBegin; 10395f80ce2aSJacob 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); 10405f80ce2aSJacob 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); 104118c3aa7eSMark pc_gamg->emax = emax; 104218c3aa7eSMark pc_gamg->emin = emin; 104318c3aa7eSMark PetscFunctionReturn(0); 104418c3aa7eSMark } 104518c3aa7eSMark 104618c3aa7eSMark /*@ 1047cab9ed1eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1048dfd5c07aSMark F. Adams 1049dfd5c07aSMark F. Adams Collective on PC 1050dfd5c07aSMark F. Adams 1051dfd5c07aSMark F. Adams Input Parameters: 10521cc46a46SBarry Smith + pc - the preconditioner context 10531cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 1054dfd5c07aSMark F. Adams 1055dfd5c07aSMark F. Adams Options Database Key: 1056147403d9SBarry Smith . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1057dfd5c07aSMark F. Adams 1058dfd5c07aSMark F. Adams Level: intermediate 1059dfd5c07aSMark F. Adams 106095452b02SPatrick Sanan Notes: 1061147403d9SBarry Smith May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1062cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1063cab9ed1eSBarry Smith 1064dfd5c07aSMark F. Adams @*/ 10651cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1066dfd5c07aSMark F. Adams { 1067dfd5c07aSMark F. Adams PetscFunctionBegin; 1068dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1069cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n)); 1070dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1071dfd5c07aSMark F. Adams } 1072dfd5c07aSMark F. Adams 10731cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1074dfd5c07aSMark F. Adams { 1075dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1076dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1077dfd5c07aSMark F. Adams 1078dfd5c07aSMark F. Adams PetscFunctionBegin; 1079dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1080dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1081dfd5c07aSMark F. Adams } 1082dfd5c07aSMark F. Adams 1083ffc955d6SMark F. Adams /*@ 1084cab9ed1eSBarry 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. 1085ffc955d6SMark F. Adams 1086ffc955d6SMark F. Adams Collective on PC 1087ffc955d6SMark F. Adams 1088ffc955d6SMark F. Adams Input Parameters: 1089cab9ed1eSBarry Smith + pc - the preconditioner context 1090cab9ed1eSBarry Smith - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not 1091ffc955d6SMark F. Adams 1092ffc955d6SMark F. Adams Options Database Key: 1093147403d9SBarry Smith . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1094ffc955d6SMark F. Adams 1095ffc955d6SMark F. Adams Level: intermediate 1096ffc955d6SMark F. Adams 1097ffc955d6SMark F. Adams @*/ 1098cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1099ffc955d6SMark F. Adams { 1100ffc955d6SMark F. Adams PetscFunctionBegin; 1101ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1102cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg)); 1103ffc955d6SMark F. Adams PetscFunctionReturn(0); 1104ffc955d6SMark F. Adams } 1105ffc955d6SMark F. Adams 1106cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1107ffc955d6SMark F. Adams { 1108ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1109ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1110ffc955d6SMark F. Adams 1111ffc955d6SMark F. Adams PetscFunctionBegin; 1112cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 1113ffc955d6SMark F. Adams PetscFunctionReturn(0); 1114ffc955d6SMark F. Adams } 1115ffc955d6SMark F. Adams 1116171cca9aSMark Adams /*@ 1117cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1118171cca9aSMark Adams 1119171cca9aSMark Adams Collective on PC 1120171cca9aSMark Adams 1121171cca9aSMark Adams Input Parameters: 1122171cca9aSMark Adams + pc - the preconditioner context 1123cf8ae1d3SMark Adams - flg - PETSC_TRUE to not force coarse grid onto one processor 1124171cca9aSMark Adams 1125171cca9aSMark Adams Options Database Key: 1126147403d9SBarry Smith . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1127171cca9aSMark Adams 1128171cca9aSMark Adams Level: intermediate 1129171cca9aSMark Adams 113039d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids() 1131171cca9aSMark Adams @*/ 1132171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1133171cca9aSMark Adams { 1134171cca9aSMark Adams PetscFunctionBegin; 1135171cca9aSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1136cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg)); 1137171cca9aSMark Adams PetscFunctionReturn(0); 1138171cca9aSMark Adams } 1139171cca9aSMark Adams 1140171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1141171cca9aSMark Adams { 1142171cca9aSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1143171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1144171cca9aSMark Adams 1145171cca9aSMark Adams PetscFunctionBegin; 1146171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 1147ffc955d6SMark F. Adams PetscFunctionReturn(0); 1148ffc955d6SMark F. Adams } 1149ffc955d6SMark F. Adams 11504ef23d27SMark F. Adams /*@ 1151ce7c7f2fSMark Adams PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1152ce7c7f2fSMark Adams 1153ce7c7f2fSMark Adams Collective on PC 1154ce7c7f2fSMark Adams 1155ce7c7f2fSMark Adams Input Parameters: 1156ce7c7f2fSMark Adams + pc - the preconditioner context 1157ce7c7f2fSMark Adams - flg - PETSC_TRUE to pin coarse grids to CPU 1158ce7c7f2fSMark Adams 1159ce7c7f2fSMark Adams Options Database Key: 1160147403d9SBarry Smith . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1161ce7c7f2fSMark Adams 1162ce7c7f2fSMark Adams Level: intermediate 1163ce7c7f2fSMark Adams 116439d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve() 1165ce7c7f2fSMark Adams @*/ 1166ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1167ce7c7f2fSMark Adams { 1168ce7c7f2fSMark Adams PetscFunctionBegin; 1169ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1170cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg)); 1171ce7c7f2fSMark Adams PetscFunctionReturn(0); 1172ce7c7f2fSMark Adams } 1173ce7c7f2fSMark Adams 1174ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1175ce7c7f2fSMark Adams { 1176ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1177ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1178ce7c7f2fSMark Adams 1179ce7c7f2fSMark Adams PetscFunctionBegin; 1180ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 1181ce7c7f2fSMark Adams PetscFunctionReturn(0); 1182ce7c7f2fSMark Adams } 1183ce7c7f2fSMark Adams 1184ce7c7f2fSMark Adams /*@ 1185147403d9SBarry Smith PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1186ce7c7f2fSMark Adams 1187ce7c7f2fSMark Adams Collective on PC 1188ce7c7f2fSMark Adams 1189ce7c7f2fSMark Adams Input Parameters: 1190ce7c7f2fSMark Adams + pc - the preconditioner context 1191ce7c7f2fSMark Adams - flg - Layout type 1192ce7c7f2fSMark Adams 1193ce7c7f2fSMark Adams Options Database Key: 1194147403d9SBarry Smith . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1195ce7c7f2fSMark Adams 1196ce7c7f2fSMark Adams Level: intermediate 1197ce7c7f2fSMark Adams 119839d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids() 1199ce7c7f2fSMark Adams @*/ 1200ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1201ce7c7f2fSMark Adams { 1202ce7c7f2fSMark Adams PetscFunctionBegin; 1203ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1204cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg)); 1205ce7c7f2fSMark Adams PetscFunctionReturn(0); 1206ce7c7f2fSMark Adams } 1207ce7c7f2fSMark Adams 1208ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1209ce7c7f2fSMark Adams { 1210ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1211ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1212ce7c7f2fSMark Adams 1213ce7c7f2fSMark Adams PetscFunctionBegin; 1214ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 1215ce7c7f2fSMark Adams PetscFunctionReturn(0); 1216ce7c7f2fSMark Adams } 1217ce7c7f2fSMark Adams 1218ce7c7f2fSMark Adams /*@ 12191cc46a46SBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 12204ef23d27SMark F. Adams 12214ef23d27SMark F. Adams Not collective on PC 12224ef23d27SMark F. Adams 12234ef23d27SMark F. Adams Input Parameters: 12241cc46a46SBarry Smith + pc - the preconditioner 12251cc46a46SBarry Smith - n - the maximum number of levels to use 12264ef23d27SMark F. Adams 12274ef23d27SMark F. Adams Options Database Key: 1228147403d9SBarry Smith . -pc_mg_levels <n> - set the maximum number of levels to allow 12294ef23d27SMark F. Adams 12304ef23d27SMark F. Adams Level: intermediate 12314ef23d27SMark F. Adams 12324ef23d27SMark F. Adams @*/ 12334ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 12344ef23d27SMark F. Adams { 12354ef23d27SMark F. Adams PetscFunctionBegin; 12364ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1237cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n)); 12384ef23d27SMark F. Adams PetscFunctionReturn(0); 12394ef23d27SMark F. Adams } 12404ef23d27SMark F. Adams 12411e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 12424ef23d27SMark F. Adams { 12434ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 12444ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12454ef23d27SMark F. Adams 12464ef23d27SMark F. Adams PetscFunctionBegin; 12479d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12484ef23d27SMark F. Adams PetscFunctionReturn(0); 12494ef23d27SMark F. Adams } 12504ef23d27SMark F. Adams 12513542efc5SMark F. Adams /*@ 12523542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12533542efc5SMark F. Adams 12543542efc5SMark F. Adams Not collective on PC 12553542efc5SMark F. Adams 12563542efc5SMark F. Adams Input Parameters: 12571cc46a46SBarry Smith + pc - the preconditioner context 1258c9567895SMark . 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 1259055c8bd0SJed Brown - n - number of threshold values provided in array 12603542efc5SMark F. Adams 12613542efc5SMark F. Adams Options Database Key: 1262147403d9SBarry Smith . -pc_gamg_threshold <threshold> - the threshold to drop edges 12633542efc5SMark F. Adams 126495452b02SPatrick Sanan Notes: 1265af3c827dSMark 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. 1266af3c827dSMark Adams Before coarsening or aggregating the graph, GAMG 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. 1267cab9ed1eSBarry Smith 1268055c8bd0SJed Brown If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening. 1269055c8bd0SJed Brown In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold(). 1270055c8bd0SJed Brown If n is greater than the total number of levels, the excess entries in threshold will not be used. 1271055c8bd0SJed Brown 12723542efc5SMark F. Adams Level: intermediate 12733542efc5SMark F. Adams 1274af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph() 12753542efc5SMark F. Adams @*/ 1276c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 12773542efc5SMark F. Adams { 12783542efc5SMark F. Adams PetscFunctionBegin; 12793542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1280055c8bd0SJed Brown if (n) PetscValidRealPointer(v,2); 1281cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n)); 12823542efc5SMark F. Adams PetscFunctionReturn(0); 12833542efc5SMark F. Adams } 12843542efc5SMark F. Adams 1285c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 12863542efc5SMark F. Adams { 1287c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1288c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1289c1eae691SMark Adams PetscInt i; 1290c1eae691SMark Adams PetscFunctionBegin; 1291055c8bd0SJed Brown for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1292055c8bd0SJed Brown for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale; 1293c1eae691SMark Adams PetscFunctionReturn(0); 1294c1eae691SMark Adams } 1295c1eae691SMark Adams 1296c1eae691SMark Adams /*@ 1297147403d9SBarry Smith PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids 1298c9567895SMark 1299c9567895SMark Collective on PC 1300c9567895SMark 1301c9567895SMark Input Parameters: 1302c9567895SMark + pc - the preconditioner context 1303c9567895SMark . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda 1304c9567895SMark - n - number of values provided in array 1305c9567895SMark 1306c9567895SMark Options Database Key: 1307147403d9SBarry Smith . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1308c9567895SMark 1309c9567895SMark Level: intermediate 1310c9567895SMark 1311c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim() 1312c9567895SMark @*/ 1313c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1314c9567895SMark { 1315c9567895SMark PetscFunctionBegin; 1316c9567895SMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1317c9567895SMark if (n) PetscValidIntPointer(v,2); 1318cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n)); 1319c9567895SMark PetscFunctionReturn(0); 1320c9567895SMark } 1321c9567895SMark 1322c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1323c9567895SMark { 1324c9567895SMark PC_MG *mg = (PC_MG*)pc->data; 1325c9567895SMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1326c9567895SMark PetscInt i; 1327c9567895SMark PetscFunctionBegin; 1328c9567895SMark for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1329c9567895SMark for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1330c9567895SMark PetscFunctionReturn(0); 1331c9567895SMark } 1332c9567895SMark 1333c9567895SMark /*@ 1334c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1335c1eae691SMark Adams 1336c1eae691SMark Adams Not collective on PC 1337c1eae691SMark Adams 1338c1eae691SMark Adams Input Parameters: 1339c1eae691SMark Adams + pc - the preconditioner context 1340147403d9SBarry Smith - scale - the threshold value reduction, usually < 1.0 1341c1eae691SMark Adams 1342c1eae691SMark Adams Options Database Key: 1343147403d9SBarry Smith . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1344c1eae691SMark Adams 1345055c8bd0SJed Brown Notes: 1346055c8bd0SJed Brown The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold(). 1347055c8bd0SJed Brown This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold(). 1348055c8bd0SJed Brown 1349c1eae691SMark Adams Level: advanced 1350c1eae691SMark Adams 1351055c8bd0SJed Brown .seealso: PCGAMGSetThreshold() 1352c1eae691SMark Adams @*/ 1353c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1354c1eae691SMark Adams { 13553542efc5SMark F. Adams PetscFunctionBegin; 1356c1eae691SMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1357cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v)); 1358c1eae691SMark Adams PetscFunctionReturn(0); 1359c1eae691SMark Adams } 1360c1eae691SMark Adams 1361c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1362c1eae691SMark Adams { 1363c1eae691SMark Adams PC_MG *mg = (PC_MG*)pc->data; 1364c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1365c1eae691SMark Adams PetscFunctionBegin; 1366c1eae691SMark Adams pc_gamg->threshold_scale = v; 13673542efc5SMark F. Adams PetscFunctionReturn(0); 13683542efc5SMark F. Adams } 13693542efc5SMark F. Adams 1370e20c40e8SBarry Smith /*@C 1371c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1372676e1743SMark F. Adams 1373676e1743SMark F. Adams Collective on PC 1374676e1743SMark F. Adams 1375676e1743SMark F. Adams Input Parameters: 1376c60c7ad4SBarry Smith + pc - the preconditioner context 1377c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1378676e1743SMark F. Adams 1379676e1743SMark F. Adams Options Database Key: 1380cab9ed1eSBarry Smith . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1381676e1743SMark F. Adams 1382676e1743SMark F. Adams Level: intermediate 1383676e1743SMark F. Adams 1384cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType 1385676e1743SMark F. Adams @*/ 138619fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1387676e1743SMark F. Adams { 1388676e1743SMark F. Adams PetscFunctionBegin; 1389676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1390cac4c232SBarry Smith PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type)); 1391676e1743SMark F. Adams PetscFunctionReturn(0); 1392676e1743SMark F. Adams } 1393676e1743SMark F. Adams 1394e20c40e8SBarry Smith /*@C 1395c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1396c60c7ad4SBarry Smith 1397c60c7ad4SBarry Smith Collective on PC 1398c60c7ad4SBarry Smith 1399c60c7ad4SBarry Smith Input Parameter: 1400c60c7ad4SBarry Smith . pc - the preconditioner context 1401c60c7ad4SBarry Smith 1402c60c7ad4SBarry Smith Output Parameter: 1403c60c7ad4SBarry Smith . type - the type of algorithm used 1404c60c7ad4SBarry Smith 1405c60c7ad4SBarry Smith Level: intermediate 1406c60c7ad4SBarry Smith 14071c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType 1408c60c7ad4SBarry Smith @*/ 1409c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1410c60c7ad4SBarry Smith { 1411c60c7ad4SBarry Smith PetscFunctionBegin; 1412c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1413cac4c232SBarry Smith PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type)); 1414c60c7ad4SBarry Smith PetscFunctionReturn(0); 1415c60c7ad4SBarry Smith } 1416c60c7ad4SBarry Smith 1417c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1418c60c7ad4SBarry Smith { 1419c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1420c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1421c60c7ad4SBarry Smith 1422c60c7ad4SBarry Smith PetscFunctionBegin; 1423c60c7ad4SBarry Smith *type = pc_gamg->type; 1424c60c7ad4SBarry Smith PetscFunctionReturn(0); 1425c60c7ad4SBarry Smith } 1426c60c7ad4SBarry Smith 14271e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1428676e1743SMark F. Adams { 14291ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 14301ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 14315f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(PC); 1432676e1743SMark F. Adams 1433676e1743SMark F. Adams PetscFunctionBegin; 1434c60c7ad4SBarry Smith pc_gamg->type = type; 14359566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(GAMGList,type,&r)); 14365f80ce2aSJacob Faibussowitsch PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 14371ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 14389566063dSJacob Faibussowitsch PetscCall((*pc_gamg->ops->destroy)(pc)); 14399566063dSJacob Faibussowitsch PetscCall(PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps))); 1440e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 14413ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 14423ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 14433ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 14443ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 14453ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 14469566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 14473ae0bb68SMark Adams pc_gamg->data_sz = 0; 14481ab5ffc9SJed Brown } 14499566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 14509566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type,&pc_gamg->gamg_type_name)); 14519566063dSJacob Faibussowitsch PetscCall((*r)(pc)); 1452676e1743SMark F. Adams PetscFunctionReturn(0); 1453676e1743SMark F. Adams } 1454676e1743SMark F. Adams 14555adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 14565adeb434SBarry Smith { 14575adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 14585adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1459e7d4b4cbSMark Adams PetscReal gc=0, oc=0; 146090db8557SMark Adams 14615adeb434SBarry Smith PetscFunctionBegin; 14629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," GAMG specific options\n")); 14639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =")); 14649566063dSJacob Faibussowitsch for (PetscInt i=0;i<mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i])); 14659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 14669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale)); 1467cab9ed1eSBarry Smith if (pc_gamg->use_aggs_in_asm) { 14689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n")); 1469cab9ed1eSBarry Smith } 1470171cca9aSMark Adams if (pc_gamg->use_parallel_coarse_grid_solver) { 14719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n")); 1472171cca9aSMark Adams } 1473ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1474ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 14759566063dSJacob Faibussowitsch /* PetscCall(PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n")); */ 1476ce7c7f2fSMark Adams } 1477ce7c7f2fSMark Adams #endif 1478ce7c7f2fSMark Adams /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */ 14799566063dSJacob Faibussowitsch /* PetscCall(PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n")); */ 1480ce7c7f2fSMark Adams /* } else { */ 14819566063dSJacob Faibussowitsch /* PetscCall(PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n")); */ 1482ce7c7f2fSMark Adams /* } */ 14835adeb434SBarry Smith if (pc_gamg->ops->view) { 14849566063dSJacob Faibussowitsch PetscCall((*pc_gamg->ops->view)(pc,viewer)); 14855adeb434SBarry Smith } 14869566063dSJacob Faibussowitsch PetscCall(PCMGGetGridComplexity(pc,&gc,&oc)); 14879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Complexity: grid = %g operator = %g\n",gc,oc)); 14885adeb434SBarry Smith PetscFunctionReturn(0); 14895adeb434SBarry Smith } 14905adeb434SBarry Smith 14914416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 14925b89ad90SMark F. Adams { 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)); 15039566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(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)); 1505bd94a7aaSJed Brown if (flag) { 15069566063dSJacob Faibussowitsch PetscCall(PCGAMGSetType(pc,tname)); 15071ab5ffc9SJed Brown } 15089566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL)); 15099566063dSJacob Faibussowitsch 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)); 15109566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL)); 15119566063dSJacob 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)); 15129566063dSJacob 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)); 15139566063dSJacob 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)); 15149566063dSJacob Faibussowitsch 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,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL)); 15159566063dSJacob 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)); 15169566063dSJacob 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)); 15179566063dSJacob 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)); 151818c3aa7eSMark n = PETSC_MG_MAXLEVELS; 15199566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag)); 152018c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1521efd3c5ceSMark Adams if (!flag) n = 1; 1522c1eae691SMark Adams i = n; 152318c3aa7eSMark do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1524c1eae691SMark Adams } 1525c9567895SMark n = PETSC_MG_MAXLEVELS; 15269566063dSJacob 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)); 1527c9567895SMark if (!flag) i = 0; 1528c9567895SMark else i = n; 1529c9567895SMark do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS); 15309566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL)); 153118c3aa7eSMark { 153218c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 153318c3aa7eSMark n = 2; 15349566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag)); 153518c3aa7eSMark if (flag) { 1536*08401ef6SPierre Jolivet PetscCheck(n == 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 15379566063dSJacob Faibussowitsch PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 153818c3aa7eSMark } 153918c3aa7eSMark } 1540b7cbab4eSMark Adams /* set options for subtype */ 15419566063dSJacob Faibussowitsch if (pc_gamg->ops->setfromoptions) PetscCall((*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc)); 154218c3aa7eSMark 15439566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 15449566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "")); 15459566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 15465b89ad90SMark F. Adams PetscFunctionReturn(0); 15475b89ad90SMark F. Adams } 15485b89ad90SMark F. Adams 15495b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 15505b89ad90SMark F. Adams /*MC 15511cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 15525b89ad90SMark F. Adams 1553280d9858SJed Brown Options Database Keys: 1554cab9ed1eSBarry Smith + -pc_gamg_type <type> - one of agg, geo, or classical 1555cab9ed1eSBarry Smith . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1556cab9ed1eSBarry Smith . -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1557cab9ed1eSBarry Smith . -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother 1558cab9ed1eSBarry Smith . -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> 1559cab9ed1eSBarry Smith equations on each process that has degrees of freedom 1560cab9ed1eSBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 15616008e27bSRichard Tran Mills . -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level 1562c1eae691SMark Adams - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1563cab9ed1eSBarry Smith 1564cab9ed1eSBarry Smith Options Database Keys for default Aggregation: 1565cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1566cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 1567cab9ed1eSBarry Smith - -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it 1568cab9ed1eSBarry Smith 1569db9745e2SBarry Smith Multigrid options: 1570db9745e2SBarry Smith + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1571db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1572db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1573cab9ed1eSBarry Smith - -pc_mg_levels <levels> - Number of levels of multigrid to use. 15745b89ad90SMark F. Adams 157595452b02SPatrick Sanan Notes: 157695452b02SPatrick Sanan In order to obtain good performance for PCGAMG for vector valued problems you must 1577db9745e2SBarry Smith Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1578db9745e2SBarry Smith Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1579db9745e2SBarry Smith See the Users Manual Chapter 4 for more details 15801cc46a46SBarry Smith 15815b89ad90SMark F. Adams Level: intermediate 1582280d9858SJed Brown 15831cc46a46SBarry Smith .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 158473f7197eSJed Brown PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig() 15855b89ad90SMark F. Adams M*/ 1586b2573a8aSBarry Smith 15878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 15885b89ad90SMark F. Adams { 15895b89ad90SMark F. Adams PC_GAMG *pc_gamg; 15905b89ad90SMark F. Adams PC_MG *mg; 15915b89ad90SMark F. Adams 15925b89ad90SMark F. Adams PetscFunctionBegin; 15931c1aac46SBarry Smith /* register AMG type */ 15949566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 15951c1aac46SBarry Smith 15965b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 15979566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); 15989566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 15995b89ad90SMark F. Adams 16005b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 16019566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&pc_gamg)); 16029566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 16035b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 16045b89ad90SMark F. Adams mg->innerctx = pc_gamg; 16055b89ad90SMark F. Adams 16069566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&pc_gamg->ops)); 16071ab5ffc9SJed Brown 16089d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 16099d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 16100a545947SLisandro Dalcin pc_gamg->data = NULL; 16115b89ad90SMark F. Adams 16129d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 16135b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 16145b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 16155b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 16165b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 16175adeb434SBarry Smith mg->view = PCView_GAMG; 16185b89ad90SMark F. Adams 16199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG)); 16209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG)); 16219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG)); 16229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG)); 16239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG)); 16249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG)); 16259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG)); 16269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG)); 16279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG)); 16289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 16299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG)); 16309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG)); 16319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG)); 16329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG)); 16339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG)); 16349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG)); 16359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG)); 16369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG)); 16379d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1638d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 16390c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1640171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1641a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1642a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1643038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 164425a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 16459566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS)); 1646c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 164718c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 16489ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 16497e6512fdSJed Brown pc_gamg->use_sa_esteig = PETSC_TRUE; 165018c3aa7eSMark pc_gamg->emin = 0; 165118c3aa7eSMark pc_gamg->emax = 0; 165218c3aa7eSMark 1653c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 16549d5b6da9SMark F. Adams 1655bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 16569566063dSJacob Faibussowitsch PetscCall(PCGAMGSetType(pc,PCGAMGAGG)); 16575b89ad90SMark F. Adams PetscFunctionReturn(0); 16585b89ad90SMark F. Adams } 16593e3471ccSMark Adams 16603e3471ccSMark Adams /*@C 16613e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 16628a690491SBarry Smith from PCInitializePackage(). 16633e3471ccSMark Adams 16643e3471ccSMark Adams Level: developer 16653e3471ccSMark Adams 16663e3471ccSMark Adams .seealso: PetscInitialize() 16673e3471ccSMark Adams @*/ 16683e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 16693e3471ccSMark Adams { 16704555aa8cSStefano Zampini PetscInt l; 16713e3471ccSMark Adams 16723e3471ccSMark Adams PetscFunctionBegin; 16733e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 16743e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 16759566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO)); 16769566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG)); 16779566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical)); 16789566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1679c1c463dbSMark Adams 1680c1c463dbSMark Adams /* general events */ 1681849bee69SMark Adams PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1682849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1683849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGFilter", PC_CLASSID, &petsc_gamg_setup_events[GAMG_FILTER])); 1684849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1685849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Square G", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SQUARE])); 1686849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1687849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1688849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1689849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1690849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1691849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1692849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1693849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG PtAP",PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1694849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Reduce",PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1695849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1696849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1697849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1698849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 16994555aa8cSStefano Zampini for (l=0;l<PETSC_MG_MAXLEVELS;l++) { 17004555aa8cSStefano Zampini char ename[32]; 17015b89ad90SMark F. Adams 17029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l)); 17039566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 17049566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l)); 17059566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 17069566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l)); 17079566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 17084555aa8cSStefano Zampini } 17094555aa8cSStefano Zampini #if defined(GAMG_STAGES) 1710849bee69SMark Adams { /* create timer stages */ 17115b89ad90SMark F. Adams char str[32]; 1712849bee69SMark Adams sprintf(str,"GAMG Level %d",0); 17139566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 17145b89ad90SMark F. Adams } 17155b89ad90SMark F. Adams #endif 17163e3471ccSMark Adams PetscFunctionReturn(0); 17173e3471ccSMark Adams } 17183e3471ccSMark Adams 17193e3471ccSMark Adams /*@C 17201c1aac46SBarry Smith PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 17211c1aac46SBarry Smith called from PetscFinalize() automatically. 17223e3471ccSMark Adams 17233e3471ccSMark Adams Level: developer 17243e3471ccSMark Adams 17253e3471ccSMark Adams .seealso: PetscFinalize() 17263e3471ccSMark Adams @*/ 17273e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 17283e3471ccSMark Adams { 17293e3471ccSMark Adams PetscFunctionBegin; 17303e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 17319566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&GAMGList)); 17323e3471ccSMark Adams PetscFunctionReturn(0); 17333e3471ccSMark Adams } 1734a36cf38bSToby Isaac 1735a36cf38bSToby Isaac /*@C 1736a36cf38bSToby Isaac PCGAMGRegister - Register a PCGAMG implementation. 1737a36cf38bSToby Isaac 1738a36cf38bSToby Isaac Input Parameters: 1739a36cf38bSToby Isaac + type - string that will be used as the name of the GAMG type. 1740a36cf38bSToby Isaac - create - function for creating the gamg context. 1741a36cf38bSToby Isaac 1742a36cf38bSToby Isaac Level: advanced 1743a36cf38bSToby Isaac 17441c1aac46SBarry Smith .seealso: PCGAMGType, PCGAMG, PCGAMGSetType() 1745a36cf38bSToby Isaac @*/ 1746a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1747a36cf38bSToby Isaac { 1748a36cf38bSToby Isaac PetscFunctionBegin; 17499566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 17509566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList,type,create)); 1751a36cf38bSToby Isaac PetscFunctionReturn(0); 1752a36cf38bSToby Isaac } 1753