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*/ 65b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */ 718c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/ 8f96513f1SMatthew G Knepley 9c9567895SMark #if defined(PETSC_HAVE_CUDA) 10c9567895SMark #include <cuda_runtime.h> 11c9567895SMark #endif 12c9567895SMark 13c9567895SMark #if defined(PETSC_HAVE_HIP) 14c9567895SMark #include <hip/hip_runtime.h> 15c9567895SMark #endif 16c9567895SMark 170cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 18*4555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3]; 19fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG; 20fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO; 210cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 220cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 230cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 240cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 25fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG; 260cbbd2e1SMark F. Adams 27b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 28*4555aa8cSStefano Zampini #if defined(GAMG_STAGES) 2918c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS]; 30b4fbaa2aSMark F. Adams #endif 31f96513f1SMatthew G Knepley 320a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL; 333e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 349d5b6da9SMark F. Adams 35d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 37d3d6bff4SMark F. Adams { 3818c3aa7eSMark PetscErrorCode ierr, level; 39d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 40d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 41d3d6bff4SMark F. Adams 42d3d6bff4SMark F. Adams PetscFunctionBegin; 4322a233eaSStefano Zampini ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 441c1aac46SBarry Smith pc_gamg->data_sz = 0; 45878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 4618c3aa7eSMark for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) { 4718c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 4818c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 4918c3aa7eSMark } 5018c3aa7eSMark pc_gamg->emin = 0; 5118c3aa7eSMark pc_gamg->emax = 0; 52a2f3521dSMark F. Adams PetscFunctionReturn(0); 53a2f3521dSMark F. Adams } 54a2f3521dSMark F. Adams 555b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 565b89ad90SMark F. Adams /* 57c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 58a147abb0SMark F. Adams of active processors. 595b89ad90SMark F. Adams 605b89ad90SMark F. Adams Input Parameter: 61a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 62a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 639d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 64c5bfad50SMark F. Adams . cr_bs - coarse block size 653530afc2SMark F. Adams In/Output Parameter: 66a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 67afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 6811e60469SMark F. Adams Output Parameter: 693530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 705b89ad90SMark F. Adams */ 715cb416c2SMark F. Adams 72171cca9aSMark 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) 735b89ad90SMark F. Adams { 74a2f3521dSMark F. Adams PetscErrorCode ierr; 759d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 76486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 77a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 783b4367a7SBarry Smith MPI_Comm comm; 79c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 803ae0bb68SMark Adams PetscInt ncrs_eq,ncrs,f_bs; 815b89ad90SMark F. Adams 825b89ad90SMark F. Adams PetscFunctionBegin; 833b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 84ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 85ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 86c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 87*4555aa8cSStefano Zampini ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);CHKERRQ(ierr); 889d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 89*4555aa8cSStefano Zampini ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);CHKERRQ(ierr); 90038e3b61SMark F. Adams 91ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 92ce7c7f2fSMark Adams 933ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 940298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 953ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 963ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 9773911c69SBarry Smith } else { 983ae0bb68SMark Adams PetscInt bs; 993ae0bb68SMark Adams ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 1003ae0bb68SMark Adams ncrs = ncrs_eq/bs; 1013ae0bb68SMark Adams } 102c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 103c9567895SMark 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. */ 104c9567895SMark #if defined(PETSC_HAVE_CUDA) 105c9567895SMark PetscShmComm pshmcomm; 106c9567895SMark PetscMPIInt locrank; 107c9567895SMark MPI_Comm loccomm; 108c9567895SMark PetscInt s_nnodes,r_nnodes, new_new_size; 109c9567895SMark cudaError_t cerr; 110c9567895SMark int devCount; 111c9567895SMark ierr = PetscShmCommGet(comm,&pshmcomm);CHKERRQ(ierr); 112c9567895SMark ierr = PetscShmCommGetMpiShmComm(pshmcomm,&loccomm);CHKERRQ(ierr); 113c9567895SMark ierr = MPI_Comm_rank(loccomm, &locrank);CHKERRQ(ierr); 114c9567895SMark s_nnodes = !locrank; 115c9567895SMark ierr = MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 116c9567895SMark if (size%r_nnodes) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes); 117c9567895SMark devCount = 0; 118c9567895SMark cerr = cudaGetDeviceCount(&devCount); 119c9567895SMark cudaGetLastError(); /* Reset the last error */ 120c9567895SMark if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 121c9567895SMark new_new_size = r_nnodes * devCount; 122c9567895SMark new_size = new_new_size; 123adf5291fSStefano Zampini ierr = PetscInfo5(pc,"Fine grid with Cuda. %D nodes. Change new active set size %d --> %d (devCount=%d #nodes=%D)\n",r_nnodes,nactive,new_size,devCount,r_nnodes);CHKERRQ(ierr); 124c9567895SMark } else { 125c9567895SMark ierr = PetscInfo(pc,"With Cuda but no device. Use heuristics.");CHKERRQ(ierr); 126c9567895SMark goto HEURISTIC; 127c9567895SMark } 128c9567895SMark #else 129c9567895SMark SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here"); 130c9567895SMark #endif 131c9567895SMark } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 132c9567895SMark if (nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level]) SETERRQ2(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]); 133c9567895SMark new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level]; 134adf5291fSStefano Zampini ierr = PetscInfo3(pc,"Manually setting reduction to %d active processes (%d/%D)\n",new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);CHKERRQ(ierr); 135c9567895SMark } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 136c9567895SMark new_size = 1; 137adf5291fSStefano Zampini ierr = PetscInfo1(pc,"Force coarsest grid reduction to %d active processes\n",new_size);CHKERRQ(ierr); 138c9567895SMark } else { 139472110cdSMark F. Adams PetscInt ncrs_eq_glob; 140c9567895SMark #if defined(PETSC_HAVE_CUDA) 141c9567895SMark HEURISTIC: 142c9567895SMark #endif 1430298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 144a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1457f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 146c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 147adf5291fSStefano Zampini ierr = PetscInfo2(pc,"Coarse grid reduction from %d to %d active processes\n",nactive,new_size);CHKERRQ(ierr); 148a2f3521dSMark F. Adams } 1492e3501ffSMark Adams if (new_size==nactive) { 150ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 151ce7c7f2fSMark Adams if (new_size < size) { 152ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 153adf5291fSStefano Zampini ierr = PetscInfo1(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr); 154ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 155b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr); 156b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr); 157ce7c7f2fSMark Adams } 158ce7c7f2fSMark Adams } 159ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1602e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 161192c0e8bSMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1; 162885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 16371959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 16471959b99SBarry Smith if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 165ce7c7f2fSMark Adams /* get new_size and rfactor */ 166ce7c7f2fSMark Adams if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 167ce7c7f2fSMark Adams /* find factor */ 168ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 169ce7c7f2fSMark Adams else { 170ce7c7f2fSMark Adams PetscReal best_fact = 0.; 171ce7c7f2fSMark Adams jj = -1; 172ce7c7f2fSMark Adams for (kk = 1 ; kk <= size ; kk++) { 173ce7c7f2fSMark Adams if (!(size%kk)) { /* a candidate */ 174ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 175ce7c7f2fSMark Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 176ce7c7f2fSMark Adams if (fact > best_fact) { 177ce7c7f2fSMark Adams best_fact = fact; jj = kk; 178ce7c7f2fSMark Adams } 179ce7c7f2fSMark Adams } 180ce7c7f2fSMark Adams } 181ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 182ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 183ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 184ce7c7f2fSMark Adams else expand_factor = rfactor; 185ce7c7f2fSMark Adams } 186ce7c7f2fSMark Adams new_size = size/rfactor; /* make new size one that is factor */ 1874cdfd227SMark if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1884cdfd227SMark *a_Amat_crs = Cmat; 189adf5291fSStefano Zampini ierr = PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr); 190ce7c7f2fSMark Adams PetscFunctionReturn(0); 191ce7c7f2fSMark Adams } 192ce7c7f2fSMark Adams } 1934cdfd227SMark ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 194a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 195785e854fSJed Brown ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 1962e3501ffSMark Adams if (pc_gamg->repart) { 197a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1985a9b9e01SMark F. Adams Mat adj; 199adf5291fSStefano Zampini ierr = PetscInfo4(pc,"Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");CHKERRQ(ierr); 200a2f3521dSMark F. Adams /* get 'adj' */ 201c5bfad50SMark F. Adams if (cr_bs == 1) { 202038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 203806fa848SBarry Smith } else { 204a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 205eb07cef2SMark F. Adams Mat tMat; 206a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 207b4fbaa2aSMark F. Adams const PetscScalar *vals; 208b4fbaa2aSMark F. Adams const PetscInt *idx; 209a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 21039d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 211d9558ea9SBarry Smith MatType mtype; 212b4fbaa2aSMark F. Adams 213e632b94dSBarry Smith ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr); 214a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 215a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 216c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 2170a545947SLisandro Dalcin ierr = MatGetRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 218c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 219c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 2200a545947SLisandro Dalcin ierr = MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 2213ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 2223ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 22358471d46SMark F. Adams } 2246876a03eSMark F. Adams 225d9558ea9SBarry Smith ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr); 2263b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 2273ae0bb68SMark Adams ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 228d9558ea9SBarry Smith ierr = MatSetType(tMat,mtype);CHKERRQ(ierr); 229a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 230a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 231e632b94dSBarry Smith ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 232eb07cef2SMark F. Adams 233a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 234c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 23522063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 236eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 237c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 238eb07cef2SMark F. Adams PetscScalar v = 1.0; 239eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 240eb07cef2SMark F. Adams } 24122063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 242eb07cef2SMark F. Adams } 243eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 244eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 245eb07cef2SMark F. Adams 246b4fbaa2aSMark F. Adams if (llev++ == -1) { 247b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2488caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 2493b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 250b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 2513bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 252b4fbaa2aSMark F. Adams } 253eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 254eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 255a2f3521dSMark F. Adams } /* create 'adj' */ 256f150b916SMark F. Adams 257a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2585a9b9e01SMark F. Adams char prefix[256]; 2595a9b9e01SMark F. Adams const char *pcpre; 260b4fbaa2aSMark F. Adams const PetscInt *is_idx; 261b4fbaa2aSMark F. Adams MatPartitioning mpart; 262a4b7d37bSMark F. Adams IS proc_is; 2632f03bc48SMark F. Adams 2643b4367a7SBarry Smith ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 2655ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 2669d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 2678caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 26859a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 26911e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 270c5df96a5SBarry Smith ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 271a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 27211e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 2735a9b9e01SMark F. Adams 2745ef31b24SMark F. Adams /* collect IS info */ 275785e854fSJed Brown ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 276a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 277a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 278c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 279ce7c7f2fSMark Adams newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ 280eb07cef2SMark F. Adams } 2815ef31b24SMark F. Adams } 282a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 283a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2845ef31b24SMark F. Adams } 2855ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2865a9b9e01SMark F. Adams 2873b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2888263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 28931cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 290ce7c7f2fSMark Adams PetscInt targetPE; 2914cdfd227SMark if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen"); 292302440fdSBarry Smith ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr); 293ce7c7f2fSMark Adams targetPE = (rank/rfactor)*expand_factor; 2943b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 295a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 296e33ef3b1SMark F. Adams 29711e60469SMark F. Adams /* 298a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 29911e60469SMark F. Adams */ 300a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 3017700e67bSMark Adams is_eq_num_prim = is_eq_num; 30211e60469SMark F. Adams /* 303a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 30411e60469SMark F. Adams */ 305c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 306c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 307a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 308ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new/cr_bs; 309a2f3521dSMark F. Adams 310a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 311885364a3SMark 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 */ 312885364a3SMark Adams { 313885364a3SMark Adams Vec src_crd, dest_crd; 314885364a3SMark 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; 315885364a3SMark Adams VecScatter vecscat; 316885364a3SMark Adams PetscScalar *array; 317885364a3SMark Adams IS isscat; 318a2f3521dSMark F. Adams /* move data (for primal equations only) */ 31922063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3203b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 3213ae0bb68SMark Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 322c0dedaeaSBarry Smith ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 32311e60469SMark F. Adams /* 3249d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 325c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 32611e60469SMark F. Adams */ 327854ce69bSBarry Smith ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr); 328a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3293ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 330c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 331a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 33211e60469SMark F. Adams } 333a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3343ae0bb68SMark Adams ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 33592a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 33611e60469SMark F. Adams /* 33711e60469SMark F. Adams Create a vector to contain the original vertex information for each element 33811e60469SMark F. Adams */ 3393ae0bb68SMark Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 3409d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3413ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3423ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3439d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 344a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 345c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 346676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 347d3d6bff4SMark F. Adams } 348038e3b61SMark F. Adams } 349eb07cef2SMark F. Adams } 350eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 351eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 35211e60469SMark F. Adams /* 35311e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 35411e60469SMark F. Adams to the correct processor 35511e60469SMark F. Adams */ 3569448b7f1SJunchao Zhang ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 35711e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 35811e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35911e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 36011e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 36111e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 36211e60469SMark F. Adams /* 36311e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 36411e60469SMark F. Adams */ 365c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 366578f55a3SPeter Brune ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 3672fa5cd67SKarl Rupp 3683ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3693ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3702fa5cd67SKarl Rupp 37111e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 3729d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3733ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3749d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 375a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 376c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 377d3d6bff4SMark F. Adams } 378038e3b61SMark F. Adams } 379038e3b61SMark F. Adams } 38011e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 38111e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 382885364a3SMark Adams } 383a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 3840cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 38511e60469SMark F. Adams /* 3867dae84e0SHong Zhang Invert for MatCreateSubMatrix 38711e60469SMark F. Adams */ 388a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 389a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 390c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 391a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 392a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 393a2f3521dSMark F. Adams } 3943cb8563fSToby Isaac if (Pcolumnperm) { 3953cb8563fSToby Isaac ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr); 3963cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3973cb8563fSToby Isaac } 398a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 3990cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 4000cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 401a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 402a2f3521dSMark F. Adams { 403a2f3521dSMark F. Adams Mat mat; 4047dae84e0SHong Zhang ierr = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 405a2f3521dSMark F. Adams *a_Amat_crs = mat; 406a2f3521dSMark F. Adams } 407038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 408a2f3521dSMark F. Adams 4090cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 41011e60469SMark F. Adams /* prolongator */ 41111e60469SMark F. Adams { 41211e60469SMark F. Adams IS findices; 413a2f3521dSMark F. Adams PetscInt Istart,Iend; 414a2f3521dSMark F. Adams Mat Pnew; 41562294041SBarry Smith 416a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 4170cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 4183b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 419c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 4207dae84e0SHong Zhang ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 42111e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 422e589036eSStefano Zampini #if defined(PETSC_HAVE_CUDA) 423e589036eSStefano Zampini ierr = MatAIJCUSPARSESetGenerateTranspose(Pnew,PETSC_TRUE);CHKERRQ(ierr); 424e589036eSStefano Zampini #endif 425c5bfad50SMark F. Adams 4260cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 4273530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 428a2f3521dSMark F. Adams 429a2f3521dSMark F. Adams /* output - repartitioned */ 430a2f3521dSMark F. Adams *a_P_inout = Pnew; 431e33ef3b1SMark F. Adams } 432a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4335b89ad90SMark F. Adams 434c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 435ce7c7f2fSMark Adams 436ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 437ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 438ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4398bca76a6SMark Adams static PetscInt llev = 2; 44039d09545SMark Adams ierr = PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);CHKERRQ(ierr); 441ce7c7f2fSMark Adams #endif 442b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr); 443b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr); 444adf5291fSStefano Zampini if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 445ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 446ce7c7f2fSMark Adams PetscMPIInt size; 447ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 448ce7c7f2fSMark Adams if (size > 1) { 449ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data; 450b470e4b4SRichard Tran Mills ierr = VecBindToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr); 451b470e4b4SRichard Tran Mills ierr = VecBindToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr); 452ce7c7f2fSMark Adams } 453ce7c7f2fSMark Adams } 454ce7c7f2fSMark Adams } 4554cdfd227SMark ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 456a2f3521dSMark F. Adams } 4575b89ad90SMark F. Adams PetscFunctionReturn(0); 4585b89ad90SMark F. Adams } 4595b89ad90SMark F. Adams 4604b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2) 4614b1575e2SStefano Zampini { 4624b1575e2SStefano Zampini PetscErrorCode ierr; 4634b1575e2SStefano Zampini const char *prefix; 4644b1575e2SStefano Zampini char addp[32]; 4654b1575e2SStefano Zampini PC_MG *mg = (PC_MG*)a_pc->data; 4664b1575e2SStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4674b1575e2SStefano Zampini 4684b1575e2SStefano Zampini PetscFunctionBegin; 4694b1575e2SStefano Zampini ierr = PCGetOptionsPrefix(a_pc,&prefix);CHKERRQ(ierr); 4704b1575e2SStefano Zampini ierr = PetscInfo1(a_pc,"Square Graph on level %D\n",pc_gamg->current_level+1);CHKERRQ(ierr); 4714b1575e2SStefano Zampini ierr = MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);CHKERRQ(ierr); 4724b1575e2SStefano Zampini ierr = MatSetOptionsPrefix(*Gmat2,prefix);CHKERRQ(ierr); 4734b1575e2SStefano Zampini ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);CHKERRQ(ierr); 4744b1575e2SStefano Zampini ierr = MatAppendOptionsPrefix(*Gmat2,addp);CHKERRQ(ierr); 475b4da3a1bSStefano Zampini if ((*Gmat2)->structurally_symmetric) { 476b4da3a1bSStefano Zampini ierr = MatProductSetType(*Gmat2,MATPRODUCT_AB);CHKERRQ(ierr); 477b4da3a1bSStefano Zampini } else { 478b4da3a1bSStefano Zampini #if defined(PETSC_HAVE_CUDA) 479b4da3a1bSStefano Zampini ierr = MatAIJCUSPARSESetGenerateTranspose(Gmat1,PETSC_TRUE);CHKERRQ(ierr); 480b4da3a1bSStefano Zampini #endif 4814b1575e2SStefano Zampini ierr = MatProductSetType(*Gmat2,MATPRODUCT_AtB);CHKERRQ(ierr); 482b4da3a1bSStefano Zampini } 4834b1575e2SStefano Zampini ierr = MatProductSetFromOptions(*Gmat2);CHKERRQ(ierr); 484*4555aa8cSStefano Zampini ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);CHKERRQ(ierr); 4854b1575e2SStefano Zampini ierr = MatProductSymbolic(*Gmat2);CHKERRQ(ierr); 486*4555aa8cSStefano Zampini ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);CHKERRQ(ierr); 487b4da3a1bSStefano Zampini ierr = MatProductClear(*Gmat2);CHKERRQ(ierr); 4884b1575e2SStefano Zampini /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 4894b1575e2SStefano Zampini (*Gmat2)->assembled = PETSC_TRUE; 4904b1575e2SStefano Zampini PetscFunctionReturn(0); 4914b1575e2SStefano Zampini } 4924b1575e2SStefano Zampini 4935b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4945b89ad90SMark F. Adams /* 4955b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4965b89ad90SMark F. Adams by setting data structures and options. 4975b89ad90SMark F. Adams 4985b89ad90SMark F. Adams Input Parameter: 4995b89ad90SMark F. Adams . pc - the preconditioner context 5005b89ad90SMark F. Adams 5015b89ad90SMark F. Adams */ 5029d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 5035b89ad90SMark F. Adams { 5045b89ad90SMark F. Adams PetscErrorCode ierr; 5059d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 5065b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 5072adcac29SMark F. Adams Mat Pmat = pc->pmat; 50818c3aa7eSMark PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS]; 5093b4367a7SBarry Smith MPI_Comm comm; 510c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 51118c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS]; 51218c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 513a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 514569f4572SMark Adams MatInfo info; 515171cca9aSMark Adams PetscBool is_last = PETSC_FALSE; 5165ef31b24SMark F. Adams 5175b89ad90SMark F. Adams PetscFunctionBegin; 5183b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 519ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 520ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 521dfd5c07aSMark F. Adams 5228abdc6daSStefano Zampini if (pc->setupcalled) { 5238abdc6daSStefano Zampini if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 524878e152fSMark F. Adams /* reset everything */ 525878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 526878e152fSMark F. Adams pc->setupcalled = 0; 527806fa848SBarry Smith } else { 52884d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 52903a628feSMark F. Adams /* just do Galerkin grids */ 53058471d46SMark F. Adams Mat B,dA,dB; 53158471d46SMark F. Adams 5329d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 533*4555aa8cSStefano Zampini PetscInt gl; 53458471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 53523ee1639SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 53658471d46SMark F. Adams /* (re)set to get dirty flag */ 53723ee1639SBarry Smith ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 53858471d46SMark F. Adams 539*4555aa8cSStefano Zampini for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) { 5408abdc6daSStefano Zampini MatReuse reuse = MAT_INITIAL_MATRIX ; 5418abdc6daSStefano Zampini 5428abdc6daSStefano Zampini /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 54323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 5448abdc6daSStefano Zampini if (B->product) { 5458abdc6daSStefano Zampini if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) { 5468abdc6daSStefano Zampini reuse = MAT_REUSE_MATRIX; 54703a628feSMark F. Adams } 5488abdc6daSStefano Zampini } 5498abdc6daSStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); } 5508abdc6daSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 5518abdc6daSStefano Zampini ierr = PetscInfo1(pc,"RAP after first solve, reuse matrix level %D\n",level);CHKERRQ(ierr); 5528abdc6daSStefano Zampini } else { 5538abdc6daSStefano Zampini ierr = PetscInfo1(pc,"RAP after first solve, new matrix level %D\n",level);CHKERRQ(ierr); 5548abdc6daSStefano Zampini } 555*4555aa8cSStefano Zampini ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);CHKERRQ(ierr); 5568abdc6daSStefano Zampini ierr = MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B);CHKERRQ(ierr); 557*4555aa8cSStefano Zampini ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);CHKERRQ(ierr); 5588abdc6daSStefano Zampini mglevels[level]->A = B; 55923ee1639SBarry Smith ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 56058471d46SMark F. Adams dB = B; 56158471d46SMark F. Adams } 5625f8cf99dSMark F. Adams } 563d5280255SMark F. Adams 564d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 56558471d46SMark F. Adams PetscFunctionReturn(0); 566eb07cef2SMark F. Adams } 567878e152fSMark F. Adams } 568f6536408SMark F. Adams 569878e152fSMark F. Adams if (!pc_gamg->data) { 570878e152fSMark F. Adams if (pc_gamg->orig_data) { 571878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 5720298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 5732fa5cd67SKarl Rupp 574878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 575878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 576878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5772fa5cd67SKarl Rupp 578785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 579878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 580806fa848SBarry Smith } else { 5811ab5ffc9SJed Brown if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 5827700e67bSMark Adams ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 5839d5b6da9SMark F. Adams } 584878e152fSMark F. Adams } 585878e152fSMark F. Adams 586878e152fSMark F. Adams /* cache original data for reuse */ 5871c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 588785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 589878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 590878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 591878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 592878e152fSMark F. Adams } 593038e3b61SMark F. Adams 594302f38e8SMark F. Adams /* get basic dims */ 595302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 596171cca9aSMark Adams ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr); 59784d3f75bSMark F. Adams 598569f4572SMark Adams ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 599569f4572SMark Adams nnz0 = info.nz_used; 600569f4572SMark Adams nnztot = info.nz_used; 601c9567895SMark ierr = PetscInfo6(pc,"level %D) N=%D, n data rows=%D, n data cols=%D, nnz/row (ave)=%d, np=%D\n",0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);CHKERRQ(ierr); 602569f4572SMark Adams 603a2f3521dSMark F. Adams /* Get A_i and R_i */ 60462294041SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 6059ab59c8bSMark Adams pc_gamg->current_level = level; 60618c3aa7eSMark if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level); 6075b89ad90SMark F. Adams level1 = level + 1; 6080cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 609*4555aa8cSStefano Zampini #if defined(GAMG_STAGES) 610a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 611b4fbaa2aSMark F. Adams #endif 612c8b0795cSMark F. Adams { /* construct prolongator */ 613725b86d8SJed Brown Mat Gmat; 6140cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 6157700e67bSMark Adams Mat Prol11; 616c8b0795cSMark F. Adams 6177700e67bSMark Adams ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 6181ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 6197700e67bSMark Adams ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 620c8b0795cSMark F. Adams 621a2f3521dSMark F. Adams /* could have failed to create new level */ 622a2f3521dSMark F. Adams if (Prol11) { 623f7df55f0SStefano Zampini const char *prefix; 624f7df55f0SStefano Zampini char addp[32]; 625f7df55f0SStefano Zampini 6269d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6270298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 628a2f3521dSMark F. Adams 629fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 630c8b0795cSMark F. Adams /* smooth */ 631fd1112cbSBarry Smith ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 632c8b0795cSMark F. Adams } 633c8b0795cSMark F. Adams 6340c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6351b18a24aSMark Adams PetscInt bs; 6361b18a24aSMark Adams ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 6370a3c815dSMark Adams ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 638ffc955d6SMark F. Adams } 639ffc955d6SMark F. Adams 640f7df55f0SStefano Zampini ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 641f7df55f0SStefano Zampini ierr = MatSetOptionsPrefix(Prol11,prefix);CHKERRQ(ierr); 642c9567895SMark ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level);CHKERRQ(ierr); 643f7df55f0SStefano Zampini ierr = MatAppendOptionsPrefix(Prol11,addp);CHKERRQ(ierr); 64491f31d3dSStefano Zampini /* Always generate the transpose with CUDA 645f7df55f0SStefano Zampini Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 646f7df55f0SStefano Zampini #if defined(PETSC_HAVE_CUDA) 647e589036eSStefano Zampini ierr = MatAIJCUSPARSESetGenerateTranspose(Prol11,PETSC_TRUE);CHKERRQ(ierr); 648f7df55f0SStefano Zampini #endif 649f7df55f0SStefano Zampini ierr = MatSetFromOptions(Prol11);CHKERRQ(ierr); 6504bde40a0SMark Adams Parr[level1] = Prol11; 6514bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 6524bde40a0SMark Adams 653a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 65441b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 655a2f3521dSMark F. Adams } /* construct prolongator scope */ 6560cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 6577f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 658171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 659569f4572SMark Adams ierr = PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr); 660*4555aa8cSStefano Zampini #if defined(GAMG_STAGES) 661a90e85d9SMark Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 662a90e85d9SMark Adams #endif 663c8b0795cSMark F. Adams break; 664c8b0795cSMark F. Adams } 6650cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 666171cca9aSMark Adams ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */ 667c9567895SMark if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?"); 668171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 6690e2909e1SMark Adams if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE; 670171cca9aSMark Adams ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr); 671a2f3521dSMark F. Adams 6720cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 673171cca9aSMark Adams ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */ 674569f4572SMark Adams ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 675569f4572SMark Adams nnztot += info.nz_used; 676c9567895SMark ierr = PetscInfo5(pc,"%D) N=%D, n data cols=%D, nnz/row (ave)=%d, %D active pes\n",level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);CHKERRQ(ierr); 677569f4572SMark Adams 678*4555aa8cSStefano Zampini #if defined(GAMG_STAGES) 679b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 680b4fbaa2aSMark F. Adams #endif 681a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6829ab59c8bSMark Adams if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) { 6839ab59c8bSMark Adams ierr = PetscInfo2(pc,"HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr); 684a90e85d9SMark Adams level++; 685a90e85d9SMark Adams break; 686a90e85d9SMark Adams } 687c8b0795cSMark F. Adams } /* levels */ 688c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 689c8b0795cSMark F. Adams 690569f4572SMark Adams ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr); 6919d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6925b89ad90SMark F. Adams fine_level = level; 6930298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 6945b89ad90SMark F. Adams 69562294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 6960ed2132dSStefano Zampini PetscErrorCode (*savesetfromoptions[PETSC_MG_MAXLEVELS])(PetscOptionItems*,KSP); 6970ed2132dSStefano Zampini 698d5280255SMark F. Adams /* set default smoothers & set operators */ 69962294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 700ffc955d6SMark F. Adams KSP smoother; 701ffc955d6SMark F. Adams PC subpc; 702a2f3521dSMark F. Adams 7039d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 704f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 705ffc955d6SMark F. Adams 706a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 707a2f3521dSMark F. Adams /* set ops */ 70823ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 709a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 710a2f3521dSMark F. Adams 711a2f3521dSMark F. Adams /* set defaults */ 7126c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 713a2f3521dSMark F. Adams 7140c3bc534SBarry Smith /* set blocks for ASM smoother that uses the 'aggregates' */ 7150c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 7162d3561bbSSatish Balay PetscInt sz; 7177a28f3e5SMark Adams IS *iss; 718a2f3521dSMark F. Adams 7192d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7207a28f3e5SMark Adams iss = ASMLocalIDsArr[level]; 7210c3bc534SBarry Smith ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr); 7220a3c815dSMark Adams ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr); 7230c3bc534SBarry Smith ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr); 7247f66b68fSMark Adams if (!sz) { 725ffc955d6SMark F. Adams IS is; 7260a3c815dSMark Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 7277a28f3e5SMark Adams ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr); 728a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 729806fa848SBarry Smith } else { 730a94c3b12SMark F. Adams PetscInt kk; 7317a28f3e5SMark Adams ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr); 732a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 7337a28f3e5SMark Adams ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr); 734a94c3b12SMark F. Adams } 7357a28f3e5SMark Adams ierr = PetscFree(iss);CHKERRQ(ierr); 736ffc955d6SMark F. Adams } 7370298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 738ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 739806fa848SBarry Smith } else { 740890ffe84SMark Adams ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 741ffc955d6SMark F. Adams } 742d5280255SMark F. Adams } 743d5280255SMark F. Adams { 744d5280255SMark F. Adams /* coarse grid */ 745d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 746d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 7470ed2132dSStefano Zampini 748d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 74923ee1639SBarry Smith ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 750cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 751d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 752d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 753d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 754d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 75571959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 75671959b99SBarry Smith if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 757d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 758d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 7599dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 7602fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 76108e36f19SMark Adams ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr); 7625b42dca8SJed Brown /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 7635b42dca8SJed Brown * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 7645b42dca8SJed Brown * view every subdomain as though they were different. */ 7655b42dca8SJed Brown ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 766d5280255SMark F. Adams } 767cf8ae1d3SMark Adams } 768d5280255SMark F. Adams 769d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 770d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 771e55864a3SBarry Smith ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); 772d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 77369aca0b8SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr); 774d5280255SMark F. Adams 77518c3aa7eSMark /* setup cheby eigen estimates from SA */ 7760ed2132dSStefano Zampini if (pc_gamg->use_sa_esteig==1) { 77718c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 77818c3aa7eSMark KSP smoother; 77918c3aa7eSMark PetscBool ischeb; 7800ed2132dSStefano Zampini 7810ed2132dSStefano Zampini savesetfromoptions[level] = NULL; 78218c3aa7eSMark ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 78318c3aa7eSMark ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr); 78418c3aa7eSMark if (ischeb) { 78518c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data; 7860ed2132dSStefano Zampini 7870ed2132dSStefano Zampini ierr = KSPSetFromOptions(smoother);CHKERRQ(ierr); /* let command line emax override using SA's eigenvalues */ 7880ed2132dSStefano Zampini if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) { 78918c3aa7eSMark PC subpc; 79018c3aa7eSMark PetscBool isjac; 79118c3aa7eSMark ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 79218c3aa7eSMark ierr = PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);CHKERRQ(ierr); 7930ed2132dSStefano Zampini if (isjac && pc_gamg->use_sa_esteig==1) { 79418c3aa7eSMark PetscReal emax,emin; 7950ed2132dSStefano Zampini 79618c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 79718c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 79818c3aa7eSMark ierr = PetscInfo4(pc,"PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",level,Aarr[level]->rmap->N,(double)emax,(double)emin);CHKERRQ(ierr); 79918c3aa7eSMark cheb->emin_computed = emin; 80018c3aa7eSMark cheb->emax_computed = emax; 80118c3aa7eSMark ierr = KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);CHKERRQ(ierr); 8020ed2132dSStefano Zampini 8030ed2132dSStefano Zampini /* We have set the eigenvalues and consumed the transformation values 8040ed2132dSStefano Zampini prevent from flagging the recomputation of the eigenvalues again in PCSetUp_MG 8050ed2132dSStefano Zampini below when setfromoptions will be called again */ 8060ed2132dSStefano Zampini savesetfromoptions[level] = smoother->ops->setfromoptions; 8070ed2132dSStefano Zampini smoother->ops->setfromoptions = NULL; 80818c3aa7eSMark } 80918c3aa7eSMark } 81018c3aa7eSMark } 81118c3aa7eSMark } 8120ed2132dSStefano Zampini } 8130ed2132dSStefano Zampini 8140ed2132dSStefano Zampini ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 8150ed2132dSStefano Zampini 8160ed2132dSStefano Zampini /* restore Chebyshev smoother for next calls */ 8170ed2132dSStefano Zampini if (pc_gamg->use_sa_esteig==1) { 8180ed2132dSStefano Zampini for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 8190ed2132dSStefano Zampini if (savesetfromoptions[level]) { 8200ed2132dSStefano Zampini KSP smoother; 8210ed2132dSStefano Zampini ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 8220ed2132dSStefano Zampini smoother->ops->setfromoptions = savesetfromoptions[level]; 8230ed2132dSStefano Zampini } 8240ed2132dSStefano Zampini } 8250ed2132dSStefano Zampini } 82618c3aa7eSMark 827d5280255SMark F. Adams /* clean up */ 828d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 829587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 830587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 8315b89ad90SMark F. Adams } 832806fa848SBarry Smith } else { 8335f8cf99dSMark F. Adams KSP smoother; 8340ed2132dSStefano Zampini 835302440fdSBarry Smith ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr); 8369d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 83723ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 8385f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 8399d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 8405f8cf99dSMark F. Adams } 8415b89ad90SMark F. Adams PetscFunctionReturn(0); 8425b89ad90SMark F. Adams } 8435b89ad90SMark F. Adams 844eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 8455b89ad90SMark F. Adams /* 8465b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 8475b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 8485b89ad90SMark F. Adams 8495b89ad90SMark F. Adams Input Parameter: 8505b89ad90SMark F. Adams . pc - the preconditioner context 8515b89ad90SMark F. Adams 8525b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 8535b89ad90SMark F. Adams */ 8545b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 8555b89ad90SMark F. Adams { 8565b89ad90SMark F. Adams PetscErrorCode ierr; 8575b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 8585b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 8595b89ad90SMark F. Adams 8605b89ad90SMark F. Adams PetscFunctionBegin; 8615b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 8629b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 8639b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 8649b8ffb57SJed Brown } 8651ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 8661ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 8675b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 8685b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 8695b89ad90SMark F. Adams PetscFunctionReturn(0); 8705b89ad90SMark F. Adams } 8715b89ad90SMark F. Adams 872676e1743SMark F. Adams /*@ 873cab9ed1eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 874676e1743SMark F. Adams 8751cc46a46SBarry Smith Logically Collective on PC 876676e1743SMark F. Adams 877676e1743SMark F. Adams Input Parameters: 8781cc46a46SBarry Smith + pc - the preconditioner context 8791cc46a46SBarry Smith - n - the number of equations 880676e1743SMark F. Adams 881676e1743SMark F. Adams 882676e1743SMark F. Adams Options Database Key: 8831cc46a46SBarry Smith . -pc_gamg_process_eq_limit <limit> 884676e1743SMark F. Adams 88595452b02SPatrick Sanan Notes: 88695452b02SPatrick 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 887cab9ed1eSBarry Smith that has degrees of freedom 888cab9ed1eSBarry Smith 889676e1743SMark F. Adams Level: intermediate 890676e1743SMark F. Adams 891c9567895SMark .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors() 892676e1743SMark F. Adams @*/ 893676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 894676e1743SMark F. Adams { 895676e1743SMark F. Adams PetscErrorCode ierr; 896676e1743SMark F. Adams 897676e1743SMark F. Adams PetscFunctionBegin; 898676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 899676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 900676e1743SMark F. Adams PetscFunctionReturn(0); 901676e1743SMark F. Adams } 902676e1743SMark F. Adams 9031e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 904676e1743SMark F. Adams { 905c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 906c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 907676e1743SMark F. Adams 908676e1743SMark F. Adams PetscFunctionBegin; 9099d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 910676e1743SMark F. Adams PetscFunctionReturn(0); 911676e1743SMark F. Adams } 912676e1743SMark F. Adams 913389730f3SMark F. Adams /*@ 914cab9ed1eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 915389730f3SMark F. Adams 916389730f3SMark F. Adams Collective on PC 917389730f3SMark F. Adams 918389730f3SMark F. Adams Input Parameters: 9191cc46a46SBarry Smith + pc - the preconditioner context 9201cc46a46SBarry Smith - n - maximum number of equations to aim for 921389730f3SMark F. Adams 922389730f3SMark F. Adams Options Database Key: 9231cc46a46SBarry Smith . -pc_gamg_coarse_eq_limit <limit> 924389730f3SMark F. Adams 92574329af1SBarry Smith Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 92674329af1SBarry Smith has less than 1000 unknowns. 92774329af1SBarry Smith 928389730f3SMark F. Adams Level: intermediate 929389730f3SMark F. Adams 930c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors() 931389730f3SMark F. Adams @*/ 932389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 933389730f3SMark F. Adams { 934389730f3SMark F. Adams PetscErrorCode ierr; 935389730f3SMark F. Adams 936389730f3SMark F. Adams PetscFunctionBegin; 937389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 938389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 939389730f3SMark F. Adams PetscFunctionReturn(0); 940389730f3SMark F. Adams } 941389730f3SMark F. Adams 9421e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 943389730f3SMark F. Adams { 944389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 945389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 946389730f3SMark F. Adams 947389730f3SMark F. Adams PetscFunctionBegin; 9489d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 949389730f3SMark F. Adams PetscFunctionReturn(0); 950389730f3SMark F. Adams } 951389730f3SMark F. Adams 952676e1743SMark F. Adams /*@ 953cab9ed1eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 954676e1743SMark F. Adams 955676e1743SMark F. Adams Collective on PC 956676e1743SMark F. Adams 957676e1743SMark F. Adams Input Parameters: 9581cc46a46SBarry Smith + pc - the preconditioner context 9591cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 960676e1743SMark F. Adams 961676e1743SMark F. Adams Options Database Key: 9621cc46a46SBarry Smith . -pc_gamg_repartition <true,false> 963676e1743SMark F. Adams 96495452b02SPatrick Sanan Notes: 96595452b02SPatrick Sanan this will generally improve the loading balancing of the work on each level 966cab9ed1eSBarry Smith 967676e1743SMark F. Adams Level: intermediate 968676e1743SMark F. Adams 969676e1743SMark F. Adams .seealso: () 970676e1743SMark F. Adams @*/ 971cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 972676e1743SMark F. Adams { 973676e1743SMark F. Adams PetscErrorCode ierr; 974676e1743SMark F. Adams 975676e1743SMark F. Adams PetscFunctionBegin; 976676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 977cab9ed1eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 978676e1743SMark F. Adams PetscFunctionReturn(0); 979676e1743SMark F. Adams } 980676e1743SMark F. Adams 981cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 982676e1743SMark F. Adams { 983c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 984c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 985676e1743SMark F. Adams 986676e1743SMark F. Adams PetscFunctionBegin; 9879d5b6da9SMark F. Adams pc_gamg->repart = n; 988676e1743SMark F. Adams PetscFunctionReturn(0); 989676e1743SMark F. Adams } 990676e1743SMark F. Adams 991dfd5c07aSMark F. Adams /*@ 99218c3aa7eSMark PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby) 99318c3aa7eSMark 99418c3aa7eSMark Collective on PC 99518c3aa7eSMark 99618c3aa7eSMark Input Parameters: 99718c3aa7eSMark + pc - the preconditioner context 99818c3aa7eSMark - n - number of its 99918c3aa7eSMark 100018c3aa7eSMark Options Database Key: 100118c3aa7eSMark . -pc_gamg_esteig_ksp_max_it <its> 100218c3aa7eSMark 100318c3aa7eSMark Notes: 100418c3aa7eSMark 100518c3aa7eSMark Level: intermediate 100618c3aa7eSMark 100718c3aa7eSMark .seealso: () 100818c3aa7eSMark @*/ 100918c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n) 101018c3aa7eSMark { 101118c3aa7eSMark PetscErrorCode ierr; 101218c3aa7eSMark 101318c3aa7eSMark PetscFunctionBegin; 101418c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 101518c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 101618c3aa7eSMark PetscFunctionReturn(0); 101718c3aa7eSMark } 101818c3aa7eSMark 101918c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n) 102018c3aa7eSMark { 102118c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 102218c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 102318c3aa7eSMark 102418c3aa7eSMark PetscFunctionBegin; 102518c3aa7eSMark pc_gamg->esteig_max_it = n; 102618c3aa7eSMark PetscFunctionReturn(0); 102718c3aa7eSMark } 102818c3aa7eSMark 102918c3aa7eSMark /*@ 103018c3aa7eSMark PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother 103118c3aa7eSMark 103218c3aa7eSMark Collective on PC 103318c3aa7eSMark 103418c3aa7eSMark Input Parameters: 103518c3aa7eSMark + pc - the preconditioner context 103618c3aa7eSMark - n - number of its 103718c3aa7eSMark 103818c3aa7eSMark Options Database Key: 103918c3aa7eSMark . -pc_gamg_use_sa_esteig <true,false> 104018c3aa7eSMark 104118c3aa7eSMark Notes: 104218c3aa7eSMark 104318c3aa7eSMark Level: intermediate 104418c3aa7eSMark 104518c3aa7eSMark .seealso: () 104618c3aa7eSMark @*/ 104718c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 104818c3aa7eSMark { 104918c3aa7eSMark PetscErrorCode ierr; 105018c3aa7eSMark 105118c3aa7eSMark PetscFunctionBegin; 105218c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 105318c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 105418c3aa7eSMark PetscFunctionReturn(0); 105518c3aa7eSMark } 105618c3aa7eSMark 10570ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 105818c3aa7eSMark { 105918c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 106018c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 106118c3aa7eSMark 106218c3aa7eSMark PetscFunctionBegin; 106318c3aa7eSMark pc_gamg->use_sa_esteig = n ? 1 : 0; 106418c3aa7eSMark PetscFunctionReturn(0); 106518c3aa7eSMark } 106618c3aa7eSMark 106718c3aa7eSMark /*@C 106818c3aa7eSMark PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby) 106918c3aa7eSMark 107018c3aa7eSMark Collective on PC 107118c3aa7eSMark 107218c3aa7eSMark Input Parameters: 107318c3aa7eSMark + pc - the preconditioner context 107418c3aa7eSMark - t - ksp type 107518c3aa7eSMark 107618c3aa7eSMark Options Database Key: 107718c3aa7eSMark . -pc_gamg_esteig_ksp_type <type> 107818c3aa7eSMark 107918c3aa7eSMark Notes: 108018c3aa7eSMark 108118c3aa7eSMark Level: intermediate 108218c3aa7eSMark 108318c3aa7eSMark .seealso: () 108418c3aa7eSMark @*/ 108518c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[]) 108618c3aa7eSMark { 108718c3aa7eSMark PetscErrorCode ierr; 108818c3aa7eSMark 108918c3aa7eSMark PetscFunctionBegin; 109018c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 109118c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));CHKERRQ(ierr); 109218c3aa7eSMark PetscFunctionReturn(0); 109318c3aa7eSMark } 109418c3aa7eSMark 109518c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[]) 109618c3aa7eSMark { 109718c3aa7eSMark PetscErrorCode ierr; 109818c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 109918c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 110018c3aa7eSMark 110118c3aa7eSMark PetscFunctionBegin; 110218c3aa7eSMark ierr = PetscStrcpy(pc_gamg->esteig_type,t);CHKERRQ(ierr); 110318c3aa7eSMark PetscFunctionReturn(0); 110418c3aa7eSMark } 110518c3aa7eSMark 110618c3aa7eSMark /*@ 110718c3aa7eSMark PCGAMGSetEigenvalues - Set eigenvalues 110818c3aa7eSMark 110918c3aa7eSMark Collective on PC 111018c3aa7eSMark 111118c3aa7eSMark Input Parameters: 111218c3aa7eSMark + pc - the preconditioner context 111318c3aa7eSMark - emax - max eigenvalue 111418c3aa7eSMark . emin - min eigenvalue 111518c3aa7eSMark 111618c3aa7eSMark Options Database Key: 111718c3aa7eSMark . -pc_gamg_eigenvalues 111818c3aa7eSMark 111918c3aa7eSMark Level: intermediate 112018c3aa7eSMark 112118c3aa7eSMark .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType() 112218c3aa7eSMark @*/ 112318c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin) 112418c3aa7eSMark { 112518c3aa7eSMark PetscErrorCode ierr; 112618c3aa7eSMark 112718c3aa7eSMark PetscFunctionBegin; 112818c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 112918c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr); 113018c3aa7eSMark PetscFunctionReturn(0); 113118c3aa7eSMark } 113218c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin) 113318c3aa7eSMark { 113418c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 113518c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 113618c3aa7eSMark 113718c3aa7eSMark PetscFunctionBegin; 113818c3aa7eSMark if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin); 113918c3aa7eSMark if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin); 114018c3aa7eSMark pc_gamg->emax = emax; 114118c3aa7eSMark pc_gamg->emin = emin; 114218c3aa7eSMark 114318c3aa7eSMark PetscFunctionReturn(0); 114418c3aa7eSMark } 114518c3aa7eSMark 114618c3aa7eSMark /*@ 1147cab9ed1eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1148dfd5c07aSMark F. Adams 1149dfd5c07aSMark F. Adams Collective on PC 1150dfd5c07aSMark F. Adams 1151dfd5c07aSMark F. Adams Input Parameters: 11521cc46a46SBarry Smith + pc - the preconditioner context 11531cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 1154dfd5c07aSMark F. Adams 1155dfd5c07aSMark F. Adams Options Database Key: 11561cc46a46SBarry Smith . -pc_gamg_reuse_interpolation <true,false> 1157dfd5c07aSMark F. Adams 1158dfd5c07aSMark F. Adams Level: intermediate 1159dfd5c07aSMark F. Adams 116095452b02SPatrick Sanan Notes: 116195452b02SPatrick Sanan this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1162cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1163cab9ed1eSBarry Smith 1164dfd5c07aSMark F. Adams .seealso: () 1165dfd5c07aSMark F. Adams @*/ 11661cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1167dfd5c07aSMark F. Adams { 1168dfd5c07aSMark F. Adams PetscErrorCode ierr; 1169dfd5c07aSMark F. Adams 1170dfd5c07aSMark F. Adams PetscFunctionBegin; 1171dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11721cc46a46SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1173dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1174dfd5c07aSMark F. Adams } 1175dfd5c07aSMark F. Adams 11761cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1177dfd5c07aSMark F. Adams { 1178dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1179dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1180dfd5c07aSMark F. Adams 1181dfd5c07aSMark F. Adams PetscFunctionBegin; 1182dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1183dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1184dfd5c07aSMark F. Adams } 1185dfd5c07aSMark F. Adams 1186ffc955d6SMark F. Adams /*@ 1187cab9ed1eSBarry 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. 1188ffc955d6SMark F. Adams 1189ffc955d6SMark F. Adams Collective on PC 1190ffc955d6SMark F. Adams 1191ffc955d6SMark F. Adams Input Parameters: 1192cab9ed1eSBarry Smith + pc - the preconditioner context 1193cab9ed1eSBarry Smith - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not 1194ffc955d6SMark F. Adams 1195ffc955d6SMark F. Adams Options Database Key: 1196cab9ed1eSBarry Smith . -pc_gamg_asm_use_agg 1197ffc955d6SMark F. Adams 1198ffc955d6SMark F. Adams Level: intermediate 1199ffc955d6SMark F. Adams 1200ffc955d6SMark F. Adams .seealso: () 1201ffc955d6SMark F. Adams @*/ 1202cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1203ffc955d6SMark F. Adams { 1204ffc955d6SMark F. Adams PetscErrorCode ierr; 1205ffc955d6SMark F. Adams 1206ffc955d6SMark F. Adams PetscFunctionBegin; 1207ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1208cab9ed1eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1209ffc955d6SMark F. Adams PetscFunctionReturn(0); 1210ffc955d6SMark F. Adams } 1211ffc955d6SMark F. Adams 1212cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1213ffc955d6SMark F. Adams { 1214ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1215ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1216ffc955d6SMark F. Adams 1217ffc955d6SMark F. Adams PetscFunctionBegin; 1218cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 1219ffc955d6SMark F. Adams PetscFunctionReturn(0); 1220ffc955d6SMark F. Adams } 1221ffc955d6SMark F. Adams 1222171cca9aSMark Adams /*@ 1223cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1224171cca9aSMark Adams 1225171cca9aSMark Adams Collective on PC 1226171cca9aSMark Adams 1227171cca9aSMark Adams Input Parameters: 1228171cca9aSMark Adams + pc - the preconditioner context 1229cf8ae1d3SMark Adams - flg - PETSC_TRUE to not force coarse grid onto one processor 1230171cca9aSMark Adams 1231171cca9aSMark Adams Options Database Key: 1232cf8ae1d3SMark Adams . -pc_gamg_use_parallel_coarse_grid_solver 1233171cca9aSMark Adams 1234171cca9aSMark Adams Level: intermediate 1235171cca9aSMark Adams 123639d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids() 1237171cca9aSMark Adams @*/ 1238171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1239171cca9aSMark Adams { 1240171cca9aSMark Adams PetscErrorCode ierr; 1241171cca9aSMark Adams 1242171cca9aSMark Adams PetscFunctionBegin; 1243171cca9aSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1244171cca9aSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1245171cca9aSMark Adams PetscFunctionReturn(0); 1246171cca9aSMark Adams } 1247171cca9aSMark Adams 1248171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1249171cca9aSMark Adams { 1250171cca9aSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1251171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1252171cca9aSMark Adams 1253171cca9aSMark Adams PetscFunctionBegin; 1254171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 1255ffc955d6SMark F. Adams PetscFunctionReturn(0); 1256ffc955d6SMark F. Adams } 1257ffc955d6SMark F. Adams 12584ef23d27SMark F. Adams /*@ 1259ce7c7f2fSMark Adams PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1260ce7c7f2fSMark Adams 1261ce7c7f2fSMark Adams Collective on PC 1262ce7c7f2fSMark Adams 1263ce7c7f2fSMark Adams Input Parameters: 1264ce7c7f2fSMark Adams + pc - the preconditioner context 1265ce7c7f2fSMark Adams - flg - PETSC_TRUE to pin coarse grids to CPU 1266ce7c7f2fSMark Adams 1267ce7c7f2fSMark Adams Options Database Key: 1268ce7c7f2fSMark Adams . -pc_gamg_cpu_pin_coarse_grids 1269ce7c7f2fSMark Adams 1270ce7c7f2fSMark Adams Level: intermediate 1271ce7c7f2fSMark Adams 127239d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve() 1273ce7c7f2fSMark Adams @*/ 1274ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1275ce7c7f2fSMark Adams { 1276ce7c7f2fSMark Adams PetscErrorCode ierr; 1277ce7c7f2fSMark Adams 1278ce7c7f2fSMark Adams PetscFunctionBegin; 1279ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1280ce7c7f2fSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1281ce7c7f2fSMark Adams PetscFunctionReturn(0); 1282ce7c7f2fSMark Adams } 1283ce7c7f2fSMark Adams 1284ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1285ce7c7f2fSMark Adams { 1286ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1287ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1288ce7c7f2fSMark Adams 1289ce7c7f2fSMark Adams PetscFunctionBegin; 1290ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 1291ce7c7f2fSMark Adams PetscFunctionReturn(0); 1292ce7c7f2fSMark Adams } 1293ce7c7f2fSMark Adams 1294ce7c7f2fSMark Adams /*@ 1295ce7c7f2fSMark Adams PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type) 1296ce7c7f2fSMark Adams 1297ce7c7f2fSMark Adams Collective on PC 1298ce7c7f2fSMark Adams 1299ce7c7f2fSMark Adams Input Parameters: 1300ce7c7f2fSMark Adams + pc - the preconditioner context 1301ce7c7f2fSMark Adams - flg - Layout type 1302ce7c7f2fSMark Adams 1303ce7c7f2fSMark Adams Options Database Key: 1304ce7c7f2fSMark Adams . -pc_gamg_coarse_grid_layout_type 1305ce7c7f2fSMark Adams 1306ce7c7f2fSMark Adams Level: intermediate 1307ce7c7f2fSMark Adams 130839d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids() 1309ce7c7f2fSMark Adams @*/ 1310ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1311ce7c7f2fSMark Adams { 1312ce7c7f2fSMark Adams PetscErrorCode ierr; 1313ce7c7f2fSMark Adams 1314ce7c7f2fSMark Adams PetscFunctionBegin; 1315ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1316ce7c7f2fSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr); 1317ce7c7f2fSMark Adams PetscFunctionReturn(0); 1318ce7c7f2fSMark Adams } 1319ce7c7f2fSMark Adams 1320ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1321ce7c7f2fSMark Adams { 1322ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1323ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1324ce7c7f2fSMark Adams 1325ce7c7f2fSMark Adams PetscFunctionBegin; 1326ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 1327ce7c7f2fSMark Adams PetscFunctionReturn(0); 1328ce7c7f2fSMark Adams } 1329ce7c7f2fSMark Adams 1330ce7c7f2fSMark Adams /*@ 13311cc46a46SBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 13324ef23d27SMark F. Adams 13334ef23d27SMark F. Adams Not collective on PC 13344ef23d27SMark F. Adams 13354ef23d27SMark F. Adams Input Parameters: 13361cc46a46SBarry Smith + pc - the preconditioner 13371cc46a46SBarry Smith - n - the maximum number of levels to use 13384ef23d27SMark F. Adams 13394ef23d27SMark F. Adams Options Database Key: 13404ef23d27SMark F. Adams . -pc_mg_levels 13414ef23d27SMark F. Adams 13424ef23d27SMark F. Adams Level: intermediate 13434ef23d27SMark F. Adams 13444ef23d27SMark F. Adams .seealso: () 13454ef23d27SMark F. Adams @*/ 13464ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 13474ef23d27SMark F. Adams { 13484ef23d27SMark F. Adams PetscErrorCode ierr; 13494ef23d27SMark F. Adams 13504ef23d27SMark F. Adams PetscFunctionBegin; 13514ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13524ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 13534ef23d27SMark F. Adams PetscFunctionReturn(0); 13544ef23d27SMark F. Adams } 13554ef23d27SMark F. Adams 13561e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 13574ef23d27SMark F. Adams { 13584ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 13594ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 13604ef23d27SMark F. Adams 13614ef23d27SMark F. Adams PetscFunctionBegin; 13629d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 13634ef23d27SMark F. Adams PetscFunctionReturn(0); 13644ef23d27SMark F. Adams } 13654ef23d27SMark F. Adams 13663542efc5SMark F. Adams /*@ 13673542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 13683542efc5SMark F. Adams 13693542efc5SMark F. Adams Not collective on PC 13703542efc5SMark F. Adams 13713542efc5SMark F. Adams Input Parameters: 13721cc46a46SBarry Smith + pc - the preconditioner context 1373c9567895SMark . 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 1374055c8bd0SJed Brown - n - number of threshold values provided in array 13753542efc5SMark F. Adams 13763542efc5SMark F. Adams Options Database Key: 13771cc46a46SBarry Smith . -pc_gamg_threshold <threshold> 13783542efc5SMark F. Adams 137995452b02SPatrick Sanan Notes: 1380af3c827dSMark 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. 1381af3c827dSMark 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. 1382cab9ed1eSBarry Smith 1383055c8bd0SJed Brown If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening. 1384055c8bd0SJed Brown In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold(). 1385055c8bd0SJed Brown If n is greater than the total number of levels, the excess entries in threshold will not be used. 1386055c8bd0SJed Brown 13873542efc5SMark F. Adams Level: intermediate 13883542efc5SMark F. Adams 1389af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph() 13903542efc5SMark F. Adams @*/ 1391c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 13923542efc5SMark F. Adams { 13933542efc5SMark F. Adams PetscErrorCode ierr; 13943542efc5SMark F. Adams 13953542efc5SMark F. Adams PetscFunctionBegin; 13963542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1397055c8bd0SJed Brown if (n) PetscValidRealPointer(v,2); 1398c1eae691SMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr); 13993542efc5SMark F. Adams PetscFunctionReturn(0); 14003542efc5SMark F. Adams } 14013542efc5SMark F. Adams 1402c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 14033542efc5SMark F. Adams { 1404c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1405c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1406c1eae691SMark Adams PetscInt i; 1407c1eae691SMark Adams PetscFunctionBegin; 1408055c8bd0SJed Brown for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1409055c8bd0SJed Brown for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale; 1410c1eae691SMark Adams PetscFunctionReturn(0); 1411c1eae691SMark Adams } 1412c1eae691SMark Adams 1413c1eae691SMark Adams /*@ 1414c9567895SMark PCGAMGSetRankReductionFactors - Set manual schedual for process reduction on coarse grids 1415c9567895SMark 1416c9567895SMark Collective on PC 1417c9567895SMark 1418c9567895SMark Input Parameters: 1419c9567895SMark + pc - the preconditioner context 1420c9567895SMark . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda 1421c9567895SMark - n - number of values provided in array 1422c9567895SMark 1423c9567895SMark Options Database Key: 1424c9567895SMark . -pc_gamg_rank_reduction_factors <factors> 1425c9567895SMark 1426c9567895SMark Level: intermediate 1427c9567895SMark 1428c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim() 1429c9567895SMark @*/ 1430c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1431c9567895SMark { 1432c9567895SMark PetscErrorCode ierr; 1433c9567895SMark 1434c9567895SMark PetscFunctionBegin; 1435c9567895SMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1436c9567895SMark if (n) PetscValidIntPointer(v,2); 1437c9567895SMark ierr = PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));CHKERRQ(ierr); 1438c9567895SMark PetscFunctionReturn(0); 1439c9567895SMark } 1440c9567895SMark 1441c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1442c9567895SMark { 1443c9567895SMark PC_MG *mg = (PC_MG*)pc->data; 1444c9567895SMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1445c9567895SMark PetscInt i; 1446c9567895SMark PetscFunctionBegin; 1447c9567895SMark for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1448c9567895SMark for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1449c9567895SMark PetscFunctionReturn(0); 1450c9567895SMark } 1451c9567895SMark 1452c9567895SMark /*@ 1453c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1454c1eae691SMark Adams 1455c1eae691SMark Adams Not collective on PC 1456c1eae691SMark Adams 1457c1eae691SMark Adams Input Parameters: 1458c1eae691SMark Adams + pc - the preconditioner context 1459c1eae691SMark Adams - scale - the threshold value reduction, ussually < 1.0 1460c1eae691SMark Adams 1461c1eae691SMark Adams Options Database Key: 1462c1eae691SMark Adams . -pc_gamg_threshold_scale <v> 1463c1eae691SMark Adams 1464055c8bd0SJed Brown Notes: 1465055c8bd0SJed Brown The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold(). 1466055c8bd0SJed Brown This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold(). 1467055c8bd0SJed Brown 1468c1eae691SMark Adams Level: advanced 1469c1eae691SMark Adams 1470055c8bd0SJed Brown .seealso: PCGAMGSetThreshold() 1471c1eae691SMark Adams @*/ 1472c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1473c1eae691SMark Adams { 1474c1eae691SMark Adams PetscErrorCode ierr; 14753542efc5SMark F. Adams 14763542efc5SMark F. Adams PetscFunctionBegin; 1477c1eae691SMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1478c1eae691SMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr); 1479c1eae691SMark Adams PetscFunctionReturn(0); 1480c1eae691SMark Adams } 1481c1eae691SMark Adams 1482c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1483c1eae691SMark Adams { 1484c1eae691SMark Adams PC_MG *mg = (PC_MG*)pc->data; 1485c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1486c1eae691SMark Adams PetscFunctionBegin; 1487c1eae691SMark Adams pc_gamg->threshold_scale = v; 14883542efc5SMark F. Adams PetscFunctionReturn(0); 14893542efc5SMark F. Adams } 14903542efc5SMark F. Adams 1491e20c40e8SBarry Smith /*@C 1492c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1493676e1743SMark F. Adams 1494676e1743SMark F. Adams Collective on PC 1495676e1743SMark F. Adams 1496676e1743SMark F. Adams Input Parameters: 1497c60c7ad4SBarry Smith + pc - the preconditioner context 1498c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1499676e1743SMark F. Adams 1500676e1743SMark F. Adams Options Database Key: 1501cab9ed1eSBarry Smith . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1502676e1743SMark F. Adams 1503676e1743SMark F. Adams Level: intermediate 1504676e1743SMark F. Adams 1505cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType 1506676e1743SMark F. Adams @*/ 150719fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1508676e1743SMark F. Adams { 1509676e1743SMark F. Adams PetscErrorCode ierr; 1510676e1743SMark F. Adams 1511676e1743SMark F. Adams PetscFunctionBegin; 1512676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1513806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1514676e1743SMark F. Adams PetscFunctionReturn(0); 1515676e1743SMark F. Adams } 1516676e1743SMark F. Adams 1517e20c40e8SBarry Smith /*@C 1518c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1519c60c7ad4SBarry Smith 1520c60c7ad4SBarry Smith Collective on PC 1521c60c7ad4SBarry Smith 1522c60c7ad4SBarry Smith Input Parameter: 1523c60c7ad4SBarry Smith . pc - the preconditioner context 1524c60c7ad4SBarry Smith 1525c60c7ad4SBarry Smith Output Parameter: 1526c60c7ad4SBarry Smith . type - the type of algorithm used 1527c60c7ad4SBarry Smith 1528c60c7ad4SBarry Smith Level: intermediate 1529c60c7ad4SBarry Smith 15301c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType 1531c60c7ad4SBarry Smith @*/ 1532c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1533c60c7ad4SBarry Smith { 1534c60c7ad4SBarry Smith PetscErrorCode ierr; 1535c60c7ad4SBarry Smith 1536c60c7ad4SBarry Smith PetscFunctionBegin; 1537c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1538c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1539c60c7ad4SBarry Smith PetscFunctionReturn(0); 1540c60c7ad4SBarry Smith } 1541c60c7ad4SBarry Smith 1542c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1543c60c7ad4SBarry Smith { 1544c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1545c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1546c60c7ad4SBarry Smith 1547c60c7ad4SBarry Smith PetscFunctionBegin; 1548c60c7ad4SBarry Smith *type = pc_gamg->type; 1549c60c7ad4SBarry Smith PetscFunctionReturn(0); 1550c60c7ad4SBarry Smith } 1551c60c7ad4SBarry Smith 15521e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1553676e1743SMark F. Adams { 15549d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 15551ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 15561ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1557676e1743SMark F. Adams 1558676e1743SMark F. Adams PetscFunctionBegin; 1559c60c7ad4SBarry Smith pc_gamg->type = type; 15601c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 15619d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 15621ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 15631ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 15641ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 1565e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 15663ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 15673ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 15683ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 15693ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 15703ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 15713ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 15723ae0bb68SMark Adams pc_gamg->data_sz = 0; 15731ab5ffc9SJed Brown } 15741ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 15751ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 15769d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1577676e1743SMark F. Adams PetscFunctionReturn(0); 1578676e1743SMark F. Adams } 1579676e1743SMark F. Adams 15809d766c59SMark Adams /* -------------------------------------------------------------------------- */ 15819d766c59SMark Adams /* 15829d766c59SMark Adams PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy 15839d766c59SMark Adams 15849d766c59SMark Adams Input Parameter: 15859d766c59SMark Adams . pc - the preconditioner context 15869d766c59SMark Adams 15879d766c59SMark Adams Output Parameter: 15889d766c59SMark Adams . gc - grid complexity = sum_i(nnz_i) / nnz_0 15899d766c59SMark Adams 15909d766c59SMark Adams Level: advanced 15919d766c59SMark Adams */ 15929d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc) 15939d766c59SMark Adams { 15949d766c59SMark Adams PetscErrorCode ierr; 15959d766c59SMark Adams PC_MG *mg = (PC_MG*)pc->data; 15969d766c59SMark Adams PC_MG_Levels **mglevels = mg->levels; 15973966268fSBarry Smith PetscInt lev; 15983966268fSBarry Smith PetscLogDouble nnz0 = 0, sgc = 0; 15999d766c59SMark Adams MatInfo info; 16003966268fSBarry Smith 16019d766c59SMark Adams PetscFunctionBegin; 1602dbf6bb8dSprj- if (!pc->setupcalled) { 1603dbf6bb8dSprj- *gc = 0; 1604dbf6bb8dSprj- PetscFunctionReturn(0); 1605dbf6bb8dSprj- } 16069d766c59SMark Adams if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels"); 16073966268fSBarry Smith for (lev=0; lev<mg->nlevels; lev++) { 160862a6e064SMark Adams Mat dB; 160962a6e064SMark Adams ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr); 161062a6e064SMark Adams ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 16113966268fSBarry Smith sgc += info.nz_used; 16129d766c59SMark Adams if (lev==mg->nlevels-1) nnz0 = info.nz_used; 16139d766c59SMark Adams } 16143966268fSBarry Smith if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0); 16153966268fSBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available"); 16169d766c59SMark Adams PetscFunctionReturn(0); 16179d766c59SMark Adams } 16189d766c59SMark Adams 16195adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 16205adeb434SBarry Smith { 1621c1eae691SMark Adams PetscErrorCode ierr,i; 16225adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 16235adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 162423b2d91dSMark Adams PetscReal gc=0; 16255adeb434SBarry Smith PetscFunctionBegin; 16265adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," GAMG specific options\n");CHKERRQ(ierr); 1627459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =");CHKERRQ(ierr); 1628b3e187dcSStefano Zampini for (i=0;i<mg->nlevels; i++) { 1629c1eae691SMark Adams ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr); 1630c1eae691SMark Adams } 1631459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1632459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr); 1633cab9ed1eSBarry Smith if (pc_gamg->use_aggs_in_asm) { 1634cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr); 1635cab9ed1eSBarry Smith } 1636171cca9aSMark Adams if (pc_gamg->use_parallel_coarse_grid_solver) { 1637171cca9aSMark Adams ierr = PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr); 1638171cca9aSMark Adams } 1639ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1640ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1641ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */ 1642ce7c7f2fSMark Adams } 1643ce7c7f2fSMark Adams #endif 1644ce7c7f2fSMark Adams /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */ 1645ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */ 1646ce7c7f2fSMark Adams /* } else { */ 1647ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */ 1648ce7c7f2fSMark Adams /* } */ 16495adeb434SBarry Smith if (pc_gamg->ops->view) { 16505adeb434SBarry Smith ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr); 16515adeb434SBarry Smith } 16529d766c59SMark Adams ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr); 16539d766c59SMark Adams ierr = PetscViewerASCIIPrintf(viewer," Complexity: grid = %g\n",gc);CHKERRQ(ierr); 16545adeb434SBarry Smith PetscFunctionReturn(0); 16555adeb434SBarry Smith } 16565adeb434SBarry Smith 16574416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 16585b89ad90SMark F. Adams { 1659676e1743SMark F. Adams PetscErrorCode ierr; 1660676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1661676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 166218c3aa7eSMark PetscBool flag,f2; 16633b4367a7SBarry Smith MPI_Comm comm; 166418c3aa7eSMark char prefix[256],tname[32]; 1665c1eae691SMark Adams PetscInt i,n; 166614a9496bSBarry Smith const char *pcpre; 16670a545947SLisandro Dalcin static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL}; 16685b89ad90SMark F. Adams PetscFunctionBegin; 16693b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1670e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr); 16711a1c1e04SBarry Smith ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1672bd94a7aaSJed Brown if (flag) { 1673bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 16741ab5ffc9SJed Brown } 167518c3aa7eSMark ierr = PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);CHKERRQ(ierr); 167618c3aa7eSMark if (flag) { 167718c3aa7eSMark ierr = PCGAMGSetEstEigKSPType(pc,tname);CHKERRQ(ierr); 167818c3aa7eSMark } 1679cab9ed1eSBarry Smith ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 168018c3aa7eSMark f2 = PETSC_TRUE; 168118c3aa7eSMark ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);CHKERRQ(ierr); 168218c3aa7eSMark if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0; 16831cc46a46SBarry Smith ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 1684a303c832SJed Brown ierr = 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);CHKERRQ(ierr); 1685cf8ae1d3SMark Adams ierr = 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);CHKERRQ(ierr); 1686ce7c7f2fSMark Adams ierr = 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);CHKERRQ(ierr); 1687a0095786SMark ierr = 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);CHKERRQ(ierr); 168894ae4db5SBarry Smith ierr = 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);CHKERRQ(ierr); 168918c3aa7eSMark ierr = PetscOptionsInt("-pc_gamg_esteig_ksp_max_it","Number of iterations of eigen estimator","PCGAMGSetEstEigKSPMaxIt",pc_gamg->esteig_max_it,&pc_gamg->esteig_max_it,NULL);CHKERRQ(ierr); 169094ae4db5SBarry Smith ierr = 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);CHKERRQ(ierr); 1691a303c832SJed Brown ierr = PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);CHKERRQ(ierr); 169218c3aa7eSMark n = PETSC_MG_MAXLEVELS; 1693c1eae691SMark Adams ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr); 169418c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1695efd3c5ceSMark Adams if (!flag) n = 1; 1696c1eae691SMark Adams i = n; 169718c3aa7eSMark do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1698c1eae691SMark Adams } 1699c9567895SMark n = PETSC_MG_MAXLEVELS; 1700c9567895SMark ierr = 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);CHKERRQ(ierr); 1701c9567895SMark if (!flag) i = 0; 1702c9567895SMark else i = n; 1703c9567895SMark do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS); 170494ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 170518c3aa7eSMark { 170618c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 170718c3aa7eSMark n = 2; 170818c3aa7eSMark ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr); 170918c3aa7eSMark if (flag) { 171018c3aa7eSMark if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 171118c3aa7eSMark ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr); 171218c3aa7eSMark } 171318c3aa7eSMark } 1714b7cbab4eSMark Adams /* set options for subtype */ 1715e55864a3SBarry Smith if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);} 171618c3aa7eSMark 171714a9496bSBarry Smith ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 171814a9496bSBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 1719676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 17205b89ad90SMark F. Adams PetscFunctionReturn(0); 17215b89ad90SMark F. Adams } 17225b89ad90SMark F. Adams 17235b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 17245b89ad90SMark F. Adams /*MC 17251cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 17265b89ad90SMark F. Adams 1727280d9858SJed Brown Options Database Keys: 1728cab9ed1eSBarry Smith + -pc_gamg_type <type> - one of agg, geo, or classical 1729cab9ed1eSBarry Smith . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1730cab9ed1eSBarry Smith . -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1731cab9ed1eSBarry 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 1732cab9ed1eSBarry 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> 1733cab9ed1eSBarry Smith equations on each process that has degrees of freedom 1734cab9ed1eSBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 17356008e27bSRichard Tran Mills . -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level 1736c1eae691SMark Adams - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1737cab9ed1eSBarry Smith 1738cab9ed1eSBarry Smith Options Database Keys for default Aggregation: 1739cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1740cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 1741cab9ed1eSBarry Smith - -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it 1742cab9ed1eSBarry Smith 1743db9745e2SBarry Smith Multigrid options: 1744db9745e2SBarry Smith + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1745db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1746db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1747cab9ed1eSBarry Smith - -pc_mg_levels <levels> - Number of levels of multigrid to use. 17485b89ad90SMark F. Adams 17491cc46a46SBarry Smith 175095452b02SPatrick Sanan Notes: 175195452b02SPatrick Sanan In order to obtain good performance for PCGAMG for vector valued problems you must 1752db9745e2SBarry Smith Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1753db9745e2SBarry Smith Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1754db9745e2SBarry Smith See the Users Manual Chapter 4 for more details 17551cc46a46SBarry Smith 17565b89ad90SMark F. Adams Level: intermediate 1757280d9858SJed Brown 17581cc46a46SBarry Smith .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 175918c3aa7eSMark PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType() 17605b89ad90SMark F. Adams M*/ 1761b2573a8aSBarry Smith 17628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 17635b89ad90SMark F. Adams { 1764c1eae691SMark Adams PetscErrorCode ierr,i; 17655b89ad90SMark F. Adams PC_GAMG *pc_gamg; 17665b89ad90SMark F. Adams PC_MG *mg; 17675b89ad90SMark F. Adams 17685b89ad90SMark F. Adams PetscFunctionBegin; 17691c1aac46SBarry Smith /* register AMG type */ 17701c1aac46SBarry Smith ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 17711c1aac46SBarry Smith 17725b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 17731c1aac46SBarry Smith ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); 17745b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 17755b89ad90SMark F. Adams 17765b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 1777b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 177869aca0b8SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr); 17795b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 17805b89ad90SMark F. Adams mg->innerctx = pc_gamg; 17815b89ad90SMark F. Adams 1782b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 17831ab5ffc9SJed Brown 17849d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 17859d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 17860a545947SLisandro Dalcin pc_gamg->data = NULL; 17875b89ad90SMark F. Adams 17889d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 17895b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 17905b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 17915b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 17925b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 17935adeb434SBarry Smith mg->view = PCView_GAMG; 17945b89ad90SMark F. Adams 179597d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 179697d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1797bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1798bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1799cab9ed1eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr); 180018c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);CHKERRQ(ierr); 180118c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);CHKERRQ(ierr); 180218c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr); 180318c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr); 18041cc46a46SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr); 1805cab9ed1eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr); 1806171cca9aSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr); 1807ce7c7f2fSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr); 1808ce7c7f2fSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr); 1809bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1810c9567895SMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG);CHKERRQ(ierr); 1811c1eae691SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr); 1812bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1813c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1814bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 18159d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1816d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 18170c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1818171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1819a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1820a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1821038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 182225a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 182318c3aa7eSMark for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.; 1824c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 182518c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 18269ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 1827d24ecf33SMark ierr = PetscStrcpy(pc_gamg->esteig_type,NULL);CHKERRQ(ierr); 182818c3aa7eSMark pc_gamg->esteig_max_it = 10; 182918c3aa7eSMark pc_gamg->use_sa_esteig = -1; 183018c3aa7eSMark pc_gamg->emin = 0; 183118c3aa7eSMark pc_gamg->emax = 0; 183218c3aa7eSMark 1833c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 18349d5b6da9SMark F. Adams 1835bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1836bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 18375b89ad90SMark F. Adams PetscFunctionReturn(0); 18385b89ad90SMark F. Adams } 18393e3471ccSMark Adams 18403e3471ccSMark Adams /*@C 18413e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 18428a690491SBarry Smith from PCInitializePackage(). 18433e3471ccSMark Adams 18443e3471ccSMark Adams Level: developer 18453e3471ccSMark Adams 18463e3471ccSMark Adams .seealso: PetscInitialize() 18473e3471ccSMark Adams @*/ 18483e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 18493e3471ccSMark Adams { 18503e3471ccSMark Adams PetscErrorCode ierr; 1851*4555aa8cSStefano Zampini PetscInt l; 18523e3471ccSMark Adams 18533e3471ccSMark Adams PetscFunctionBegin; 18543e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 18553e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 18563e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 18573e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 18588e6d0c30SPeter Brune ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 18593e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1860c1c463dbSMark Adams 1861c1c463dbSMark Adams /* general events */ 1862fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr); 1863fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr); 1864fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1865fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1866c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1867c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1868fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr); 1869c1c463dbSMark Adams 18705b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 18715b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 18725b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 18735b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 18745b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 18755b89ad90SMark F. Adams ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 18765b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 18775b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1878bb235841SBarry Smith ierr = PetscLogEventRegister(" search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 18795b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 18805b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 18815b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 18825b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 18835b89ad90SMark F. Adams ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 18845b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 18855b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 18865b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 1887*4555aa8cSStefano Zampini for (l=0;l<PETSC_MG_MAXLEVELS;l++) { 1888*4555aa8cSStefano Zampini char ename[32]; 18895b89ad90SMark F. Adams 1890*4555aa8cSStefano Zampini ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l);CHKERRQ(ierr); 1891*4555aa8cSStefano Zampini ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]);CHKERRQ(ierr); 1892*4555aa8cSStefano Zampini ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l);CHKERRQ(ierr); 1893*4555aa8cSStefano Zampini ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]);CHKERRQ(ierr); 1894*4555aa8cSStefano Zampini ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l);CHKERRQ(ierr); 1895*4555aa8cSStefano Zampini ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]);CHKERRQ(ierr); 1896*4555aa8cSStefano Zampini } 18975b89ad90SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 18985b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 18995b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 19005b89ad90SMark F. Adams /* create timer stages */ 1901*4555aa8cSStefano Zampini #if defined(GAMG_STAGES) 19025b89ad90SMark F. Adams { 19035b89ad90SMark F. Adams char str[32]; 19045b89ad90SMark F. Adams PetscInt lidx; 19055b89ad90SMark F. Adams sprintf(str,"MG Level %d (finest)",0); 19065b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 19075b89ad90SMark F. Adams for (lidx=1; lidx<9; lidx++) { 1908c9567895SMark sprintf(str,"MG Level %d",(int)lidx); 19095b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 19105b89ad90SMark F. Adams } 19115b89ad90SMark F. Adams } 19125b89ad90SMark F. Adams #endif 19133e3471ccSMark Adams PetscFunctionReturn(0); 19143e3471ccSMark Adams } 19153e3471ccSMark Adams 19163e3471ccSMark Adams /*@C 19171c1aac46SBarry Smith PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 19181c1aac46SBarry Smith called from PetscFinalize() automatically. 19193e3471ccSMark Adams 19203e3471ccSMark Adams Level: developer 19213e3471ccSMark Adams 19223e3471ccSMark Adams .seealso: PetscFinalize() 19233e3471ccSMark Adams @*/ 19243e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 19253e3471ccSMark Adams { 19263e3471ccSMark Adams PetscErrorCode ierr; 19273e3471ccSMark Adams 19283e3471ccSMark Adams PetscFunctionBegin; 19293e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 19303e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 19313e3471ccSMark Adams PetscFunctionReturn(0); 19323e3471ccSMark Adams } 1933a36cf38bSToby Isaac 1934a36cf38bSToby Isaac /*@C 1935a36cf38bSToby Isaac PCGAMGRegister - Register a PCGAMG implementation. 1936a36cf38bSToby Isaac 1937a36cf38bSToby Isaac Input Parameters: 1938a36cf38bSToby Isaac + type - string that will be used as the name of the GAMG type. 1939a36cf38bSToby Isaac - create - function for creating the gamg context. 1940a36cf38bSToby Isaac 1941a36cf38bSToby Isaac Level: advanced 1942a36cf38bSToby Isaac 19431c1aac46SBarry Smith .seealso: PCGAMGType, PCGAMG, PCGAMGSetType() 1944a36cf38bSToby Isaac @*/ 1945a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1946a36cf38bSToby Isaac { 1947a36cf38bSToby Isaac PetscErrorCode ierr; 1948a36cf38bSToby Isaac 1949a36cf38bSToby Isaac PetscFunctionBegin; 1950a36cf38bSToby Isaac ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1951a36cf38bSToby Isaac ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr); 1952a36cf38bSToby Isaac PetscFunctionReturn(0); 1953a36cf38bSToby Isaac } 1954a36cf38bSToby Isaac 1955