15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4b45d2f2cSJed Brown #include "petsc-private/matimpl.h" 5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h> 7*5b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */ 8f96513f1SMatthew G Knepley 90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 11b4fbaa2aSMark F. Adams #endif 120cbbd2e1SMark F. Adams 130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG; 150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO; 160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 200cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG; 21a2f3521dSMark F. Adams PetscLogEvent PC_GAMGKKTProl_AGG; 220cbbd2e1SMark F. Adams #endif 230cbbd2e1SMark F. Adams 24b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 25b4fbaa2aSMark F. Adams 26b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 270cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 28b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 29b4fbaa2aSMark F. Adams #endif 30f96513f1SMatthew G Knepley 31140e18c1SBarry Smith static PetscFunctionList GAMGList = 0; 323e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 339d5b6da9SMark F. Adams 34d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 35d3d6bff4SMark F. Adams #undef __FUNCT__ 36d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 37d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 38d3d6bff4SMark F. Adams { 39d3d6bff4SMark F. Adams PetscErrorCode ierr; 40d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 41d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 42d3d6bff4SMark F. Adams 43d3d6bff4SMark F. Adams PetscFunctionBegin; 44a2f3521dSMark F. Adams if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */ 45ce94432eSBarry Smith PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__); 469d5b6da9SMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 4758471d46SMark F. Adams } 480298fd71SBarry Smith pc_gamg->data = NULL; pc_gamg->data_sz = 0; 49878e152fSMark F. Adams 50878e152fSMark F. Adams if (pc_gamg->orig_data) { 51878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 52878e152fSMark F. Adams } 53a2f3521dSMark F. Adams PetscFunctionReturn(0); 54a2f3521dSMark F. Adams } 55a2f3521dSMark F. Adams 56a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */ 57a2f3521dSMark F. Adams typedef struct { 58a2f3521dSMark F. Adams Mat A11,A21,A12,Amat; 59a2f3521dSMark F. Adams IS prim_is,constr_is; 60a2f3521dSMark F. Adams } GAMGKKTMat; 61a2f3521dSMark F. Adams 62a2f3521dSMark F. Adams #undef __FUNCT__ 63a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate" 64a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out) 65a2f3521dSMark F. Adams { 66a2f3521dSMark F. Adams PetscFunctionBegin; 67a2f3521dSMark F. Adams out->Amat = A; 68a2f3521dSMark F. Adams if (iskkt) { 69a2f3521dSMark F. Adams PetscErrorCode ierr; 70a2f3521dSMark F. Adams IS is_constraint, is_prime; 71a2f3521dSMark F. Adams PetscInt nmin,nmax; 72a2f3521dSMark F. Adams 73a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr); 74a2f3521dSMark F. Adams ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr); 75a2f3521dSMark F. Adams ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr); 762fa5cd67SKarl Rupp 77a2f3521dSMark F. Adams out->prim_is = is_prime; 78a2f3521dSMark F. Adams out->constr_is = is_constraint; 79a2f3521dSMark F. Adams 80a2f3521dSMark F. Adams ierr = MatGetSubMatrix(A, is_prime, is_prime, MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr); 81a2f3521dSMark F. Adams ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr); 82a2f3521dSMark F. Adams ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr); 83806fa848SBarry Smith } else { 84a2f3521dSMark F. Adams out->A11 = A; 850298fd71SBarry Smith out->A21 = NULL; 860298fd71SBarry Smith out->A12 = NULL; 870298fd71SBarry Smith out->prim_is = NULL; 880298fd71SBarry Smith out->constr_is = NULL; 89a2f3521dSMark F. Adams } 90a2f3521dSMark F. Adams PetscFunctionReturn(0); 91a2f3521dSMark F. Adams } 92a2f3521dSMark F. Adams 93a2f3521dSMark F. Adams #undef __FUNCT__ 94a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy" 95a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat) 96a2f3521dSMark F. Adams { 97a2f3521dSMark F. Adams PetscErrorCode ierr; 98a2f3521dSMark F. Adams 99a2f3521dSMark F. Adams PetscFunctionBegin; 100a2f3521dSMark F. Adams if (mat->A11 && mat->A11 != mat->Amat) { 101a2f3521dSMark F. Adams ierr = MatDestroy(&mat->A11);CHKERRQ(ierr); 102a2f3521dSMark F. Adams } 103a2f3521dSMark F. Adams ierr = MatDestroy(&mat->A21);CHKERRQ(ierr); 104a2f3521dSMark F. Adams ierr = MatDestroy(&mat->A12);CHKERRQ(ierr); 105a2f3521dSMark F. Adams 106a2f3521dSMark F. Adams ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr); 107a2f3521dSMark F. Adams ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr); 108d3d6bff4SMark F. Adams PetscFunctionReturn(0); 109d3d6bff4SMark F. Adams } 110d3d6bff4SMark F. Adams 1115b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1125b89ad90SMark F. Adams /* 113a147abb0SMark F. Adams createLevel: create coarse op with RAP. repartition and/or reduce number 114a147abb0SMark F. Adams of active processors. 1155b89ad90SMark F. Adams 1165b89ad90SMark F. Adams Input Parameter: 117a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 118a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 1199d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 120c5bfad50SMark F. Adams . cr_bs - coarse block size 1219d5b6da9SMark F. Adams . isLast - 122a2f3521dSMark F. Adams . stokes - 1233530afc2SMark F. Adams In/Output Parameter: 124a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 125afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 12611e60469SMark F. Adams Output Parameter: 1273530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 1285b89ad90SMark F. Adams */ 1295cb416c2SMark F. Adams 1305b89ad90SMark F. Adams #undef __FUNCT__ 1318263b398SMark F. Adams #define __FUNCT__ "createLevel" 1322fa5cd67SKarl Rupp static PetscErrorCode createLevel(const PC pc,const Mat Amat_fine,const PetscInt cr_bs,const PetscBool isLast,const PetscBool stokes,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc) 1335b89ad90SMark F. Adams { 134a2f3521dSMark F. Adams PetscErrorCode ierr; 1359d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 136486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1379d5b6da9SMark F. Adams const PetscBool repart = pc_gamg->repart; 1389d5b6da9SMark F. Adams const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 139a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 1403b4367a7SBarry Smith MPI_Comm comm; 141c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 142c5bfad50SMark F. Adams PetscInt ncrs_eq,ncrs_prim,f_bs; 1435b89ad90SMark F. Adams 1445b89ad90SMark F. Adams PetscFunctionBegin; 1453b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 1463b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1473b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 148c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 14911e60469SMark F. Adams /* RAP */ 1509d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 151038e3b61SMark F. Adams 152a2f3521dSMark F. Adams /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/ 153a2f3521dSMark F. Adams ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 15471959b99SBarry Smith if (pc_gamg->data_sz % (pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_sz %D not divisible by (pc_gamg->data_cell_cols %D *pc_gamg->data_cell_rows %D)",pc_gamg->data_sz,pc_gamg->data_cell_cols,pc_gamg->data_cell_rows); 1550298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 156a2f3521dSMark F. Adams 157c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 158a2f3521dSMark F. Adams { 159472110cdSMark F. Adams PetscInt ncrs_eq_glob; 1600298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 161c5df96a5SBarry Smith new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 162c5df96a5SBarry Smith if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1; 163c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 164c5df96a5SBarry Smith if (isLast) new_size = 1; 165a2f3521dSMark F. Adams } 166f852f58cSMark F. Adams 1672fa5cd67SKarl Rupp if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 1682fa5cd67SKarl Rupp else { 169a2f3521dSMark F. 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; 170a2f3521dSMark F. Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old; 171a2f3521dSMark F. Adams IS is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices; 1725a9b9e01SMark F. Adams VecScatter vecscat; 17322063be5SMark F. Adams PetscScalar *array; 17422063be5SMark F. Adams Vec src_crd, dest_crd; 175e33ef3b1SMark F. Adams 17671959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 17771959b99SBarry Smith if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 1780cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1790cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 180b4fbaa2aSMark F. Adams #endif 181a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 182c5df96a5SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr); 183a2f3521dSMark F. Adams if (repart && !stokes) { 184a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1855a9b9e01SMark F. Adams Mat adj; 1865a9b9e01SMark F. Adams 187a2f3521dSMark F. Adams if (pc_gamg->verbose>0) { 1883b4367a7SBarry Smith if (pc_gamg->verbose==1) PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq); 189a2f3521dSMark F. Adams else { 190a2f3521dSMark F. Adams PetscInt n; 1913b4367a7SBarry Smith ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr); 1923b4367a7SBarry Smith PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n); 193a2f3521dSMark F. Adams } 194a2f3521dSMark F. Adams } 1955a9b9e01SMark F. Adams 196a2f3521dSMark F. Adams /* get 'adj' */ 197c5bfad50SMark F. Adams if (cr_bs == 1) { 198038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 199806fa848SBarry Smith } else { 200a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 201eb07cef2SMark F. Adams Mat tMat; 202a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 203b4fbaa2aSMark F. Adams const PetscScalar *vals; 204b4fbaa2aSMark F. Adams const PetscInt *idx; 205a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 2069057884aSMark F. Adams static PetscInt llev = 0; 207b4fbaa2aSMark F. Adams 208a2f3521dSMark F. Adams ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr); 209a2f3521dSMark F. Adams ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr); 210a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 211a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 212c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 21358471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 214c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 215c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 21658471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 217a2f3521dSMark F. Adams if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim; 218c5bfad50SMark F. Adams if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim; 21958471d46SMark F. Adams } 2206876a03eSMark F. Adams 2213b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 222806fa848SBarry Smith ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 223a2f3521dSMark F. Adams ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr); 224a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 225a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 22658471d46SMark F. Adams ierr = PetscFree(d_nnz);CHKERRQ(ierr); 2275f8cf99dSMark F. Adams ierr = PetscFree(o_nnz);CHKERRQ(ierr); 228eb07cef2SMark F. Adams 229a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 230c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 23122063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 232eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 233c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 234eb07cef2SMark F. Adams PetscScalar v = 1.0; 235eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 236eb07cef2SMark F. Adams } 23722063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 238eb07cef2SMark F. Adams } 239eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 240eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 241eb07cef2SMark F. Adams 242b4fbaa2aSMark F. Adams if (llev++ == -1) { 243b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2448caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 2453b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 246b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 2473bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 248b4fbaa2aSMark F. Adams } 249b4fbaa2aSMark F. Adams 250eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 251eb07cef2SMark F. Adams 252eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 253a2f3521dSMark F. Adams } /* create 'adj' */ 254f150b916SMark F. Adams 255a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2565a9b9e01SMark F. Adams char prefix[256]; 2575a9b9e01SMark F. Adams const char *pcpre; 258b4fbaa2aSMark F. Adams const PetscInt *is_idx; 259b4fbaa2aSMark F. Adams MatPartitioning mpart; 260a4b7d37bSMark F. Adams IS proc_is; 261a2f3521dSMark F. Adams PetscInt targetPE; 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 */ 274a2f3521dSMark F. Adams ierr = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr); 275a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 276a2f3521dSMark F. Adams targetPE = 1; /* bring to "front" of machine */ 277c5df96a5SBarry Smith /*targetPE = size/new_size;*/ /* spread partitioning across machine */ 278a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 279c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 280a2f3521dSMark F. Adams newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 281eb07cef2SMark F. Adams } 2825ef31b24SMark F. Adams } 283a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 284a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2855ef31b24SMark F. Adams } 2865ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2875a9b9e01SMark F. Adams 2883b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2898263b398SMark F. Adams if (newproc_idx != 0) { 2908263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 2915ef31b24SMark F. Adams } 292806fa848SBarry Smith } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 293a2f3521dSMark F. Adams 294a2f3521dSMark F. Adams PetscInt rfactor,targetPE; 2955a9b9e01SMark F. Adams /* find factor */ 296c5df96a5SBarry Smith if (new_size == 1) rfactor = size; /* easy */ 2975a9b9e01SMark F. Adams else { 2985a9b9e01SMark F. Adams PetscReal best_fact = 0.; 2995a9b9e01SMark F. Adams jj = -1; 300c5df96a5SBarry Smith for (kk = 1 ; kk <= size ; kk++) { 301c5df96a5SBarry Smith if (size%kk==0) { /* a candidate */ 302c5df96a5SBarry Smith PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 3035a9b9e01SMark F. Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 3045a9b9e01SMark F. Adams if (fact > best_fact) { 3055a9b9e01SMark F. Adams best_fact = fact; jj = kk; 3065a9b9e01SMark F. Adams } 3075a9b9e01SMark F. Adams } 3085a9b9e01SMark F. Adams } 3095a9b9e01SMark F. Adams if (jj != -1) rfactor = jj; 310a2f3521dSMark F. Adams else rfactor = 1; /* does this happen .. a prime */ 3115a9b9e01SMark F. Adams } 312c5df96a5SBarry Smith new_size = size/rfactor; 3135a9b9e01SMark F. Adams 314c5df96a5SBarry Smith if (new_size==nactive) { 315a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 3165a9b9e01SMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 317a2f3521dSMark F. Adams if (pc_gamg->verbose>0) { 3183b4367a7SBarry Smith PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq); 319a2f3521dSMark F. Adams } 3205a9b9e01SMark F. Adams PetscFunctionReturn(0); 3215a9b9e01SMark F. Adams } 3225a9b9e01SMark F. Adams 3233b4367a7SBarry Smith if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq); 324c5df96a5SBarry Smith targetPE = rank/rfactor; 3253b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 3265a9b9e01SMark F. Adams 327a2f3521dSMark F. Adams if (stokes) { 3283b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr); 3295a9b9e01SMark F. Adams } 330a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 331e33ef3b1SMark F. Adams 33211e60469SMark F. Adams /* 333a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 33411e60469SMark F. Adams */ 335a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 336a2f3521dSMark F. Adams if (stokes) { 337a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr); 338806fa848SBarry Smith } else is_eq_num_prim = is_eq_num; 33911e60469SMark F. Adams /* 340a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 34111e60469SMark F. Adams */ 342c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 343c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 344a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 345a2f3521dSMark F. Adams if (stokes) { 346c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc_prim, size, counts);CHKERRQ(ierr); 347a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr); 348c5df96a5SBarry Smith ncrs_prim_new = counts[rank]/cr_bs; /* nodes */ 349806fa848SBarry Smith } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */ 350a2f3521dSMark F. Adams 351a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 3520cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3530cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 354b4fbaa2aSMark F. Adams #endif 355a2f3521dSMark F. Adams 356a2f3521dSMark F. Adams /* move data (for primal equations only) */ 35722063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3583b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 359a2f3521dSMark F. Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr); 360470e26f8SMark F. Adams ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */ 36111e60469SMark F. Adams /* 3629d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 363c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 36411e60469SMark F. Adams */ 365a2f3521dSMark F. Adams ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr); 366a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 367a2f3521dSMark F. Adams for (ii=0,jj=0; ii<ncrs_prim; ii++) { 368c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 369a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 37011e60469SMark F. Adams } 371a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3723b4367a7SBarry Smith ierr = ISCreateGeneral(comm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 37392a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 37411e60469SMark F. Adams /* 37511e60469SMark F. Adams Create a vector to contain the original vertex information for each element 37611e60469SMark F. Adams */ 377a2f3521dSMark F. Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr); 3789d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 379a2f3521dSMark F. Adams const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows; 380a2f3521dSMark F. Adams for (ii=0; ii<ncrs_prim; ii++) { 3819d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 382a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 383c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 384676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 385d3d6bff4SMark F. Adams } 386038e3b61SMark F. Adams } 387eb07cef2SMark F. Adams } 388eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 389eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 39011e60469SMark F. Adams /* 39111e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 39211e60469SMark F. Adams to the correct processor 39311e60469SMark F. Adams */ 3940298fd71SBarry Smith ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 39511e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 39611e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 39711e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 39811e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 39911e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 40011e60469SMark F. Adams /* 40111e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 40211e60469SMark F. Adams */ 403c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 404a2f3521dSMark F. Adams ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr); 4052fa5cd67SKarl Rupp 406a2f3521dSMark F. Adams pc_gamg->data_sz = node_data_sz*ncrs_prim_new; 407a2f3521dSMark F. Adams strideNew = ncrs_prim_new*ndata_rows; 4082fa5cd67SKarl Rupp 40911e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 4109d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 411a2f3521dSMark F. Adams for (ii=0; ii<ncrs_prim_new; ii++) { 4129d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 413a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 414c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 415d3d6bff4SMark F. Adams } 416038e3b61SMark F. Adams } 417038e3b61SMark F. Adams } 41811e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 41911e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 420a2f3521dSMark F. Adams 421a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 4220cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4230cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 424ed3f9983SMark F. Adams #endif 425a2f3521dSMark F. Adams 42611e60469SMark F. Adams /* 42711e60469SMark F. Adams Invert for MatGetSubMatrix 42811e60469SMark F. Adams */ 429a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 430a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 431c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 432a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 433a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 434a2f3521dSMark F. Adams } 435a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 4360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4370cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 4380cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 439ed3f9983SMark F. Adams #endif 440a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 441a2f3521dSMark F. Adams { 442a2f3521dSMark F. Adams Mat mat; 443806fa848SBarry Smith ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 444a2f3521dSMark F. Adams *a_Amat_crs = mat; 445c5bfad50SMark F. Adams 446c5bfad50SMark F. Adams if (!PETSC_TRUE) { 447c5bfad50SMark F. Adams PetscInt cbs, rbs; 448c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr); 449c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 450c5bfad50SMark F. Adams ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr); 451c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr); 452c5bfad50SMark F. Adams } 453a2f3521dSMark F. Adams } 454038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 455a2f3521dSMark F. Adams 4560cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4570cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 458ed3f9983SMark F. Adams #endif 45911e60469SMark F. Adams /* prolongator */ 46011e60469SMark F. Adams { 46111e60469SMark F. Adams IS findices; 462a2f3521dSMark F. Adams PetscInt Istart,Iend; 463a2f3521dSMark F. Adams Mat Pnew; 464a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 4650cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4660cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 467ed3f9983SMark F. Adams #endif 4683b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 469c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 470806fa848SBarry Smith ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 47111e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 472c5bfad50SMark F. Adams 473c5bfad50SMark F. Adams if (!PETSC_TRUE) { 474c5bfad50SMark F. Adams PetscInt cbs, rbs; 475c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr); 476c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 477c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr); 478c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 479c5bfad50SMark F. Adams } 4800cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4810cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 482ed3f9983SMark F. Adams #endif 4833530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 484a2f3521dSMark F. Adams 485a2f3521dSMark F. Adams /* output - repartitioned */ 486a2f3521dSMark F. Adams *a_P_inout = Pnew; 487e33ef3b1SMark F. Adams } 488a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4895b89ad90SMark F. Adams 490c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 491a2f3521dSMark F. Adams } 4925a9b9e01SMark F. Adams 493a2f3521dSMark F. Adams /* outout matrix data */ 494c8b0795cSMark F. Adams if (!PETSC_TRUE) { 495c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 496c8b0795cSMark F. Adams if (llev==0) { 497c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 4983b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 499c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 500c8b0795cSMark F. Adams ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr); 501c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 502c8b0795cSMark F. Adams } 503c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 5043b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 505c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 506c8b0795cSMark F. Adams ierr = MatView(Cmat, viewer);CHKERRQ(ierr); 507c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 508c8b0795cSMark F. Adams } 5095b89ad90SMark F. Adams PetscFunctionReturn(0); 5105b89ad90SMark F. Adams } 5115b89ad90SMark F. Adams 5125b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 5135b89ad90SMark F. Adams /* 5145b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 5155b89ad90SMark F. Adams by setting data structures and options. 5165b89ad90SMark F. Adams 5175b89ad90SMark F. Adams Input Parameter: 5185b89ad90SMark F. Adams . pc - the preconditioner context 5195b89ad90SMark F. Adams 5205b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 5215b89ad90SMark F. Adams 5225b89ad90SMark F. Adams Notes: 5235b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 5245b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 5255b89ad90SMark F. Adams */ 5265b89ad90SMark F. Adams #undef __FUNCT__ 5275b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 5289d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 5295b89ad90SMark F. Adams { 5305b89ad90SMark F. Adams PetscErrorCode ierr; 5319d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 5325b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 5332adcac29SMark F. Adams Mat Pmat = pc->pmat; 534a2f3521dSMark F. Adams PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 5353b4367a7SBarry Smith MPI_Comm comm; 536c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 537587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 538c8b0795cSMark F. Adams PetscReal emaxs[GAMG_MAXLEVELS]; 539e696ed0bSMark F. Adams IS *ASMLocalIDsArr[GAMG_MAXLEVELS]; 540a2f3521dSMark F. Adams GAMGKKTMat kktMatsArr[GAMG_MAXLEVELS]; 541a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 542737a81a9SMark F. Adams MatInfo info; 5435169fedaSMark F. Adams PetscBool stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol); 5445ef31b24SMark F. Adams 5455b89ad90SMark F. Adams PetscFunctionBegin; 5463b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5473b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 5483b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 549dfd5c07aSMark F. Adams 5503b4367a7SBarry Smith if (pc_gamg->verbose>2) PetscPrintf(comm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled); 551dfd5c07aSMark F. Adams 55284d3f75bSMark F. Adams if (pc_gamg->setup_count++ > 0) { 553878e152fSMark F. Adams if (redo_mesh_setup) { 554878e152fSMark F. Adams /* reset everything */ 555878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 556878e152fSMark F. Adams pc->setupcalled = 0; 557806fa848SBarry Smith } else { 55884d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 55903a628feSMark F. Adams /* just do Galerkin grids */ 56058471d46SMark F. Adams Mat B,dA,dB; 56158471d46SMark F. Adams 56271959b99SBarry Smith if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 5639d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 56458471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5650298fd71SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,NULL);CHKERRQ(ierr); 56658471d46SMark F. Adams /* (re)set to get dirty flag */ 5679d5b6da9SMark F. Adams ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 56858471d46SMark F. Adams 5692fb0b348SMark F. Adams for (level=pc_gamg->Nlevels-2; level>=0; level--) { 57003a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 5712fb0b348SMark F. Adams if (pc_gamg->setup_count==2 && (pc_gamg->repart || level==0)) { 57203a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 573084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 5742fa5cd67SKarl Rupp 57503a628feSMark F. Adams mglevels[level]->A = B; 576806fa848SBarry Smith } else { 5770298fd71SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B,NULL);CHKERRQ(ierr); 57858471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 57903a628feSMark F. Adams } 58058471d46SMark F. Adams ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 58158471d46SMark F. Adams dB = B; 58258471d46SMark F. Adams } 5835f8cf99dSMark F. Adams } 584d5280255SMark F. Adams 585d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 586d5280255SMark F. Adams 587d5280255SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 5881ac9965eSMark F. Adams ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr); 58958471d46SMark F. Adams PetscFunctionReturn(0); 590eb07cef2SMark F. Adams } 591878e152fSMark F. Adams } 592f6536408SMark F. Adams 5930298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 594eb07cef2SMark F. Adams 5959d1b15c3SMark F. Adams ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr); 5969d1b15c3SMark F. Adams 597878e152fSMark F. Adams if (!pc_gamg->data) { 598878e152fSMark F. Adams if (pc_gamg->orig_data) { 599878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 6000298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 6012fa5cd67SKarl Rupp 602878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 603878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 604878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 6052fa5cd67SKarl Rupp 606878e152fSMark F. Adams ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr); 607878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 608806fa848SBarry Smith } else { 6091ab5ffc9SJed Brown if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 6103b4367a7SBarry Smith if (stokes) SETERRQ(comm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems"); 6111ab5ffc9SJed Brown ierr = pc_gamg->ops->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr); 6129d5b6da9SMark F. Adams } 613878e152fSMark F. Adams } 614878e152fSMark F. Adams 615878e152fSMark F. Adams /* cache original data for reuse */ 616878e152fSMark F. Adams if (!pc_gamg->orig_data && redo_mesh_setup) { 617878e152fSMark F. Adams ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr); 618878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 619878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 620878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 621878e152fSMark F. Adams } 622038e3b61SMark F. Adams 623302f38e8SMark F. Adams /* get basic dims */ 6242fa5cd67SKarl Rupp if (stokes) bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */ 6252fa5cd67SKarl Rupp else { 626302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 627302f38e8SMark F. Adams } 628a2f3521dSMark F. Adams 629a2f3521dSMark F. Adams ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr); 630c8b0795cSMark F. Adams if (pc_gamg->verbose) { 63184f9421dSMark F. Adams PetscInt NN = M; 63284f9421dSMark F. Adams if (pc_gamg->verbose==1) { 63384f9421dSMark F. Adams ierr = MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr); 6343bf036e2SBarry Smith ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr); 635806fa848SBarry Smith } else { 636806fa848SBarry Smith ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 63784f9421dSMark F. Adams } 638b2a4f308SMark F. Adams nnz0 = info.nz_used; 639b2a4f308SMark F. Adams nnztot = info.nz_used; 6403b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 641c5df96a5SBarry Smith rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 642c5df96a5SBarry Smith (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr); 643c8b0795cSMark F. Adams } 64484d3f75bSMark F. Adams 645a2f3521dSMark F. Adams /* Get A_i and R_i */ 646c5df96a5SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */ 647c5df96a5SBarry Smith level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (size==1 || nactivepe>1); */ 6480205a208SMark F. Adams level++) { 6495b89ad90SMark F. Adams level1 = level + 1; 6500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6510cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 652a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 653a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 654b4fbaa2aSMark F. Adams #endif 655a2f3521dSMark F. Adams #endif 656a2f3521dSMark F. Adams /* deal with Stokes, get sub matrices */ 6579d1b15c3SMark F. Adams if (level > 0) { 658a2f3521dSMark F. Adams ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr); 6599d1b15c3SMark F. Adams } 660c8b0795cSMark F. Adams { /* construct prolongator */ 661725b86d8SJed Brown Mat Gmat; 6620cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 663a2f3521dSMark F. Adams Mat Prol11,Prol22; 664c8b0795cSMark F. Adams 6651ab5ffc9SJed Brown ierr = pc_gamg->ops->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr); 6661ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 6671ab5ffc9SJed Brown ierr = pc_gamg->ops->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr); 668c8b0795cSMark F. Adams 669a2f3521dSMark F. Adams /* could have failed to create new level */ 670a2f3521dSMark F. Adams if (Prol11) { 6719d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6720298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 673a2f3521dSMark F. Adams 674a2f3521dSMark F. Adams if (stokes) { 6751ab5ffc9SJed Brown if (!pc_gamg->ops->formkktprol) SETERRQ(comm,PETSC_ERR_USER,"Stokes not supportd by AMG method."); 676a2f3521dSMark F. Adams /* R A12 == (T = A21 P)'; G = T' T; coarsen G; form plain agg with G */ 6771ab5ffc9SJed Brown ierr = pc_gamg->ops->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr); 678c8b0795cSMark F. Adams } 679c8b0795cSMark F. Adams 6801ab5ffc9SJed Brown if (pc_gamg->ops->optprol) { 681c8b0795cSMark F. Adams /* smooth */ 6821ab5ffc9SJed Brown ierr = pc_gamg->ops->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr); 683c8b0795cSMark F. Adams } 684c8b0795cSMark F. Adams 685a2f3521dSMark F. Adams if (stokes) { 686da80777bSKarl Rupp IS is_row[2]; 687da80777bSKarl Rupp Mat a[4]; 688da80777bSKarl Rupp 689da80777bSKarl Rupp is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is; 6900298fd71SBarry Smith a[0] = Prol11; a[1] = NULL; a[2] = NULL; a[3] = Prol22; 6913b4367a7SBarry Smith ierr = MatCreateNest(comm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr); 6922fa5cd67SKarl Rupp } else Parr[level1] = Prol11; 6930298fd71SBarry Smith } else Parr[level1] = NULL; 694ffc955d6SMark F. Adams 695ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 696806fa848SBarry Smith ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 697ffc955d6SMark F. Adams } 698ffc955d6SMark F. Adams 699a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 70041b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 701a2f3521dSMark F. Adams } /* construct prolongator scope */ 7020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7030cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 704c8b0795cSMark F. Adams #endif 7059d5b6da9SMark F. Adams /* cache eigen estimate */ 7069d5b6da9SMark F. Adams if (pc_gamg->emax_id != -1) { 7079d5b6da9SMark F. Adams PetscBool flag; 708806fa848SBarry Smith ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr); 7099d5b6da9SMark F. Adams if (!flag) emaxs[level] = -1.; 710806fa848SBarry Smith } else emaxs[level] = -1.; 7112adcac29SMark F. Adams if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 712c8b0795cSMark F. Adams if (!Parr[level1]) { 713806fa848SBarry Smith if (pc_gamg->verbose) { 7143b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr); 715806fa848SBarry Smith } 716c8b0795cSMark F. Adams break; 717c8b0795cSMark F. Adams } 7180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7190cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 720b4fbaa2aSMark F. Adams #endif 721a2f3521dSMark F. Adams 722a2f3521dSMark F. Adams ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2), 723806fa848SBarry Smith stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr); 724a2f3521dSMark F. Adams 7250cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7260cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 727b4fbaa2aSMark F. Adams #endif 728a2f3521dSMark F. Adams ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr); 729a2f3521dSMark F. Adams 730a2f3521dSMark F. Adams if (pc_gamg->verbose > 0) { 7310cbbd2e1SMark F. Adams PetscInt NN = M; 7320cbbd2e1SMark F. Adams if (pc_gamg->verbose==1) { 733a2f3521dSMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr); 7343bf036e2SBarry Smith ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr); 735806fa848SBarry Smith } else { 736806fa848SBarry Smith ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 7370cbbd2e1SMark F. Adams } 738a2f3521dSMark F. Adams 7390cbbd2e1SMark F. Adams nnztot += info.nz_used; 7403b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 741c5df96a5SBarry Smith rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols, 742806fa848SBarry Smith (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr); 743c8b0795cSMark F. Adams } 744a2f3521dSMark F. Adams 745a2f3521dSMark F. Adams /* stop if one node -- could pull back for singular problems */ 746c8b0795cSMark F. Adams if (M/pc_gamg->data_cell_cols < 2) { 74781708dccSMark F. Adams level++; 74881708dccSMark F. Adams break; 74981708dccSMark F. Adams } 7500cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 751b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 752b4fbaa2aSMark F. Adams #endif 753c8b0795cSMark F. Adams } /* levels */ 754c8b0795cSMark F. Adams 755c8b0795cSMark F. Adams if (pc_gamg->data) { 756c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 7570298fd71SBarry Smith pc_gamg->data = NULL; 7585b89ad90SMark F. Adams } 759c8b0795cSMark F. Adams 7603b4367a7SBarry Smith if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 7619d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 7625b89ad90SMark F. Adams fine_level = level; 7630298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 7645b89ad90SMark F. Adams 76584d3f75bSMark F. Adams /* simple setup */ 76684d3f75bSMark F. Adams if (!PETSC_TRUE) { 76784d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 76884d3f75bSMark F. Adams for (lidx=0,level=pc_gamg->Nlevels-1; 76984d3f75bSMark F. Adams lidx<fine_level; 77084d3f75bSMark F. Adams lidx++, level--) { 77184d3f75bSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr); 77284d3f75bSMark F. Adams ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 77384d3f75bSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 77484d3f75bSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 77584d3f75bSMark F. Adams } 77684d3f75bSMark F. Adams ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 77784d3f75bSMark F. Adams 77884d3f75bSMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 779806fa848SBarry Smith } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 780d5280255SMark F. Adams /* set default smoothers & set operators */ 7819d5b6da9SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 782587fa25dSMark F. Adams lidx <= fine_level; 783587fa25dSMark F. Adams lidx++, level--) { 784ffc955d6SMark F. Adams KSP smoother; 785ffc955d6SMark F. Adams PC subpc; 786a2f3521dSMark F. Adams 7879d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 788f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 789ffc955d6SMark F. Adams 790a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 791a2f3521dSMark F. Adams /* set ops */ 792a2f3521dSMark F. Adams ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 793a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 794a2f3521dSMark F. Adams 795a2f3521dSMark F. Adams /* create field split PC, get subsmoother */ 796a2f3521dSMark F. Adams if (stokes) { 797a2f3521dSMark F. Adams KSP *ksps; 798a2f3521dSMark F. Adams PetscInt nn; 799a2f3521dSMark F. Adams ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr); 800a2f3521dSMark F. Adams ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr); 801a2f3521dSMark F. Adams ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr); 802a2f3521dSMark F. Adams smoother = ksps[0]; 803a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 804a2f3521dSMark F. Adams ierr = PetscFree(ksps);CHKERRQ(ierr); 805a2f3521dSMark F. Adams } 806a2f3521dSMark F. Adams ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr); 807a2f3521dSMark F. Adams 808a2f3521dSMark F. Adams /* set defaults */ 8096c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 810a2f3521dSMark F. Adams 811ffc955d6SMark F. Adams /* override defaults and command line args (!) */ 812ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 8132d3561bbSSatish Balay PetscInt sz; 8142d3561bbSSatish Balay IS *is; 815a2f3521dSMark F. Adams 8162d3561bbSSatish Balay sz = nASMBlocksArr[level]; 8172d3561bbSSatish Balay is = ASMLocalIDsArr[level]; 818ffc955d6SMark F. Adams ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr); 819ffc955d6SMark F. Adams if (sz==0) { 820ffc955d6SMark F. Adams IS is; 821ffc955d6SMark F. Adams PetscInt my0,kk; 822ffc955d6SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr); 823ffc955d6SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 8240298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr); 825a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 826806fa848SBarry Smith } else { 827a94c3b12SMark F. Adams PetscInt kk; 8280298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr); 829a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 830a94c3b12SMark F. Adams ierr = ISDestroy(&is[kk]);CHKERRQ(ierr); 831a94c3b12SMark F. Adams } 832ffc955d6SMark F. Adams ierr = PetscFree(is);CHKERRQ(ierr); 833ffc955d6SMark F. Adams } 834ffc955d6SMark F. Adams ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr); 835ffc955d6SMark F. Adams 8360298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 837ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 838ffc955d6SMark F. Adams ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr); 839806fa848SBarry Smith } else { 8409d5b6da9SMark F. Adams ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr); 841ffc955d6SMark F. Adams } 842d5280255SMark F. Adams } 843d5280255SMark F. Adams { 844d5280255SMark F. Adams /* coarse grid */ 845d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 846d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 847d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 848d5280255SMark F. Adams ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 849d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 850d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 851d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 852d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 85371959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 85471959b99SBarry Smith if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 855d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 856d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 8579dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 8582fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 859*5b42dca8SJed Brown /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 860*5b42dca8SJed Brown * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 861*5b42dca8SJed Brown * view every subdomain as though they were different. */ 862*5b42dca8SJed Brown ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 863d5280255SMark F. Adams } 864d5280255SMark F. Adams 865d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 866d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 867d5280255SMark F. Adams ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); 868d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 8693b4367a7SBarry Smith if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 870d5280255SMark F. Adams 871d5280255SMark F. Adams /* create cheby smoothers */ 872d5280255SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 873d5280255SMark F. Adams lidx <= fine_level; 874d5280255SMark F. Adams lidx++, level--) { 875d5280255SMark F. Adams KSP smoother; 876ffc955d6SMark F. Adams PetscBool flag; 877d5280255SMark F. Adams PC subpc; 878d5280255SMark F. Adams 879ffc955d6SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 880a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 881a2f3521dSMark F. Adams 882a2f3521dSMark F. Adams /* create field split PC, get subsmoother */ 883a2f3521dSMark F. Adams if (stokes) { 884a2f3521dSMark F. Adams KSP *ksps; 885a2f3521dSMark F. Adams PetscInt nn; 886a2f3521dSMark F. Adams ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr); 887a2f3521dSMark F. Adams smoother = ksps[0]; 888a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 889a2f3521dSMark F. Adams ierr = PetscFree(ksps);CHKERRQ(ierr); 890a2f3521dSMark F. Adams } 891ffc955d6SMark F. Adams 892ffc955d6SMark F. Adams /* do my own cheby */ 8936c9de887SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr); 894ffc955d6SMark F. Adams if (flag) { 895ffc955d6SMark F. Adams PetscReal emax, emin; 896251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr); 897ffc955d6SMark F. Adams if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 898587fa25dSMark F. Adams else { /* eigen estimate 'emax' */ 899e696ed0bSMark F. Adams KSP eksp; 900e696ed0bSMark F. Adams Mat Lmat = Aarr[level]; 901b2a4f308SMark F. Adams Vec bb, xx; 902038e3b61SMark F. Adams 9035745e0f5SMark F. Adams ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr); 9045745e0f5SMark F. Adams ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr); 905fc4362bfSMark F. Adams { 906fc4362bfSMark F. Adams PetscRandom rctx; 9073b4367a7SBarry Smith ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr); 908fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 909fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 910fc4362bfSMark F. Adams ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 9115b89ad90SMark F. Adams } 912a2f3521dSMark F. Adams 913e696ed0bSMark F. Adams /* zeroing out BC rows -- needed for crazy matrices */ 914e696ed0bSMark F. Adams { 915e696ed0bSMark F. Adams PetscInt Istart,Iend,ncols,jj,Ii; 916e696ed0bSMark F. Adams PetscScalar zero = 0.0; 917e696ed0bSMark F. Adams ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr); 918e696ed0bSMark F. Adams for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 919e696ed0bSMark F. Adams ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 920e696ed0bSMark F. Adams if (ncols <= 1) { 921e696ed0bSMark F. Adams ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 922a94c3b12SMark F. Adams } 923e696ed0bSMark F. Adams ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 924a94c3b12SMark F. Adams } 925a94c3b12SMark F. Adams ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 926a94c3b12SMark F. Adams ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 927a94c3b12SMark F. Adams } 928e696ed0bSMark F. Adams 9293b4367a7SBarry Smith ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr); 930806fa848SBarry Smith ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr); 931fc4362bfSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 9321a166f3bSJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 9331a166f3bSJed Brown ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 934f6536408SMark F. Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 935f6536408SMark F. Adams 936f6536408SMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 937f6536408SMark F. Adams ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 938fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 9395a9b9e01SMark F. Adams 940d3d0db20SJed Brown /* set PC type to be same as smoother */ 941ffc955d6SMark F. Adams ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr); 942b2a4f308SMark F. Adams 9435a9b9e01SMark F. Adams /* solve - keep stuff out of logging */ 9445a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 9455a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 946fc4362bfSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 9475a9b9e01SMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 9485a9b9e01SMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 9495a9b9e01SMark F. Adams 950fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 9515a9b9e01SMark F. Adams 952fc4362bfSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 953fc4362bfSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 954fc4362bfSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 955f6536408SMark F. Adams 956ffc955d6SMark F. Adams if (pc_gamg->verbose > 0) { 957a94c3b12SMark F. Adams PetscInt N1, tt; 958a94c3b12SMark F. Adams ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr); 9593b4367a7SBarry Smith PetscPrintf(comm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 960f6536408SMark F. Adams } 961fc4362bfSMark F. Adams } 962038e3b61SMark F. Adams { 963c5bfad50SMark F. Adams PetscInt N1, N0; 9640298fd71SBarry Smith ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr); 9650298fd71SBarry Smith ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr); 966f6536408SMark F. Adams /* heuristic - is this crap? */ 967b4ec6429SMark F. Adams /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */ 9685e7c91beSJed Brown emin = emax * pc_gamg->eigtarget[0]; 9695e7c91beSJed Brown emax *= pc_gamg->eigtarget[1]; 970038e3b61SMark F. Adams } 9716c9de887SHong Zhang ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); 972ffc955d6SMark F. Adams } /* setup checby flag */ 973ffc955d6SMark F. Adams } /* non-coarse levels */ 974737a81a9SMark F. Adams 975d5280255SMark F. Adams /* clean up */ 976d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 977587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 978587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 9795b89ad90SMark F. Adams } 9800cbbd2e1SMark F. Adams 9810cbbd2e1SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 9820cbbd2e1SMark F. Adams 9831ac9965eSMark F. Adams if (PETSC_TRUE) { 98458471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 9859d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 98658471d46SMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 98758471d46SMark F. Adams } 988806fa848SBarry Smith } else { 9895f8cf99dSMark F. Adams KSP smoother; 9903b4367a7SBarry Smith if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__); 9919d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 9925f8cf99dSMark F. Adams ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 9935f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 9949d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 9955f8cf99dSMark F. Adams } 9965b89ad90SMark F. Adams PetscFunctionReturn(0); 9975b89ad90SMark F. Adams } 9985b89ad90SMark F. Adams 999eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 10005b89ad90SMark F. Adams /* 10015b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 10025b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 10035b89ad90SMark F. Adams 10045b89ad90SMark F. Adams Input Parameter: 10055b89ad90SMark F. Adams . pc - the preconditioner context 10065b89ad90SMark F. Adams 10075b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 10085b89ad90SMark F. Adams */ 10095b89ad90SMark F. Adams #undef __FUNCT__ 10105b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 10115b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 10125b89ad90SMark F. Adams { 10135b89ad90SMark F. Adams PetscErrorCode ierr; 10145b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 10155b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 10165b89ad90SMark F. Adams 10175b89ad90SMark F. Adams PetscFunctionBegin; 10185b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 10199b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 10209b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 10219b8ffb57SJed Brown } 10221ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 10231ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 10245b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 10255b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 10265b89ad90SMark F. Adams PetscFunctionReturn(0); 10275b89ad90SMark F. Adams } 10285b89ad90SMark F. Adams 1029676e1743SMark F. Adams 1030676e1743SMark F. Adams #undef __FUNCT__ 1031676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 1032676e1743SMark F. Adams /*@ 1033676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 1034676e1743SMark F. Adams processor reduction. 1035676e1743SMark F. Adams 1036676e1743SMark F. Adams Not Collective on PC 1037676e1743SMark F. Adams 1038676e1743SMark F. Adams Input Parameters: 1039676e1743SMark F. Adams . pc - the preconditioner context 1040676e1743SMark F. Adams 1041676e1743SMark F. Adams 1042676e1743SMark F. Adams Options Database Key: 1043676e1743SMark F. Adams . -pc_gamg_process_eq_limit 1044676e1743SMark F. Adams 1045676e1743SMark F. Adams Level: intermediate 1046676e1743SMark F. Adams 1047676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1048676e1743SMark F. Adams 1049676e1743SMark F. Adams .seealso: () 1050676e1743SMark F. Adams @*/ 1051676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 1052676e1743SMark F. Adams { 1053676e1743SMark F. Adams PetscErrorCode ierr; 1054676e1743SMark F. Adams 1055676e1743SMark F. Adams PetscFunctionBegin; 1056676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1057676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1058676e1743SMark F. Adams PetscFunctionReturn(0); 1059676e1743SMark F. Adams } 1060676e1743SMark F. Adams 1061676e1743SMark F. Adams #undef __FUNCT__ 1062676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 10631e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 1064676e1743SMark F. Adams { 1065c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1066c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1067676e1743SMark F. Adams 1068676e1743SMark F. Adams PetscFunctionBegin; 10699d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 1070676e1743SMark F. Adams PetscFunctionReturn(0); 1071676e1743SMark F. Adams } 1072676e1743SMark F. Adams 1073676e1743SMark F. Adams #undef __FUNCT__ 1074389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 1075389730f3SMark F. Adams /*@ 1076389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 1077389730f3SMark F. Adams 1078389730f3SMark F. Adams Collective on PC 1079389730f3SMark F. Adams 1080389730f3SMark F. Adams Input Parameters: 1081389730f3SMark F. Adams . pc - the preconditioner context 1082389730f3SMark F. Adams 1083389730f3SMark F. Adams 1084389730f3SMark F. Adams Options Database Key: 1085389730f3SMark F. Adams . -pc_gamg_coarse_eq_limit 1086389730f3SMark F. Adams 1087389730f3SMark F. Adams Level: intermediate 1088389730f3SMark F. Adams 1089389730f3SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1090389730f3SMark F. Adams 1091389730f3SMark F. Adams .seealso: () 1092389730f3SMark F. Adams @*/ 1093389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 1094389730f3SMark F. Adams { 1095389730f3SMark F. Adams PetscErrorCode ierr; 1096389730f3SMark F. Adams 1097389730f3SMark F. Adams PetscFunctionBegin; 1098389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1099389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1100389730f3SMark F. Adams PetscFunctionReturn(0); 1101389730f3SMark F. Adams } 1102389730f3SMark F. Adams 1103389730f3SMark F. Adams #undef __FUNCT__ 1104389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 11051e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1106389730f3SMark F. Adams { 1107389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1108389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1109389730f3SMark F. Adams 1110389730f3SMark F. Adams PetscFunctionBegin; 11119d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 1112389730f3SMark F. Adams PetscFunctionReturn(0); 1113389730f3SMark F. Adams } 1114389730f3SMark F. Adams 1115389730f3SMark F. Adams #undef __FUNCT__ 11168263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 1117676e1743SMark F. Adams /*@ 11188263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 1119676e1743SMark F. Adams 1120676e1743SMark F. Adams Collective on PC 1121676e1743SMark F. Adams 1122676e1743SMark F. Adams Input Parameters: 1123676e1743SMark F. Adams . pc - the preconditioner context 1124676e1743SMark F. Adams 1125676e1743SMark F. Adams 1126676e1743SMark F. Adams Options Database Key: 11278263b398SMark F. Adams . -pc_gamg_repartition 1128676e1743SMark F. Adams 1129676e1743SMark F. Adams Level: intermediate 1130676e1743SMark F. Adams 1131676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1132676e1743SMark F. Adams 1133676e1743SMark F. Adams .seealso: () 1134676e1743SMark F. Adams @*/ 11358263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1136676e1743SMark F. Adams { 1137676e1743SMark F. Adams PetscErrorCode ierr; 1138676e1743SMark F. Adams 1139676e1743SMark F. Adams PetscFunctionBegin; 1140676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11418263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1142676e1743SMark F. Adams PetscFunctionReturn(0); 1143676e1743SMark F. Adams } 1144676e1743SMark F. Adams 1145676e1743SMark F. Adams #undef __FUNCT__ 11468263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 11471e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 1148676e1743SMark F. Adams { 1149c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1150c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1151676e1743SMark F. Adams 1152676e1743SMark F. Adams PetscFunctionBegin; 11539d5b6da9SMark F. Adams pc_gamg->repart = n; 1154676e1743SMark F. Adams PetscFunctionReturn(0); 1155676e1743SMark F. Adams } 1156676e1743SMark F. Adams 1157676e1743SMark F. Adams #undef __FUNCT__ 1158dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl" 1159dfd5c07aSMark F. Adams /*@ 1160dfd5c07aSMark F. Adams PCGAMGSetReuseProl - Reuse prlongation 1161dfd5c07aSMark F. Adams 1162dfd5c07aSMark F. Adams Collective on PC 1163dfd5c07aSMark F. Adams 1164dfd5c07aSMark F. Adams Input Parameters: 1165dfd5c07aSMark F. Adams . pc - the preconditioner context 1166dfd5c07aSMark F. Adams 1167dfd5c07aSMark F. Adams 1168dfd5c07aSMark F. Adams Options Database Key: 1169dfd5c07aSMark F. Adams . -pc_gamg_reuse_interpolation 1170dfd5c07aSMark F. Adams 1171dfd5c07aSMark F. Adams Level: intermediate 1172dfd5c07aSMark F. Adams 1173dfd5c07aSMark F. Adams Concepts: Unstructured multrigrid preconditioner 1174dfd5c07aSMark F. Adams 1175dfd5c07aSMark F. Adams .seealso: () 1176dfd5c07aSMark F. Adams @*/ 1177dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n) 1178dfd5c07aSMark F. Adams { 1179dfd5c07aSMark F. Adams PetscErrorCode ierr; 1180dfd5c07aSMark F. Adams 1181dfd5c07aSMark F. Adams PetscFunctionBegin; 1182dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1183dfd5c07aSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1184dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1185dfd5c07aSMark F. Adams } 1186dfd5c07aSMark F. Adams 1187dfd5c07aSMark F. Adams #undef __FUNCT__ 1188dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG" 11891e6b0712SBarry Smith static PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n) 1190dfd5c07aSMark F. Adams { 1191dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1192dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1193dfd5c07aSMark F. Adams 1194dfd5c07aSMark F. Adams PetscFunctionBegin; 1195dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1196dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1197dfd5c07aSMark F. Adams } 1198dfd5c07aSMark F. Adams 1199dfd5c07aSMark F. Adams #undef __FUNCT__ 1200ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs" 1201ffc955d6SMark F. Adams /*@ 1202ffc955d6SMark F. Adams PCGAMGSetUseASMAggs - 1203ffc955d6SMark F. Adams 1204ffc955d6SMark F. Adams Collective on PC 1205ffc955d6SMark F. Adams 1206ffc955d6SMark F. Adams Input Parameters: 1207ffc955d6SMark F. Adams . pc - the preconditioner context 1208ffc955d6SMark F. Adams 1209ffc955d6SMark F. Adams 1210ffc955d6SMark F. Adams Options Database Key: 1211ffc955d6SMark F. Adams . -pc_gamg_use_agg_gasm 1212ffc955d6SMark F. Adams 1213ffc955d6SMark F. Adams Level: intermediate 1214ffc955d6SMark F. Adams 1215ffc955d6SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1216ffc955d6SMark F. Adams 1217ffc955d6SMark F. Adams .seealso: () 1218ffc955d6SMark F. Adams @*/ 1219ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1220ffc955d6SMark F. Adams { 1221ffc955d6SMark F. Adams PetscErrorCode ierr; 1222ffc955d6SMark F. Adams 1223ffc955d6SMark F. Adams PetscFunctionBegin; 1224ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1225ffc955d6SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1226ffc955d6SMark F. Adams PetscFunctionReturn(0); 1227ffc955d6SMark F. Adams } 1228ffc955d6SMark F. Adams 1229ffc955d6SMark F. Adams #undef __FUNCT__ 1230ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 12311e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1232ffc955d6SMark F. Adams { 1233ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1234ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1235ffc955d6SMark F. Adams 1236ffc955d6SMark F. Adams PetscFunctionBegin; 1237ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = n; 1238ffc955d6SMark F. Adams PetscFunctionReturn(0); 1239ffc955d6SMark F. Adams } 1240ffc955d6SMark F. Adams 1241ffc955d6SMark F. Adams #undef __FUNCT__ 12424ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 12434ef23d27SMark F. Adams /*@ 12444ef23d27SMark F. Adams PCGAMGSetNlevels - 12454ef23d27SMark F. Adams 12464ef23d27SMark F. Adams Not collective on PC 12474ef23d27SMark F. Adams 12484ef23d27SMark F. Adams Input Parameters: 12494ef23d27SMark F. Adams . pc - the preconditioner context 12504ef23d27SMark F. Adams 12514ef23d27SMark F. Adams 12524ef23d27SMark F. Adams Options Database Key: 12534ef23d27SMark F. Adams . -pc_mg_levels 12544ef23d27SMark F. Adams 12554ef23d27SMark F. Adams Level: intermediate 12564ef23d27SMark F. Adams 12574ef23d27SMark F. Adams Concepts: Unstructured multrigrid preconditioner 12584ef23d27SMark F. Adams 12594ef23d27SMark F. Adams .seealso: () 12604ef23d27SMark F. Adams @*/ 12614ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 12624ef23d27SMark F. Adams { 12634ef23d27SMark F. Adams PetscErrorCode ierr; 12644ef23d27SMark F. Adams 12654ef23d27SMark F. Adams PetscFunctionBegin; 12664ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12674ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 12684ef23d27SMark F. Adams PetscFunctionReturn(0); 12694ef23d27SMark F. Adams } 12704ef23d27SMark F. Adams 12714ef23d27SMark F. Adams #undef __FUNCT__ 12724ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 12731e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 12744ef23d27SMark F. Adams { 12754ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 12764ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12774ef23d27SMark F. Adams 12784ef23d27SMark F. Adams PetscFunctionBegin; 12799d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12804ef23d27SMark F. Adams PetscFunctionReturn(0); 12814ef23d27SMark F. Adams } 12824ef23d27SMark F. Adams 12834ef23d27SMark F. Adams #undef __FUNCT__ 12843542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 12853542efc5SMark F. Adams /*@ 12863542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12873542efc5SMark F. Adams 12883542efc5SMark F. Adams Not collective on PC 12893542efc5SMark F. Adams 12903542efc5SMark F. Adams Input Parameters: 12913542efc5SMark F. Adams . pc - the preconditioner context 12923542efc5SMark F. Adams 12933542efc5SMark F. Adams 12943542efc5SMark F. Adams Options Database Key: 12953542efc5SMark F. Adams . -pc_gamg_threshold 12963542efc5SMark F. Adams 12973542efc5SMark F. Adams Level: intermediate 12983542efc5SMark F. Adams 12993542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 13003542efc5SMark F. Adams 13013542efc5SMark F. Adams .seealso: () 13023542efc5SMark F. Adams @*/ 13033542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 13043542efc5SMark F. Adams { 13053542efc5SMark F. Adams PetscErrorCode ierr; 13063542efc5SMark F. Adams 13073542efc5SMark F. Adams PetscFunctionBegin; 13083542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13093542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 13103542efc5SMark F. Adams PetscFunctionReturn(0); 13113542efc5SMark F. Adams } 13123542efc5SMark F. Adams 13133542efc5SMark F. Adams #undef __FUNCT__ 13143542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 13151e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 13163542efc5SMark F. Adams { 1317c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1318c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 13193542efc5SMark F. Adams 13203542efc5SMark F. Adams PetscFunctionBegin; 13219d5b6da9SMark F. Adams pc_gamg->threshold = n; 13223542efc5SMark F. Adams PetscFunctionReturn(0); 13233542efc5SMark F. Adams } 13243542efc5SMark F. Adams 13253542efc5SMark F. Adams #undef __FUNCT__ 13269d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1327676e1743SMark F. Adams /*@ 13289d5b6da9SMark F. Adams PCGAMGSetType - Set solution method - calls sub create method 1329676e1743SMark F. Adams 1330676e1743SMark F. Adams Collective on PC 1331676e1743SMark F. Adams 1332676e1743SMark F. Adams Input Parameters: 1333676e1743SMark F. Adams . pc - the preconditioner context 1334676e1743SMark F. Adams 1335676e1743SMark F. Adams 1336676e1743SMark F. Adams Options Database Key: 13373542efc5SMark F. Adams . -pc_gamg_type 1338676e1743SMark F. Adams 1339676e1743SMark F. Adams Level: intermediate 1340676e1743SMark F. Adams 1341676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1342676e1743SMark F. Adams 1343676e1743SMark F. Adams .seealso: () 1344676e1743SMark F. Adams @*/ 134519fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1346676e1743SMark F. Adams { 1347676e1743SMark F. Adams PetscErrorCode ierr; 1348676e1743SMark F. Adams 1349676e1743SMark F. Adams PetscFunctionBegin; 1350676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1351806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1352676e1743SMark F. Adams PetscFunctionReturn(0); 1353676e1743SMark F. Adams } 1354676e1743SMark F. Adams 1355676e1743SMark F. Adams #undef __FUNCT__ 13569d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 13571e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1358676e1743SMark F. Adams { 13599d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 13601ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 13611ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1362676e1743SMark F. Adams 1363676e1743SMark F. Adams PetscFunctionBegin; 13641c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 13659d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 13661ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 13671ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 13681ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 13691ab5ffc9SJed Brown } 13701ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 13711ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 13729d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1373676e1743SMark F. Adams PetscFunctionReturn(0); 1374676e1743SMark F. Adams } 1375676e1743SMark F. Adams 13765b89ad90SMark F. Adams #undef __FUNCT__ 13775b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 13785b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 13795b89ad90SMark F. Adams { 1380676e1743SMark F. Adams PetscErrorCode ierr; 1381676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1382676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1383676e1743SMark F. Adams PetscBool flag; 13845e7c91beSJed Brown PetscInt two = 2; 13853b4367a7SBarry Smith MPI_Comm comm; 13865b89ad90SMark F. Adams 13875b89ad90SMark F. Adams PetscFunctionBegin; 13883b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1389676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr); 1390676e1743SMark F. Adams { 1391b7cbab4eSMark Adams /* -pc_gamg_type */ 1392b7cbab4eSMark Adams { 1393b7cbab4eSMark Adams char tname[256] = PCGAMGAGG; 13941ab5ffc9SJed Brown const char *deftype = pc_gamg->gamg_type_name ? pc_gamg->gamg_type_name : tname; 13951ab5ffc9SJed Brown ierr = PetscOptionsList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1396b7cbab4eSMark Adams /* call PCCreateGAMG_XYZ */ 13971ab5ffc9SJed Brown if (flag || !pc_gamg->gamg_type_name) { 13981ab5ffc9SJed Brown ierr = PCGAMGSetType(pc, flag ? tname : deftype);CHKERRQ(ierr); 13991ab5ffc9SJed Brown } 1400b7cbab4eSMark Adams } 140175b74bdaSMark F. Adams /* -pc_gamg_verbose */ 14029d5b6da9SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 14039d5b6da9SMark F. Adams "none", pc_gamg->verbose, 14040298fd71SBarry Smith &pc_gamg->verbose, NULL);CHKERRQ(ierr); 14058263b398SMark F. Adams /* -pc_gamg_repartition */ 14068263b398SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_repartition", 14078263b398SMark F. Adams "Repartion coarse grids (false)", 14088263b398SMark F. Adams "PCGAMGRepartitioning", 14099d5b6da9SMark F. Adams pc_gamg->repart, 14109d5b6da9SMark F. Adams &pc_gamg->repart, 1411806fa848SBarry Smith &flag);CHKERRQ(ierr); 1412dfd5c07aSMark F. Adams /* -pc_gamg_reuse_interpolation */ 1413dfd5c07aSMark F. Adams ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation", 1414dfd5c07aSMark F. Adams "Reuse prolongation operator (true)", 1415dfd5c07aSMark F. Adams "PCGAMGReuseProl", 1416dfd5c07aSMark F. Adams pc_gamg->reuse_prol, 1417dfd5c07aSMark F. Adams &pc_gamg->reuse_prol, 1418dfd5c07aSMark F. Adams &flag);CHKERRQ(ierr); 1419ffc955d6SMark F. Adams /* -pc_gamg_use_agg_gasm */ 1420ffc955d6SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm", 1421ffc955d6SMark F. Adams "Use aggregation agragates for GASM smoother (false)", 1422ffc955d6SMark F. Adams "PCGAMGUseASMAggs", 1423ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm, 1424ffc955d6SMark F. Adams &pc_gamg->use_aggs_in_gasm, 1425806fa848SBarry Smith &flag);CHKERRQ(ierr); 1426c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 1427676e1743SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1428676e1743SMark F. Adams "Limit (goal) on number of equations per process on coarse grids", 1429676e1743SMark F. Adams "PCGAMGSetProcEqLim", 14309d5b6da9SMark F. Adams pc_gamg->min_eq_proc, 14319d5b6da9SMark F. Adams &pc_gamg->min_eq_proc, 1432806fa848SBarry Smith &flag);CHKERRQ(ierr); 1433389730f3SMark F. Adams /* -pc_gamg_coarse_eq_limit */ 1434389730f3SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1435389730f3SMark F. Adams "Limit on number of equations for the coarse grid", 1436389730f3SMark F. Adams "PCGAMGSetCoarseEqLim", 14379d5b6da9SMark F. Adams pc_gamg->coarse_eq_limit, 14389d5b6da9SMark F. Adams &pc_gamg->coarse_eq_limit, 1439806fa848SBarry Smith &flag);CHKERRQ(ierr); 1440c20e4228SMark F. Adams /* -pc_gamg_threshold */ 14413542efc5SMark F. Adams ierr = PetscOptionsReal("-pc_gamg_threshold", 14423542efc5SMark F. Adams "Relative threshold to use for dropping edges in aggregation graph", 14433542efc5SMark F. Adams "PCGAMGSetThreshold", 14449d5b6da9SMark F. Adams pc_gamg->threshold, 14459d5b6da9SMark F. Adams &pc_gamg->threshold, 1446806fa848SBarry Smith &flag);CHKERRQ(ierr); 1447806fa848SBarry Smith if (flag && pc_gamg->verbose) { 14483b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr); 1449806fa848SBarry Smith } 1450b7cbab4eSMark Adams /* -pc_gamg_eigtarget */ 14510298fd71SBarry Smith ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr); 14524ef23d27SMark F. Adams ierr = PetscOptionsInt("-pc_mg_levels", 14534ef23d27SMark F. Adams "Set number of MG levels", 14544ef23d27SMark F. Adams "PCGAMGSetNlevels", 14559d5b6da9SMark F. Adams pc_gamg->Nlevels, 14569d5b6da9SMark F. Adams &pc_gamg->Nlevels, 1457806fa848SBarry Smith &flag);CHKERRQ(ierr); 1458b7cbab4eSMark Adams 1459b7cbab4eSMark Adams /* set options for subtype */ 14601ab5ffc9SJed Brown if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(pc);CHKERRQ(ierr);} 1461676e1743SMark F. Adams } 1462676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 14635b89ad90SMark F. Adams PetscFunctionReturn(0); 14645b89ad90SMark F. Adams } 14655b89ad90SMark F. Adams 14665b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 14675b89ad90SMark F. Adams /*MC 1468856830bfSJed Brown PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 14699d5b6da9SMark F. Adams - This is the entry point to GAMG, registered in pcregis.c 14705b89ad90SMark F. Adams 1471280d9858SJed Brown Options Database Keys: 14725b89ad90SMark F. Adams Multigrid options(inherited) 1473280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1474280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1475280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 14768c1c2452SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 14775b89ad90SMark F. Adams 14785b89ad90SMark F. Adams Level: intermediate 1479280d9858SJed Brown 14805b89ad90SMark F. Adams Concepts: multigrid 14815b89ad90SMark F. Adams 14825b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1483280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 14845b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 14855b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 14865b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 14875b89ad90SMark F. Adams M*/ 1488b2573a8aSBarry Smith 14895b89ad90SMark F. Adams #undef __FUNCT__ 14905b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 14918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 14925b89ad90SMark F. Adams { 14935b89ad90SMark F. Adams PetscErrorCode ierr; 14945b89ad90SMark F. Adams PC_GAMG *pc_gamg; 14955b89ad90SMark F. Adams PC_MG *mg; 14960cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 14979d5b6da9SMark F. Adams static long count = 0; 14985ee9c036SSatish Balay #endif 14995b89ad90SMark F. Adams 15005b89ad90SMark F. Adams PetscFunctionBegin; 15015b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 15025b89ad90SMark F. Adams ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 15035b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 15045b89ad90SMark F. Adams 15055b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 15065b89ad90SMark F. Adams ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr); 15075b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 1508ce4cda84SJed Brown mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 15095b89ad90SMark F. Adams mg->innerctx = pc_gamg; 15105b89ad90SMark F. Adams 15111ab5ffc9SJed Brown ierr = PetscNewLog(pc,struct _PCGAMGOps,&pc_gamg->ops);CHKERRQ(ierr); 15121ab5ffc9SJed Brown 15139d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 15149d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 15159d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 15169d5b6da9SMark F. Adams pc_gamg->data = 0; 15175b89ad90SMark F. Adams 15189d5b6da9SMark F. Adams /* register AMG type */ 15193e3471ccSMark Adams #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 15203e3471ccSMark Adams ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 15213e3471ccSMark Adams #endif 15229d5b6da9SMark F. Adams 15239d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 15245b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 15255b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 15265b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 15275b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 15285b89ad90SMark F. Adams 1529bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1530bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1531bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 1532bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseProl_C",PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr); 1533bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1534bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1535bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1536bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 15379d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1538d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 1539ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1540038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 1541c8b0795cSMark F. Adams pc_gamg->coarse_eq_limit = 800; 1542d3042614SMark Adams pc_gamg->threshold = 0.; 15439d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 15449d5b6da9SMark F. Adams pc_gamg->verbose = 0; 15459d5b6da9SMark F. Adams pc_gamg->emax_id = -1; 15465e7c91beSJed Brown pc_gamg->eigtarget[0] = 0.05; 15475e7c91beSJed Brown pc_gamg->eigtarget[1] = 1.05; 15489d5b6da9SMark F. Adams 15490cbbd2e1SMark F. Adams /* private events */ 15500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1551785cba28SMark F. Adams if (count++ == 0) { 1552806fa848SBarry Smith ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 1553806fa848SBarry Smith ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 15540cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 15550cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 15560cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 1557806fa848SBarry Smith ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 1558806fa848SBarry Smith ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 1559806fa848SBarry Smith ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1560806fa848SBarry Smith ierr = PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 1561806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 1562806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 1563806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 1564806fa848SBarry Smith ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 1565806fa848SBarry Smith ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 1566806fa848SBarry Smith ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 1567806fa848SBarry Smith ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 1568806fa848SBarry Smith ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 1569f852f58cSMark F. Adams 15700cbbd2e1SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 15710cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 15720cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1573b4fbaa2aSMark F. Adams /* create timer stages */ 1574b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1575b4fbaa2aSMark F. Adams { 1576b4fbaa2aSMark F. Adams char str[32]; 1577b4fbaa2aSMark F. Adams PetscInt lidx; 1578806fa848SBarry Smith sprintf(str,"MG Level %d (finest)",0); 1579806fa848SBarry Smith ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 1580b4fbaa2aSMark F. Adams for (lidx=1; lidx<9; lidx++) { 1581b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1582806fa848SBarry Smith ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 1583b4fbaa2aSMark F. Adams } 1584b4fbaa2aSMark F. Adams } 1585b4fbaa2aSMark F. Adams #endif 1586b4fbaa2aSMark F. Adams } 1587b4fbaa2aSMark F. Adams #endif 15880cbbd2e1SMark F. Adams /* general events */ 15890cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 1590806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr); 1591806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr); 1592806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1593806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1594806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1595806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1596806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr); 1597806fa848SBarry Smith ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr); 15980cbbd2e1SMark F. Adams #endif 15990cbbd2e1SMark F. Adams 16005b89ad90SMark F. Adams PetscFunctionReturn(0); 16015b89ad90SMark F. Adams } 16023e3471ccSMark Adams 16033e3471ccSMark Adams #undef __FUNCT__ 16043e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage" 16053e3471ccSMark Adams /*@C 16063e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 16073e3471ccSMark Adams from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG() 16083e3471ccSMark Adams when using static libraries. 16093e3471ccSMark Adams 16103e3471ccSMark Adams Level: developer 16113e3471ccSMark Adams 16123e3471ccSMark Adams .keywords: PC, PCGAMG, initialize, package 16133e3471ccSMark Adams .seealso: PetscInitialize() 16143e3471ccSMark Adams @*/ 16153e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 16163e3471ccSMark Adams { 16173e3471ccSMark Adams PetscErrorCode ierr; 16183e3471ccSMark Adams 16193e3471ccSMark Adams PetscFunctionBegin; 16203e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 16213e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 16223e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 16233e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 16243e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 16253e3471ccSMark Adams PetscFunctionReturn(0); 16263e3471ccSMark Adams } 16273e3471ccSMark Adams 16283e3471ccSMark Adams #undef __FUNCT__ 16293e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage" 16303e3471ccSMark Adams /*@C 16313e3471ccSMark Adams PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is 16323e3471ccSMark Adams called from PetscFinalize(). 16333e3471ccSMark Adams 16343e3471ccSMark Adams Level: developer 16353e3471ccSMark Adams 16363e3471ccSMark Adams .keywords: Petsc, destroy, package 16373e3471ccSMark Adams .seealso: PetscFinalize() 16383e3471ccSMark Adams @*/ 16393e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 16403e3471ccSMark Adams { 16413e3471ccSMark Adams PetscErrorCode ierr; 16423e3471ccSMark Adams 16433e3471ccSMark Adams PetscFunctionBegin; 16443e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 16453e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 16463e3471ccSMark Adams PetscFunctionReturn(0); 16473e3471ccSMark Adams } 1648