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> 7f96513f1SMatthew G Knepley 80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 10b4fbaa2aSMark F. Adams #endif 110cbbd2e1SMark F. Adams 120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 130cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG; 140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO; 150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG; 20a2f3521dSMark F. Adams PetscLogEvent PC_GAMGKKTProl_AGG; 210cbbd2e1SMark F. Adams #endif 220cbbd2e1SMark F. Adams 23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 24b4fbaa2aSMark F. Adams 25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 28b4fbaa2aSMark F. Adams #endif 29f96513f1SMatthew G Knepley 309d5b6da9SMark F. Adams static PetscFList GAMGList = 0; 319d5b6da9SMark F. Adams 32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 33d3d6bff4SMark F. Adams #undef __FUNCT__ 34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 36d3d6bff4SMark F. Adams { 37d3d6bff4SMark F. Adams PetscErrorCode ierr; 38d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 39d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 40d3d6bff4SMark F. Adams 41d3d6bff4SMark F. Adams PetscFunctionBegin; 42a2f3521dSMark F. Adams if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */ 43878e152fSMark F. Adams PetscPrintf(((PetscObject)pc)->comm,"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__); 449d5b6da9SMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 4558471d46SMark F. Adams } 46a2f3521dSMark F. Adams pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0; 47878e152fSMark F. Adams 48878e152fSMark F. Adams if (pc_gamg->orig_data) { 49878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 50878e152fSMark F. Adams } 51878e152fSMark F. Adams 52a2f3521dSMark F. Adams PetscFunctionReturn(0); 53a2f3521dSMark F. Adams } 54a2f3521dSMark F. Adams 55a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */ 56a2f3521dSMark F. Adams typedef struct{ 57a2f3521dSMark F. Adams Mat A11,A21,A12,Amat; 58a2f3521dSMark F. Adams IS prim_is,constr_is; 59a2f3521dSMark F. Adams }GAMGKKTMat; 60a2f3521dSMark F. Adams 61a2f3521dSMark F. Adams #undef __FUNCT__ 62a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate" 63a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out) 64a2f3521dSMark F. Adams { 65a2f3521dSMark F. Adams PetscFunctionBegin; 66a2f3521dSMark F. Adams out->Amat = A; 67a2f3521dSMark F. Adams if (iskkt){ 68a2f3521dSMark F. Adams PetscErrorCode ierr; 69a2f3521dSMark F. Adams IS is_constraint, is_prime; 70a2f3521dSMark F. Adams PetscInt nmin,nmax; 71a2f3521dSMark F. Adams 72a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr); 73a2f3521dSMark F. Adams ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr); 74a2f3521dSMark F. Adams ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr); 75a2f3521dSMark F. Adams out->prim_is = is_prime; 76a2f3521dSMark F. Adams out->constr_is = is_constraint; 77a2f3521dSMark F. Adams 78a2f3521dSMark F. Adams ierr = MatGetSubMatrix(A, is_prime, is_prime, MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr); 79a2f3521dSMark F. Adams ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr); 80a2f3521dSMark F. Adams ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr); 81806fa848SBarry Smith } else { 82a2f3521dSMark F. Adams out->A11 = A; 83a2f3521dSMark F. Adams out->A21 = PETSC_NULL; 84a2f3521dSMark F. Adams out->A12 = PETSC_NULL; 85a2f3521dSMark F. Adams out->prim_is = PETSC_NULL; 86a2f3521dSMark F. Adams out->constr_is = PETSC_NULL; 87a2f3521dSMark F. Adams } 88a2f3521dSMark F. Adams PetscFunctionReturn(0); 89a2f3521dSMark F. Adams } 90a2f3521dSMark F. Adams 91a2f3521dSMark F. Adams #undef __FUNCT__ 92a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy" 93a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat) 94a2f3521dSMark F. Adams { 95a2f3521dSMark F. Adams PetscErrorCode ierr; 96a2f3521dSMark F. Adams 97a2f3521dSMark F. Adams PetscFunctionBegin; 98a2f3521dSMark F. Adams if (mat->A11 && mat->A11 != mat->Amat) { 99a2f3521dSMark F. Adams ierr = MatDestroy(&mat->A11);CHKERRQ(ierr); 100a2f3521dSMark F. Adams } 101a2f3521dSMark F. Adams ierr = MatDestroy(&mat->A21);CHKERRQ(ierr); 102a2f3521dSMark F. Adams ierr = MatDestroy(&mat->A12);CHKERRQ(ierr); 103a2f3521dSMark F. Adams 104a2f3521dSMark F. Adams ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr); 105a2f3521dSMark F. Adams ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr); 106a2f3521dSMark F. Adams 107d3d6bff4SMark F. Adams PetscFunctionReturn(0); 108d3d6bff4SMark F. Adams } 109d3d6bff4SMark F. Adams 1105b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 1115b89ad90SMark F. Adams /* 112a147abb0SMark F. Adams createLevel: create coarse op with RAP. repartition and/or reduce number 113a147abb0SMark F. Adams of active processors. 1145b89ad90SMark F. Adams 1155b89ad90SMark F. Adams Input Parameter: 116a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 117a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 1189d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 119c5bfad50SMark F. Adams . cr_bs - coarse block size 1209d5b6da9SMark F. Adams . isLast - 121a2f3521dSMark F. Adams . stokes - 1223530afc2SMark F. Adams In/Output Parameter: 123a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 124afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 12511e60469SMark F. Adams Output Parameter: 1263530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 1275b89ad90SMark F. Adams */ 1285cb416c2SMark F. Adams 1295b89ad90SMark F. Adams #undef __FUNCT__ 1308263b398SMark F. Adams #define __FUNCT__ "createLevel" 1310cbbd2e1SMark F. Adams static PetscErrorCode createLevel(const PC pc, 1329d5b6da9SMark F. Adams const Mat Amat_fine, 133c5bfad50SMark F. Adams const PetscInt cr_bs, 1349d5b6da9SMark F. Adams const PetscBool isLast, 135a2f3521dSMark F. Adams const PetscBool stokes, 1363530afc2SMark F. Adams Mat *a_P_inout, 137a2f3521dSMark F. Adams Mat *a_Amat_crs, 138a2f3521dSMark F. Adams PetscMPIInt *a_nactive_proc 13911e60469SMark F. Adams ) 1405b89ad90SMark F. Adams { 141a2f3521dSMark F. Adams PetscErrorCode ierr; 1429d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 143486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1449d5b6da9SMark F. Adams const PetscBool repart = pc_gamg->repart; 1459d5b6da9SMark F. Adams const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 146a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 1479d5b6da9SMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat_fine)->comm; 1485a9b9e01SMark F. Adams PetscMPIInt mype,npe,new_npe,nactive=*a_nactive_proc; 149c5bfad50SMark F. Adams PetscInt ncrs_eq,ncrs_prim,f_bs; 1505b89ad90SMark F. Adams 1515b89ad90SMark F. Adams PetscFunctionBegin; 15211e60469SMark F. Adams ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr); 15311e60469SMark F. Adams ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr); 154c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 15511e60469SMark F. Adams /* RAP */ 1569d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 157038e3b61SMark F. Adams 158a2f3521dSMark F. Adams /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/ 159a2f3521dSMark F. Adams ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 160a2f3521dSMark F. Adams assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0); 161a2f3521dSMark F. Adams ierr = MatGetLocalSize(Cmat, &ncrs_eq, PETSC_NULL);CHKERRQ(ierr); 162a2f3521dSMark F. Adams 163a2f3521dSMark F. Adams /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */ 164a2f3521dSMark F. Adams { 165472110cdSMark F. Adams PetscInt ncrs_eq_glob; 166a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &ncrs_eq_glob, PETSC_NULL);CHKERRQ(ierr); 167472110cdSMark F. Adams new_npe = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 168472110cdSMark F. Adams if (new_npe == 0 || ncrs_eq_glob < coarse_max) new_npe = 1; 1695a9b9e01SMark F. Adams else if (new_npe >= nactive) new_npe = nactive; /* no change, rare */ 1709d5b6da9SMark F. Adams if (isLast) new_npe = 1; 171a2f3521dSMark F. Adams } 172f852f58cSMark F. Adams 1735a9b9e01SMark F. Adams if (!repart && new_npe==nactive) { 174a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 175806fa848SBarry Smith } else { 176a2f3521dSMark 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; 177a2f3521dSMark F. Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old; 178a2f3521dSMark F. Adams IS is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices; 1795a9b9e01SMark F. Adams VecScatter vecscat; 18022063be5SMark F. Adams PetscScalar *array; 18122063be5SMark F. Adams Vec src_crd, dest_crd; 182e33ef3b1SMark F. Adams 183c5bfad50SMark F. Adams nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0); 1840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1850cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 186b4fbaa2aSMark F. Adams #endif 187a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 188a2f3521dSMark F. Adams ierr = PetscMalloc(npe*sizeof(PetscInt), &counts);CHKERRQ(ierr); 189a2f3521dSMark F. Adams if (repart && !stokes) { 190a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1915a9b9e01SMark F. Adams Mat adj; 1925a9b9e01SMark F. Adams 193a2f3521dSMark F. Adams if (pc_gamg->verbose>0) { 194a2f3521dSMark F. Adams if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq); 195a2f3521dSMark F. Adams else { 196a2f3521dSMark F. Adams PetscInt n; 197a2f3521dSMark F. Adams ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm);CHKERRQ(ierr); 198a2f3521dSMark F. Adams PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n); 199a2f3521dSMark F. Adams } 200a2f3521dSMark F. Adams } 2015a9b9e01SMark F. Adams 202a2f3521dSMark F. Adams /* get 'adj' */ 203c5bfad50SMark F. Adams if (cr_bs == 1) { 204038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 205806fa848SBarry Smith } else{ 206a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 207eb07cef2SMark F. Adams Mat tMat; 208a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 209b4fbaa2aSMark F. Adams const PetscScalar *vals; 210b4fbaa2aSMark F. Adams const PetscInt *idx; 211a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 2129057884aSMark F. Adams static PetscInt llev = 0; 213b4fbaa2aSMark F. Adams 214a2f3521dSMark F. Adams ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr); 215a2f3521dSMark F. Adams ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr); 216a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 217a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 218c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cr_bs, jj++) { 21958471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 220c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 221c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 22258471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 223a2f3521dSMark F. Adams if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim; 224c5bfad50SMark F. Adams if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim; 22558471d46SMark F. Adams } 2266876a03eSMark F. Adams 227a2f3521dSMark F. Adams ierr = MatCreate(wcomm, &tMat);CHKERRQ(ierr); 228806fa848SBarry Smith ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 229a2f3521dSMark F. Adams ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr); 230a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 231a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 23258471d46SMark F. Adams ierr = PetscFree(d_nnz);CHKERRQ(ierr); 2335f8cf99dSMark F. Adams ierr = PetscFree(o_nnz);CHKERRQ(ierr); 234eb07cef2SMark F. Adams 235a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 236c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 23722063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 238eb07cef2SMark F. Adams for (jj = 0 ; jj < ncols ; jj++){ 239c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 240eb07cef2SMark F. Adams PetscScalar v = 1.0; 241eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 242eb07cef2SMark F. Adams } 24322063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 244eb07cef2SMark F. Adams } 245eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 246eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 247eb07cef2SMark F. Adams 248b4fbaa2aSMark F. Adams if (llev++ == -1) { 249b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2508caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 251b4fbaa2aSMark F. Adams PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 252b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 2533bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 254b4fbaa2aSMark F. Adams } 255b4fbaa2aSMark F. Adams 256eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 257eb07cef2SMark F. Adams 258eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 259a2f3521dSMark F. Adams } /* create 'adj' */ 260f150b916SMark F. Adams 261a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2625a9b9e01SMark F. Adams char prefix[256]; 2635a9b9e01SMark F. Adams const char *pcpre; 264b4fbaa2aSMark F. Adams const PetscInt *is_idx; 265b4fbaa2aSMark F. Adams MatPartitioning mpart; 266a4b7d37bSMark F. Adams IS proc_is; 267a2f3521dSMark F. Adams PetscInt targetPE; 2682f03bc48SMark F. Adams 2695a9b9e01SMark F. Adams ierr = MatPartitioningCreate(wcomm, &mpart);CHKERRQ(ierr); 2705ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 2719d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 2728caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr); 27359a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 27411e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 2755ef31b24SMark F. Adams ierr = MatPartitioningSetNParts(mpart, new_npe);CHKERRQ(ierr); 276a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 27711e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 2785a9b9e01SMark F. Adams 2795ef31b24SMark F. Adams /* collect IS info */ 280a2f3521dSMark F. Adams ierr = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr); 281a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 282a2f3521dSMark F. Adams targetPE = 1; /* bring to "front" of machine */ 283a2f3521dSMark F. Adams /*targetPE = npe/new_npe;*/ /* spread partitioning across machine */ 284a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++){ 285c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++){ 286a2f3521dSMark F. Adams newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 287eb07cef2SMark F. Adams } 2885ef31b24SMark F. Adams } 289a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 290a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2915ef31b24SMark F. Adams } 2925ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2935a9b9e01SMark F. Adams 294806fa848SBarry Smith ierr = ISCreateGeneral(wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2958263b398SMark F. Adams if (newproc_idx != 0) { 2968263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 2975ef31b24SMark F. Adams } 298806fa848SBarry Smith } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 299a2f3521dSMark F. Adams 300a2f3521dSMark F. Adams PetscInt rfactor,targetPE; 3015a9b9e01SMark F. Adams /* find factor */ 3025a9b9e01SMark F. Adams if (new_npe == 1) rfactor = npe; /* easy */ 3035a9b9e01SMark F. Adams else { 3045a9b9e01SMark F. Adams PetscReal best_fact = 0.; 3055a9b9e01SMark F. Adams jj = -1; 3065a9b9e01SMark F. Adams for (kk = 1 ; kk <= npe ; kk++){ 3075a9b9e01SMark F. Adams if (npe%kk==0) { /* a candidate */ 3085a9b9e01SMark F. Adams PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe; 3095a9b9e01SMark F. Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 3105a9b9e01SMark F. Adams if (fact > best_fact) { 3115a9b9e01SMark F. Adams best_fact = fact; jj = kk; 3125a9b9e01SMark F. Adams } 3135a9b9e01SMark F. Adams } 3145a9b9e01SMark F. Adams } 3155a9b9e01SMark F. Adams if (jj != -1) rfactor = jj; 316a2f3521dSMark F. Adams else rfactor = 1; /* does this happen .. a prime */ 3175a9b9e01SMark F. Adams } 3185a9b9e01SMark F. Adams new_npe = npe/rfactor; 3195a9b9e01SMark F. Adams 3205a9b9e01SMark F. Adams if (new_npe==nactive) { 321a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 3225a9b9e01SMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 323a2f3521dSMark F. Adams if (pc_gamg->verbose>0){ 324a2f3521dSMark F. Adams PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq); 325a2f3521dSMark F. Adams } 3265a9b9e01SMark F. Adams PetscFunctionReturn(0); 3275a9b9e01SMark F. Adams } 3285a9b9e01SMark F. Adams 329a2f3521dSMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq); 3305a9b9e01SMark F. Adams targetPE = mype/rfactor; 331a2f3521dSMark F. Adams ierr = ISCreateStride(wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 3325a9b9e01SMark F. Adams 333a2f3521dSMark F. Adams if (stokes) { 334c5bfad50SMark F. Adams ierr = ISCreateStride(wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr); 3355a9b9e01SMark F. Adams } 336a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 337e33ef3b1SMark F. Adams 33811e60469SMark F. Adams /* 339a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 34011e60469SMark F. Adams */ 341a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 342a2f3521dSMark F. Adams if (stokes) { 343a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr); 344806fa848SBarry Smith } else is_eq_num_prim = is_eq_num; 34511e60469SMark F. Adams /* 346a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 34711e60469SMark F. Adams */ 348a2f3521dSMark F. Adams ierr = ISPartitioningCount(is_eq_newproc, npe, counts);CHKERRQ(ierr); 349a2f3521dSMark F. Adams ncrs_eq_new = counts[mype]; 350a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 351a2f3521dSMark F. Adams if (stokes) { 352a2f3521dSMark F. Adams ierr = ISPartitioningCount(is_eq_newproc_prim, npe, counts);CHKERRQ(ierr); 353a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr); 354c5bfad50SMark F. Adams ncrs_prim_new = counts[mype]/cr_bs; /* nodes */ 355806fa848SBarry Smith } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */ 356a2f3521dSMark F. Adams 357a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 3580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3590cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 360b4fbaa2aSMark F. Adams #endif 361a2f3521dSMark F. Adams 362a2f3521dSMark F. Adams /* move data (for primal equations only) */ 36322063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3643bf036e2SBarry Smith ierr = VecCreate(wcomm, &dest_crd);CHKERRQ(ierr); 365a2f3521dSMark F. Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr); 366470e26f8SMark F. Adams ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */ 36711e60469SMark F. Adams /* 3689d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 369c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 37011e60469SMark F. Adams */ 371a2f3521dSMark F. Adams ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr); 372a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 373a2f3521dSMark F. Adams for (ii=0,jj=0; ii<ncrs_prim ; ii++) { 374c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 375a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 37611e60469SMark F. Adams } 377a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 378806fa848SBarry Smith ierr = ISCreateGeneral(wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 37992a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 38011e60469SMark F. Adams /* 38111e60469SMark F. Adams Create a vector to contain the original vertex information for each element 38211e60469SMark F. Adams */ 383a2f3521dSMark F. Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr); 3849d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols ; jj++) { 385a2f3521dSMark F. Adams const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows; 386a2f3521dSMark F. Adams for (ii=0 ; ii<ncrs_prim ; ii++) { 3879d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows ; kk++) { 388a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 389c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 390676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 391d3d6bff4SMark F. Adams } 392038e3b61SMark F. Adams } 393eb07cef2SMark F. Adams } 394eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 395eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 39611e60469SMark F. Adams /* 39711e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 39811e60469SMark F. Adams to the correct processor 39911e60469SMark F. Adams */ 400806fa848SBarry Smith ierr = VecScatterCreate(src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 40111e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 40211e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40311e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40411e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 40511e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 40611e60469SMark F. Adams /* 40711e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 40811e60469SMark F. Adams */ 409c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 410a2f3521dSMark F. Adams ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr); 411a2f3521dSMark F. Adams pc_gamg->data_sz = node_data_sz*ncrs_prim_new; 412a2f3521dSMark F. Adams strideNew = ncrs_prim_new*ndata_rows; 41311e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 4149d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols ; jj++) { 415a2f3521dSMark F. Adams for (ii=0 ; ii<ncrs_prim_new ; ii++) { 4169d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows ; kk++) { 417a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 418c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 419d3d6bff4SMark F. Adams } 420038e3b61SMark F. Adams } 421038e3b61SMark F. Adams } 42211e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 42311e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 424a2f3521dSMark F. Adams 425a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 4260cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4270cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 428ed3f9983SMark F. Adams #endif 429a2f3521dSMark F. Adams 43011e60469SMark F. Adams /* 43111e60469SMark F. Adams Invert for MatGetSubMatrix 43211e60469SMark F. Adams */ 433a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 434a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 435c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 436a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 437a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 438a2f3521dSMark F. Adams } 439a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 4400cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4410cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 4420cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 443ed3f9983SMark F. Adams #endif 444a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 445a2f3521dSMark F. Adams { 446a2f3521dSMark F. Adams Mat mat; 447806fa848SBarry Smith ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 448a2f3521dSMark F. Adams *a_Amat_crs = mat; 449c5bfad50SMark F. Adams 450c5bfad50SMark F. Adams if (!PETSC_TRUE){ 451c5bfad50SMark F. Adams PetscInt cbs, rbs; 452c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr); 453806fa848SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 454c5bfad50SMark F. Adams ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr); 455806fa848SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",mype,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr); 456c5bfad50SMark F. Adams } 457a2f3521dSMark F. Adams } 458038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 459a2f3521dSMark F. Adams 4600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4610cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 462ed3f9983SMark F. Adams #endif 46311e60469SMark F. Adams /* prolongator */ 46411e60469SMark F. Adams { 46511e60469SMark F. Adams IS findices; 466a2f3521dSMark F. Adams PetscInt Istart,Iend; 467a2f3521dSMark F. Adams Mat Pnew; 468a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 4690cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4700cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 471ed3f9983SMark F. Adams #endif 47211e60469SMark F. Adams ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 473c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 474806fa848SBarry Smith ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 47511e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 476c5bfad50SMark F. Adams 477c5bfad50SMark F. Adams if (!PETSC_TRUE){ 478c5bfad50SMark F. Adams PetscInt cbs, rbs; 479c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr); 480806fa848SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 481c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr); 482806fa848SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 483c5bfad50SMark F. Adams } 4840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4850cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 486ed3f9983SMark F. Adams #endif 4873530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 488a2f3521dSMark F. Adams 489a2f3521dSMark F. Adams /* output - repartitioned */ 490a2f3521dSMark F. Adams *a_P_inout = Pnew; 491e33ef3b1SMark F. Adams } 492a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4935b89ad90SMark F. Adams 4945a9b9e01SMark F. Adams *a_nactive_proc = new_npe; /* output */ 495a2f3521dSMark F. Adams } 4965a9b9e01SMark F. Adams 497a2f3521dSMark F. Adams /* outout matrix data */ 498c8b0795cSMark F. Adams if (!PETSC_TRUE) { 499c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 500c8b0795cSMark F. Adams if (llev==0) { 501c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 502c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 503c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 504c8b0795cSMark F. Adams ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr); 505c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 506c8b0795cSMark F. Adams } 507c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 508c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 509c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 510c8b0795cSMark F. Adams ierr = MatView(Cmat, viewer);CHKERRQ(ierr); 511c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 512c8b0795cSMark F. Adams } 513c8b0795cSMark F. Adams 5145b89ad90SMark F. Adams PetscFunctionReturn(0); 5155b89ad90SMark F. Adams } 5165b89ad90SMark F. Adams 5175b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 5185b89ad90SMark F. Adams /* 5195b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 5205b89ad90SMark F. Adams by setting data structures and options. 5215b89ad90SMark F. Adams 5225b89ad90SMark F. Adams Input Parameter: 5235b89ad90SMark F. Adams . pc - the preconditioner context 5245b89ad90SMark F. Adams 5255b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 5265b89ad90SMark F. Adams 5275b89ad90SMark F. Adams Notes: 5285b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 5295b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 5305b89ad90SMark F. Adams */ 5315b89ad90SMark F. Adams #undef __FUNCT__ 5325b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 5339d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 5345b89ad90SMark F. Adams { 5355b89ad90SMark F. Adams PetscErrorCode ierr; 5369d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 5375b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 5382adcac29SMark F. Adams Mat Pmat = pc->pmat; 539a2f3521dSMark F. Adams PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 5409d5b6da9SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 5413530afc2SMark F. Adams PetscMPIInt mype,npe,nactivepe; 542587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 543c8b0795cSMark F. Adams PetscReal emaxs[GAMG_MAXLEVELS]; 544e696ed0bSMark F. Adams IS *ASMLocalIDsArr[GAMG_MAXLEVELS]; 545a2f3521dSMark F. Adams GAMGKKTMat kktMatsArr[GAMG_MAXLEVELS]; 546a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 547737a81a9SMark F. Adams MatInfo info; 548*dfd5c07aSMark F. Adams PetscBool stokes = PETSC_FALSE, redo_mesh_setup = !pc_gamg->reuse_prol; 5495ef31b24SMark F. Adams 5505b89ad90SMark F. Adams PetscFunctionBegin; 55184d3f75bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 55284d3f75bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 553*dfd5c07aSMark F. Adams 554a2f3521dSMark F. Adams if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled); 555*dfd5c07aSMark F. Adams 55684d3f75bSMark F. Adams if (pc_gamg->setup_count++ > 0) { 557878e152fSMark F. Adams if (redo_mesh_setup) { 558878e152fSMark F. Adams /* reset everything */ 559878e152fSMark F. Adams ierr = PCReset_MG(pc); CHKERRQ(ierr); 560878e152fSMark F. Adams pc->setupcalled = 0; 561806fa848SBarry Smith } else { 56284d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 56303a628feSMark F. Adams /* just do Galerkin grids */ 56458471d46SMark F. Adams Mat B,dA,dB; 565d5280255SMark F. Adams assert(pc->setupcalled); 56658471d46SMark F. Adams 5679d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 56858471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5699d5b6da9SMark F. Adams ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 57058471d46SMark F. Adams /* (re)set to get dirty flag */ 5719d5b6da9SMark F. Adams ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 57258471d46SMark F. Adams 5739d5b6da9SMark F. Adams for (level=pc_gamg->Nlevels-2; level>-1; level--) { 57403a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 57584d3f75bSMark F. Adams if (pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) { 57603a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 577084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 57803a628feSMark F. Adams mglevels[level]->A = B; 579806fa848SBarry Smith } else { 580084a8fe3SJed Brown ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 58158471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 58203a628feSMark F. Adams } 58358471d46SMark F. Adams ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 58458471d46SMark F. Adams dB = B; 58558471d46SMark F. Adams } 5865f8cf99dSMark F. Adams } 587d5280255SMark F. Adams 588d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 589d5280255SMark F. Adams 590d5280255SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 591*dfd5c07aSMark F. Adams /* ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr); */ 592d5280255SMark F. Adams 59358471d46SMark F. Adams PetscFunctionReturn(0); 594eb07cef2SMark F. Adams } 595878e152fSMark F. Adams } 596f6536408SMark F. Adams 597302f38e8SMark F. Adams ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 598eb07cef2SMark F. Adams 5999d1b15c3SMark F. Adams ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr); 6009d1b15c3SMark F. Adams 601878e152fSMark F. Adams if (!pc_gamg->data) { 602878e152fSMark F. Adams if (pc_gamg->orig_data) { 603878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 604878e152fSMark F. Adams ierr = MatGetLocalSize(Pmat, &qq, PETSC_NULL);CHKERRQ(ierr); 605878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 606878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 607878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 608878e152fSMark F. Adams ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr); 609878e152fSMark F. Adams for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 610806fa848SBarry Smith } else { 6119d5b6da9SMark F. Adams if (!pc_gamg->createdefaultdata){ 612c666adbfSMark F. Adams SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 613302f38e8SMark F. Adams } 614302f38e8SMark F. Adams if (stokes) { 615c666adbfSMark F. Adams SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems"); 616eb07cef2SMark F. Adams } 6179d1b15c3SMark F. Adams ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr); 6189d5b6da9SMark F. Adams } 619878e152fSMark F. Adams } 620878e152fSMark F. Adams 621878e152fSMark F. Adams /* cache original data for reuse */ 622878e152fSMark F. Adams if (!pc_gamg->orig_data && redo_mesh_setup) { 623878e152fSMark F. Adams ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr); 624878e152fSMark F. Adams for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 625878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 626878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 627878e152fSMark F. Adams } 628038e3b61SMark F. Adams 629302f38e8SMark F. Adams /* get basic dims */ 630302f38e8SMark F. Adams if (stokes) { 631a2f3521dSMark F. Adams bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */ 632806fa848SBarry Smith } else { 633302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 634302f38e8SMark F. Adams } 635a2f3521dSMark F. Adams 636a2f3521dSMark F. Adams ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr); 637c8b0795cSMark F. Adams if (pc_gamg->verbose) { 63884f9421dSMark F. Adams PetscInt NN = M; 63984f9421dSMark F. Adams if (pc_gamg->verbose==1) { 64084f9421dSMark F. Adams ierr = MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr); 6413bf036e2SBarry Smith ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr); 642806fa848SBarry Smith } else { 643806fa848SBarry Smith ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 64484f9421dSMark F. Adams } 645b2a4f308SMark F. Adams nnz0 = info.nz_used; 646b2a4f308SMark F. Adams nnztot = info.nz_used; 647806fa848SBarry Smith ierr = PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 648a2f3521dSMark F. Adams mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 649806fa848SBarry Smith (int)(nnz0/(PetscReal)NN),npe);CHKERRQ(ierr); 650c8b0795cSMark F. Adams } 65184d3f75bSMark F. Adams 652a2f3521dSMark F. Adams /* Get A_i and R_i */ 6538f4b7eb5SMark F. Adams for (level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */ 6549d5b6da9SMark F. Adams level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 6550205a208SMark F. Adams level++){ 6565b89ad90SMark F. Adams level1 = level + 1; 6570cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6580cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 659a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 660a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 661b4fbaa2aSMark F. Adams #endif 662a2f3521dSMark F. Adams #endif 663a2f3521dSMark F. Adams /* deal with Stokes, get sub matrices */ 6649d1b15c3SMark F. Adams if (level > 0) { 665a2f3521dSMark F. Adams ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr); 6669d1b15c3SMark F. Adams } 667c8b0795cSMark F. Adams { /* construct prolongator */ 668725b86d8SJed Brown Mat Gmat; 6690cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 670a2f3521dSMark F. Adams Mat Prol11,Prol22; 671c8b0795cSMark F. Adams 672a2f3521dSMark F. Adams ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr); 6730cbbd2e1SMark F. Adams ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 674a2f3521dSMark F. Adams ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr); 675c8b0795cSMark F. Adams 676a2f3521dSMark F. Adams /* could have failed to create new level */ 677a2f3521dSMark F. Adams if (Prol11){ 6789d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6799d1b15c3SMark F. Adams ierr = MatGetBlockSizes(Prol11, PETSC_NULL, &bs);CHKERRQ(ierr); 680a2f3521dSMark F. Adams 681a2f3521dSMark F. Adams if (stokes) { 682a2f3521dSMark F. Adams if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method."); 683a2f3521dSMark F. Adams /* R A12 == (T = A21 P)'; G = T' T; coarsen G; form plain agg with G */ 684a2f3521dSMark F. Adams ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr); 685c8b0795cSMark F. Adams } 686c8b0795cSMark F. Adams 687a2f3521dSMark F. Adams if (pc_gamg->optprol){ 688c8b0795cSMark F. Adams /* smooth */ 689a2f3521dSMark F. Adams ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr); 690c8b0795cSMark F. Adams } 691c8b0795cSMark F. Adams 692a2f3521dSMark F. Adams if (stokes) { 693a2f3521dSMark F. Adams IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is}; 694a2f3521dSMark F. Adams Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 }; 695a2f3521dSMark F. Adams ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr); 696806fa848SBarry Smith } else { 697a2f3521dSMark F. Adams Parr[level1] = Prol11; 698a2f3521dSMark F. Adams } 699806fa848SBarry Smith } else Parr[level1] = PETSC_NULL; 700ffc955d6SMark F. Adams 701ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 702806fa848SBarry Smith ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 703ffc955d6SMark F. Adams } 704ffc955d6SMark F. Adams 705a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 70641b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 707a2f3521dSMark F. Adams } /* construct prolongator scope */ 7080cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7090cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 710c8b0795cSMark F. Adams #endif 7119d5b6da9SMark F. Adams /* cache eigen estimate */ 7129d5b6da9SMark F. Adams if (pc_gamg->emax_id != -1){ 7139d5b6da9SMark F. Adams PetscBool flag; 714806fa848SBarry Smith ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr); 7159d5b6da9SMark F. Adams if (!flag) emaxs[level] = -1.; 716806fa848SBarry Smith } else emaxs[level] = -1.; 7172adcac29SMark F. Adams if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 718c8b0795cSMark F. Adams if (!Parr[level1]) { 719806fa848SBarry Smith if (pc_gamg->verbose) { 7203bf036e2SBarry Smith ierr = PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);CHKERRQ(ierr); 721806fa848SBarry Smith } 722c8b0795cSMark F. Adams break; 723c8b0795cSMark F. Adams } 7240cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7250cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 726b4fbaa2aSMark F. Adams #endif 727a2f3521dSMark F. Adams 728a2f3521dSMark F. Adams ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2), 729806fa848SBarry Smith stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr); 730a2f3521dSMark F. Adams 7310cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7320cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 733b4fbaa2aSMark F. Adams #endif 734a2f3521dSMark F. Adams ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr); 735a2f3521dSMark F. Adams 736a2f3521dSMark F. Adams if (pc_gamg->verbose > 0){ 7370cbbd2e1SMark F. Adams PetscInt NN = M; 7380cbbd2e1SMark F. Adams if (pc_gamg->verbose==1) { 739a2f3521dSMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr); 7403bf036e2SBarry Smith ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr); 741806fa848SBarry Smith } else { 742806fa848SBarry Smith ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 7430cbbd2e1SMark F. Adams } 744a2f3521dSMark F. Adams 7450cbbd2e1SMark F. Adams nnztot += info.nz_used; 746806fa848SBarry Smith ierr = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 747a2f3521dSMark F. Adams mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols, 748806fa848SBarry Smith (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr); 749c8b0795cSMark F. Adams } 750a2f3521dSMark F. Adams 751a2f3521dSMark F. Adams /* stop if one node -- could pull back for singular problems */ 752c8b0795cSMark F. Adams if (M/pc_gamg->data_cell_cols < 2) { 75381708dccSMark F. Adams level++; 75481708dccSMark F. Adams break; 75581708dccSMark F. Adams } 7560cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 757b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 758b4fbaa2aSMark F. Adams #endif 759c8b0795cSMark F. Adams } /* levels */ 760c8b0795cSMark F. Adams 761c8b0795cSMark F. Adams if (pc_gamg->data) { 762c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 763a2f3521dSMark F. Adams pc_gamg->data = PETSC_NULL; 7645b89ad90SMark F. Adams } 765c8b0795cSMark F. Adams 7660cbbd2e1SMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 7679d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 7685b89ad90SMark F. Adams fine_level = level; 7699d5b6da9SMark F. Adams ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr); 7705b89ad90SMark F. Adams 77184d3f75bSMark F. Adams /* simple setup */ 77284d3f75bSMark F. Adams if (!PETSC_TRUE){ 77384d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 77484d3f75bSMark F. Adams for (lidx=0,level=pc_gamg->Nlevels-1; 77584d3f75bSMark F. Adams lidx<fine_level; 77684d3f75bSMark F. Adams lidx++, level--){ 77784d3f75bSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr); 77884d3f75bSMark F. Adams ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 77984d3f75bSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 78084d3f75bSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 78184d3f75bSMark F. Adams } 78284d3f75bSMark F. Adams ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 78384d3f75bSMark F. Adams 78484d3f75bSMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 785806fa848SBarry Smith } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 786d5280255SMark F. Adams /* set default smoothers & set operators */ 7879d5b6da9SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 788587fa25dSMark F. Adams lidx <= fine_level; 789587fa25dSMark F. Adams lidx++, level--) { 790ffc955d6SMark F. Adams KSP smoother; 791ffc955d6SMark F. Adams PC subpc; 792a2f3521dSMark F. Adams 7939d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 794f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 795ffc955d6SMark F. Adams 796a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 797a2f3521dSMark F. Adams /* set ops */ 798a2f3521dSMark F. Adams ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 799a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 800a2f3521dSMark F. Adams 801a2f3521dSMark F. Adams /* create field split PC, get subsmoother */ 802a2f3521dSMark F. Adams if (stokes) { 803a2f3521dSMark F. Adams KSP *ksps; 804a2f3521dSMark F. Adams PetscInt nn; 805a2f3521dSMark F. Adams ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr); 806a2f3521dSMark F. Adams ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr); 807a2f3521dSMark F. Adams ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr); 808a2f3521dSMark F. Adams smoother = ksps[0]; 809a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 810a2f3521dSMark F. Adams ierr = PetscFree(ksps);CHKERRQ(ierr); 811a2f3521dSMark F. Adams } 812a2f3521dSMark F. Adams ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr); 813a2f3521dSMark F. Adams 814a2f3521dSMark F. Adams /* set defaults */ 8156c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 816a2f3521dSMark F. Adams 817ffc955d6SMark F. Adams /* override defaults and command line args (!) */ 818ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 8192d3561bbSSatish Balay PetscInt sz; 8202d3561bbSSatish Balay IS *is; 821a2f3521dSMark F. Adams 8222d3561bbSSatish Balay sz = nASMBlocksArr[level]; 8232d3561bbSSatish Balay is = ASMLocalIDsArr[level]; 824ffc955d6SMark F. Adams ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr); 825ffc955d6SMark F. Adams if (sz==0){ 826ffc955d6SMark F. Adams IS is; 827ffc955d6SMark F. Adams PetscInt my0,kk; 828ffc955d6SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr); 829ffc955d6SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 83006b43e7eSDmitry Karpeev ierr = PCGASMSetSubdomains(subpc, 1, &is, PETSC_NULL);CHKERRQ(ierr); 831a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 832806fa848SBarry Smith } else { 833a94c3b12SMark F. Adams PetscInt kk; 83406b43e7eSDmitry Karpeev ierr = PCGASMSetSubdomains(subpc, sz, is, PETSC_NULL);CHKERRQ(ierr); 835a94c3b12SMark F. Adams for (kk=0;kk<sz;kk++){ 836a94c3b12SMark F. Adams ierr = ISDestroy(&is[kk]);CHKERRQ(ierr); 837a94c3b12SMark F. Adams } 838ffc955d6SMark F. Adams ierr = PetscFree(is);CHKERRQ(ierr); 839ffc955d6SMark F. Adams } 840ffc955d6SMark F. Adams ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr); 841ffc955d6SMark F. Adams 842ffc955d6SMark F. Adams ASMLocalIDsArr[level] = PETSC_NULL; 843ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 844ffc955d6SMark F. Adams ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr); 845806fa848SBarry Smith } else { 8469d5b6da9SMark F. Adams ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr); 847ffc955d6SMark F. Adams } 848d5280255SMark F. Adams } 849d5280255SMark F. Adams { 850d5280255SMark F. Adams /* coarse grid */ 851d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 852d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 853d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 854d5280255SMark F. Adams ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 855d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 856d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 857d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 858d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 859d5280255SMark F. Adams ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); assert(ii==1); 860d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 861d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 862d5280255SMark F. Adams } 863d5280255SMark F. Adams 864d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 865d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 866d5280255SMark F. Adams ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); 867d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 868d5280255SMark F. Adams if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 869d5280255SMark F. Adams 870d5280255SMark F. Adams /* create cheby smoothers */ 871d5280255SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 872d5280255SMark F. Adams lidx <= fine_level; 873d5280255SMark F. Adams lidx++, level--) { 874d5280255SMark F. Adams KSP smoother; 875ffc955d6SMark F. Adams PetscBool flag; 876d5280255SMark F. Adams PC subpc; 877d5280255SMark F. Adams 878ffc955d6SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 879a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 880a2f3521dSMark F. Adams 881a2f3521dSMark F. Adams /* create field split PC, get subsmoother */ 882a2f3521dSMark F. Adams if (stokes) { 883a2f3521dSMark F. Adams KSP *ksps; 884a2f3521dSMark F. Adams PetscInt nn; 885a2f3521dSMark F. Adams ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr); 886a2f3521dSMark F. Adams smoother = ksps[0]; 887a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 888a2f3521dSMark F. Adams ierr = PetscFree(ksps);CHKERRQ(ierr); 889a2f3521dSMark F. Adams } 890ffc955d6SMark F. Adams 891ffc955d6SMark F. Adams /* do my own cheby */ 8926c9de887SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr); 893ffc955d6SMark F. Adams if (flag) { 894ffc955d6SMark F. Adams PetscReal emax, emin; 895251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr); 896ffc955d6SMark F. Adams if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 897587fa25dSMark F. Adams else{ /* eigen estimate 'emax' */ 898e696ed0bSMark F. Adams KSP eksp; 899e696ed0bSMark F. Adams Mat Lmat = Aarr[level]; 900b2a4f308SMark F. Adams Vec bb, xx; 901038e3b61SMark F. Adams 9025745e0f5SMark F. Adams ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr); 9035745e0f5SMark F. Adams ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr); 904fc4362bfSMark F. Adams { 905fc4362bfSMark F. Adams PetscRandom rctx; 906fc4362bfSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 907fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 908fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 909fc4362bfSMark F. Adams ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 9105b89ad90SMark F. Adams } 911a2f3521dSMark F. Adams 912e696ed0bSMark F. Adams /* zeroing out BC rows -- needed for crazy matrices */ 913e696ed0bSMark F. Adams { 914e696ed0bSMark F. Adams PetscInt Istart,Iend,ncols,jj,Ii; 915e696ed0bSMark F. Adams PetscScalar zero = 0.0; 916e696ed0bSMark F. Adams ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr); 917e696ed0bSMark F. Adams for (Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++) { 918e696ed0bSMark F. Adams ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 919e696ed0bSMark F. Adams if(ncols <= 1) { 920e696ed0bSMark F. Adams ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 921a94c3b12SMark F. Adams } 922e696ed0bSMark F. Adams ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 923a94c3b12SMark F. Adams } 924a94c3b12SMark F. Adams ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 925a94c3b12SMark F. Adams ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 926a94c3b12SMark F. Adams } 927e696ed0bSMark F. Adams 928fc4362bfSMark F. Adams ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr); 929806fa848SBarry Smith ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr); 930fc4362bfSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 9311a166f3bSJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 9321a166f3bSJed Brown ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 933f6536408SMark F. Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 934f6536408SMark F. Adams 935f6536408SMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 936f6536408SMark F. Adams ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 937fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 9385a9b9e01SMark F. Adams 939d3d0db20SJed Brown /* set PC type to be same as smoother */ 940ffc955d6SMark F. Adams ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr); 941b2a4f308SMark F. Adams 9425a9b9e01SMark F. Adams /* solve - keep stuff out of logging */ 9435a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 9445a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 945fc4362bfSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 9465a9b9e01SMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 9475a9b9e01SMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 9485a9b9e01SMark F. Adams 949fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 9505a9b9e01SMark F. Adams 951fc4362bfSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 952fc4362bfSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 953fc4362bfSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 954f6536408SMark F. Adams 955ffc955d6SMark F. Adams if (pc_gamg->verbose > 0) { 956a94c3b12SMark F. Adams PetscInt N1, tt; 957a94c3b12SMark F. Adams ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr); 958a94c3b12SMark F. Adams PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 959f6536408SMark F. Adams } 960fc4362bfSMark F. Adams } 961038e3b61SMark F. Adams { 962c5bfad50SMark F. Adams PetscInt N1, N0; 963c5bfad50SMark F. Adams ierr = MatGetSize(Aarr[level], &N1, PETSC_NULL);CHKERRQ(ierr); 964c5bfad50SMark F. Adams ierr = MatGetSize(Aarr[level+1], &N0, PETSC_NULL);CHKERRQ(ierr); 965f6536408SMark F. Adams /* heuristic - is this crap? */ 966b4ec6429SMark F. Adams /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */ 9675e7c91beSJed Brown emin = emax * pc_gamg->eigtarget[0]; 9685e7c91beSJed Brown emax *= pc_gamg->eigtarget[1]; 969038e3b61SMark F. Adams } 9706c9de887SHong Zhang ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); 971ffc955d6SMark F. Adams } /* setup checby flag */ 972ffc955d6SMark F. Adams } /* non-coarse levels */ 973737a81a9SMark F. Adams 974d5280255SMark F. Adams /* clean up */ 975d5280255SMark F. Adams for (level=1;level<pc_gamg->Nlevels;level++){ 976587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 977587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 9785b89ad90SMark F. Adams } 9790cbbd2e1SMark F. Adams 9800cbbd2e1SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 9810cbbd2e1SMark F. Adams 982a94c3b12SMark F. Adams if (PETSC_FALSE){ 98358471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 9849d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 98558471d46SMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 98658471d46SMark F. Adams } 987806fa848SBarry Smith } else { 9885f8cf99dSMark F. Adams KSP smoother; 989b43b03e9SMark F. Adams if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__); 9909d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 9915f8cf99dSMark F. Adams ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 9925f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 9939d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 9945f8cf99dSMark F. Adams } 99584d3f75bSMark 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); 10195b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 10205b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 10215b89ad90SMark F. Adams PetscFunctionReturn(0); 10225b89ad90SMark F. Adams } 10235b89ad90SMark F. Adams 1024676e1743SMark F. Adams 1025676e1743SMark F. Adams #undef __FUNCT__ 1026676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 1027676e1743SMark F. Adams /*@ 1028676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 1029676e1743SMark F. Adams processor reduction. 1030676e1743SMark F. Adams 1031676e1743SMark F. Adams Not Collective on PC 1032676e1743SMark F. Adams 1033676e1743SMark F. Adams Input Parameters: 1034676e1743SMark F. Adams . pc - the preconditioner context 1035676e1743SMark F. Adams 1036676e1743SMark F. Adams 1037676e1743SMark F. Adams Options Database Key: 1038676e1743SMark F. Adams . -pc_gamg_process_eq_limit 1039676e1743SMark F. Adams 1040676e1743SMark F. Adams Level: intermediate 1041676e1743SMark F. Adams 1042676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1043676e1743SMark F. Adams 1044676e1743SMark F. Adams .seealso: () 1045676e1743SMark F. Adams @*/ 1046676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 1047676e1743SMark F. Adams { 1048676e1743SMark F. Adams PetscErrorCode ierr; 1049676e1743SMark F. Adams 1050676e1743SMark F. Adams PetscFunctionBegin; 1051676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1052676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1053676e1743SMark F. Adams PetscFunctionReturn(0); 1054676e1743SMark F. Adams } 1055676e1743SMark F. Adams 1056676e1743SMark F. Adams EXTERN_C_BEGIN 1057676e1743SMark F. Adams #undef __FUNCT__ 1058676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 1059676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 1060676e1743SMark F. Adams { 1061c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1062c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1063676e1743SMark F. Adams 1064676e1743SMark F. Adams PetscFunctionBegin; 10659d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 1066676e1743SMark F. Adams PetscFunctionReturn(0); 1067676e1743SMark F. Adams } 1068676e1743SMark F. Adams EXTERN_C_END 1069676e1743SMark F. Adams 1070676e1743SMark F. Adams #undef __FUNCT__ 1071389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 1072389730f3SMark F. Adams /*@ 1073389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 1074389730f3SMark F. Adams 1075389730f3SMark F. Adams Collective on PC 1076389730f3SMark F. Adams 1077389730f3SMark F. Adams Input Parameters: 1078389730f3SMark F. Adams . pc - the preconditioner context 1079389730f3SMark F. Adams 1080389730f3SMark F. Adams 1081389730f3SMark F. Adams Options Database Key: 1082389730f3SMark F. Adams . -pc_gamg_coarse_eq_limit 1083389730f3SMark F. Adams 1084389730f3SMark F. Adams Level: intermediate 1085389730f3SMark F. Adams 1086389730f3SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1087389730f3SMark F. Adams 1088389730f3SMark F. Adams .seealso: () 1089389730f3SMark F. Adams @*/ 1090389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 1091389730f3SMark F. Adams { 1092389730f3SMark F. Adams PetscErrorCode ierr; 1093389730f3SMark F. Adams 1094389730f3SMark F. Adams PetscFunctionBegin; 1095389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1096389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1097389730f3SMark F. Adams PetscFunctionReturn(0); 1098389730f3SMark F. Adams } 1099389730f3SMark F. Adams 1100389730f3SMark F. Adams EXTERN_C_BEGIN 1101389730f3SMark F. Adams #undef __FUNCT__ 1102389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 1103389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1104389730f3SMark F. Adams { 1105389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1106389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1107389730f3SMark F. Adams 1108389730f3SMark F. Adams PetscFunctionBegin; 11099d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 1110389730f3SMark F. Adams PetscFunctionReturn(0); 1111389730f3SMark F. Adams } 1112389730f3SMark F. Adams EXTERN_C_END 1113389730f3SMark F. Adams 1114389730f3SMark F. Adams #undef __FUNCT__ 11158263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 1116676e1743SMark F. Adams /*@ 11178263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 1118676e1743SMark F. Adams 1119676e1743SMark F. Adams Collective on PC 1120676e1743SMark F. Adams 1121676e1743SMark F. Adams Input Parameters: 1122676e1743SMark F. Adams . pc - the preconditioner context 1123676e1743SMark F. Adams 1124676e1743SMark F. Adams 1125676e1743SMark F. Adams Options Database Key: 11268263b398SMark F. Adams . -pc_gamg_repartition 1127676e1743SMark F. Adams 1128676e1743SMark F. Adams Level: intermediate 1129676e1743SMark F. Adams 1130676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1131676e1743SMark F. Adams 1132676e1743SMark F. Adams .seealso: () 1133676e1743SMark F. Adams @*/ 11348263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1135676e1743SMark F. Adams { 1136676e1743SMark F. Adams PetscErrorCode ierr; 1137676e1743SMark F. Adams 1138676e1743SMark F. Adams PetscFunctionBegin; 1139676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11408263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1141676e1743SMark F. Adams PetscFunctionReturn(0); 1142676e1743SMark F. Adams } 1143676e1743SMark F. Adams 1144676e1743SMark F. Adams EXTERN_C_BEGIN 1145676e1743SMark F. Adams #undef __FUNCT__ 11468263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 11478263b398SMark F. Adams 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 EXTERN_C_END 1157676e1743SMark F. Adams 1158676e1743SMark F. Adams #undef __FUNCT__ 1159*dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl" 1160*dfd5c07aSMark F. Adams /*@ 1161*dfd5c07aSMark F. Adams PCGAMGSetReuseProl - Reuse prlongation 1162*dfd5c07aSMark F. Adams 1163*dfd5c07aSMark F. Adams Collective on PC 1164*dfd5c07aSMark F. Adams 1165*dfd5c07aSMark F. Adams Input Parameters: 1166*dfd5c07aSMark F. Adams . pc - the preconditioner context 1167*dfd5c07aSMark F. Adams 1168*dfd5c07aSMark F. Adams 1169*dfd5c07aSMark F. Adams Options Database Key: 1170*dfd5c07aSMark F. Adams . -pc_gamg_reuse_interpolation 1171*dfd5c07aSMark F. Adams 1172*dfd5c07aSMark F. Adams Level: intermediate 1173*dfd5c07aSMark F. Adams 1174*dfd5c07aSMark F. Adams Concepts: Unstructured multrigrid preconditioner 1175*dfd5c07aSMark F. Adams 1176*dfd5c07aSMark F. Adams .seealso: () 1177*dfd5c07aSMark F. Adams @*/ 1178*dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n) 1179*dfd5c07aSMark F. Adams { 1180*dfd5c07aSMark F. Adams PetscErrorCode ierr; 1181*dfd5c07aSMark F. Adams 1182*dfd5c07aSMark F. Adams PetscFunctionBegin; 1183*dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1184*dfd5c07aSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1185*dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1186*dfd5c07aSMark F. Adams } 1187*dfd5c07aSMark F. Adams 1188*dfd5c07aSMark F. Adams EXTERN_C_BEGIN 1189*dfd5c07aSMark F. Adams #undef __FUNCT__ 1190*dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG" 1191*dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n) 1192*dfd5c07aSMark F. Adams { 1193*dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1194*dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1195*dfd5c07aSMark F. Adams 1196*dfd5c07aSMark F. Adams PetscFunctionBegin; 1197*dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1198*dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1199*dfd5c07aSMark F. Adams } 1200*dfd5c07aSMark F. Adams EXTERN_C_END 1201*dfd5c07aSMark F. Adams 1202*dfd5c07aSMark F. Adams #undef __FUNCT__ 1203ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs" 1204ffc955d6SMark F. Adams /*@ 1205ffc955d6SMark F. Adams PCGAMGSetUseASMAggs - 1206ffc955d6SMark F. Adams 1207ffc955d6SMark F. Adams Collective on PC 1208ffc955d6SMark F. Adams 1209ffc955d6SMark F. Adams Input Parameters: 1210ffc955d6SMark F. Adams . pc - the preconditioner context 1211ffc955d6SMark F. Adams 1212ffc955d6SMark F. Adams 1213ffc955d6SMark F. Adams Options Database Key: 1214ffc955d6SMark F. Adams . -pc_gamg_use_agg_gasm 1215ffc955d6SMark F. Adams 1216ffc955d6SMark F. Adams Level: intermediate 1217ffc955d6SMark F. Adams 1218ffc955d6SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1219ffc955d6SMark F. Adams 1220ffc955d6SMark F. Adams .seealso: () 1221ffc955d6SMark F. Adams @*/ 1222ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1223ffc955d6SMark F. Adams { 1224ffc955d6SMark F. Adams PetscErrorCode ierr; 1225ffc955d6SMark F. Adams 1226ffc955d6SMark F. Adams PetscFunctionBegin; 1227ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1228ffc955d6SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1229ffc955d6SMark F. Adams PetscFunctionReturn(0); 1230ffc955d6SMark F. Adams } 1231ffc955d6SMark F. Adams 1232ffc955d6SMark F. Adams EXTERN_C_BEGIN 1233ffc955d6SMark F. Adams #undef __FUNCT__ 1234ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 1235ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1236ffc955d6SMark F. Adams { 1237ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1238ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1239ffc955d6SMark F. Adams 1240ffc955d6SMark F. Adams PetscFunctionBegin; 1241ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = n; 1242ffc955d6SMark F. Adams PetscFunctionReturn(0); 1243ffc955d6SMark F. Adams } 1244ffc955d6SMark F. Adams EXTERN_C_END 1245ffc955d6SMark F. Adams 1246ffc955d6SMark F. Adams #undef __FUNCT__ 12474ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 12484ef23d27SMark F. Adams /*@ 12494ef23d27SMark F. Adams PCGAMGSetNlevels - 12504ef23d27SMark F. Adams 12514ef23d27SMark F. Adams Not collective on PC 12524ef23d27SMark F. Adams 12534ef23d27SMark F. Adams Input Parameters: 12544ef23d27SMark F. Adams . pc - the preconditioner context 12554ef23d27SMark F. Adams 12564ef23d27SMark F. Adams 12574ef23d27SMark F. Adams Options Database Key: 12584ef23d27SMark F. Adams . -pc_mg_levels 12594ef23d27SMark F. Adams 12604ef23d27SMark F. Adams Level: intermediate 12614ef23d27SMark F. Adams 12624ef23d27SMark F. Adams Concepts: Unstructured multrigrid preconditioner 12634ef23d27SMark F. Adams 12644ef23d27SMark F. Adams .seealso: () 12654ef23d27SMark F. Adams @*/ 12664ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 12674ef23d27SMark F. Adams { 12684ef23d27SMark F. Adams PetscErrorCode ierr; 12694ef23d27SMark F. Adams 12704ef23d27SMark F. Adams PetscFunctionBegin; 12714ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12724ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 12734ef23d27SMark F. Adams PetscFunctionReturn(0); 12744ef23d27SMark F. Adams } 12754ef23d27SMark F. Adams 12764ef23d27SMark F. Adams EXTERN_C_BEGIN 12774ef23d27SMark F. Adams #undef __FUNCT__ 12784ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 12794ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 12804ef23d27SMark F. Adams { 12814ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 12824ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12834ef23d27SMark F. Adams 12844ef23d27SMark F. Adams PetscFunctionBegin; 12859d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12864ef23d27SMark F. Adams PetscFunctionReturn(0); 12874ef23d27SMark F. Adams } 12884ef23d27SMark F. Adams EXTERN_C_END 12894ef23d27SMark F. Adams 12904ef23d27SMark F. Adams #undef __FUNCT__ 12913542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 12923542efc5SMark F. Adams /*@ 12933542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12943542efc5SMark F. Adams 12953542efc5SMark F. Adams Not collective on PC 12963542efc5SMark F. Adams 12973542efc5SMark F. Adams Input Parameters: 12983542efc5SMark F. Adams . pc - the preconditioner context 12993542efc5SMark F. Adams 13003542efc5SMark F. Adams 13013542efc5SMark F. Adams Options Database Key: 13023542efc5SMark F. Adams . -pc_gamg_threshold 13033542efc5SMark F. Adams 13043542efc5SMark F. Adams Level: intermediate 13053542efc5SMark F. Adams 13063542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 13073542efc5SMark F. Adams 13083542efc5SMark F. Adams .seealso: () 13093542efc5SMark F. Adams @*/ 13103542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 13113542efc5SMark F. Adams { 13123542efc5SMark F. Adams PetscErrorCode ierr; 13133542efc5SMark F. Adams 13143542efc5SMark F. Adams PetscFunctionBegin; 13153542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13163542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 13173542efc5SMark F. Adams PetscFunctionReturn(0); 13183542efc5SMark F. Adams } 13193542efc5SMark F. Adams 13203542efc5SMark F. Adams EXTERN_C_BEGIN 13213542efc5SMark F. Adams #undef __FUNCT__ 13223542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 13233542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 13243542efc5SMark F. Adams { 1325c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1326c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 13273542efc5SMark F. Adams 13283542efc5SMark F. Adams PetscFunctionBegin; 13299d5b6da9SMark F. Adams pc_gamg->threshold = n; 13303542efc5SMark F. Adams PetscFunctionReturn(0); 13313542efc5SMark F. Adams } 13323542efc5SMark F. Adams EXTERN_C_END 13333542efc5SMark F. Adams 13343542efc5SMark F. Adams #undef __FUNCT__ 13359d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1336676e1743SMark F. Adams /*@ 13379d5b6da9SMark F. Adams PCGAMGSetType - Set solution method - calls sub create method 1338676e1743SMark F. Adams 1339676e1743SMark F. Adams Collective on PC 1340676e1743SMark F. Adams 1341676e1743SMark F. Adams Input Parameters: 1342676e1743SMark F. Adams . pc - the preconditioner context 1343676e1743SMark F. Adams 1344676e1743SMark F. Adams 1345676e1743SMark F. Adams Options Database Key: 13463542efc5SMark F. Adams . -pc_gamg_type 1347676e1743SMark F. Adams 1348676e1743SMark F. Adams Level: intermediate 1349676e1743SMark F. Adams 1350676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1351676e1743SMark F. Adams 1352676e1743SMark F. Adams .seealso: () 1353676e1743SMark F. Adams @*/ 135419fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1355676e1743SMark F. Adams { 1356676e1743SMark F. Adams PetscErrorCode ierr; 1357676e1743SMark F. Adams 1358676e1743SMark F. Adams PetscFunctionBegin; 1359676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1360806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1361676e1743SMark F. Adams PetscFunctionReturn(0); 1362676e1743SMark F. Adams } 1363676e1743SMark F. Adams 1364676e1743SMark F. Adams EXTERN_C_BEGIN 1365676e1743SMark F. Adams #undef __FUNCT__ 13669d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 136719fd82e9SBarry Smith PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1368676e1743SMark F. Adams { 13699d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 1370676e1743SMark F. Adams 1371676e1743SMark F. Adams PetscFunctionBegin; 1372806fa848SBarry Smith ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr); 13739d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 13749d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1375676e1743SMark F. Adams PetscFunctionReturn(0); 1376676e1743SMark F. Adams } 1377676e1743SMark F. Adams EXTERN_C_END 1378676e1743SMark F. Adams 13795b89ad90SMark F. Adams #undef __FUNCT__ 13805b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 13815b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 13825b89ad90SMark F. Adams { 1383676e1743SMark F. Adams PetscErrorCode ierr; 1384676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1385676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1386676e1743SMark F. Adams PetscBool flag; 13875e7c91beSJed Brown PetscInt two = 2; 1388b43b03e9SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 13895b89ad90SMark F. Adams 13905b89ad90SMark F. Adams PetscFunctionBegin; 1391676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1392676e1743SMark F. Adams { 139375b74bdaSMark F. Adams /* -pc_gamg_verbose */ 13949d5b6da9SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 13959d5b6da9SMark F. Adams "none", pc_gamg->verbose, 1396806fa848SBarry Smith &pc_gamg->verbose, PETSC_NULL);CHKERRQ(ierr); 13978263b398SMark F. Adams /* -pc_gamg_repartition */ 13988263b398SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_repartition", 13998263b398SMark F. Adams "Repartion coarse grids (false)", 14008263b398SMark F. Adams "PCGAMGRepartitioning", 14019d5b6da9SMark F. Adams pc_gamg->repart, 14029d5b6da9SMark F. Adams &pc_gamg->repart, 1403806fa848SBarry Smith &flag);CHKERRQ(ierr); 1404*dfd5c07aSMark F. Adams /* -pc_gamg_reuse_interpolation */ 1405*dfd5c07aSMark F. Adams ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation", 1406*dfd5c07aSMark F. Adams "Reuse prolongation operator (true)", 1407*dfd5c07aSMark F. Adams "PCGAMGReuseProl", 1408*dfd5c07aSMark F. Adams pc_gamg->reuse_prol, 1409*dfd5c07aSMark F. Adams &pc_gamg->reuse_prol, 1410*dfd5c07aSMark F. Adams &flag);CHKERRQ(ierr); 1411ffc955d6SMark F. Adams /* -pc_gamg_use_agg_gasm */ 1412ffc955d6SMark F. Adams ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm", 1413ffc955d6SMark F. Adams "Use aggregation agragates for GASM smoother (false)", 1414ffc955d6SMark F. Adams "PCGAMGUseASMAggs", 1415ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm, 1416ffc955d6SMark F. Adams &pc_gamg->use_aggs_in_gasm, 1417806fa848SBarry Smith &flag);CHKERRQ(ierr); 1418c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 1419676e1743SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1420676e1743SMark F. Adams "Limit (goal) on number of equations per process on coarse grids", 1421676e1743SMark F. Adams "PCGAMGSetProcEqLim", 14229d5b6da9SMark F. Adams pc_gamg->min_eq_proc, 14239d5b6da9SMark F. Adams &pc_gamg->min_eq_proc, 1424806fa848SBarry Smith &flag); CHKERRQ(ierr); 1425389730f3SMark F. Adams /* -pc_gamg_coarse_eq_limit */ 1426389730f3SMark F. Adams ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1427389730f3SMark F. Adams "Limit on number of equations for the coarse grid", 1428389730f3SMark F. Adams "PCGAMGSetCoarseEqLim", 14299d5b6da9SMark F. Adams pc_gamg->coarse_eq_limit, 14309d5b6da9SMark F. Adams &pc_gamg->coarse_eq_limit, 1431806fa848SBarry Smith &flag);CHKERRQ(ierr); 1432c20e4228SMark F. Adams /* -pc_gamg_threshold */ 14333542efc5SMark F. Adams ierr = PetscOptionsReal("-pc_gamg_threshold", 14343542efc5SMark F. Adams "Relative threshold to use for dropping edges in aggregation graph", 14353542efc5SMark F. Adams "PCGAMGSetThreshold", 14369d5b6da9SMark F. Adams pc_gamg->threshold, 14379d5b6da9SMark F. Adams &pc_gamg->threshold, 1438806fa848SBarry Smith &flag);CHKERRQ(ierr); 1439806fa848SBarry Smith if (flag && pc_gamg->verbose) { 1440806fa848SBarry Smith ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr); 1441806fa848SBarry Smith } 14424ef23d27SMark F. Adams 14435e7c91beSJed Brown ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr); 14444ef23d27SMark F. Adams ierr = PetscOptionsInt("-pc_mg_levels", 14454ef23d27SMark F. Adams "Set number of MG levels", 14464ef23d27SMark F. Adams "PCGAMGSetNlevels", 14479d5b6da9SMark F. Adams pc_gamg->Nlevels, 14489d5b6da9SMark F. Adams &pc_gamg->Nlevels, 1449806fa848SBarry Smith &flag);CHKERRQ(ierr); 1450676e1743SMark F. Adams } 1451676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1452676e1743SMark F. Adams 14535b89ad90SMark F. Adams PetscFunctionReturn(0); 14545b89ad90SMark F. Adams } 14555b89ad90SMark F. Adams 14565b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 14575b89ad90SMark F. Adams /*MC 1458856830bfSJed Brown PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 14599d5b6da9SMark F. Adams - This is the entry point to GAMG, registered in pcregis.c 14605b89ad90SMark F. Adams 1461280d9858SJed Brown Options Database Keys: 14625b89ad90SMark F. Adams Multigrid options(inherited) 1463280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1464280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1465280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 14668c1c2452SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 14675b89ad90SMark F. Adams 14685b89ad90SMark F. Adams Level: intermediate 1469280d9858SJed Brown 14705b89ad90SMark F. Adams Concepts: multigrid 14715b89ad90SMark F. Adams 14725b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1473280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 14745b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 14755b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 14765b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 14775b89ad90SMark F. Adams M*/ 14785b89ad90SMark F. Adams EXTERN_C_BEGIN 14795b89ad90SMark F. Adams #undef __FUNCT__ 14805b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 14815b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 14825b89ad90SMark F. Adams { 14835b89ad90SMark F. Adams PetscErrorCode ierr; 14845b89ad90SMark F. Adams PC_GAMG *pc_gamg; 14855b89ad90SMark F. Adams PC_MG *mg; 14860cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 14879d5b6da9SMark F. Adams static long count = 0; 14885ee9c036SSatish Balay #endif 14895b89ad90SMark F. Adams 14905b89ad90SMark F. Adams PetscFunctionBegin; 149160a6d8f8SMark F. Adams 14925b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 14935b89ad90SMark F. Adams ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 14945b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 14955b89ad90SMark F. Adams 14965b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 14975b89ad90SMark F. Adams ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr); 14985b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 1499ce4cda84SJed Brown mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 15005b89ad90SMark F. Adams mg->innerctx = pc_gamg; 15015b89ad90SMark F. Adams 15029d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 15039d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 15049d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 15059d5b6da9SMark F. Adams pc_gamg->data = 0; 15065b89ad90SMark F. Adams 15079d5b6da9SMark F. Adams /* register AMG type */ 15089d5b6da9SMark F. Adams if (!GAMGList){ 15099d5b6da9SMark F. Adams ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 15109d5b6da9SMark F. Adams ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 15119d5b6da9SMark F. Adams } 15129d5b6da9SMark F. Adams 15139d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 15145b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 15155b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 15165b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 15175b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 15185b89ad90SMark F. Adams 15195b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1520676e1743SMark F. Adams "PCGAMGSetProcEqLim_C", 1521676e1743SMark F. Adams "PCGAMGSetProcEqLim_GAMG", 1522806fa848SBarry Smith PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1523676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1524389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_C", 1525389730f3SMark F. Adams "PCGAMGSetCoarseEqLim_GAMG", 1526806fa848SBarry Smith PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1527389730f3SMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 15288263b398SMark F. Adams "PCGAMGSetRepartitioning_C", 15298263b398SMark F. Adams "PCGAMGSetRepartitioning_GAMG", 1530806fa848SBarry Smith PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 1531676e1743SMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1532*dfd5c07aSMark F. Adams "PCGAMGSetReuseProl_C", 1533*dfd5c07aSMark F. Adams "PCGAMGSetReuseProl_GAMG", 1534*dfd5c07aSMark F. Adams PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr); 1535*dfd5c07aSMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1536ffc955d6SMark F. Adams "PCGAMGSetUseASMAggs_C", 1537ffc955d6SMark F. Adams "PCGAMGSetUseASMAggs_GAMG", 1538806fa848SBarry Smith PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1539ffc955d6SMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 15403542efc5SMark F. Adams "PCGAMGSetThreshold_C", 15413542efc5SMark F. Adams "PCGAMGSetThreshold_GAMG", 1542806fa848SBarry Smith PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 15433542efc5SMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 15449d5b6da9SMark F. Adams "PCGAMGSetType_C", 15459d5b6da9SMark F. Adams "PCGAMGSetType_GAMG", 1546806fa848SBarry Smith PCGAMGSetType_GAMG);CHKERRQ(ierr); 15479d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1548*dfd5c07aSMark F. Adams pc_gamg->reuse_prol = PETSC_TRUE; 1549ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 15509d5b6da9SMark F. Adams pc_gamg->min_eq_proc = 100; 1551c8b0795cSMark F. Adams pc_gamg->coarse_eq_limit = 800; 1552a2f3521dSMark F. Adams pc_gamg->threshold = 0.001; 15539d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 15549d5b6da9SMark F. Adams pc_gamg->verbose = 0; 15559d5b6da9SMark F. Adams pc_gamg->emax_id = -1; 15565e7c91beSJed Brown pc_gamg->eigtarget[0] = 0.05; 15575e7c91beSJed Brown pc_gamg->eigtarget[1] = 1.05; 15589d5b6da9SMark F. Adams 15590cbbd2e1SMark F. Adams /* private events */ 15600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1561785cba28SMark F. Adams if (count++ == 0) { 1562806fa848SBarry Smith ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 1563806fa848SBarry Smith ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 15640cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 15650cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 15660cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 1567806fa848SBarry Smith ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 1568806fa848SBarry Smith ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 1569806fa848SBarry Smith ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1570806fa848SBarry Smith ierr = PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 1571806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 1572806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 1573806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 1574806fa848SBarry Smith ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 1575806fa848SBarry Smith ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 1576806fa848SBarry Smith ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 1577806fa848SBarry Smith ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 1578806fa848SBarry Smith ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 1579f852f58cSMark F. Adams 15800cbbd2e1SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 15810cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 15820cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1583b4fbaa2aSMark F. Adams /* create timer stages */ 1584b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1585b4fbaa2aSMark F. Adams { 1586b4fbaa2aSMark F. Adams char str[32]; 1587b4fbaa2aSMark F. Adams PetscInt lidx; 1588806fa848SBarry Smith sprintf(str,"MG Level %d (finest)",0); 1589806fa848SBarry Smith ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 1590b4fbaa2aSMark F. Adams for (lidx=1;lidx<9;lidx++){ 1591b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1592806fa848SBarry Smith ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 1593b4fbaa2aSMark F. Adams } 1594b4fbaa2aSMark F. Adams } 1595b4fbaa2aSMark F. Adams #endif 1596b4fbaa2aSMark F. Adams } 1597b4fbaa2aSMark F. Adams #endif 15980cbbd2e1SMark F. Adams /* general events */ 15990cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 1600806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr); 1601806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr); 1602806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1603806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1604806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1605806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1606806fa848SBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr); 1607806fa848SBarry Smith ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr); 16080cbbd2e1SMark F. Adams #endif 16090cbbd2e1SMark F. Adams 161060a6d8f8SMark F. Adams /* instantiate derived type */ 161160a6d8f8SMark F. Adams ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr); 161260a6d8f8SMark F. Adams { 161360a6d8f8SMark F. Adams char tname[256] = GAMGAGG; 1614806fa848SBarry Smith ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), PETSC_NULL);CHKERRQ(ierr); 161560a6d8f8SMark F. Adams ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr); 161660a6d8f8SMark F. Adams } 161760a6d8f8SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 16185b89ad90SMark F. Adams PetscFunctionReturn(0); 16195b89ad90SMark F. Adams } 16205b89ad90SMark F. Adams EXTERN_C_END 1621