15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4af0996ceSBarry Smith #include <petsc/private/matimpl.h> 5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 618c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/ 7f96513f1SMatthew G Knepley 8c9567895SMark #if defined(PETSC_HAVE_CUDA) 9c9567895SMark #include <cuda_runtime.h> 10c9567895SMark #endif 11c9567895SMark 12c9567895SMark #if defined(PETSC_HAVE_HIP) 13c9567895SMark #include <hip/hip_runtime.h> 14c9567895SMark #endif 15c9567895SMark 160cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 174555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3]; 18fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG; 19fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO; 200cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 210cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 220cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 230cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 24fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG; 250cbbd2e1SMark F. Adams 26b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 274555aa8cSStefano Zampini #if defined(GAMG_STAGES) 2818c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS]; 29b4fbaa2aSMark F. Adams #endif 30f96513f1SMatthew G Knepley 310a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL; 323e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 339d5b6da9SMark F. Adams 34d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 36d3d6bff4SMark F. Adams { 3718c3aa7eSMark PetscErrorCode ierr, level; 38d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 39d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 40d3d6bff4SMark F. Adams 41d3d6bff4SMark F. Adams PetscFunctionBegin; 4222a233eaSStefano Zampini ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 431c1aac46SBarry Smith pc_gamg->data_sz = 0; 44878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 4518c3aa7eSMark for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) { 4618c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 4718c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 4818c3aa7eSMark } 4918c3aa7eSMark pc_gamg->emin = 0; 5018c3aa7eSMark pc_gamg->emax = 0; 51a2f3521dSMark F. Adams PetscFunctionReturn(0); 52a2f3521dSMark F. Adams } 53a2f3521dSMark F. Adams 545b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 555b89ad90SMark F. Adams /* 56c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 57a147abb0SMark F. Adams of active processors. 585b89ad90SMark F. Adams 595b89ad90SMark F. Adams Input Parameter: 60a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 61a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 629d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 63c5bfad50SMark F. Adams . cr_bs - coarse block size 643530afc2SMark F. Adams In/Output Parameter: 65a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 66afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 6711e60469SMark F. Adams Output Parameter: 683530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 695b89ad90SMark F. Adams */ 705cb416c2SMark F. Adams 71171cca9aSMark 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) 725b89ad90SMark F. Adams { 73a2f3521dSMark F. Adams PetscErrorCode ierr; 749d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 75486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 76a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 773b4367a7SBarry Smith MPI_Comm comm; 78c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 793ae0bb68SMark Adams PetscInt ncrs_eq,ncrs,f_bs; 805b89ad90SMark F. Adams 815b89ad90SMark F. Adams PetscFunctionBegin; 823b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 83ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 84ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 85c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 864555aa8cSStefano Zampini ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);CHKERRQ(ierr); 879d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 884555aa8cSStefano Zampini ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);CHKERRQ(ierr); 89038e3b61SMark F. Adams 90ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 91ce7c7f2fSMark Adams 923ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 930298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 943ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 953ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 9673911c69SBarry Smith } else { 973ae0bb68SMark Adams PetscInt bs; 983ae0bb68SMark Adams ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 993ae0bb68SMark Adams ncrs = ncrs_eq/bs; 1003ae0bb68SMark Adams } 101c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 102c9567895SMark 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. */ 103c9567895SMark #if defined(PETSC_HAVE_CUDA) 104c9567895SMark PetscShmComm pshmcomm; 105c9567895SMark PetscMPIInt locrank; 106c9567895SMark MPI_Comm loccomm; 107c9567895SMark PetscInt s_nnodes,r_nnodes, new_new_size; 108c9567895SMark cudaError_t cerr; 109c9567895SMark int devCount; 110c9567895SMark ierr = PetscShmCommGet(comm,&pshmcomm);CHKERRQ(ierr); 111c9567895SMark ierr = PetscShmCommGetMpiShmComm(pshmcomm,&loccomm);CHKERRQ(ierr); 11255b25c41SPierre Jolivet ierr = MPI_Comm_rank(loccomm, &locrank);CHKERRMPI(ierr); 113c9567895SMark s_nnodes = !locrank; 11455b25c41SPierre Jolivet ierr = MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 1152c71b3e2SJacob Faibussowitsch PetscCheckFalse(size%r_nnodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes); 116c9567895SMark devCount = 0; 117c9567895SMark cerr = cudaGetDeviceCount(&devCount); 118c9567895SMark cudaGetLastError(); /* Reset the last error */ 119c9567895SMark if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 120c9567895SMark new_new_size = r_nnodes * devCount; 121c9567895SMark new_size = new_new_size; 1227d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 123c9567895SMark } else { 124c9567895SMark ierr = PetscInfo(pc,"With Cuda but no device. Use heuristics.");CHKERRQ(ierr); 125c9567895SMark goto HEURISTIC; 126c9567895SMark } 127c9567895SMark #else 128c9567895SMark SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here"); 129c9567895SMark #endif 130c9567895SMark } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 1312c71b3e2SJacob Faibussowitsch PetscCheckFalse(nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level],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]); 132c9567895SMark new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level]; 1337d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 134c9567895SMark } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 135c9567895SMark new_size = 1; 1367d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Force coarsest grid reduction to %d active processes\n",new_size);CHKERRQ(ierr); 137c9567895SMark } else { 138472110cdSMark F. Adams PetscInt ncrs_eq_glob; 139c9567895SMark #if defined(PETSC_HAVE_CUDA) 140c9567895SMark HEURISTIC: 141c9567895SMark #endif 1420298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 143a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1447f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 145c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 1467d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Coarse grid reduction from %d to %d active processes\n",nactive,new_size);CHKERRQ(ierr); 147a2f3521dSMark F. Adams } 1482e3501ffSMark Adams if (new_size==nactive) { 149ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 150ce7c7f2fSMark Adams if (new_size < size) { 151ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 1527d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr); 153ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 154b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr); 155b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr); 156ce7c7f2fSMark Adams } 157ce7c7f2fSMark Adams } 158ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1592e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 160192c0e8bSMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1; 161885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 16271959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 1632c71b3e2SJacob Faibussowitsch PetscCheckFalse(ncrs_eq % cr_bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 164ce7c7f2fSMark Adams /* get new_size and rfactor */ 165ce7c7f2fSMark Adams if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 166ce7c7f2fSMark Adams /* find factor */ 167ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 168ce7c7f2fSMark Adams else { 169ce7c7f2fSMark Adams PetscReal best_fact = 0.; 170ce7c7f2fSMark Adams jj = -1; 171ce7c7f2fSMark Adams for (kk = 1 ; kk <= size ; kk++) { 172ce7c7f2fSMark Adams if (!(size%kk)) { /* a candidate */ 173ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 174ce7c7f2fSMark Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 175ce7c7f2fSMark Adams if (fact > best_fact) { 176ce7c7f2fSMark Adams best_fact = fact; jj = kk; 177ce7c7f2fSMark Adams } 178ce7c7f2fSMark Adams } 179ce7c7f2fSMark Adams } 180ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 181ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 182ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 183ce7c7f2fSMark Adams else expand_factor = rfactor; 184ce7c7f2fSMark Adams } 185ce7c7f2fSMark Adams new_size = size/rfactor; /* make new size one that is factor */ 1864cdfd227SMark if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1874cdfd227SMark *a_Amat_crs = Cmat; 1887d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr); 189ce7c7f2fSMark Adams PetscFunctionReturn(0); 190ce7c7f2fSMark Adams } 191ce7c7f2fSMark Adams } 1924cdfd227SMark ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 193a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 194785e854fSJed Brown ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 1952e3501ffSMark Adams if (pc_gamg->repart) { 196a5b23f4aSJose E. Roman /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 1975a9b9e01SMark F. Adams Mat adj; 1987d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 199a2f3521dSMark F. Adams /* get 'adj' */ 200c5bfad50SMark F. Adams if (cr_bs == 1) { 201038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 202806fa848SBarry Smith } else { 203a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 204eb07cef2SMark F. Adams Mat tMat; 205a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 206b4fbaa2aSMark F. Adams const PetscScalar *vals; 207b4fbaa2aSMark F. Adams const PetscInt *idx; 208a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 20939d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 210d9558ea9SBarry Smith MatType mtype; 211b4fbaa2aSMark F. Adams 212e632b94dSBarry Smith ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr); 213a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 214a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 215c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 2160a545947SLisandro Dalcin ierr = MatGetRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 217c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 218c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 2190a545947SLisandro Dalcin ierr = MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 2203ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 2213ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 22258471d46SMark F. Adams } 2236876a03eSMark F. Adams 224d9558ea9SBarry Smith ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr); 2253b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 2263ae0bb68SMark Adams ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 227d9558ea9SBarry Smith ierr = MatSetType(tMat,mtype);CHKERRQ(ierr); 228a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 229a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 230e632b94dSBarry Smith ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 231eb07cef2SMark F. Adams 232a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 233c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 23422063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 235eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 236c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 237eb07cef2SMark F. Adams PetscScalar v = 1.0; 238eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 239eb07cef2SMark F. Adams } 24022063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 241eb07cef2SMark F. Adams } 242eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 243eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 244eb07cef2SMark F. Adams 245b4fbaa2aSMark F. Adams if (llev++ == -1) { 246b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2478caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 2483b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 249b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 2503bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 251b4fbaa2aSMark F. Adams } 252eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 253eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 254a2f3521dSMark F. Adams } /* create 'adj' */ 255f150b916SMark F. Adams 256a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2575a9b9e01SMark F. Adams char prefix[256]; 2585a9b9e01SMark F. Adams const char *pcpre; 259b4fbaa2aSMark F. Adams const PetscInt *is_idx; 260b4fbaa2aSMark F. Adams MatPartitioning mpart; 261a4b7d37bSMark F. Adams IS proc_is; 2622f03bc48SMark F. Adams 2633b4367a7SBarry Smith ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 2645ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 2659d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 2668caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 26759a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 26811e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 269c5df96a5SBarry Smith ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 270a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 27111e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 2725a9b9e01SMark F. Adams 2735ef31b24SMark F. Adams /* collect IS info */ 274785e854fSJed Brown ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 275a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 276a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 277c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 278ce7c7f2fSMark Adams newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ 279eb07cef2SMark F. Adams } 2805ef31b24SMark F. Adams } 281a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 282a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2835ef31b24SMark F. Adams } 2845ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2855a9b9e01SMark F. Adams 2863b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2878263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 28831cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 289ce7c7f2fSMark Adams PetscInt targetPE; 2902c71b3e2SJacob Faibussowitsch PetscCheckFalse(new_size==nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen"); 2917d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr); 292ce7c7f2fSMark Adams targetPE = (rank/rfactor)*expand_factor; 2933b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 294a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 295e33ef3b1SMark F. Adams 29611e60469SMark F. Adams /* 297a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 29811e60469SMark F. Adams */ 299a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 3007700e67bSMark Adams is_eq_num_prim = is_eq_num; 30111e60469SMark F. Adams /* 302a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 30311e60469SMark F. Adams */ 304c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 305c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 306a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 307ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new/cr_bs; 308a2f3521dSMark F. Adams 309a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 310885364a3SMark 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 */ 311885364a3SMark Adams { 312885364a3SMark Adams Vec src_crd, dest_crd; 313885364a3SMark 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; 314885364a3SMark Adams VecScatter vecscat; 315885364a3SMark Adams PetscScalar *array; 316885364a3SMark Adams IS isscat; 317a2f3521dSMark F. Adams /* move data (for primal equations only) */ 31822063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3193b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 3203ae0bb68SMark Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 321c0dedaeaSBarry Smith ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 32211e60469SMark F. Adams /* 3239d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 324c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 32511e60469SMark F. Adams */ 326854ce69bSBarry Smith ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr); 327a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3283ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 329c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 330a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 33111e60469SMark F. Adams } 332a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3333ae0bb68SMark Adams ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 33492a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 33511e60469SMark F. Adams /* 33611e60469SMark F. Adams Create a vector to contain the original vertex information for each element 33711e60469SMark F. Adams */ 3383ae0bb68SMark Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 3399d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3403ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3413ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3429d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 343a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 344c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 345676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 346d3d6bff4SMark F. Adams } 347038e3b61SMark F. Adams } 348eb07cef2SMark F. Adams } 349eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 350eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 35111e60469SMark F. Adams /* 35211e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 35311e60469SMark F. Adams to the correct processor 35411e60469SMark F. Adams */ 3559448b7f1SJunchao Zhang ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 35611e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 35711e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35811e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35911e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 36011e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 36111e60469SMark F. Adams /* 36211e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 36311e60469SMark F. Adams */ 364c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 365578f55a3SPeter Brune ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 3662fa5cd67SKarl Rupp 3673ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3683ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3692fa5cd67SKarl Rupp 37011e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 3719d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3723ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3739d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 374a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 375c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 376d3d6bff4SMark F. Adams } 377038e3b61SMark F. Adams } 378038e3b61SMark F. Adams } 37911e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 38011e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 381885364a3SMark Adams } 382a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 3830cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 38411e60469SMark F. Adams /* 3857dae84e0SHong Zhang Invert for MatCreateSubMatrix 38611e60469SMark F. Adams */ 387a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 388a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 389c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 390a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 391a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 392a2f3521dSMark F. Adams } 3933cb8563fSToby Isaac if (Pcolumnperm) { 3943cb8563fSToby Isaac ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr); 3953cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3963cb8563fSToby Isaac } 397a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 3980cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 3990cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 400a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 401a2f3521dSMark F. Adams { 402a2f3521dSMark F. Adams Mat mat; 4037dae84e0SHong Zhang ierr = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 404a2f3521dSMark F. Adams *a_Amat_crs = mat; 405a2f3521dSMark F. Adams } 406038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 407a2f3521dSMark F. Adams 4080cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 40911e60469SMark F. Adams /* prolongator */ 41011e60469SMark F. Adams { 41111e60469SMark F. Adams IS findices; 412a2f3521dSMark F. Adams PetscInt Istart,Iend; 413a2f3521dSMark F. Adams Mat Pnew; 41462294041SBarry Smith 415a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 4160cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 4173b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 418c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 4197dae84e0SHong Zhang ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 42011e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 4211a2c6b5cSJunchao Zhang ierr = MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr); 422c5bfad50SMark F. Adams 4230cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 4243530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 425a2f3521dSMark F. Adams 426a2f3521dSMark F. Adams /* output - repartitioned */ 427a2f3521dSMark F. Adams *a_P_inout = Pnew; 428e33ef3b1SMark F. Adams } 429a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4305b89ad90SMark F. Adams 431c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 432ce7c7f2fSMark Adams 433ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 434ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 435ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4368bca76a6SMark Adams static PetscInt llev = 2; 4377d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Pinning level %D to the CPU\n",llev++);CHKERRQ(ierr); 438ce7c7f2fSMark Adams #endif 439b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr); 440b470e4b4SRichard Tran Mills ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr); 441adf5291fSStefano Zampini if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 442ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 443ce7c7f2fSMark Adams PetscMPIInt size; 444ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 445ce7c7f2fSMark Adams if (size > 1) { 446ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data; 447b470e4b4SRichard Tran Mills ierr = VecBindToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr); 448b470e4b4SRichard Tran Mills ierr = VecBindToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr); 449ce7c7f2fSMark Adams } 450ce7c7f2fSMark Adams } 451ce7c7f2fSMark Adams } 4524cdfd227SMark ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 453a2f3521dSMark F. Adams } 4545b89ad90SMark F. Adams PetscFunctionReturn(0); 4555b89ad90SMark F. Adams } 4565b89ad90SMark F. Adams 4574b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2) 4584b1575e2SStefano Zampini { 4594b1575e2SStefano Zampini PetscErrorCode ierr; 4604b1575e2SStefano Zampini const char *prefix; 4614b1575e2SStefano Zampini char addp[32]; 4624b1575e2SStefano Zampini PC_MG *mg = (PC_MG*)a_pc->data; 4634b1575e2SStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4644b1575e2SStefano Zampini 4654b1575e2SStefano Zampini PetscFunctionBegin; 4664b1575e2SStefano Zampini ierr = PCGetOptionsPrefix(a_pc,&prefix);CHKERRQ(ierr); 4677d3de750SJacob Faibussowitsch ierr = PetscInfo(a_pc,"Square Graph on level %D\n",pc_gamg->current_level+1);CHKERRQ(ierr); 4684b1575e2SStefano Zampini ierr = MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);CHKERRQ(ierr); 4694b1575e2SStefano Zampini ierr = MatSetOptionsPrefix(*Gmat2,prefix);CHKERRQ(ierr); 4704b1575e2SStefano Zampini ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);CHKERRQ(ierr); 4714b1575e2SStefano Zampini ierr = MatAppendOptionsPrefix(*Gmat2,addp);CHKERRQ(ierr); 472b4da3a1bSStefano Zampini if ((*Gmat2)->structurally_symmetric) { 473b4da3a1bSStefano Zampini ierr = MatProductSetType(*Gmat2,MATPRODUCT_AB);CHKERRQ(ierr); 474b4da3a1bSStefano Zampini } else { 4751a2c6b5cSJunchao Zhang ierr = MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr); 4764b1575e2SStefano Zampini ierr = MatProductSetType(*Gmat2,MATPRODUCT_AtB);CHKERRQ(ierr); 477b4da3a1bSStefano Zampini } 4784b1575e2SStefano Zampini ierr = MatProductSetFromOptions(*Gmat2);CHKERRQ(ierr); 4794555aa8cSStefano Zampini ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);CHKERRQ(ierr); 4804b1575e2SStefano Zampini ierr = MatProductSymbolic(*Gmat2);CHKERRQ(ierr); 4814555aa8cSStefano Zampini ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);CHKERRQ(ierr); 482b4da3a1bSStefano Zampini ierr = MatProductClear(*Gmat2);CHKERRQ(ierr); 4834b1575e2SStefano Zampini /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 4844b1575e2SStefano Zampini (*Gmat2)->assembled = PETSC_TRUE; 4854b1575e2SStefano Zampini PetscFunctionReturn(0); 4864b1575e2SStefano Zampini } 4874b1575e2SStefano Zampini 4885b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4895b89ad90SMark F. Adams /* 4905b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4915b89ad90SMark F. Adams by setting data structures and options. 4925b89ad90SMark F. Adams 4935b89ad90SMark F. Adams Input Parameter: 4945b89ad90SMark F. Adams . pc - the preconditioner context 4955b89ad90SMark F. Adams 4965b89ad90SMark F. Adams */ 4979d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 4985b89ad90SMark F. Adams { 4995b89ad90SMark F. Adams PetscErrorCode ierr; 5009d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 5015b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 5022adcac29SMark F. Adams Mat Pmat = pc->pmat; 50318c3aa7eSMark PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS]; 5043b4367a7SBarry Smith MPI_Comm comm; 505c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 50618c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS]; 50718c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 508a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 509569f4572SMark Adams MatInfo info; 510171cca9aSMark Adams PetscBool is_last = PETSC_FALSE; 5115ef31b24SMark F. Adams 5125b89ad90SMark F. Adams PetscFunctionBegin; 5133b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 514ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 515ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 516dfd5c07aSMark F. Adams 5178abdc6daSStefano Zampini if (pc->setupcalled) { 5188abdc6daSStefano Zampini if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 519878e152fSMark F. Adams /* reset everything */ 520878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 521878e152fSMark F. Adams pc->setupcalled = 0; 522806fa848SBarry Smith } else { 52384d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 52403a628feSMark F. Adams /* just do Galerkin grids */ 52558471d46SMark F. Adams Mat B,dA,dB; 52658471d46SMark F. Adams 5279d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 5284555aa8cSStefano Zampini PetscInt gl; 52958471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 53023ee1639SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 53158471d46SMark F. Adams /* (re)set to get dirty flag */ 53223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 53358471d46SMark F. Adams 5344555aa8cSStefano Zampini for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) { 5358abdc6daSStefano Zampini MatReuse reuse = MAT_INITIAL_MATRIX ; 5368abdc6daSStefano Zampini 5378abdc6daSStefano Zampini /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 53823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 5398abdc6daSStefano Zampini if (B->product) { 5408abdc6daSStefano Zampini if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) { 5418abdc6daSStefano Zampini reuse = MAT_REUSE_MATRIX; 54203a628feSMark F. Adams } 5438abdc6daSStefano Zampini } 5448abdc6daSStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); } 5458abdc6daSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 5467d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"RAP after first solve, reuse matrix level %D\n",level);CHKERRQ(ierr); 5478abdc6daSStefano Zampini } else { 5487d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"RAP after first solve, new matrix level %D\n",level);CHKERRQ(ierr); 5498abdc6daSStefano Zampini } 5504555aa8cSStefano Zampini ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);CHKERRQ(ierr); 5518abdc6daSStefano Zampini ierr = MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B);CHKERRQ(ierr); 5524555aa8cSStefano Zampini ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);CHKERRQ(ierr); 55363b77682SMark Adams if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 55423ee1639SBarry Smith ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 55558471d46SMark F. Adams dB = B; 55658471d46SMark F. Adams } 5575f8cf99dSMark F. Adams } 558d5280255SMark F. Adams 559d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 56058471d46SMark F. Adams PetscFunctionReturn(0); 561eb07cef2SMark F. Adams } 562878e152fSMark F. Adams } 563f6536408SMark F. Adams 564878e152fSMark F. Adams if (!pc_gamg->data) { 565878e152fSMark F. Adams if (pc_gamg->orig_data) { 566878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 5670298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 5682fa5cd67SKarl Rupp 569878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 570878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 571878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5722fa5cd67SKarl Rupp 573785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 574878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 575806fa848SBarry Smith } else { 5762c71b3e2SJacob Faibussowitsch PetscCheckFalse(!pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 5777700e67bSMark Adams ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 5789d5b6da9SMark F. Adams } 579878e152fSMark F. Adams } 580878e152fSMark F. Adams 581878e152fSMark F. Adams /* cache original data for reuse */ 5821c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 583785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 584878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 585878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 586878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 587878e152fSMark F. Adams } 588038e3b61SMark F. Adams 589302f38e8SMark F. Adams /* get basic dims */ 590302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 591171cca9aSMark Adams ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr); 59284d3f75bSMark F. Adams 593569f4572SMark Adams ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 594569f4572SMark Adams nnz0 = info.nz_used; 595569f4572SMark Adams nnztot = info.nz_used; 5967d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 597569f4572SMark Adams 598a2f3521dSMark F. Adams /* Get A_i and R_i */ 59962294041SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 6009ab59c8bSMark Adams pc_gamg->current_level = level; 6012c71b3e2SJacob Faibussowitsch PetscCheckFalse(level >= PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level); 6025b89ad90SMark F. Adams level1 = level + 1; 6030cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 6044555aa8cSStefano Zampini #if defined(GAMG_STAGES) 605a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 606b4fbaa2aSMark F. Adams #endif 607c8b0795cSMark F. Adams { /* construct prolongator */ 608725b86d8SJed Brown Mat Gmat; 6090cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 6107700e67bSMark Adams Mat Prol11; 611c8b0795cSMark F. Adams 6127700e67bSMark Adams ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 6131ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 6147700e67bSMark Adams ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 615c8b0795cSMark F. Adams 616a2f3521dSMark F. Adams /* could have failed to create new level */ 617a2f3521dSMark F. Adams if (Prol11) { 618f7df55f0SStefano Zampini const char *prefix; 619f7df55f0SStefano Zampini char addp[32]; 620f7df55f0SStefano Zampini 6219d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6220298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 623a2f3521dSMark F. Adams 624fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 625c8b0795cSMark F. Adams /* smooth */ 626fd1112cbSBarry Smith ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 627c8b0795cSMark F. Adams } 628c8b0795cSMark F. Adams 6290c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6301b18a24aSMark Adams PetscInt bs; 6311b18a24aSMark Adams ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 6320a3c815dSMark Adams ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 633ffc955d6SMark F. Adams } 634ffc955d6SMark F. Adams 635f7df55f0SStefano Zampini ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 636f7df55f0SStefano Zampini ierr = MatSetOptionsPrefix(Prol11,prefix);CHKERRQ(ierr); 637c9567895SMark ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level);CHKERRQ(ierr); 638f7df55f0SStefano Zampini ierr = MatAppendOptionsPrefix(Prol11,addp);CHKERRQ(ierr); 63991f31d3dSStefano Zampini /* Always generate the transpose with CUDA 640f7df55f0SStefano Zampini Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 6411a2c6b5cSJunchao Zhang ierr = MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr); 642f7df55f0SStefano Zampini ierr = MatSetFromOptions(Prol11);CHKERRQ(ierr); 6434bde40a0SMark Adams Parr[level1] = Prol11; 6444bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 6454bde40a0SMark Adams 646a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 64741b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 648a2f3521dSMark F. Adams } /* construct prolongator scope */ 6490cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 6507f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 651171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 6527d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr); 6534555aa8cSStefano Zampini #if defined(GAMG_STAGES) 654a90e85d9SMark Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 655a90e85d9SMark Adams #endif 656c8b0795cSMark F. Adams break; 657c8b0795cSMark F. Adams } 6580cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 659171cca9aSMark Adams ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */ 6602c71b3e2SJacob Faibussowitsch PetscCheckFalse(is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?"); 661171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 6620e2909e1SMark Adams if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE; 663171cca9aSMark Adams ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr); 664a2f3521dSMark F. Adams 6650cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 666171cca9aSMark Adams ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */ 667569f4572SMark Adams ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 668569f4572SMark Adams nnztot += info.nz_used; 6697d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 670569f4572SMark Adams 6714555aa8cSStefano Zampini #if defined(GAMG_STAGES) 672b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 673b4fbaa2aSMark F. Adams #endif 674a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6759ab59c8bSMark Adams if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) { 6767d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr); 677a90e85d9SMark Adams level++; 678a90e85d9SMark Adams break; 679a90e85d9SMark Adams } 680c8b0795cSMark F. Adams } /* levels */ 681c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 682c8b0795cSMark F. Adams 6837d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr); 6849d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6855b89ad90SMark F. Adams fine_level = level; 6860298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 6875b89ad90SMark F. Adams 68862294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 6890ed2132dSStefano Zampini PetscErrorCode (*savesetfromoptions[PETSC_MG_MAXLEVELS])(PetscOptionItems*,KSP); 6900ed2132dSStefano Zampini 691d5280255SMark F. Adams /* set default smoothers & set operators */ 69262294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 693ffc955d6SMark F. Adams KSP smoother; 694ffc955d6SMark F. Adams PC subpc; 695a2f3521dSMark F. Adams 6969d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 697f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 698ffc955d6SMark F. Adams 699a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 700a2f3521dSMark F. Adams /* set ops */ 70123ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 702a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 703a2f3521dSMark F. Adams 704a2f3521dSMark F. Adams /* set defaults */ 7056c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 706a2f3521dSMark F. Adams 7070c3bc534SBarry Smith /* set blocks for ASM smoother that uses the 'aggregates' */ 7080c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 7092d3561bbSSatish Balay PetscInt sz; 7107a28f3e5SMark Adams IS *iss; 711a2f3521dSMark F. Adams 7122d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7137a28f3e5SMark Adams iss = ASMLocalIDsArr[level]; 7140c3bc534SBarry Smith ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr); 7150a3c815dSMark Adams ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr); 7160c3bc534SBarry Smith ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr); 7177f66b68fSMark Adams if (!sz) { 718ffc955d6SMark F. Adams IS is; 7190a3c815dSMark Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 7207a28f3e5SMark Adams ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr); 721a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 722806fa848SBarry Smith } else { 723a94c3b12SMark F. Adams PetscInt kk; 7247a28f3e5SMark Adams ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr); 725a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 7267a28f3e5SMark Adams ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr); 727a94c3b12SMark F. Adams } 7287a28f3e5SMark Adams ierr = PetscFree(iss);CHKERRQ(ierr); 729ffc955d6SMark F. Adams } 7300298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 731ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 732806fa848SBarry Smith } else { 733890ffe84SMark Adams ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 734ffc955d6SMark F. Adams } 735d5280255SMark F. Adams } 736d5280255SMark F. Adams { 737d5280255SMark F. Adams /* coarse grid */ 738d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 739d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 7400ed2132dSStefano Zampini 741d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 74223ee1639SBarry Smith ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 743cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 744d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 745d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 746d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 747d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 74871959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 7492c71b3e2SJacob Faibussowitsch PetscCheckFalse(ii != 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 750d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 751d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 7529dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 7532fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 75408e36f19SMark Adams ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr); 755d5280255SMark F. Adams } 756cf8ae1d3SMark Adams } 757d5280255SMark F. Adams 758d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 759d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 760e55864a3SBarry Smith ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); 761d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 76269aca0b8SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr); 763d5280255SMark F. Adams 76418c3aa7eSMark /* setup cheby eigen estimates from SA */ 7650ed2132dSStefano Zampini if (pc_gamg->use_sa_esteig==1) { 76618c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 76718c3aa7eSMark KSP smoother; 76818c3aa7eSMark PetscBool ischeb; 7690ed2132dSStefano Zampini 7700ed2132dSStefano Zampini savesetfromoptions[level] = NULL; 77118c3aa7eSMark ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 77218c3aa7eSMark ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr); 77318c3aa7eSMark if (ischeb) { 77418c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data; 7750ed2132dSStefano Zampini 7760ed2132dSStefano Zampini ierr = KSPSetFromOptions(smoother);CHKERRQ(ierr); /* let command line emax override using SA's eigenvalues */ 7770ed2132dSStefano Zampini if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) { 77818c3aa7eSMark PC subpc; 77918c3aa7eSMark PetscBool isjac; 78018c3aa7eSMark ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 78118c3aa7eSMark ierr = PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);CHKERRQ(ierr); 7820ed2132dSStefano Zampini if (isjac && pc_gamg->use_sa_esteig==1) { 78318c3aa7eSMark PetscReal emax,emin; 7840ed2132dSStefano Zampini 78518c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 78618c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 7877d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 78818c3aa7eSMark cheb->emin_computed = emin; 78918c3aa7eSMark cheb->emax_computed = emax; 79018c3aa7eSMark ierr = KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);CHKERRQ(ierr); 7910ed2132dSStefano Zampini 7920ed2132dSStefano Zampini /* We have set the eigenvalues and consumed the transformation values 7930ed2132dSStefano Zampini prevent from flagging the recomputation of the eigenvalues again in PCSetUp_MG 7940ed2132dSStefano Zampini below when setfromoptions will be called again */ 7950ed2132dSStefano Zampini savesetfromoptions[level] = smoother->ops->setfromoptions; 7960ed2132dSStefano Zampini smoother->ops->setfromoptions = NULL; 79718c3aa7eSMark } 79818c3aa7eSMark } 79918c3aa7eSMark } 80018c3aa7eSMark } 8010ed2132dSStefano Zampini } 8020ed2132dSStefano Zampini 8030ed2132dSStefano Zampini ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 8040ed2132dSStefano Zampini 8050ed2132dSStefano Zampini /* restore Chebyshev smoother for next calls */ 8060ed2132dSStefano Zampini if (pc_gamg->use_sa_esteig==1) { 8070ed2132dSStefano Zampini for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 8080ed2132dSStefano Zampini if (savesetfromoptions[level]) { 8090ed2132dSStefano Zampini KSP smoother; 8100ed2132dSStefano Zampini ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 8110ed2132dSStefano Zampini smoother->ops->setfromoptions = savesetfromoptions[level]; 8120ed2132dSStefano Zampini } 8130ed2132dSStefano Zampini } 8140ed2132dSStefano Zampini } 81518c3aa7eSMark 816d5280255SMark F. Adams /* clean up */ 817d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 818587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 819587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 8205b89ad90SMark F. Adams } 821806fa848SBarry Smith } else { 8225f8cf99dSMark F. Adams KSP smoother; 8230ed2132dSStefano Zampini 824302440fdSBarry Smith ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr); 8259d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 82623ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 8275f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 8289d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 8295f8cf99dSMark F. Adams } 8305b89ad90SMark F. Adams PetscFunctionReturn(0); 8315b89ad90SMark F. Adams } 8325b89ad90SMark F. Adams 833eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 8345b89ad90SMark F. Adams /* 8355b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 8365b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 8375b89ad90SMark F. Adams 8385b89ad90SMark F. Adams Input Parameter: 8395b89ad90SMark F. Adams . pc - the preconditioner context 8405b89ad90SMark F. Adams 8415b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 8425b89ad90SMark F. Adams */ 8435b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 8445b89ad90SMark F. Adams { 8455b89ad90SMark F. Adams PetscErrorCode ierr; 8465b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 8475b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 8485b89ad90SMark F. Adams 8495b89ad90SMark F. Adams PetscFunctionBegin; 8505b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 8519b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 8529b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 8539b8ffb57SJed Brown } 8541ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 8551ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 8565b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 8575b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 8585b89ad90SMark F. Adams PetscFunctionReturn(0); 8595b89ad90SMark F. Adams } 8605b89ad90SMark F. Adams 861676e1743SMark F. Adams /*@ 862cab9ed1eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 863676e1743SMark F. Adams 8641cc46a46SBarry Smith Logically Collective on PC 865676e1743SMark F. Adams 866676e1743SMark F. Adams Input Parameters: 8671cc46a46SBarry Smith + pc - the preconditioner context 8681cc46a46SBarry Smith - n - the number of equations 869676e1743SMark F. Adams 870676e1743SMark F. Adams Options Database Key: 8711cc46a46SBarry Smith . -pc_gamg_process_eq_limit <limit> 872676e1743SMark F. Adams 87395452b02SPatrick Sanan Notes: 87495452b02SPatrick 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 875cab9ed1eSBarry Smith that has degrees of freedom 876cab9ed1eSBarry Smith 877676e1743SMark F. Adams Level: intermediate 878676e1743SMark F. Adams 879c9567895SMark .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors() 880676e1743SMark F. Adams @*/ 881676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 882676e1743SMark F. Adams { 883676e1743SMark F. Adams PetscErrorCode ierr; 884676e1743SMark F. Adams 885676e1743SMark F. Adams PetscFunctionBegin; 886676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 887676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 888676e1743SMark F. Adams PetscFunctionReturn(0); 889676e1743SMark F. Adams } 890676e1743SMark F. Adams 8911e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 892676e1743SMark F. Adams { 893c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 894c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 895676e1743SMark F. Adams 896676e1743SMark F. Adams PetscFunctionBegin; 8979d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 898676e1743SMark F. Adams PetscFunctionReturn(0); 899676e1743SMark F. Adams } 900676e1743SMark F. Adams 901389730f3SMark F. Adams /*@ 902cab9ed1eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 903389730f3SMark F. Adams 904389730f3SMark F. Adams Collective on PC 905389730f3SMark F. Adams 906389730f3SMark F. Adams Input Parameters: 9071cc46a46SBarry Smith + pc - the preconditioner context 9081cc46a46SBarry Smith - n - maximum number of equations to aim for 909389730f3SMark F. Adams 910389730f3SMark F. Adams Options Database Key: 9111cc46a46SBarry Smith . -pc_gamg_coarse_eq_limit <limit> 912389730f3SMark F. Adams 91374329af1SBarry Smith Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 91474329af1SBarry Smith has less than 1000 unknowns. 91574329af1SBarry Smith 916389730f3SMark F. Adams Level: intermediate 917389730f3SMark F. Adams 918c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors() 919389730f3SMark F. Adams @*/ 920389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 921389730f3SMark F. Adams { 922389730f3SMark F. Adams PetscErrorCode ierr; 923389730f3SMark F. Adams 924389730f3SMark F. Adams PetscFunctionBegin; 925389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 926389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 927389730f3SMark F. Adams PetscFunctionReturn(0); 928389730f3SMark F. Adams } 929389730f3SMark F. Adams 9301e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 931389730f3SMark F. Adams { 932389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 933389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 934389730f3SMark F. Adams 935389730f3SMark F. Adams PetscFunctionBegin; 9369d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 937389730f3SMark F. Adams PetscFunctionReturn(0); 938389730f3SMark F. Adams } 939389730f3SMark F. Adams 940676e1743SMark F. Adams /*@ 941cab9ed1eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 942676e1743SMark F. Adams 943676e1743SMark F. Adams Collective on PC 944676e1743SMark F. Adams 945676e1743SMark F. Adams Input Parameters: 9461cc46a46SBarry Smith + pc - the preconditioner context 9471cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 948676e1743SMark F. Adams 949676e1743SMark F. Adams Options Database Key: 9501cc46a46SBarry Smith . -pc_gamg_repartition <true,false> 951676e1743SMark F. Adams 95295452b02SPatrick Sanan Notes: 95395452b02SPatrick Sanan this will generally improve the loading balancing of the work on each level 954cab9ed1eSBarry Smith 955676e1743SMark F. Adams Level: intermediate 956676e1743SMark F. Adams 957676e1743SMark F. Adams .seealso: () 958676e1743SMark F. Adams @*/ 959cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 960676e1743SMark F. Adams { 961676e1743SMark F. Adams PetscErrorCode ierr; 962676e1743SMark F. Adams 963676e1743SMark F. Adams PetscFunctionBegin; 964676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 965cab9ed1eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 966676e1743SMark F. Adams PetscFunctionReturn(0); 967676e1743SMark F. Adams } 968676e1743SMark F. Adams 969cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 970676e1743SMark F. Adams { 971c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 972c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 973676e1743SMark F. Adams 974676e1743SMark F. Adams PetscFunctionBegin; 9759d5b6da9SMark F. Adams pc_gamg->repart = n; 976676e1743SMark F. Adams PetscFunctionReturn(0); 977676e1743SMark F. Adams } 978676e1743SMark F. Adams 979dfd5c07aSMark F. Adams /*@ 98018c3aa7eSMark PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby) 98118c3aa7eSMark 98218c3aa7eSMark Collective on PC 98318c3aa7eSMark 98418c3aa7eSMark Input Parameters: 98518c3aa7eSMark + pc - the preconditioner context 98618c3aa7eSMark - n - number of its 98718c3aa7eSMark 98818c3aa7eSMark Options Database Key: 98918c3aa7eSMark . -pc_gamg_esteig_ksp_max_it <its> 99018c3aa7eSMark 99118c3aa7eSMark Notes: 99218c3aa7eSMark 99318c3aa7eSMark Level: intermediate 99418c3aa7eSMark 99518c3aa7eSMark .seealso: () 99618c3aa7eSMark @*/ 99718c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n) 99818c3aa7eSMark { 99918c3aa7eSMark PetscErrorCode ierr; 100018c3aa7eSMark 100118c3aa7eSMark PetscFunctionBegin; 100218c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 100318c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 100418c3aa7eSMark PetscFunctionReturn(0); 100518c3aa7eSMark } 100618c3aa7eSMark 100718c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n) 100818c3aa7eSMark { 100918c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 101018c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 101118c3aa7eSMark 101218c3aa7eSMark PetscFunctionBegin; 101318c3aa7eSMark pc_gamg->esteig_max_it = n; 101418c3aa7eSMark PetscFunctionReturn(0); 101518c3aa7eSMark } 101618c3aa7eSMark 101718c3aa7eSMark /*@ 101818c3aa7eSMark PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother 101918c3aa7eSMark 102018c3aa7eSMark Collective on PC 102118c3aa7eSMark 102218c3aa7eSMark Input Parameters: 102318c3aa7eSMark + pc - the preconditioner context 102418c3aa7eSMark - n - number of its 102518c3aa7eSMark 102618c3aa7eSMark Options Database Key: 102718c3aa7eSMark . -pc_gamg_use_sa_esteig <true,false> 102818c3aa7eSMark 102918c3aa7eSMark Notes: 103018c3aa7eSMark 103118c3aa7eSMark Level: intermediate 103218c3aa7eSMark 103318c3aa7eSMark .seealso: () 103418c3aa7eSMark @*/ 103518c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 103618c3aa7eSMark { 103718c3aa7eSMark PetscErrorCode ierr; 103818c3aa7eSMark 103918c3aa7eSMark PetscFunctionBegin; 104018c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 104118c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 104218c3aa7eSMark PetscFunctionReturn(0); 104318c3aa7eSMark } 104418c3aa7eSMark 10450ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 104618c3aa7eSMark { 104718c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 104818c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 104918c3aa7eSMark 105018c3aa7eSMark PetscFunctionBegin; 105118c3aa7eSMark pc_gamg->use_sa_esteig = n ? 1 : 0; 105218c3aa7eSMark PetscFunctionReturn(0); 105318c3aa7eSMark } 105418c3aa7eSMark 105518c3aa7eSMark /*@C 105618c3aa7eSMark PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby) 105718c3aa7eSMark 105818c3aa7eSMark Collective on PC 105918c3aa7eSMark 106018c3aa7eSMark Input Parameters: 106118c3aa7eSMark + pc - the preconditioner context 106218c3aa7eSMark - t - ksp type 106318c3aa7eSMark 106418c3aa7eSMark Options Database Key: 106518c3aa7eSMark . -pc_gamg_esteig_ksp_type <type> 106618c3aa7eSMark 106718c3aa7eSMark Notes: 106818c3aa7eSMark 106918c3aa7eSMark Level: intermediate 107018c3aa7eSMark 107118c3aa7eSMark .seealso: () 107218c3aa7eSMark @*/ 107318c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[]) 107418c3aa7eSMark { 107518c3aa7eSMark PetscErrorCode ierr; 107618c3aa7eSMark 107718c3aa7eSMark PetscFunctionBegin; 107818c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 107918c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));CHKERRQ(ierr); 108018c3aa7eSMark PetscFunctionReturn(0); 108118c3aa7eSMark } 108218c3aa7eSMark 108318c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[]) 108418c3aa7eSMark { 108518c3aa7eSMark PetscErrorCode ierr; 108618c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 108718c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 108818c3aa7eSMark 108918c3aa7eSMark PetscFunctionBegin; 109018c3aa7eSMark ierr = PetscStrcpy(pc_gamg->esteig_type,t);CHKERRQ(ierr); 109118c3aa7eSMark PetscFunctionReturn(0); 109218c3aa7eSMark } 109318c3aa7eSMark 109418c3aa7eSMark /*@ 109518c3aa7eSMark PCGAMGSetEigenvalues - Set eigenvalues 109618c3aa7eSMark 109718c3aa7eSMark Collective on PC 109818c3aa7eSMark 109918c3aa7eSMark Input Parameters: 110018c3aa7eSMark + pc - the preconditioner context 110118c3aa7eSMark - emax - max eigenvalue 110218c3aa7eSMark . emin - min eigenvalue 110318c3aa7eSMark 110418c3aa7eSMark Options Database Key: 110518c3aa7eSMark . -pc_gamg_eigenvalues 110618c3aa7eSMark 110718c3aa7eSMark Level: intermediate 110818c3aa7eSMark 110918c3aa7eSMark .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType() 111018c3aa7eSMark @*/ 111118c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin) 111218c3aa7eSMark { 111318c3aa7eSMark PetscErrorCode ierr; 111418c3aa7eSMark 111518c3aa7eSMark PetscFunctionBegin; 111618c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 111718c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr); 111818c3aa7eSMark PetscFunctionReturn(0); 111918c3aa7eSMark } 112041ffd417SStefano Zampini 112118c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin) 112218c3aa7eSMark { 112318c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 112418c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 112518c3aa7eSMark 112618c3aa7eSMark PetscFunctionBegin; 11272c71b3e2SJacob Faibussowitsch PetscCheckFalse(emax <= emin,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin); 11282c71b3e2SJacob Faibussowitsch PetscCheckFalse(emax*emin <= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin); 112918c3aa7eSMark pc_gamg->emax = emax; 113018c3aa7eSMark pc_gamg->emin = emin; 113118c3aa7eSMark 113218c3aa7eSMark PetscFunctionReturn(0); 113318c3aa7eSMark } 113418c3aa7eSMark 113518c3aa7eSMark /*@ 1136cab9ed1eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1137dfd5c07aSMark F. Adams 1138dfd5c07aSMark F. Adams Collective on PC 1139dfd5c07aSMark F. Adams 1140dfd5c07aSMark F. Adams Input Parameters: 11411cc46a46SBarry Smith + pc - the preconditioner context 11421cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 1143dfd5c07aSMark F. Adams 1144dfd5c07aSMark F. Adams Options Database Key: 11451cc46a46SBarry Smith . -pc_gamg_reuse_interpolation <true,false> 1146dfd5c07aSMark F. Adams 1147dfd5c07aSMark F. Adams Level: intermediate 1148dfd5c07aSMark F. Adams 114995452b02SPatrick Sanan Notes: 115095452b02SPatrick Sanan this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1151cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1152cab9ed1eSBarry Smith 1153dfd5c07aSMark F. Adams .seealso: () 1154dfd5c07aSMark F. Adams @*/ 11551cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1156dfd5c07aSMark F. Adams { 1157dfd5c07aSMark F. Adams PetscErrorCode ierr; 1158dfd5c07aSMark F. Adams 1159dfd5c07aSMark F. Adams PetscFunctionBegin; 1160dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11611cc46a46SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1162dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1163dfd5c07aSMark F. Adams } 1164dfd5c07aSMark F. Adams 11651cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1166dfd5c07aSMark F. Adams { 1167dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1168dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1169dfd5c07aSMark F. Adams 1170dfd5c07aSMark F. Adams PetscFunctionBegin; 1171dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1172dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1173dfd5c07aSMark F. Adams } 1174dfd5c07aSMark F. Adams 1175ffc955d6SMark F. Adams /*@ 1176cab9ed1eSBarry 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. 1177ffc955d6SMark F. Adams 1178ffc955d6SMark F. Adams Collective on PC 1179ffc955d6SMark F. Adams 1180ffc955d6SMark F. Adams Input Parameters: 1181cab9ed1eSBarry Smith + pc - the preconditioner context 1182cab9ed1eSBarry Smith - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not 1183ffc955d6SMark F. Adams 1184ffc955d6SMark F. Adams Options Database Key: 1185cab9ed1eSBarry Smith . -pc_gamg_asm_use_agg 1186ffc955d6SMark F. Adams 1187ffc955d6SMark F. Adams Level: intermediate 1188ffc955d6SMark F. Adams 1189ffc955d6SMark F. Adams .seealso: () 1190ffc955d6SMark F. Adams @*/ 1191cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1192ffc955d6SMark F. Adams { 1193ffc955d6SMark F. Adams PetscErrorCode ierr; 1194ffc955d6SMark F. Adams 1195ffc955d6SMark F. Adams PetscFunctionBegin; 1196ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1197cab9ed1eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1198ffc955d6SMark F. Adams PetscFunctionReturn(0); 1199ffc955d6SMark F. Adams } 1200ffc955d6SMark F. Adams 1201cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1202ffc955d6SMark F. Adams { 1203ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1204ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1205ffc955d6SMark F. Adams 1206ffc955d6SMark F. Adams PetscFunctionBegin; 1207cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 1208ffc955d6SMark F. Adams PetscFunctionReturn(0); 1209ffc955d6SMark F. Adams } 1210ffc955d6SMark F. Adams 1211171cca9aSMark Adams /*@ 1212cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1213171cca9aSMark Adams 1214171cca9aSMark Adams Collective on PC 1215171cca9aSMark Adams 1216171cca9aSMark Adams Input Parameters: 1217171cca9aSMark Adams + pc - the preconditioner context 1218cf8ae1d3SMark Adams - flg - PETSC_TRUE to not force coarse grid onto one processor 1219171cca9aSMark Adams 1220171cca9aSMark Adams Options Database Key: 1221cf8ae1d3SMark Adams . -pc_gamg_use_parallel_coarse_grid_solver 1222171cca9aSMark Adams 1223171cca9aSMark Adams Level: intermediate 1224171cca9aSMark Adams 122539d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids() 1226171cca9aSMark Adams @*/ 1227171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1228171cca9aSMark Adams { 1229171cca9aSMark Adams PetscErrorCode ierr; 1230171cca9aSMark Adams 1231171cca9aSMark Adams PetscFunctionBegin; 1232171cca9aSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1233171cca9aSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1234171cca9aSMark Adams PetscFunctionReturn(0); 1235171cca9aSMark Adams } 1236171cca9aSMark Adams 1237171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1238171cca9aSMark Adams { 1239171cca9aSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1240171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1241171cca9aSMark Adams 1242171cca9aSMark Adams PetscFunctionBegin; 1243171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 1244ffc955d6SMark F. Adams PetscFunctionReturn(0); 1245ffc955d6SMark F. Adams } 1246ffc955d6SMark F. Adams 12474ef23d27SMark F. Adams /*@ 1248ce7c7f2fSMark Adams PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1249ce7c7f2fSMark Adams 1250ce7c7f2fSMark Adams Collective on PC 1251ce7c7f2fSMark Adams 1252ce7c7f2fSMark Adams Input Parameters: 1253ce7c7f2fSMark Adams + pc - the preconditioner context 1254ce7c7f2fSMark Adams - flg - PETSC_TRUE to pin coarse grids to CPU 1255ce7c7f2fSMark Adams 1256ce7c7f2fSMark Adams Options Database Key: 1257ce7c7f2fSMark Adams . -pc_gamg_cpu_pin_coarse_grids 1258ce7c7f2fSMark Adams 1259ce7c7f2fSMark Adams Level: intermediate 1260ce7c7f2fSMark Adams 126139d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve() 1262ce7c7f2fSMark Adams @*/ 1263ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1264ce7c7f2fSMark Adams { 1265ce7c7f2fSMark Adams PetscErrorCode ierr; 1266ce7c7f2fSMark Adams 1267ce7c7f2fSMark Adams PetscFunctionBegin; 1268ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1269ce7c7f2fSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1270ce7c7f2fSMark Adams PetscFunctionReturn(0); 1271ce7c7f2fSMark Adams } 1272ce7c7f2fSMark Adams 1273ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1274ce7c7f2fSMark Adams { 1275ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1276ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1277ce7c7f2fSMark Adams 1278ce7c7f2fSMark Adams PetscFunctionBegin; 1279ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 1280ce7c7f2fSMark Adams PetscFunctionReturn(0); 1281ce7c7f2fSMark Adams } 1282ce7c7f2fSMark Adams 1283ce7c7f2fSMark Adams /*@ 1284ce7c7f2fSMark Adams PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type) 1285ce7c7f2fSMark Adams 1286ce7c7f2fSMark Adams Collective on PC 1287ce7c7f2fSMark Adams 1288ce7c7f2fSMark Adams Input Parameters: 1289ce7c7f2fSMark Adams + pc - the preconditioner context 1290ce7c7f2fSMark Adams - flg - Layout type 1291ce7c7f2fSMark Adams 1292ce7c7f2fSMark Adams Options Database Key: 1293ce7c7f2fSMark Adams . -pc_gamg_coarse_grid_layout_type 1294ce7c7f2fSMark Adams 1295ce7c7f2fSMark Adams Level: intermediate 1296ce7c7f2fSMark Adams 129739d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids() 1298ce7c7f2fSMark Adams @*/ 1299ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1300ce7c7f2fSMark Adams { 1301ce7c7f2fSMark Adams PetscErrorCode ierr; 1302ce7c7f2fSMark Adams 1303ce7c7f2fSMark Adams PetscFunctionBegin; 1304ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1305ce7c7f2fSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr); 1306ce7c7f2fSMark Adams PetscFunctionReturn(0); 1307ce7c7f2fSMark Adams } 1308ce7c7f2fSMark Adams 1309ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1310ce7c7f2fSMark Adams { 1311ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1312ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1313ce7c7f2fSMark Adams 1314ce7c7f2fSMark Adams PetscFunctionBegin; 1315ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 1316ce7c7f2fSMark Adams PetscFunctionReturn(0); 1317ce7c7f2fSMark Adams } 1318ce7c7f2fSMark Adams 1319ce7c7f2fSMark Adams /*@ 13201cc46a46SBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 13214ef23d27SMark F. Adams 13224ef23d27SMark F. Adams Not collective on PC 13234ef23d27SMark F. Adams 13244ef23d27SMark F. Adams Input Parameters: 13251cc46a46SBarry Smith + pc - the preconditioner 13261cc46a46SBarry Smith - n - the maximum number of levels to use 13274ef23d27SMark F. Adams 13284ef23d27SMark F. Adams Options Database Key: 13294ef23d27SMark F. Adams . -pc_mg_levels 13304ef23d27SMark F. Adams 13314ef23d27SMark F. Adams Level: intermediate 13324ef23d27SMark F. Adams 13334ef23d27SMark F. Adams .seealso: () 13344ef23d27SMark F. Adams @*/ 13354ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 13364ef23d27SMark F. Adams { 13374ef23d27SMark F. Adams PetscErrorCode ierr; 13384ef23d27SMark F. Adams 13394ef23d27SMark F. Adams PetscFunctionBegin; 13404ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13414ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 13424ef23d27SMark F. Adams PetscFunctionReturn(0); 13434ef23d27SMark F. Adams } 13444ef23d27SMark F. Adams 13451e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 13464ef23d27SMark F. Adams { 13474ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 13484ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 13494ef23d27SMark F. Adams 13504ef23d27SMark F. Adams PetscFunctionBegin; 13519d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 13524ef23d27SMark F. Adams PetscFunctionReturn(0); 13534ef23d27SMark F. Adams } 13544ef23d27SMark F. Adams 13553542efc5SMark F. Adams /*@ 13563542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 13573542efc5SMark F. Adams 13583542efc5SMark F. Adams Not collective on PC 13593542efc5SMark F. Adams 13603542efc5SMark F. Adams Input Parameters: 13611cc46a46SBarry Smith + pc - the preconditioner context 1362c9567895SMark . 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 1363055c8bd0SJed Brown - n - number of threshold values provided in array 13643542efc5SMark F. Adams 13653542efc5SMark F. Adams Options Database Key: 13661cc46a46SBarry Smith . -pc_gamg_threshold <threshold> 13673542efc5SMark F. Adams 136895452b02SPatrick Sanan Notes: 1369af3c827dSMark 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. 1370af3c827dSMark 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. 1371cab9ed1eSBarry Smith 1372055c8bd0SJed Brown If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening. 1373055c8bd0SJed Brown In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold(). 1374055c8bd0SJed Brown If n is greater than the total number of levels, the excess entries in threshold will not be used. 1375055c8bd0SJed Brown 13763542efc5SMark F. Adams Level: intermediate 13773542efc5SMark F. Adams 1378af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph() 13793542efc5SMark F. Adams @*/ 1380c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 13813542efc5SMark F. Adams { 13823542efc5SMark F. Adams PetscErrorCode ierr; 13833542efc5SMark F. Adams 13843542efc5SMark F. Adams PetscFunctionBegin; 13853542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1386055c8bd0SJed Brown if (n) PetscValidRealPointer(v,2); 1387c1eae691SMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr); 13883542efc5SMark F. Adams PetscFunctionReturn(0); 13893542efc5SMark F. Adams } 13903542efc5SMark F. Adams 1391c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 13923542efc5SMark F. Adams { 1393c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1394c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1395c1eae691SMark Adams PetscInt i; 1396c1eae691SMark Adams PetscFunctionBegin; 1397055c8bd0SJed Brown for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1398055c8bd0SJed Brown for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale; 1399c1eae691SMark Adams PetscFunctionReturn(0); 1400c1eae691SMark Adams } 1401c1eae691SMark Adams 1402c1eae691SMark Adams /*@ 1403c9567895SMark PCGAMGSetRankReductionFactors - Set manual schedual for process reduction on coarse grids 1404c9567895SMark 1405c9567895SMark Collective on PC 1406c9567895SMark 1407c9567895SMark Input Parameters: 1408c9567895SMark + pc - the preconditioner context 1409c9567895SMark . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda 1410c9567895SMark - n - number of values provided in array 1411c9567895SMark 1412c9567895SMark Options Database Key: 1413c9567895SMark . -pc_gamg_rank_reduction_factors <factors> 1414c9567895SMark 1415c9567895SMark Level: intermediate 1416c9567895SMark 1417c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim() 1418c9567895SMark @*/ 1419c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1420c9567895SMark { 1421c9567895SMark PetscErrorCode ierr; 1422c9567895SMark 1423c9567895SMark PetscFunctionBegin; 1424c9567895SMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1425c9567895SMark if (n) PetscValidIntPointer(v,2); 1426c9567895SMark ierr = PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));CHKERRQ(ierr); 1427c9567895SMark PetscFunctionReturn(0); 1428c9567895SMark } 1429c9567895SMark 1430c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1431c9567895SMark { 1432c9567895SMark PC_MG *mg = (PC_MG*)pc->data; 1433c9567895SMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1434c9567895SMark PetscInt i; 1435c9567895SMark PetscFunctionBegin; 1436c9567895SMark for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1437c9567895SMark for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1438c9567895SMark PetscFunctionReturn(0); 1439c9567895SMark } 1440c9567895SMark 1441c9567895SMark /*@ 1442c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1443c1eae691SMark Adams 1444c1eae691SMark Adams Not collective on PC 1445c1eae691SMark Adams 1446c1eae691SMark Adams Input Parameters: 1447c1eae691SMark Adams + pc - the preconditioner context 1448c1eae691SMark Adams - scale - the threshold value reduction, ussually < 1.0 1449c1eae691SMark Adams 1450c1eae691SMark Adams Options Database Key: 1451c1eae691SMark Adams . -pc_gamg_threshold_scale <v> 1452c1eae691SMark Adams 1453055c8bd0SJed Brown Notes: 1454055c8bd0SJed Brown The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold(). 1455055c8bd0SJed Brown This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold(). 1456055c8bd0SJed Brown 1457c1eae691SMark Adams Level: advanced 1458c1eae691SMark Adams 1459055c8bd0SJed Brown .seealso: PCGAMGSetThreshold() 1460c1eae691SMark Adams @*/ 1461c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1462c1eae691SMark Adams { 1463c1eae691SMark Adams PetscErrorCode ierr; 14643542efc5SMark F. Adams 14653542efc5SMark F. Adams PetscFunctionBegin; 1466c1eae691SMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1467c1eae691SMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr); 1468c1eae691SMark Adams PetscFunctionReturn(0); 1469c1eae691SMark Adams } 1470c1eae691SMark Adams 1471c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1472c1eae691SMark Adams { 1473c1eae691SMark Adams PC_MG *mg = (PC_MG*)pc->data; 1474c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1475c1eae691SMark Adams PetscFunctionBegin; 1476c1eae691SMark Adams pc_gamg->threshold_scale = v; 14773542efc5SMark F. Adams PetscFunctionReturn(0); 14783542efc5SMark F. Adams } 14793542efc5SMark F. Adams 1480e20c40e8SBarry Smith /*@C 1481c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1482676e1743SMark F. Adams 1483676e1743SMark F. Adams Collective on PC 1484676e1743SMark F. Adams 1485676e1743SMark F. Adams Input Parameters: 1486c60c7ad4SBarry Smith + pc - the preconditioner context 1487c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1488676e1743SMark F. Adams 1489676e1743SMark F. Adams Options Database Key: 1490cab9ed1eSBarry Smith . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1491676e1743SMark F. Adams 1492676e1743SMark F. Adams Level: intermediate 1493676e1743SMark F. Adams 1494cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType 1495676e1743SMark F. Adams @*/ 149619fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1497676e1743SMark F. Adams { 1498676e1743SMark F. Adams PetscErrorCode ierr; 1499676e1743SMark F. Adams 1500676e1743SMark F. Adams PetscFunctionBegin; 1501676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1502806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1503676e1743SMark F. Adams PetscFunctionReturn(0); 1504676e1743SMark F. Adams } 1505676e1743SMark F. Adams 1506e20c40e8SBarry Smith /*@C 1507c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1508c60c7ad4SBarry Smith 1509c60c7ad4SBarry Smith Collective on PC 1510c60c7ad4SBarry Smith 1511c60c7ad4SBarry Smith Input Parameter: 1512c60c7ad4SBarry Smith . pc - the preconditioner context 1513c60c7ad4SBarry Smith 1514c60c7ad4SBarry Smith Output Parameter: 1515c60c7ad4SBarry Smith . type - the type of algorithm used 1516c60c7ad4SBarry Smith 1517c60c7ad4SBarry Smith Level: intermediate 1518c60c7ad4SBarry Smith 15191c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType 1520c60c7ad4SBarry Smith @*/ 1521c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1522c60c7ad4SBarry Smith { 1523c60c7ad4SBarry Smith PetscErrorCode ierr; 1524c60c7ad4SBarry Smith 1525c60c7ad4SBarry Smith PetscFunctionBegin; 1526c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1527c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1528c60c7ad4SBarry Smith PetscFunctionReturn(0); 1529c60c7ad4SBarry Smith } 1530c60c7ad4SBarry Smith 1531c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1532c60c7ad4SBarry Smith { 1533c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1534c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1535c60c7ad4SBarry Smith 1536c60c7ad4SBarry Smith PetscFunctionBegin; 1537c60c7ad4SBarry Smith *type = pc_gamg->type; 1538c60c7ad4SBarry Smith PetscFunctionReturn(0); 1539c60c7ad4SBarry Smith } 1540c60c7ad4SBarry Smith 15411e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1542676e1743SMark F. Adams { 15439d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 15441ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 15451ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1546676e1743SMark F. Adams 1547676e1743SMark F. Adams PetscFunctionBegin; 1548c60c7ad4SBarry Smith pc_gamg->type = type; 15491c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 15502c71b3e2SJacob Faibussowitsch PetscCheckFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 15511ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 15521ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 15531ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 1554e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 15553ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 15563ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 15573ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 15583ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 15593ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 15603ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 15613ae0bb68SMark Adams pc_gamg->data_sz = 0; 15621ab5ffc9SJed Brown } 15631ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 15641ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 15659d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1566676e1743SMark F. Adams PetscFunctionReturn(0); 1567676e1743SMark F. Adams } 1568676e1743SMark F. Adams 15695adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 15705adeb434SBarry Smith { 1571c1eae691SMark Adams PetscErrorCode ierr,i; 15725adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 15735adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1574*e7d4b4cbSMark Adams PetscReal gc=0, oc=0; 15755adeb434SBarry Smith PetscFunctionBegin; 15765adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," GAMG specific options\n");CHKERRQ(ierr); 1577459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =");CHKERRQ(ierr); 1578b3e187dcSStefano Zampini for (i=0;i<mg->nlevels; i++) { 1579c1eae691SMark Adams ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr); 1580c1eae691SMark Adams } 1581459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1582459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr); 1583cab9ed1eSBarry Smith if (pc_gamg->use_aggs_in_asm) { 1584cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr); 1585cab9ed1eSBarry Smith } 1586171cca9aSMark Adams if (pc_gamg->use_parallel_coarse_grid_solver) { 1587171cca9aSMark Adams ierr = PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr); 1588171cca9aSMark Adams } 1589ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1590ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1591ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */ 1592ce7c7f2fSMark Adams } 1593ce7c7f2fSMark Adams #endif 1594ce7c7f2fSMark Adams /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */ 1595ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */ 1596ce7c7f2fSMark Adams /* } else { */ 1597ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */ 1598ce7c7f2fSMark Adams /* } */ 15995adeb434SBarry Smith if (pc_gamg->ops->view) { 16005adeb434SBarry Smith ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr); 16015adeb434SBarry Smith } 1602*e7d4b4cbSMark Adams ierr = PCMGGetGridComplexity(pc,&gc,&oc);CHKERRQ(ierr); 1603*e7d4b4cbSMark Adams ierr = PetscViewerASCIIPrintf(viewer," Complexity: grid = %g operator = %g\n",gc,oc);CHKERRQ(ierr); 16045adeb434SBarry Smith PetscFunctionReturn(0); 16055adeb434SBarry Smith } 16065adeb434SBarry Smith 16074416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 16085b89ad90SMark F. Adams { 1609676e1743SMark F. Adams PetscErrorCode ierr; 1610676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1611676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 161218c3aa7eSMark PetscBool flag,f2; 16133b4367a7SBarry Smith MPI_Comm comm; 161418c3aa7eSMark char prefix[256],tname[32]; 1615c1eae691SMark Adams PetscInt i,n; 161614a9496bSBarry Smith const char *pcpre; 16170a545947SLisandro Dalcin static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL}; 16185b89ad90SMark F. Adams PetscFunctionBegin; 16193b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1620e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr); 16211a1c1e04SBarry Smith ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1622bd94a7aaSJed Brown if (flag) { 1623bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 16241ab5ffc9SJed Brown } 162518c3aa7eSMark ierr = PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);CHKERRQ(ierr); 162618c3aa7eSMark if (flag) { 162718c3aa7eSMark ierr = PCGAMGSetEstEigKSPType(pc,tname);CHKERRQ(ierr); 162818c3aa7eSMark } 1629cab9ed1eSBarry Smith ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 163018c3aa7eSMark f2 = PETSC_TRUE; 163118c3aa7eSMark ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);CHKERRQ(ierr); 163218c3aa7eSMark if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0; 16331cc46a46SBarry Smith ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 1634a303c832SJed 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); 1635cf8ae1d3SMark 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); 1636ce7c7f2fSMark 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); 1637a0095786SMark 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); 163894ae4db5SBarry 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); 163918c3aa7eSMark 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); 164094ae4db5SBarry 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); 1641a303c832SJed 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); 164218c3aa7eSMark n = PETSC_MG_MAXLEVELS; 1643c1eae691SMark Adams ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr); 164418c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1645efd3c5ceSMark Adams if (!flag) n = 1; 1646c1eae691SMark Adams i = n; 164718c3aa7eSMark do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1648c1eae691SMark Adams } 1649c9567895SMark n = PETSC_MG_MAXLEVELS; 1650c9567895SMark 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); 1651c9567895SMark if (!flag) i = 0; 1652c9567895SMark else i = n; 1653c9567895SMark do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS); 165494ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 165518c3aa7eSMark { 165618c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 165718c3aa7eSMark n = 2; 165818c3aa7eSMark ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr); 165918c3aa7eSMark if (flag) { 16602c71b3e2SJacob Faibussowitsch PetscCheckFalse(n != 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 166118c3aa7eSMark ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr); 166218c3aa7eSMark } 166318c3aa7eSMark } 1664b7cbab4eSMark Adams /* set options for subtype */ 1665e55864a3SBarry Smith if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);} 166618c3aa7eSMark 166714a9496bSBarry Smith ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 166814a9496bSBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 1669676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 16705b89ad90SMark F. Adams PetscFunctionReturn(0); 16715b89ad90SMark F. Adams } 16725b89ad90SMark F. Adams 16735b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 16745b89ad90SMark F. Adams /*MC 16751cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 16765b89ad90SMark F. Adams 1677280d9858SJed Brown Options Database Keys: 1678cab9ed1eSBarry Smith + -pc_gamg_type <type> - one of agg, geo, or classical 1679cab9ed1eSBarry Smith . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1680cab9ed1eSBarry Smith . -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1681cab9ed1eSBarry 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 1682cab9ed1eSBarry 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> 1683cab9ed1eSBarry Smith equations on each process that has degrees of freedom 1684cab9ed1eSBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 16856008e27bSRichard Tran Mills . -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level 1686c1eae691SMark Adams - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1687cab9ed1eSBarry Smith 1688cab9ed1eSBarry Smith Options Database Keys for default Aggregation: 1689cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1690cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 1691cab9ed1eSBarry Smith - -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it 1692cab9ed1eSBarry Smith 1693db9745e2SBarry Smith Multigrid options: 1694db9745e2SBarry Smith + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1695db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1696db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1697cab9ed1eSBarry Smith - -pc_mg_levels <levels> - Number of levels of multigrid to use. 16985b89ad90SMark F. Adams 169995452b02SPatrick Sanan Notes: 170095452b02SPatrick Sanan In order to obtain good performance for PCGAMG for vector valued problems you must 1701db9745e2SBarry Smith Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1702db9745e2SBarry Smith Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1703db9745e2SBarry Smith See the Users Manual Chapter 4 for more details 17041cc46a46SBarry Smith 17055b89ad90SMark F. Adams Level: intermediate 1706280d9858SJed Brown 17071cc46a46SBarry Smith .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 170818c3aa7eSMark PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType() 17095b89ad90SMark F. Adams M*/ 1710b2573a8aSBarry Smith 17118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 17125b89ad90SMark F. Adams { 1713c1eae691SMark Adams PetscErrorCode ierr,i; 17145b89ad90SMark F. Adams PC_GAMG *pc_gamg; 17155b89ad90SMark F. Adams PC_MG *mg; 17165b89ad90SMark F. Adams 17175b89ad90SMark F. Adams PetscFunctionBegin; 17181c1aac46SBarry Smith /* register AMG type */ 17191c1aac46SBarry Smith ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 17201c1aac46SBarry Smith 17215b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 17221c1aac46SBarry Smith ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); 17235b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 17245b89ad90SMark F. Adams 17255b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 1726b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 172769aca0b8SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr); 17285b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 17295b89ad90SMark F. Adams mg->innerctx = pc_gamg; 17305b89ad90SMark F. Adams 1731b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 17321ab5ffc9SJed Brown 17339d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 17349d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 17350a545947SLisandro Dalcin pc_gamg->data = NULL; 17365b89ad90SMark F. Adams 17379d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 17385b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 17395b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 17405b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 17415b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 17425adeb434SBarry Smith mg->view = PCView_GAMG; 17435b89ad90SMark F. Adams 174497d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 174597d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1746bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1747bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1748cab9ed1eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr); 174918c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);CHKERRQ(ierr); 175018c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);CHKERRQ(ierr); 175118c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr); 175218c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr); 17531cc46a46SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr); 1754cab9ed1eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr); 1755171cca9aSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr); 1756ce7c7f2fSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr); 1757ce7c7f2fSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr); 1758bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1759c9567895SMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG);CHKERRQ(ierr); 1760c1eae691SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr); 1761bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1762c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1763bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 17649d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1765d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 17660c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1767171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1768a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1769a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1770038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 177125a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 177218c3aa7eSMark for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.; 1773c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 177418c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 17759ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 1776d24ecf33SMark ierr = PetscStrcpy(pc_gamg->esteig_type,NULL);CHKERRQ(ierr); 177718c3aa7eSMark pc_gamg->esteig_max_it = 10; 177818c3aa7eSMark pc_gamg->use_sa_esteig = -1; 177918c3aa7eSMark pc_gamg->emin = 0; 178018c3aa7eSMark pc_gamg->emax = 0; 178118c3aa7eSMark 1782c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 17839d5b6da9SMark F. Adams 1784bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1785bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 17865b89ad90SMark F. Adams PetscFunctionReturn(0); 17875b89ad90SMark F. Adams } 17883e3471ccSMark Adams 17893e3471ccSMark Adams /*@C 17903e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 17918a690491SBarry Smith from PCInitializePackage(). 17923e3471ccSMark Adams 17933e3471ccSMark Adams Level: developer 17943e3471ccSMark Adams 17953e3471ccSMark Adams .seealso: PetscInitialize() 17963e3471ccSMark Adams @*/ 17973e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 17983e3471ccSMark Adams { 17993e3471ccSMark Adams PetscErrorCode ierr; 18004555aa8cSStefano Zampini PetscInt l; 18013e3471ccSMark Adams 18023e3471ccSMark Adams PetscFunctionBegin; 18033e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 18043e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 18053e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 18063e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 18078e6d0c30SPeter Brune ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 18083e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1809c1c463dbSMark Adams 1810c1c463dbSMark Adams /* general events */ 1811fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr); 1812fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr); 1813fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1814fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1815c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1816c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1817fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr); 1818c1c463dbSMark Adams 18195b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 18205b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 18215b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 18225b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 18235b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 18245b89ad90SMark F. Adams ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 18255b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 18265b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1827bb235841SBarry Smith ierr = PetscLogEventRegister(" search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 18285b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 18295b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 18305b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 18315b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 18325b89ad90SMark F. Adams ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 18335b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 18345b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 18355b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 18364555aa8cSStefano Zampini for (l=0;l<PETSC_MG_MAXLEVELS;l++) { 18374555aa8cSStefano Zampini char ename[32]; 18385b89ad90SMark F. Adams 18394555aa8cSStefano Zampini ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l);CHKERRQ(ierr); 18404555aa8cSStefano Zampini ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]);CHKERRQ(ierr); 18414555aa8cSStefano Zampini ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l);CHKERRQ(ierr); 18424555aa8cSStefano Zampini ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]);CHKERRQ(ierr); 18434555aa8cSStefano Zampini ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l);CHKERRQ(ierr); 18444555aa8cSStefano Zampini ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]);CHKERRQ(ierr); 18454555aa8cSStefano Zampini } 18465b89ad90SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 18475b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 18485b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 18495b89ad90SMark F. Adams /* create timer stages */ 18504555aa8cSStefano Zampini #if defined(GAMG_STAGES) 18515b89ad90SMark F. Adams { 18525b89ad90SMark F. Adams char str[32]; 18535b89ad90SMark F. Adams PetscInt lidx; 18545b89ad90SMark F. Adams sprintf(str,"MG Level %d (finest)",0); 18555b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 18565b89ad90SMark F. Adams for (lidx=1; lidx<9; lidx++) { 1857c9567895SMark sprintf(str,"MG Level %d",(int)lidx); 18585b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 18595b89ad90SMark F. Adams } 18605b89ad90SMark F. Adams } 18615b89ad90SMark F. Adams #endif 18623e3471ccSMark Adams PetscFunctionReturn(0); 18633e3471ccSMark Adams } 18643e3471ccSMark Adams 18653e3471ccSMark Adams /*@C 18661c1aac46SBarry Smith PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 18671c1aac46SBarry Smith called from PetscFinalize() automatically. 18683e3471ccSMark Adams 18693e3471ccSMark Adams Level: developer 18703e3471ccSMark Adams 18713e3471ccSMark Adams .seealso: PetscFinalize() 18723e3471ccSMark Adams @*/ 18733e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 18743e3471ccSMark Adams { 18753e3471ccSMark Adams PetscErrorCode ierr; 18763e3471ccSMark Adams 18773e3471ccSMark Adams PetscFunctionBegin; 18783e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 18793e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 18803e3471ccSMark Adams PetscFunctionReturn(0); 18813e3471ccSMark Adams } 1882a36cf38bSToby Isaac 1883a36cf38bSToby Isaac /*@C 1884a36cf38bSToby Isaac PCGAMGRegister - Register a PCGAMG implementation. 1885a36cf38bSToby Isaac 1886a36cf38bSToby Isaac Input Parameters: 1887a36cf38bSToby Isaac + type - string that will be used as the name of the GAMG type. 1888a36cf38bSToby Isaac - create - function for creating the gamg context. 1889a36cf38bSToby Isaac 1890a36cf38bSToby Isaac Level: advanced 1891a36cf38bSToby Isaac 18921c1aac46SBarry Smith .seealso: PCGAMGType, PCGAMG, PCGAMGSetType() 1893a36cf38bSToby Isaac @*/ 1894a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1895a36cf38bSToby Isaac { 1896a36cf38bSToby Isaac PetscErrorCode ierr; 1897a36cf38bSToby Isaac 1898a36cf38bSToby Isaac PetscFunctionBegin; 1899a36cf38bSToby Isaac ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1900a36cf38bSToby Isaac ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr); 1901a36cf38bSToby Isaac PetscFunctionReturn(0); 1902a36cf38bSToby Isaac } 1903a36cf38bSToby Isaac 1904