15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4aaa7dc30SBarry Smith #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> 75b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */ 8f96513f1SMatthew G Knepley 90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 11b4fbaa2aSMark F. Adams #endif 120cbbd2e1SMark F. Adams 130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG; 150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO; 160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 200cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG; 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 30140e18c1SBarry Smith static PetscFunctionList GAMGList = 0; 313e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 329d5b6da9SMark F. Adams 33d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 34d3d6bff4SMark F. Adams #undef __FUNCT__ 35d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 37d3d6bff4SMark F. Adams { 38d3d6bff4SMark F. Adams PetscErrorCode ierr; 39d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 40d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 41d3d6bff4SMark F. Adams 42d3d6bff4SMark F. Adams PetscFunctionBegin; 43a2f3521dSMark F. Adams if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */ 44ce94432eSBarry Smith PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__); 459d5b6da9SMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 4658471d46SMark F. Adams } 470298fd71SBarry Smith pc_gamg->data = NULL; pc_gamg->data_sz = 0; 48878e152fSMark F. Adams 49878e152fSMark F. Adams if (pc_gamg->orig_data) { 50878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 51878e152fSMark F. Adams } 52a2f3521dSMark F. Adams PetscFunctionReturn(0); 53a2f3521dSMark F. Adams } 54a2f3521dSMark F. Adams 555b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 565b89ad90SMark F. Adams /* 57a147abb0SMark F. Adams createLevel: create coarse op with RAP. repartition and/or reduce number 58a147abb0SMark F. Adams of active processors. 595b89ad90SMark F. Adams 605b89ad90SMark F. Adams Input Parameter: 61a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 62a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 639d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 64c5bfad50SMark F. Adams . cr_bs - coarse block size 659d5b6da9SMark F. Adams . isLast - 663530afc2SMark F. Adams In/Output Parameter: 67a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 68afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 6911e60469SMark F. Adams Output Parameter: 703530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 715b89ad90SMark F. Adams */ 725cb416c2SMark F. Adams 735b89ad90SMark F. Adams #undef __FUNCT__ 748263b398SMark F. Adams #define __FUNCT__ "createLevel" 757700e67bSMark Adams static PetscErrorCode createLevel(const PC pc,const Mat Amat_fine,const PetscInt cr_bs,const PetscBool isLast, 767700e67bSMark Adams Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc) 775b89ad90SMark F. Adams { 78a2f3521dSMark F. Adams PetscErrorCode ierr; 799d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 80486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 819d5b6da9SMark F. Adams const PetscBool repart = pc_gamg->repart; 829d5b6da9SMark F. Adams const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 83a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 843b4367a7SBarry Smith MPI_Comm comm; 85c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 863ae0bb68SMark Adams PetscInt ncrs_eq,ncrs,f_bs; 875b89ad90SMark F. Adams 885b89ad90SMark F. Adams PetscFunctionBegin; 893b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 903b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 913b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 92c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 9311e60469SMark F. Adams /* RAP */ 949d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 95038e3b61SMark F. Adams 963ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 970298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 983ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 993ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 1003ae0bb68SMark Adams } 1013ae0bb68SMark Adams else { 1023ae0bb68SMark Adams PetscInt bs; 1033ae0bb68SMark Adams ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 1043ae0bb68SMark Adams ncrs = ncrs_eq/bs; 1053ae0bb68SMark Adams } 106a2f3521dSMark F. Adams 107c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 108a2f3521dSMark F. Adams { 109472110cdSMark F. Adams PetscInt ncrs_eq_glob; 1100298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 111c5df96a5SBarry Smith new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 112c5df96a5SBarry Smith if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1; 113c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 114c5df96a5SBarry Smith if (isLast) new_size = 1; 115a2f3521dSMark F. Adams } 116f852f58cSMark F. Adams 1172fa5cd67SKarl Rupp if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 1182fa5cd67SKarl Rupp else { 1193ae0bb68SMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old; 120885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 121e33ef3b1SMark F. Adams 12271959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 12371959b99SBarry Smith if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 1240cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1250cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 126b4fbaa2aSMark F. Adams #endif 127a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 128785e854fSJed Brown ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 1297700e67bSMark Adams if (repart) { 130a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1315a9b9e01SMark F. Adams Mat adj; 1325a9b9e01SMark F. Adams 133a2f3521dSMark F. Adams if (pc_gamg->verbose>0) { 1343b4367a7SBarry Smith if (pc_gamg->verbose==1) PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq); 135a2f3521dSMark F. Adams else { 136a2f3521dSMark F. Adams PetscInt n; 1373b4367a7SBarry Smith ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr); 1383b4367a7SBarry Smith PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n); 139a2f3521dSMark F. Adams } 140a2f3521dSMark F. Adams } 1415a9b9e01SMark F. Adams 142a2f3521dSMark F. Adams /* get 'adj' */ 143c5bfad50SMark F. Adams if (cr_bs == 1) { 144038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 145806fa848SBarry Smith } else { 146a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 147eb07cef2SMark F. Adams Mat tMat; 148a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 149b4fbaa2aSMark F. Adams const PetscScalar *vals; 150b4fbaa2aSMark F. Adams const PetscInt *idx; 151a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 1529057884aSMark F. Adams static PetscInt llev = 0; 153b4fbaa2aSMark F. Adams 154578f55a3SPeter Brune ierr = PetscMalloc1(ncrs, &d_nnz);CHKERRQ(ierr); 155578f55a3SPeter Brune ierr = PetscMalloc1(ncrs, &o_nnz);CHKERRQ(ierr); 156a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 157a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 158c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 15958471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 160c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 161c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 16258471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 1633ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 1643ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 16558471d46SMark F. Adams } 1666876a03eSMark F. Adams 1673b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 1683ae0bb68SMark Adams ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 169a2f3521dSMark F. Adams ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr); 170a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 171a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 17258471d46SMark F. Adams ierr = PetscFree(d_nnz);CHKERRQ(ierr); 1735f8cf99dSMark F. Adams ierr = PetscFree(o_nnz);CHKERRQ(ierr); 174eb07cef2SMark F. Adams 175a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 176c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 17722063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 178eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 179c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 180eb07cef2SMark F. Adams PetscScalar v = 1.0; 181eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 182eb07cef2SMark F. Adams } 18322063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 184eb07cef2SMark F. Adams } 185eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 186eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 187eb07cef2SMark F. Adams 188b4fbaa2aSMark F. Adams if (llev++ == -1) { 189b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 1908caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 1913b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 192b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 1933bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 194b4fbaa2aSMark F. Adams } 195b4fbaa2aSMark F. Adams 196eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 197eb07cef2SMark F. Adams 198eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 199a2f3521dSMark F. Adams } /* create 'adj' */ 200f150b916SMark F. Adams 201a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2025a9b9e01SMark F. Adams char prefix[256]; 2035a9b9e01SMark F. Adams const char *pcpre; 204b4fbaa2aSMark F. Adams const PetscInt *is_idx; 205b4fbaa2aSMark F. Adams MatPartitioning mpart; 206a4b7d37bSMark F. Adams IS proc_is; 207a2f3521dSMark F. Adams PetscInt targetPE; 2082f03bc48SMark F. Adams 2093b4367a7SBarry Smith ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 2105ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 2119d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 2128caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 21359a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 21411e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 215c5df96a5SBarry Smith ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 216a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 21711e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 2185a9b9e01SMark F. Adams 2195ef31b24SMark F. Adams /* collect IS info */ 220785e854fSJed Brown ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 221a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 222a2f3521dSMark F. Adams targetPE = 1; /* bring to "front" of machine */ 223c5df96a5SBarry Smith /*targetPE = size/new_size;*/ /* spread partitioning across machine */ 224a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 225c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 226a2f3521dSMark F. Adams newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 227eb07cef2SMark F. Adams } 2285ef31b24SMark F. Adams } 229a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 230a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2315ef31b24SMark F. Adams } 2325ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2335a9b9e01SMark F. Adams 2343b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2358263b398SMark F. Adams if (newproc_idx != 0) { 2368263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 2375ef31b24SMark F. Adams } 238806fa848SBarry Smith } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 239a2f3521dSMark F. Adams 240a2f3521dSMark F. Adams PetscInt rfactor,targetPE; 2415a9b9e01SMark F. Adams /* find factor */ 242c5df96a5SBarry Smith if (new_size == 1) rfactor = size; /* easy */ 2435a9b9e01SMark F. Adams else { 2445a9b9e01SMark F. Adams PetscReal best_fact = 0.; 2455a9b9e01SMark F. Adams jj = -1; 246c5df96a5SBarry Smith for (kk = 1 ; kk <= size ; kk++) { 247c5df96a5SBarry Smith if (size%kk==0) { /* a candidate */ 248c5df96a5SBarry Smith PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 2495a9b9e01SMark F. Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 2505a9b9e01SMark F. Adams if (fact > best_fact) { 2515a9b9e01SMark F. Adams best_fact = fact; jj = kk; 2525a9b9e01SMark F. Adams } 2535a9b9e01SMark F. Adams } 2545a9b9e01SMark F. Adams } 2555a9b9e01SMark F. Adams if (jj != -1) rfactor = jj; 256a2f3521dSMark F. Adams else rfactor = 1; /* does this happen .. a prime */ 2575a9b9e01SMark F. Adams } 258c5df96a5SBarry Smith new_size = size/rfactor; 2595a9b9e01SMark F. Adams 260c5df96a5SBarry Smith if (new_size==nactive) { 261a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 2625a9b9e01SMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 263a2f3521dSMark F. Adams if (pc_gamg->verbose>0) { 2643b4367a7SBarry Smith PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq); 265a2f3521dSMark F. Adams } 2665a9b9e01SMark F. Adams PetscFunctionReturn(0); 2675a9b9e01SMark F. Adams } 2685a9b9e01SMark F. Adams 2693b4367a7SBarry Smith if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq); 270c5df96a5SBarry Smith targetPE = rank/rfactor; 2713b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 272a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 273e33ef3b1SMark F. Adams 27411e60469SMark F. Adams /* 275a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 27611e60469SMark F. Adams */ 277a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 2787700e67bSMark Adams is_eq_num_prim = is_eq_num; 27911e60469SMark F. Adams /* 280a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 28111e60469SMark F. Adams */ 282c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 283c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 284a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 2853ae0bb68SMark Adams ncrs_new = ncrs_eq_new/cr_bs; /* eqs */ 286a2f3521dSMark F. Adams 287a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 2880cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 2890cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 290b4fbaa2aSMark F. Adams #endif 291885364a3SMark Adams /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */ 292885364a3SMark Adams { 293885364a3SMark Adams Vec src_crd, dest_crd; 294885364a3SMark 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; 295885364a3SMark Adams VecScatter vecscat; 296885364a3SMark Adams PetscScalar *array; 297885364a3SMark Adams IS isscat; 298a2f3521dSMark F. Adams 299a2f3521dSMark F. Adams /* move data (for primal equations only) */ 30022063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3013b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 3023ae0bb68SMark Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 303c0dedaeaSBarry Smith ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 30411e60469SMark F. Adams /* 3059d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 306c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 30711e60469SMark F. Adams */ 308578f55a3SPeter Brune ierr = PetscMalloc1((ncrs*node_data_sz), &tidx);CHKERRQ(ierr); 309a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3103ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 311c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 312a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 31311e60469SMark F. Adams } 314a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 3153ae0bb68SMark Adams ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 31692a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 31711e60469SMark F. Adams /* 31811e60469SMark F. Adams Create a vector to contain the original vertex information for each element 31911e60469SMark F. Adams */ 3203ae0bb68SMark Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 3219d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3223ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3233ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3249d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 325a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 326c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 327676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 328d3d6bff4SMark F. Adams } 329038e3b61SMark F. Adams } 330eb07cef2SMark F. Adams } 331eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 332eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 33311e60469SMark F. Adams /* 33411e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 33511e60469SMark F. Adams to the correct processor 33611e60469SMark F. Adams */ 3370298fd71SBarry Smith ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 33811e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 33911e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34011e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34111e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 34211e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 34311e60469SMark F. Adams /* 34411e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 34511e60469SMark F. Adams */ 346c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 347578f55a3SPeter Brune ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 3482fa5cd67SKarl Rupp 3493ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3503ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3512fa5cd67SKarl Rupp 35211e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 3539d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3543ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3559d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 356a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 357c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 358d3d6bff4SMark F. Adams } 359038e3b61SMark F. Adams } 360038e3b61SMark F. Adams } 36111e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 36211e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 363885364a3SMark Adams } 364a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 3650cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3660cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 367ed3f9983SMark F. Adams #endif 368a2f3521dSMark F. Adams 36911e60469SMark F. Adams /* 37011e60469SMark F. Adams Invert for MatGetSubMatrix 37111e60469SMark F. Adams */ 372a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 373a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 374c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 375a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 376a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 377a2f3521dSMark F. Adams } 378a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 3790cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3800cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 3810cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 382ed3f9983SMark F. Adams #endif 383a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 384a2f3521dSMark F. Adams { 385a2f3521dSMark F. Adams Mat mat; 386806fa848SBarry Smith ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 387a2f3521dSMark F. Adams *a_Amat_crs = mat; 388c5bfad50SMark F. Adams 389c5bfad50SMark F. Adams if (!PETSC_TRUE) { 390c5bfad50SMark F. Adams PetscInt cbs, rbs; 391c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr); 392c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 393c5bfad50SMark F. Adams ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr); 394c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr); 395c5bfad50SMark F. Adams } 396a2f3521dSMark F. Adams } 397038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 398a2f3521dSMark F. Adams 3990cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4000cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 401ed3f9983SMark F. Adams #endif 40211e60469SMark F. Adams /* prolongator */ 40311e60469SMark F. Adams { 40411e60469SMark F. Adams IS findices; 405a2f3521dSMark F. Adams PetscInt Istart,Iend; 406a2f3521dSMark F. Adams Mat Pnew; 407a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 4080cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4090cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 410ed3f9983SMark F. Adams #endif 4113b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 412c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 413806fa848SBarry Smith ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 41411e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 415c5bfad50SMark F. Adams 416c5bfad50SMark F. Adams if (!PETSC_TRUE) { 417c5bfad50SMark F. Adams PetscInt cbs, rbs; 418c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr); 419c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 420c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr); 421c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 422c5bfad50SMark F. Adams } 4230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4240cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 425ed3f9983SMark F. Adams #endif 4263530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 427a2f3521dSMark F. Adams 428a2f3521dSMark F. Adams /* output - repartitioned */ 429a2f3521dSMark F. Adams *a_P_inout = Pnew; 430e33ef3b1SMark F. Adams } 431a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4325b89ad90SMark F. Adams 433c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 434a2f3521dSMark F. Adams } 4355a9b9e01SMark F. Adams 436a2f3521dSMark F. Adams /* outout matrix data */ 437c8b0795cSMark F. Adams if (!PETSC_TRUE) { 438c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 439c8b0795cSMark F. Adams if (llev==0) { 440c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 4413b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 442c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 443c8b0795cSMark F. Adams ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr); 444c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 445c8b0795cSMark F. Adams } 446c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 4473b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 448c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 449c8b0795cSMark F. Adams ierr = MatView(Cmat, viewer);CHKERRQ(ierr); 450c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 451c8b0795cSMark F. Adams } 4525b89ad90SMark F. Adams PetscFunctionReturn(0); 4535b89ad90SMark F. Adams } 4545b89ad90SMark F. Adams 4555b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4565b89ad90SMark F. Adams /* 4575b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4585b89ad90SMark F. Adams by setting data structures and options. 4595b89ad90SMark F. Adams 4605b89ad90SMark F. Adams Input Parameter: 4615b89ad90SMark F. Adams . pc - the preconditioner context 4625b89ad90SMark F. Adams 4635b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 4645b89ad90SMark F. Adams 4655b89ad90SMark F. Adams Notes: 4665b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 4675b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 4685b89ad90SMark F. Adams */ 4695b89ad90SMark F. Adams #undef __FUNCT__ 4705b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 4719d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 4725b89ad90SMark F. Adams { 4735b89ad90SMark F. Adams PetscErrorCode ierr; 4749d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 4755b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4762adcac29SMark F. Adams Mat Pmat = pc->pmat; 477a2f3521dSMark F. Adams PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 4783b4367a7SBarry Smith MPI_Comm comm; 479c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 480587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 481c8b0795cSMark F. Adams PetscReal emaxs[GAMG_MAXLEVELS]; 482e696ed0bSMark F. Adams IS *ASMLocalIDsArr[GAMG_MAXLEVELS]; 483a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 484737a81a9SMark F. Adams MatInfo info; 4857700e67bSMark Adams PetscBool redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol); 4865ef31b24SMark F. Adams 4875b89ad90SMark F. Adams PetscFunctionBegin; 4883b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 4893b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 4903b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 491dfd5c07aSMark F. Adams 4923b4367a7SBarry Smith if (pc_gamg->verbose>2) PetscPrintf(comm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled); 493dfd5c07aSMark F. Adams 49484d3f75bSMark F. Adams if (pc_gamg->setup_count++ > 0) { 495878e152fSMark F. Adams if (redo_mesh_setup) { 496878e152fSMark F. Adams /* reset everything */ 497878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 498878e152fSMark F. Adams pc->setupcalled = 0; 499806fa848SBarry Smith } else { 50084d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 50103a628feSMark F. Adams /* just do Galerkin grids */ 50258471d46SMark F. Adams Mat B,dA,dB; 50358471d46SMark F. Adams 50471959b99SBarry Smith if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 5059d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 50658471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 50723ee1639SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 50858471d46SMark F. Adams /* (re)set to get dirty flag */ 50923ee1639SBarry Smith ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 51058471d46SMark F. Adams 5112fb0b348SMark F. Adams for (level=pc_gamg->Nlevels-2; level>=0; level--) { 51203a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 5130a97e771SToby Isaac if (pc_gamg->setup_count==2) { 51403a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 515084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 5162fa5cd67SKarl Rupp 51703a628feSMark F. Adams mglevels[level]->A = B; 518806fa848SBarry Smith } else { 51923ee1639SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 52058471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 52103a628feSMark F. Adams } 52223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 52358471d46SMark F. Adams dB = B; 52458471d46SMark F. Adams } 5255f8cf99dSMark F. Adams } 526d5280255SMark F. Adams 527d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 528d5280255SMark F. Adams 529d5280255SMark F. Adams /* PCSetUp_MG seems to insists on setting this to GMRES */ 5301ac9965eSMark F. Adams ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr); 53158471d46SMark F. Adams PetscFunctionReturn(0); 532eb07cef2SMark F. Adams } 533878e152fSMark F. Adams } 534f6536408SMark F. Adams 535878e152fSMark F. Adams if (!pc_gamg->data) { 536878e152fSMark F. Adams if (pc_gamg->orig_data) { 537878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 5380298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 5392fa5cd67SKarl Rupp 540878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 541878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 542878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5432fa5cd67SKarl Rupp 544785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 545878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 546806fa848SBarry Smith } else { 5471ab5ffc9SJed Brown if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 5487700e67bSMark Adams ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 5499d5b6da9SMark F. Adams } 550878e152fSMark F. Adams } 551878e152fSMark F. Adams 552878e152fSMark F. Adams /* cache original data for reuse */ 553878e152fSMark F. Adams if (!pc_gamg->orig_data && redo_mesh_setup) { 554785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 555878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 556878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 557878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 558878e152fSMark F. Adams } 559038e3b61SMark F. Adams 560302f38e8SMark F. Adams /* get basic dims */ 561302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 562a2f3521dSMark F. Adams 563a2f3521dSMark F. Adams ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr); 564c8b0795cSMark F. Adams if (pc_gamg->verbose) { 56584f9421dSMark F. Adams PetscInt NN = M; 56684f9421dSMark F. Adams if (pc_gamg->verbose==1) { 56784f9421dSMark F. Adams ierr = MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr); 5683bf036e2SBarry Smith ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr); 569806fa848SBarry Smith } else { 570806fa848SBarry Smith ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 57184f9421dSMark F. Adams } 572b2a4f308SMark F. Adams nnz0 = info.nz_used; 573b2a4f308SMark F. Adams nnztot = info.nz_used; 5743b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 575c5df96a5SBarry Smith rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 576c5df96a5SBarry Smith (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr); 577c8b0795cSMark F. Adams } 57884d3f75bSMark F. Adams 579a2f3521dSMark F. Adams /* Get A_i and R_i */ 580c5df96a5SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */ 581c5df96a5SBarry Smith level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (size==1 || nactivepe>1); */ 5820205a208SMark F. Adams level++) { 5835b89ad90SMark F. Adams level1 = level + 1; 5840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5850cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 586a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 587a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 588b4fbaa2aSMark F. Adams #endif 589a2f3521dSMark F. Adams #endif 590c8b0795cSMark F. Adams { /* construct prolongator */ 591725b86d8SJed Brown Mat Gmat; 5920cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 5937700e67bSMark Adams Mat Prol11; 594c8b0795cSMark F. Adams 5957700e67bSMark Adams ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 5961ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 5977700e67bSMark Adams ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 598c8b0795cSMark F. Adams 599a2f3521dSMark F. Adams /* could have failed to create new level */ 600a2f3521dSMark F. Adams if (Prol11) { 6019d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6020298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 603a2f3521dSMark F. Adams 6041ab5ffc9SJed Brown if (pc_gamg->ops->optprol) { 605c8b0795cSMark F. Adams /* smooth */ 6067700e67bSMark Adams ierr = pc_gamg->ops->optprol(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 607c8b0795cSMark F. Adams } 608c8b0795cSMark F. Adams 6097700e67bSMark Adams Parr[level1] = Prol11; 6100298fd71SBarry Smith } else Parr[level1] = NULL; 611ffc955d6SMark F. Adams 612ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 6131b18a24aSMark Adams PetscInt bs; 6141b18a24aSMark Adams ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 615806fa848SBarry Smith ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 616ffc955d6SMark F. Adams } 617ffc955d6SMark F. Adams 618a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 61941b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 620a2f3521dSMark F. Adams } /* construct prolongator scope */ 6210cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6220cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 623c8b0795cSMark F. Adams #endif 6249d5b6da9SMark F. Adams /* cache eigen estimate */ 6259d5b6da9SMark F. Adams if (pc_gamg->emax_id != -1) { 6269d5b6da9SMark F. Adams PetscBool flag; 6277700e67bSMark Adams ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr); 6289d5b6da9SMark F. Adams if (!flag) emaxs[level] = -1.; 629806fa848SBarry Smith } else emaxs[level] = -1.; 6302adcac29SMark F. Adams if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 631c8b0795cSMark F. Adams if (!Parr[level1]) { 632806fa848SBarry Smith if (pc_gamg->verbose) { 6333b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr); 634806fa848SBarry Smith } 635c8b0795cSMark F. Adams break; 636c8b0795cSMark F. Adams } 6370cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6380cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 639b4fbaa2aSMark F. Adams #endif 640a2f3521dSMark F. Adams 641a2f3521dSMark F. Adams ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2), 6427700e67bSMark Adams &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr); 643a2f3521dSMark F. Adams 6440cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6450cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 646b4fbaa2aSMark F. Adams #endif 647a2f3521dSMark F. Adams ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr); 648a2f3521dSMark F. Adams 649a2f3521dSMark F. Adams if (pc_gamg->verbose > 0) { 6500cbbd2e1SMark F. Adams PetscInt NN = M; 6510cbbd2e1SMark F. Adams if (pc_gamg->verbose==1) { 652a2f3521dSMark F. Adams ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr); 6533bf036e2SBarry Smith ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr); 654806fa848SBarry Smith } else { 655806fa848SBarry Smith ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 6560cbbd2e1SMark F. Adams } 657a2f3521dSMark F. Adams 6580cbbd2e1SMark F. Adams nnztot += info.nz_used; 6593b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 660c5df96a5SBarry Smith rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols, 661806fa848SBarry Smith (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr); 662c8b0795cSMark F. Adams } 663a2f3521dSMark F. Adams 664a2f3521dSMark F. Adams /* stop if one node -- could pull back for singular problems */ 6653ae0bb68SMark Adams if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M < 2)) { 66681708dccSMark F. Adams level++; 66781708dccSMark F. Adams break; 66881708dccSMark F. Adams } 6690cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 670b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 671b4fbaa2aSMark F. Adams #endif 672c8b0795cSMark F. Adams } /* levels */ 673c8b0795cSMark F. Adams 674c8b0795cSMark F. Adams if (pc_gamg->data) { 675c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 6760298fd71SBarry Smith pc_gamg->data = NULL; 6775b89ad90SMark F. Adams } 678c8b0795cSMark F. Adams 6793b4367a7SBarry Smith if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 6809d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6815b89ad90SMark F. Adams fine_level = level; 6820298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 6835b89ad90SMark F. Adams 68484d3f75bSMark F. Adams /* simple setup */ 68584d3f75bSMark F. Adams if (!PETSC_TRUE) { 68684d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 68784d3f75bSMark F. Adams for (lidx=0,level=pc_gamg->Nlevels-1; 68884d3f75bSMark F. Adams lidx<fine_level; 68984d3f75bSMark F. Adams lidx++, level--) { 69084d3f75bSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr); 69123ee1639SBarry Smith ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr); 69284d3f75bSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 69384d3f75bSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 69484d3f75bSMark F. Adams } 69523ee1639SBarry Smith ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr); 69684d3f75bSMark F. Adams 69784d3f75bSMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 698806fa848SBarry Smith } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 699d5280255SMark F. Adams /* set default smoothers & set operators */ 7009d5b6da9SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 701587fa25dSMark F. Adams lidx <= fine_level; 702587fa25dSMark F. Adams lidx++, level--) { 703ffc955d6SMark F. Adams KSP smoother; 704ffc955d6SMark F. Adams PC subpc; 705a2f3521dSMark F. Adams 7069d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 707f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 708ffc955d6SMark F. Adams 709a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 710a2f3521dSMark F. Adams /* set ops */ 71123ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 712a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 713a2f3521dSMark F. Adams 714a2f3521dSMark F. Adams /* set defaults */ 7156c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 716a2f3521dSMark F. Adams 7171b18a24aSMark Adams /* set blocks for GASM smoother that uses the 'aggregates' */ 718ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 7192d3561bbSSatish Balay PetscInt sz; 7202d3561bbSSatish Balay IS *is; 721a2f3521dSMark F. Adams 7222d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7232d3561bbSSatish Balay is = ASMLocalIDsArr[level]; 724ffc955d6SMark F. Adams ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr); 7251b18a24aSMark Adams ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr); 726ffc955d6SMark F. Adams if (sz==0) { 727ffc955d6SMark F. Adams IS is; 728ffc955d6SMark F. Adams PetscInt my0,kk; 729ffc955d6SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr); 730ffc955d6SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 7310298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr); 732a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 733806fa848SBarry Smith } else { 734a94c3b12SMark F. Adams PetscInt kk; 7350298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr); 736a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 737a94c3b12SMark F. Adams ierr = ISDestroy(&is[kk]);CHKERRQ(ierr); 738a94c3b12SMark F. Adams } 739ffc955d6SMark F. Adams ierr = PetscFree(is);CHKERRQ(ierr); 740ffc955d6SMark F. Adams } 7410298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 742ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 743ffc955d6SMark F. Adams ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr); 744806fa848SBarry Smith } else { 745890ffe84SMark Adams ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 746ffc955d6SMark F. Adams } 747d5280255SMark F. Adams } 748d5280255SMark F. Adams { 749d5280255SMark F. Adams /* coarse grid */ 750d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 751d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 752d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 75323ee1639SBarry Smith ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 754d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 755d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 756d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 757d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 75871959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 75971959b99SBarry Smith if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 760d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 761d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 7629dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 7632fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 7645b42dca8SJed Brown /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 7655b42dca8SJed Brown * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 7665b42dca8SJed Brown * view every subdomain as though they were different. */ 7675b42dca8SJed Brown ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 768d5280255SMark F. Adams } 769d5280255SMark F. Adams 770d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 771d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 772d5280255SMark F. Adams ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); 773d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 7743b4367a7SBarry Smith if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 775d5280255SMark F. Adams 776d5280255SMark F. Adams /* create cheby smoothers */ 777d5280255SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 778d5280255SMark F. Adams lidx <= fine_level; 779d5280255SMark F. Adams lidx++, level--) { 780d5280255SMark F. Adams KSP smoother; 781890ffe84SMark Adams PetscBool flag,flag2; 782d5280255SMark F. Adams PC subpc; 783d5280255SMark F. Adams 784ffc955d6SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 785a2f3521dSMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 786a2f3521dSMark F. Adams 787ffc955d6SMark F. Adams /* do my own cheby */ 7886c9de887SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr); 789ffc955d6SMark F. Adams if (flag) { 790ffc955d6SMark F. Adams PetscReal emax, emin; 791251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr); 792890ffe84SMark Adams ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr); 793890ffe84SMark Adams if ((flag||flag2) && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC but lets acccept SOR because it is close and safe (always lower) */ 794890ffe84SMark Adams else { /* eigen estimate 'emax' -- this is done in cheby */ 795e696ed0bSMark F. Adams KSP eksp; 796e696ed0bSMark F. Adams Mat Lmat = Aarr[level]; 797b2a4f308SMark F. Adams Vec bb, xx; 798038e3b61SMark F. Adams 7992a7a6963SBarry Smith ierr = MatCreateVecs(Lmat, &bb, 0);CHKERRQ(ierr); 8002a7a6963SBarry Smith ierr = MatCreateVecs(Lmat, &xx, 0);CHKERRQ(ierr); 801fc4362bfSMark F. Adams { 802fc4362bfSMark F. Adams PetscRandom rctx; 8033b4367a7SBarry Smith ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr); 804fc4362bfSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 805fc4362bfSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 806fc4362bfSMark F. Adams ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 8075b89ad90SMark F. Adams } 808a2f3521dSMark F. Adams 809e696ed0bSMark F. Adams /* zeroing out BC rows -- needed for crazy matrices */ 810e696ed0bSMark F. Adams { 811e696ed0bSMark F. Adams PetscInt Istart,Iend,ncols,jj,Ii; 812e696ed0bSMark F. Adams PetscScalar zero = 0.0; 813e696ed0bSMark F. Adams ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr); 814e696ed0bSMark F. Adams for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 815e696ed0bSMark F. Adams ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 816e696ed0bSMark F. Adams if (ncols <= 1) { 817e696ed0bSMark F. Adams ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 818a94c3b12SMark F. Adams } 819e696ed0bSMark F. Adams ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 820a94c3b12SMark F. Adams } 821a94c3b12SMark F. Adams ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 822a94c3b12SMark F. Adams ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 823a94c3b12SMark F. Adams } 824e696ed0bSMark F. Adams 8253b4367a7SBarry Smith ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr); 826806fa848SBarry Smith ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr); 827fc4362bfSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 8281a166f3bSJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 8291a166f3bSJed Brown ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 830f6536408SMark F. Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 831f6536408SMark F. Adams 832f6536408SMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 83323ee1639SBarry Smith ierr = KSPSetOperators(eksp, Lmat, Lmat);CHKERRQ(ierr); 834fc4362bfSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 8355a9b9e01SMark F. Adams 836d3d0db20SJed Brown /* set PC type to be same as smoother */ 837ffc955d6SMark F. Adams ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr); 838b2a4f308SMark F. Adams 8395a9b9e01SMark F. Adams /* solve - keep stuff out of logging */ 8405a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 8415a9b9e01SMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 842fc4362bfSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 8435a9b9e01SMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 8445a9b9e01SMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 8455a9b9e01SMark F. Adams 846fc4362bfSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 8475a9b9e01SMark F. Adams 848fc4362bfSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 849fc4362bfSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 850fc4362bfSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 851f6536408SMark F. Adams 852ffc955d6SMark F. Adams if (pc_gamg->verbose > 0) { 853a94c3b12SMark F. Adams PetscInt N1, tt; 854a94c3b12SMark F. Adams ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr); 8553b4367a7SBarry Smith PetscPrintf(comm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 856f6536408SMark F. Adams } 857fc4362bfSMark F. Adams } 858038e3b61SMark F. Adams { 859c5bfad50SMark F. Adams PetscInt N1, N0; 8600298fd71SBarry Smith ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr); 8610298fd71SBarry Smith ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr); 862f6536408SMark F. Adams /* heuristic - is this crap? */ 863b4ec6429SMark F. Adams /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */ 8645e7c91beSJed Brown emin = emax * pc_gamg->eigtarget[0]; 8655e7c91beSJed Brown emax *= pc_gamg->eigtarget[1]; 866038e3b61SMark F. Adams } 8676c9de887SHong Zhang ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); 868ffc955d6SMark F. Adams } /* setup checby flag */ 869ffc955d6SMark F. Adams } /* non-coarse levels */ 870737a81a9SMark F. Adams 871d5280255SMark F. Adams /* clean up */ 872d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 873587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 874587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 8755b89ad90SMark F. Adams } 8760cbbd2e1SMark F. Adams 8770cbbd2e1SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 8780cbbd2e1SMark F. Adams 8791ac9965eSMark F. Adams if (PETSC_TRUE) { 88058471d46SMark F. Adams KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 8819d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 88258471d46SMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 88358471d46SMark F. Adams } 884806fa848SBarry Smith } else { 8855f8cf99dSMark F. Adams KSP smoother; 8863b4367a7SBarry Smith if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__); 8879d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 88823ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 8895f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 8909d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 8915f8cf99dSMark F. Adams } 8925b89ad90SMark F. Adams PetscFunctionReturn(0); 8935b89ad90SMark F. Adams } 8945b89ad90SMark F. Adams 895eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 8965b89ad90SMark F. Adams /* 8975b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 8985b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 8995b89ad90SMark F. Adams 9005b89ad90SMark F. Adams Input Parameter: 9015b89ad90SMark F. Adams . pc - the preconditioner context 9025b89ad90SMark F. Adams 9035b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 9045b89ad90SMark F. Adams */ 9055b89ad90SMark F. Adams #undef __FUNCT__ 9065b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 9075b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 9085b89ad90SMark F. Adams { 9095b89ad90SMark F. Adams PetscErrorCode ierr; 9105b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9115b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 9125b89ad90SMark F. Adams 9135b89ad90SMark F. Adams PetscFunctionBegin; 9145b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 9159b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 9169b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 9179b8ffb57SJed Brown } 9181ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 9191ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 9205b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 9215b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 9225b89ad90SMark F. Adams PetscFunctionReturn(0); 9235b89ad90SMark F. Adams } 9245b89ad90SMark F. Adams 925676e1743SMark F. Adams 926676e1743SMark F. Adams #undef __FUNCT__ 927676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 928676e1743SMark F. Adams /*@ 929676e1743SMark F. Adams PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 930676e1743SMark F. Adams processor reduction. 931676e1743SMark F. Adams 932676e1743SMark F. Adams Not Collective on PC 933676e1743SMark F. Adams 934676e1743SMark F. Adams Input Parameters: 935676e1743SMark F. Adams . pc - the preconditioner context 936676e1743SMark F. Adams 937676e1743SMark F. Adams 938676e1743SMark F. Adams Options Database Key: 939676e1743SMark F. Adams . -pc_gamg_process_eq_limit 940676e1743SMark F. Adams 941676e1743SMark F. Adams Level: intermediate 942676e1743SMark F. Adams 943676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 944676e1743SMark F. Adams 945676e1743SMark F. Adams .seealso: () 946676e1743SMark F. Adams @*/ 947676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 948676e1743SMark F. Adams { 949676e1743SMark F. Adams PetscErrorCode ierr; 950676e1743SMark F. Adams 951676e1743SMark F. Adams PetscFunctionBegin; 952676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 953676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 954676e1743SMark F. Adams PetscFunctionReturn(0); 955676e1743SMark F. Adams } 956676e1743SMark F. Adams 957676e1743SMark F. Adams #undef __FUNCT__ 958676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 9591e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 960676e1743SMark F. Adams { 961c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 962c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 963676e1743SMark F. Adams 964676e1743SMark F. Adams PetscFunctionBegin; 9659d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 966676e1743SMark F. Adams PetscFunctionReturn(0); 967676e1743SMark F. Adams } 968676e1743SMark F. Adams 969676e1743SMark F. Adams #undef __FUNCT__ 970389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 971389730f3SMark F. Adams /*@ 972389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 973389730f3SMark F. Adams 974389730f3SMark F. Adams Collective on PC 975389730f3SMark F. Adams 976389730f3SMark F. Adams Input Parameters: 977389730f3SMark F. Adams . pc - the preconditioner context 978389730f3SMark F. Adams 979389730f3SMark F. Adams 980389730f3SMark F. Adams Options Database Key: 981389730f3SMark F. Adams . -pc_gamg_coarse_eq_limit 982389730f3SMark F. Adams 983389730f3SMark F. Adams Level: intermediate 984389730f3SMark F. Adams 985389730f3SMark F. Adams Concepts: Unstructured multrigrid preconditioner 986389730f3SMark F. Adams 987389730f3SMark F. Adams .seealso: () 988389730f3SMark F. Adams @*/ 989389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 990389730f3SMark F. Adams { 991389730f3SMark F. Adams PetscErrorCode ierr; 992389730f3SMark F. Adams 993389730f3SMark F. Adams PetscFunctionBegin; 994389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 995389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 996389730f3SMark F. Adams PetscFunctionReturn(0); 997389730f3SMark F. Adams } 998389730f3SMark F. Adams 999389730f3SMark F. Adams #undef __FUNCT__ 1000389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 10011e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1002389730f3SMark F. Adams { 1003389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1004389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1005389730f3SMark F. Adams 1006389730f3SMark F. Adams PetscFunctionBegin; 10079d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 1008389730f3SMark F. Adams PetscFunctionReturn(0); 1009389730f3SMark F. Adams } 1010389730f3SMark F. Adams 1011389730f3SMark F. Adams #undef __FUNCT__ 10128263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 1013676e1743SMark F. Adams /*@ 10148263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 1015676e1743SMark F. Adams 1016676e1743SMark F. Adams Collective on PC 1017676e1743SMark F. Adams 1018676e1743SMark F. Adams Input Parameters: 1019676e1743SMark F. Adams . pc - the preconditioner context 1020676e1743SMark F. Adams 1021676e1743SMark F. Adams 1022676e1743SMark F. Adams Options Database Key: 10238263b398SMark F. Adams . -pc_gamg_repartition 1024676e1743SMark F. Adams 1025676e1743SMark F. Adams Level: intermediate 1026676e1743SMark F. Adams 1027676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1028676e1743SMark F. Adams 1029676e1743SMark F. Adams .seealso: () 1030676e1743SMark F. Adams @*/ 10318263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1032676e1743SMark F. Adams { 1033676e1743SMark F. Adams PetscErrorCode ierr; 1034676e1743SMark F. Adams 1035676e1743SMark F. Adams PetscFunctionBegin; 1036676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10378263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1038676e1743SMark F. Adams PetscFunctionReturn(0); 1039676e1743SMark F. Adams } 1040676e1743SMark F. Adams 1041676e1743SMark F. Adams #undef __FUNCT__ 10428263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 10431e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 1044676e1743SMark F. Adams { 1045c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1046c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1047676e1743SMark F. Adams 1048676e1743SMark F. Adams PetscFunctionBegin; 10499d5b6da9SMark F. Adams pc_gamg->repart = n; 1050676e1743SMark F. Adams PetscFunctionReturn(0); 1051676e1743SMark F. Adams } 1052676e1743SMark F. Adams 1053676e1743SMark F. Adams #undef __FUNCT__ 1054dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl" 1055dfd5c07aSMark F. Adams /*@ 1056dfd5c07aSMark F. Adams PCGAMGSetReuseProl - Reuse prlongation 1057dfd5c07aSMark F. Adams 1058dfd5c07aSMark F. Adams Collective on PC 1059dfd5c07aSMark F. Adams 1060dfd5c07aSMark F. Adams Input Parameters: 1061dfd5c07aSMark F. Adams . pc - the preconditioner context 1062dfd5c07aSMark F. Adams 1063dfd5c07aSMark F. Adams 1064dfd5c07aSMark F. Adams Options Database Key: 1065dfd5c07aSMark F. Adams . -pc_gamg_reuse_interpolation 1066dfd5c07aSMark F. Adams 1067dfd5c07aSMark F. Adams Level: intermediate 1068dfd5c07aSMark F. Adams 1069dfd5c07aSMark F. Adams Concepts: Unstructured multrigrid preconditioner 1070dfd5c07aSMark F. Adams 1071dfd5c07aSMark F. Adams .seealso: () 1072dfd5c07aSMark F. Adams @*/ 1073dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n) 1074dfd5c07aSMark F. Adams { 1075dfd5c07aSMark F. Adams PetscErrorCode ierr; 1076dfd5c07aSMark F. Adams 1077dfd5c07aSMark F. Adams PetscFunctionBegin; 1078dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1079dfd5c07aSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1080dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1081dfd5c07aSMark F. Adams } 1082dfd5c07aSMark F. Adams 1083dfd5c07aSMark F. Adams #undef __FUNCT__ 1084dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG" 10851e6b0712SBarry Smith static PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n) 1086dfd5c07aSMark F. Adams { 1087dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1088dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1089dfd5c07aSMark F. Adams 1090dfd5c07aSMark F. Adams PetscFunctionBegin; 1091dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1092dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1093dfd5c07aSMark F. Adams } 1094dfd5c07aSMark F. Adams 1095dfd5c07aSMark F. Adams #undef __FUNCT__ 1096ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs" 1097ffc955d6SMark F. Adams /*@ 1098ffc955d6SMark F. Adams PCGAMGSetUseASMAggs - 1099ffc955d6SMark F. Adams 1100ffc955d6SMark F. Adams Collective on PC 1101ffc955d6SMark F. Adams 1102ffc955d6SMark F. Adams Input Parameters: 1103ffc955d6SMark F. Adams . pc - the preconditioner context 1104ffc955d6SMark F. Adams 1105ffc955d6SMark F. Adams 1106ffc955d6SMark F. Adams Options Database Key: 1107ffc955d6SMark F. Adams . -pc_gamg_use_agg_gasm 1108ffc955d6SMark F. Adams 1109ffc955d6SMark F. Adams Level: intermediate 1110ffc955d6SMark F. Adams 1111ffc955d6SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1112ffc955d6SMark F. Adams 1113ffc955d6SMark F. Adams .seealso: () 1114ffc955d6SMark F. Adams @*/ 1115ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1116ffc955d6SMark F. Adams { 1117ffc955d6SMark F. Adams PetscErrorCode ierr; 1118ffc955d6SMark F. Adams 1119ffc955d6SMark F. Adams PetscFunctionBegin; 1120ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1121ffc955d6SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1122ffc955d6SMark F. Adams PetscFunctionReturn(0); 1123ffc955d6SMark F. Adams } 1124ffc955d6SMark F. Adams 1125ffc955d6SMark F. Adams #undef __FUNCT__ 1126ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 11271e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1128ffc955d6SMark F. Adams { 1129ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1130ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1131ffc955d6SMark F. Adams 1132ffc955d6SMark F. Adams PetscFunctionBegin; 1133ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = n; 1134ffc955d6SMark F. Adams PetscFunctionReturn(0); 1135ffc955d6SMark F. Adams } 1136ffc955d6SMark F. Adams 1137ffc955d6SMark F. Adams #undef __FUNCT__ 11384ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 11394ef23d27SMark F. Adams /*@ 11404ef23d27SMark F. Adams PCGAMGSetNlevels - 11414ef23d27SMark F. Adams 11424ef23d27SMark F. Adams Not collective on PC 11434ef23d27SMark F. Adams 11444ef23d27SMark F. Adams Input Parameters: 11454ef23d27SMark F. Adams . pc - the preconditioner context 11464ef23d27SMark F. Adams 11474ef23d27SMark F. Adams 11484ef23d27SMark F. Adams Options Database Key: 11494ef23d27SMark F. Adams . -pc_mg_levels 11504ef23d27SMark F. Adams 11514ef23d27SMark F. Adams Level: intermediate 11524ef23d27SMark F. Adams 11534ef23d27SMark F. Adams Concepts: Unstructured multrigrid preconditioner 11544ef23d27SMark F. Adams 11554ef23d27SMark F. Adams .seealso: () 11564ef23d27SMark F. Adams @*/ 11574ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 11584ef23d27SMark F. Adams { 11594ef23d27SMark F. Adams PetscErrorCode ierr; 11604ef23d27SMark F. Adams 11614ef23d27SMark F. Adams PetscFunctionBegin; 11624ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11634ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 11644ef23d27SMark F. Adams PetscFunctionReturn(0); 11654ef23d27SMark F. Adams } 11664ef23d27SMark F. Adams 11674ef23d27SMark F. Adams #undef __FUNCT__ 11684ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 11691e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 11704ef23d27SMark F. Adams { 11714ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 11724ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 11734ef23d27SMark F. Adams 11744ef23d27SMark F. Adams PetscFunctionBegin; 11759d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 11764ef23d27SMark F. Adams PetscFunctionReturn(0); 11774ef23d27SMark F. Adams } 11784ef23d27SMark F. Adams 11794ef23d27SMark F. Adams #undef __FUNCT__ 11803542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 11813542efc5SMark F. Adams /*@ 11823542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 11833542efc5SMark F. Adams 11843542efc5SMark F. Adams Not collective on PC 11853542efc5SMark F. Adams 11863542efc5SMark F. Adams Input Parameters: 11873542efc5SMark F. Adams . pc - the preconditioner context 11883542efc5SMark F. Adams 11893542efc5SMark F. Adams 11903542efc5SMark F. Adams Options Database Key: 11913542efc5SMark F. Adams . -pc_gamg_threshold 11923542efc5SMark F. Adams 11933542efc5SMark F. Adams Level: intermediate 11943542efc5SMark F. Adams 11953542efc5SMark F. Adams Concepts: Unstructured multrigrid preconditioner 11963542efc5SMark F. Adams 11973542efc5SMark F. Adams .seealso: () 11983542efc5SMark F. Adams @*/ 11993542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 12003542efc5SMark F. Adams { 12013542efc5SMark F. Adams PetscErrorCode ierr; 12023542efc5SMark F. Adams 12033542efc5SMark F. Adams PetscFunctionBegin; 12043542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12053542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 12063542efc5SMark F. Adams PetscFunctionReturn(0); 12073542efc5SMark F. Adams } 12083542efc5SMark F. Adams 12093542efc5SMark F. Adams #undef __FUNCT__ 12103542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 12111e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 12123542efc5SMark F. Adams { 1213c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1214c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12153542efc5SMark F. Adams 12163542efc5SMark F. Adams PetscFunctionBegin; 12179d5b6da9SMark F. Adams pc_gamg->threshold = n; 12183542efc5SMark F. Adams PetscFunctionReturn(0); 12193542efc5SMark F. Adams } 12203542efc5SMark F. Adams 12213542efc5SMark F. Adams #undef __FUNCT__ 12229d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1223676e1743SMark F. Adams /*@ 1224*c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1225676e1743SMark F. Adams 1226676e1743SMark F. Adams Collective on PC 1227676e1743SMark F. Adams 1228676e1743SMark F. Adams Input Parameters: 1229*c60c7ad4SBarry Smith + pc - the preconditioner context 1230*c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1231676e1743SMark F. Adams 1232676e1743SMark F. Adams Options Database Key: 1233*c60c7ad4SBarry Smith . -pc_gamg_type <agg,geo,classical> 1234676e1743SMark F. Adams 1235676e1743SMark F. Adams Level: intermediate 1236676e1743SMark F. Adams 1237676e1743SMark F. Adams Concepts: Unstructured multrigrid preconditioner 1238676e1743SMark F. Adams 1239*c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG 1240676e1743SMark F. Adams @*/ 124119fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1242676e1743SMark F. Adams { 1243676e1743SMark F. Adams PetscErrorCode ierr; 1244676e1743SMark F. Adams 1245676e1743SMark F. Adams PetscFunctionBegin; 1246676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1247806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1248676e1743SMark F. Adams PetscFunctionReturn(0); 1249676e1743SMark F. Adams } 1250676e1743SMark F. Adams 1251676e1743SMark F. Adams #undef __FUNCT__ 1252*c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType" 1253*c60c7ad4SBarry Smith /*@ 1254*c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1255*c60c7ad4SBarry Smith 1256*c60c7ad4SBarry Smith Collective on PC 1257*c60c7ad4SBarry Smith 1258*c60c7ad4SBarry Smith Input Parameter: 1259*c60c7ad4SBarry Smith . pc - the preconditioner context 1260*c60c7ad4SBarry Smith 1261*c60c7ad4SBarry Smith Output Parameter: 1262*c60c7ad4SBarry Smith . type - the type of algorithm used 1263*c60c7ad4SBarry Smith 1264*c60c7ad4SBarry Smith Level: intermediate 1265*c60c7ad4SBarry Smith 1266*c60c7ad4SBarry Smith Concepts: Unstructured multrigrid preconditioner 1267*c60c7ad4SBarry Smith 1268*c60c7ad4SBarry Smith .seealso: PCGAMGSetType() 1269*c60c7ad4SBarry Smith @*/ 1270*c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1271*c60c7ad4SBarry Smith { 1272*c60c7ad4SBarry Smith PetscErrorCode ierr; 1273*c60c7ad4SBarry Smith 1274*c60c7ad4SBarry Smith PetscFunctionBegin; 1275*c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1276*c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1277*c60c7ad4SBarry Smith PetscFunctionReturn(0); 1278*c60c7ad4SBarry Smith } 1279*c60c7ad4SBarry Smith 1280*c60c7ad4SBarry Smith #undef __FUNCT__ 1281*c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG" 1282*c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1283*c60c7ad4SBarry Smith { 1284*c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1285*c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1286*c60c7ad4SBarry Smith 1287*c60c7ad4SBarry Smith PetscFunctionBegin; 1288*c60c7ad4SBarry Smith *type = pc_gamg->type; 1289*c60c7ad4SBarry Smith PetscFunctionReturn(0); 1290*c60c7ad4SBarry Smith } 1291*c60c7ad4SBarry Smith 1292*c60c7ad4SBarry Smith #undef __FUNCT__ 12939d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 12941e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1295676e1743SMark F. Adams { 12969d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 12971ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 12981ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1299676e1743SMark F. Adams 1300676e1743SMark F. Adams PetscFunctionBegin; 1301*c60c7ad4SBarry Smith pc_gamg->type = type; 13021c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 13039d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 13041ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 13053ae0bb68SMark Adams /* there was something here - kill it */ 13061ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 13071ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 13083ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 13093ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 13103ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 13113ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 13123ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 13133ae0bb68SMark Adams if (pc_gamg->data_sz) { 13143ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 13153ae0bb68SMark Adams pc_gamg->data_sz = 0; 1316*c60c7ad4SBarry Smith } else if (pc_gamg->data) { 13173ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); /* can this happen ? */ 13183ae0bb68SMark Adams } 13191ab5ffc9SJed Brown } 13201ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 13211ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 13229d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1323676e1743SMark F. Adams PetscFunctionReturn(0); 1324676e1743SMark F. Adams } 1325676e1743SMark F. Adams 13265b89ad90SMark F. Adams #undef __FUNCT__ 13275b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 13285b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 13295b89ad90SMark F. Adams { 1330676e1743SMark F. Adams PetscErrorCode ierr; 1331676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1332676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1333676e1743SMark F. Adams PetscBool flag; 13345e7c91beSJed Brown PetscInt two = 2; 13353b4367a7SBarry Smith MPI_Comm comm; 13365b89ad90SMark F. Adams 13375b89ad90SMark F. Adams PetscFunctionBegin; 13383b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1339676e1743SMark F. Adams ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr); 1340676e1743SMark F. Adams { 1341b7cbab4eSMark Adams /* -pc_gamg_type */ 1342b7cbab4eSMark Adams { 1343bd94a7aaSJed Brown char tname[256]; 13441a1c1e04SBarry Smith ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1345bd94a7aaSJed Brown if (flag) { 1346bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 13471ab5ffc9SJed Brown } 1348b7cbab4eSMark Adams } 134975b74bdaSMark F. Adams /* -pc_gamg_verbose */ 135094ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none", pc_gamg->verbose,&pc_gamg->verbose, NULL);CHKERRQ(ierr); 13518263b398SMark F. Adams /* -pc_gamg_repartition */ 135294ae4db5SBarry Smith ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 1353dfd5c07aSMark F. Adams /* -pc_gamg_reuse_interpolation */ 135494ae4db5SBarry Smith ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseProl",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 1355ffc955d6SMark F. Adams /* -pc_gamg_use_agg_gasm */ 135694ae4db5SBarry Smith ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);CHKERRQ(ierr); 1357c20e4228SMark F. Adams /* -pc_gamg_process_eq_limit */ 135894ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL);CHKERRQ(ierr); 1359389730f3SMark F. Adams /* -pc_gamg_coarse_eq_limit */ 136094ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL);CHKERRQ(ierr); 1361c20e4228SMark F. Adams /* -pc_gamg_threshold */ 136294ae4db5SBarry Smith ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr); 1363806fa848SBarry Smith if (flag && pc_gamg->verbose) { 13643b4367a7SBarry Smith ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr); 1365806fa848SBarry Smith } 1366b7cbab4eSMark Adams /* -pc_gamg_eigtarget */ 13670298fd71SBarry Smith ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr); 136894ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 1369b7cbab4eSMark Adams 1370b7cbab4eSMark Adams /* set options for subtype */ 13711ab5ffc9SJed Brown if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(pc);CHKERRQ(ierr);} 1372676e1743SMark F. Adams } 1373676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 13745b89ad90SMark F. Adams PetscFunctionReturn(0); 13755b89ad90SMark F. Adams } 13765b89ad90SMark F. Adams 13775b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 13785b89ad90SMark F. Adams /*MC 1379856830bfSJed Brown PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 13809d5b6da9SMark F. Adams - This is the entry point to GAMG, registered in pcregis.c 13815b89ad90SMark F. Adams 1382280d9858SJed Brown Options Database Keys: 13835b89ad90SMark F. Adams Multigrid options(inherited) 1384280d9858SJed Brown + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1385280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1386280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 13878c1c2452SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 13885b89ad90SMark F. Adams 13895b89ad90SMark F. Adams Level: intermediate 1390280d9858SJed Brown 13915b89ad90SMark F. Adams Concepts: multigrid 13925b89ad90SMark F. Adams 13935b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1394280d9858SJed Brown PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 13955b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 13965b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 13975b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 13985b89ad90SMark F. Adams M*/ 1399b2573a8aSBarry Smith 14005b89ad90SMark F. Adams #undef __FUNCT__ 14015b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 14028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 14035b89ad90SMark F. Adams { 14045b89ad90SMark F. Adams PetscErrorCode ierr; 14055b89ad90SMark F. Adams PC_GAMG *pc_gamg; 14065b89ad90SMark F. Adams PC_MG *mg; 14070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 14089d5b6da9SMark F. Adams static long count = 0; 14095ee9c036SSatish Balay #endif 14105b89ad90SMark F. Adams 14115b89ad90SMark F. Adams PetscFunctionBegin; 14125b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 14135b89ad90SMark F. Adams ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 14145b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 14155b89ad90SMark F. Adams 14165b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 1417b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 14185b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 1419ce4cda84SJed Brown mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 14205b89ad90SMark F. Adams mg->innerctx = pc_gamg; 14215b89ad90SMark F. Adams 1422b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 14231ab5ffc9SJed Brown 14249d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 14259d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 14269d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 14279d5b6da9SMark F. Adams pc_gamg->data = 0; 14285b89ad90SMark F. Adams 14299d5b6da9SMark F. Adams /* register AMG type */ 14303e3471ccSMark Adams ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 14319d5b6da9SMark F. Adams 14329d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 14335b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 14345b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 14355b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 14365b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 14375b89ad90SMark F. Adams 1438bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1439bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1440bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 1441bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseProl_C",PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr); 1442bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1443bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1444bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1445*c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1446bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 14479d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1448d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 1449ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1450038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 1451c8b0795cSMark F. Adams pc_gamg->coarse_eq_limit = 800; 1452d3042614SMark Adams pc_gamg->threshold = 0.; 14539d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 14549d5b6da9SMark F. Adams pc_gamg->verbose = 0; 14559d5b6da9SMark F. Adams pc_gamg->emax_id = -1; 14565e7c91beSJed Brown pc_gamg->eigtarget[0] = 0.05; 14575e7c91beSJed Brown pc_gamg->eigtarget[1] = 1.05; 14589d5b6da9SMark F. Adams 14590cbbd2e1SMark F. Adams /* private events */ 14600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1461785cba28SMark F. Adams if (count++ == 0) { 1462806fa848SBarry Smith ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 1463806fa848SBarry Smith ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 14640cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 14650cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 14660cbbd2e1SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 1467806fa848SBarry Smith ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 1468806fa848SBarry Smith ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 1469806fa848SBarry Smith ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1470806fa848SBarry Smith ierr = PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 1471806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 1472806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 1473806fa848SBarry Smith ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 1474806fa848SBarry Smith ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 1475806fa848SBarry Smith ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 1476806fa848SBarry Smith ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 1477806fa848SBarry Smith ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 1478806fa848SBarry Smith ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 1479f852f58cSMark F. Adams 14800cbbd2e1SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 14810cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 14820cbbd2e1SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1483b4fbaa2aSMark F. Adams /* create timer stages */ 1484b4fbaa2aSMark F. Adams #if defined GAMG_STAGES 1485b4fbaa2aSMark F. Adams { 1486b4fbaa2aSMark F. Adams char str[32]; 1487b4fbaa2aSMark F. Adams PetscInt lidx; 1488806fa848SBarry Smith sprintf(str,"MG Level %d (finest)",0); 1489806fa848SBarry Smith ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 1490b4fbaa2aSMark F. Adams for (lidx=1; lidx<9; lidx++) { 1491b4fbaa2aSMark F. Adams sprintf(str,"MG Level %d",lidx); 1492806fa848SBarry Smith ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 1493b4fbaa2aSMark F. Adams } 1494b4fbaa2aSMark F. Adams } 1495b4fbaa2aSMark F. Adams #endif 1496b4fbaa2aSMark F. Adams } 1497b4fbaa2aSMark F. Adams #endif 1498bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1499bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 15005b89ad90SMark F. Adams PetscFunctionReturn(0); 15015b89ad90SMark F. Adams } 15023e3471ccSMark Adams 15033e3471ccSMark Adams #undef __FUNCT__ 15043e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage" 15053e3471ccSMark Adams /*@C 15063e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 15073e3471ccSMark Adams from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG() 15083e3471ccSMark Adams when using static libraries. 15093e3471ccSMark Adams 15103e3471ccSMark Adams Level: developer 15113e3471ccSMark Adams 15123e3471ccSMark Adams .keywords: PC, PCGAMG, initialize, package 15133e3471ccSMark Adams .seealso: PetscInitialize() 15143e3471ccSMark Adams @*/ 15153e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 15163e3471ccSMark Adams { 15173e3471ccSMark Adams PetscErrorCode ierr; 15183e3471ccSMark Adams 15193e3471ccSMark Adams PetscFunctionBegin; 15203e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 15213e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 15223e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 15233e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 15248e6d0c30SPeter Brune ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 15253e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1526c1c463dbSMark Adams 1527c1c463dbSMark Adams /* general events */ 1528c1c463dbSMark Adams #if defined PETSC_USE_LOG 1529c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr); 1530c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr); 1531c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1532c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1533c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1534c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1535c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr); 1536c1c463dbSMark Adams #endif 1537c1c463dbSMark Adams 15383e3471ccSMark Adams PetscFunctionReturn(0); 15393e3471ccSMark Adams } 15403e3471ccSMark Adams 15413e3471ccSMark Adams #undef __FUNCT__ 15423e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage" 15433e3471ccSMark Adams /*@C 15443e3471ccSMark Adams PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is 15453e3471ccSMark Adams called from PetscFinalize(). 15463e3471ccSMark Adams 15473e3471ccSMark Adams Level: developer 15483e3471ccSMark Adams 15493e3471ccSMark Adams .keywords: Petsc, destroy, package 15503e3471ccSMark Adams .seealso: PetscFinalize() 15513e3471ccSMark Adams @*/ 15523e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 15533e3471ccSMark Adams { 15543e3471ccSMark Adams PetscErrorCode ierr; 15553e3471ccSMark Adams 15563e3471ccSMark Adams PetscFunctionBegin; 15573e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 15583e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 15593e3471ccSMark Adams PetscFunctionReturn(0); 15603e3471ccSMark Adams } 1561