1df826632SBarry Smith #define PETSCKSP_DLL 2df826632SBarry Smith 3df826632SBarry Smith /* 4df826632SBarry Smith This file defines a "solve the problem redistributely on each subgroup of processor" preconditioner. 5df826632SBarry Smith */ 6df826632SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 7df826632SBarry Smith #include "petscksp.h" 8df826632SBarry Smith 9df826632SBarry Smith typedef struct { 10df826632SBarry Smith KSP ksp; 11df826632SBarry Smith Vec x,b; 12df826632SBarry Smith VecScatter scatter; 13911f9fe8SBarry Smith IS is; 14181dd334SBarry Smith PetscInt dcnt,*drows; /* these are the local rows that have only diagonal entry */ 15181dd334SBarry Smith PetscScalar *diag; 16181dd334SBarry Smith Vec work; 17df826632SBarry Smith } PC_Redistribute; 18df826632SBarry Smith 19df826632SBarry Smith #undef __FUNCT__ 20df826632SBarry Smith #define __FUNCT__ "PCView_Redistribute" 21df826632SBarry Smith static PetscErrorCode PCView_Redistribute(PC pc,PetscViewer viewer) 22df826632SBarry Smith { 23df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 24df826632SBarry Smith PetscErrorCode ierr; 25ace3abfcSBarry Smith PetscBool iascii,isstring; 26181dd334SBarry Smith PetscInt ncnt,N; 27df826632SBarry Smith 28df826632SBarry Smith PetscFunctionBegin; 292692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 302692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 31df826632SBarry Smith if (iascii) { 32181dd334SBarry Smith ierr = MPI_Allreduce(&red->dcnt,&ncnt,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 33181dd334SBarry Smith ierr = MatGetSize(pc->pmat,&N,PETSC_NULL);CHKERRQ(ierr); 34572f72c7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Number rows eliminated %D Percentage rows eliminated %g\n",ncnt,100.0*((PetscReal)ncnt)/((PetscReal)N));CHKERRQ(ierr); 35df826632SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Redistribute preconditioner: \n");CHKERRQ(ierr); 36df826632SBarry Smith ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr); 37df826632SBarry Smith } else if (isstring) { 38df826632SBarry Smith ierr = PetscViewerStringSPrintf(viewer," Redistribute preconditioner");CHKERRQ(ierr); 39df826632SBarry Smith ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr); 40df826632SBarry Smith } else { 4165e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PC redistribute",((PetscObject)viewer)->type_name); 42df826632SBarry Smith } 43df826632SBarry Smith PetscFunctionReturn(0); 44df826632SBarry Smith } 45df826632SBarry Smith 46df826632SBarry Smith #undef __FUNCT__ 47df826632SBarry Smith #define __FUNCT__ "PCSetUp_Redistribute" 48df826632SBarry Smith static PetscErrorCode PCSetUp_Redistribute(PC pc) 49df826632SBarry Smith { 50df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 51df826632SBarry Smith PetscErrorCode ierr; 52df826632SBarry Smith MPI_Comm comm; 53181dd334SBarry Smith PetscInt rstart,rend,i,nz,cnt,*rows,ncnt,dcnt,*drows; 5426283091SBarry Smith PetscLayout map,nmap; 555a586d82SBarry Smith PetscMPIInt size,imdex,tag,n; 56e8dd6687SHong Zhang PetscInt *source = PETSC_NULL; 57df826632SBarry Smith PetscMPIInt *nprocs = PETSC_NULL,nrecvs; 58df826632SBarry Smith PetscInt j,nsends; 59df826632SBarry Smith PetscInt *owner = PETSC_NULL,*starts = PETSC_NULL,count,slen; 60e8dd6687SHong Zhang PetscInt *rvalues,*svalues,recvtotal; 61df826632SBarry Smith PetscMPIInt *onodes1,*olengths1; 62df826632SBarry Smith MPI_Request *send_waits = PETSC_NULL,*recv_waits = PETSC_NULL; 63df826632SBarry Smith MPI_Status recv_status,*send_status; 64181dd334SBarry Smith Vec tvec,diag; 65ca320bd4SBarry Smith Mat tmat; 66181dd334SBarry Smith const PetscScalar *d; 67df826632SBarry Smith 68df826632SBarry Smith PetscFunctionBegin; 69dc9360f3SBarry Smith if (pc->setupcalled) { 70dc9360f3SBarry Smith ierr = KSPGetOperators(red->ksp,PETSC_NULL,&tmat,PETSC_NULL);CHKERRQ(ierr); 71dc9360f3SBarry Smith ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_REUSE_MATRIX,&tmat);CHKERRQ(ierr); 721d805cfdSBarry Smith ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 73dc9360f3SBarry Smith } else { 74b862ddfaSBarry Smith PetscInt NN; 75b862ddfaSBarry Smith 76df826632SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 77df826632SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 78361c1e09SMatthew Knepley ierr = PetscObjectGetNewTag((PetscObject)pc,&tag);CHKERRQ(ierr); 79df826632SBarry Smith 80ca320bd4SBarry Smith /* count non-diagonal rows on process */ 81df826632SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&rstart,&rend);CHKERRQ(ierr); 82ca320bd4SBarry Smith cnt = 0; 83df826632SBarry Smith for (i=rstart; i<rend; i++) { 84df826632SBarry Smith ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 85df826632SBarry Smith if (nz > 1) cnt++; 86911f9fe8SBarry Smith ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 87df826632SBarry Smith } 88df826632SBarry Smith ierr = PetscMalloc(cnt*sizeof(PetscInt),&rows);CHKERRQ(ierr); 89181dd334SBarry Smith ierr = PetscMalloc((rend - rstart - cnt)*sizeof(PetscInt),&drows);CHKERRQ(ierr); 90ca320bd4SBarry Smith 91ca320bd4SBarry Smith /* list non-diagonal rows on process */ 92181dd334SBarry Smith cnt = 0; dcnt = 0; 93df826632SBarry Smith for (i=rstart; i<rend; i++) { 94df826632SBarry Smith ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 95df826632SBarry Smith if (nz > 1) rows[cnt++] = i; 96181dd334SBarry Smith else drows[dcnt++] = i - rstart; 97911f9fe8SBarry Smith ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 98df826632SBarry Smith } 99ca320bd4SBarry Smith 10026283091SBarry Smith /* create PetscLayout for non-diagonal rows on each process */ 10126283091SBarry Smith ierr = PetscLayoutCreate(comm,&map);CHKERRQ(ierr); 10226283091SBarry Smith ierr = PetscLayoutSetLocalSize(map,cnt);CHKERRQ(ierr); 10326283091SBarry Smith ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 10426283091SBarry Smith ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 105df826632SBarry Smith rstart = map->rstart; 106df826632SBarry Smith rend = map->rend; 107df826632SBarry Smith 10826283091SBarry Smith /* create PetscLayout for load-balanced non-diagonal rows on each process */ 10926283091SBarry Smith ierr = PetscLayoutCreate(comm,&nmap);CHKERRQ(ierr); 110df826632SBarry Smith ierr = MPI_Allreduce(&cnt,&ncnt,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 11126283091SBarry Smith ierr = PetscLayoutSetSize(nmap,ncnt);CHKERRQ(ierr); 11226283091SBarry Smith ierr = PetscLayoutSetBlockSize(nmap,1);CHKERRQ(ierr); 11326283091SBarry Smith ierr = PetscLayoutSetUp(nmap);CHKERRQ(ierr); 114df826632SBarry Smith 115b862ddfaSBarry Smith ierr = MatGetSize(pc->pmat,&NN,PETSC_NULL);CHKERRQ(ierr); 116b862ddfaSBarry Smith ierr = PetscInfo2(pc,"Number of diagonal rows eliminated %d, percentage eliminated %g\n",NN-ncnt,((PetscReal)(NN-ncnt))/((PetscReal)(NN))); 117ca320bd4SBarry Smith /* 118ca320bd4SBarry Smith this code is taken from VecScatterCreate_PtoS() 119ca320bd4SBarry Smith Determines what rows need to be moved where to 120ca320bd4SBarry Smith load balance the non-diagonal rows 121ca320bd4SBarry Smith */ 122df826632SBarry Smith /* count number of contributors to each processor */ 123df826632SBarry Smith ierr = PetscMalloc2(size,PetscMPIInt,&nprocs,cnt,PetscInt,&owner);CHKERRQ(ierr); 124df826632SBarry Smith ierr = PetscMemzero(nprocs,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 125df826632SBarry Smith j = 0; 126df826632SBarry Smith nsends = 0; 127df826632SBarry Smith for (i=rstart; i<rend; i++) { 128df826632SBarry Smith if (i < nmap->range[j]) j = 0; 129df826632SBarry Smith for (; j<size; j++) { 130df826632SBarry Smith if (i < nmap->range[j+1]) { 131df826632SBarry Smith if (!nprocs[j]++) nsends++; 132ca320bd4SBarry Smith owner[i-rstart] = j; 133df826632SBarry Smith break; 134df826632SBarry Smith } 135df826632SBarry Smith } 136df826632SBarry Smith } 137df826632SBarry Smith /* inform other processors of number of messages and max length*/ 138df826632SBarry Smith ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,nprocs,&nrecvs);CHKERRQ(ierr); 139df826632SBarry Smith ierr = PetscGatherMessageLengths(comm,nsends,nrecvs,nprocs,&onodes1,&olengths1);CHKERRQ(ierr); 140df826632SBarry Smith ierr = PetscSortMPIIntWithArray(nrecvs,onodes1,olengths1);CHKERRQ(ierr); 141df826632SBarry Smith recvtotal = 0; for (i=0; i<nrecvs; i++) recvtotal += olengths1[i]; 142df826632SBarry Smith 143df826632SBarry Smith /* post receives: rvalues - rows I will own; count - nu */ 144df826632SBarry Smith ierr = PetscMalloc3(recvtotal,PetscInt,&rvalues,nrecvs,PetscInt,&source,nrecvs,MPI_Request,&recv_waits);CHKERRQ(ierr); 145df826632SBarry Smith count = 0; 146df826632SBarry Smith for (i=0; i<nrecvs; i++) { 147df826632SBarry Smith ierr = MPI_Irecv((rvalues+count),olengths1[i],MPIU_INT,onodes1[i],tag,comm,recv_waits+i);CHKERRQ(ierr); 148df826632SBarry Smith count += olengths1[i]; 149df826632SBarry Smith } 150df826632SBarry Smith 151df826632SBarry Smith /* do sends: 152df826632SBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 153df826632SBarry Smith the ith processor 154df826632SBarry Smith */ 155911f9fe8SBarry Smith ierr = PetscMalloc3(cnt,PetscInt,&svalues,nsends,MPI_Request,&send_waits,size,PetscInt,&starts);CHKERRQ(ierr); 156df826632SBarry Smith starts[0] = 0; 157df826632SBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 158df826632SBarry Smith for (i=0; i<cnt; i++) { 159df826632SBarry Smith svalues[starts[owner[i]]++] = rows[i]; 160df826632SBarry Smith } 161181dd334SBarry Smith for (i=0; i<cnt; i++) rows[i] = rows[i] - rstart; 162181dd334SBarry Smith red->drows = drows; 163181dd334SBarry Smith red->dcnt = dcnt; 164ca320bd4SBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 165181dd334SBarry Smith 166df826632SBarry Smith starts[0] = 0; 167911f9fe8SBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 168df826632SBarry Smith count = 0; 169df826632SBarry Smith for (i=0; i<size; i++) { 170df826632SBarry Smith if (nprocs[i]) { 171df826632SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 172df826632SBarry Smith } 173df826632SBarry Smith } 174df826632SBarry Smith 175df826632SBarry Smith /* wait on receives */ 176df826632SBarry Smith count = nrecvs; 177df826632SBarry Smith slen = 0; 178df826632SBarry Smith while (count) { 179df826632SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 180df826632SBarry Smith /* unpack receives into our local space */ 181df826632SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 182df826632SBarry Smith slen += n; 183df826632SBarry Smith count--; 184df826632SBarry Smith } 185e32f2f54SBarry Smith if (slen != recvtotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Total message lengths %D not expected %D",slen,recvtotal); 186df826632SBarry Smith 18770b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,slen,rvalues,PETSC_COPY_VALUES,&red->is);CHKERRQ(ierr); 188911f9fe8SBarry Smith 189ca320bd4SBarry Smith /* free up all work space */ 190df826632SBarry Smith ierr = PetscFree(olengths1);CHKERRQ(ierr); 191df826632SBarry Smith ierr = PetscFree(onodes1);CHKERRQ(ierr); 192df826632SBarry Smith ierr = PetscFree3(rvalues,source,recv_waits);CHKERRQ(ierr); 193ca320bd4SBarry Smith ierr = PetscFree2(nprocs,owner);CHKERRQ(ierr); 194ca320bd4SBarry Smith if (nsends) { /* wait on sends */ 195df826632SBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 196df826632SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 197df826632SBarry Smith ierr = PetscFree(send_status);CHKERRQ(ierr); 198df826632SBarry Smith } 199df826632SBarry Smith ierr = PetscFree3(svalues,send_waits,starts);CHKERRQ(ierr); 20026283091SBarry Smith ierr = PetscLayoutDestroy(map);CHKERRQ(ierr); 20126283091SBarry Smith ierr = PetscLayoutDestroy(nmap);CHKERRQ(ierr); 202df826632SBarry Smith 203ca320bd4SBarry Smith ierr = VecCreateMPI(comm,slen,PETSC_DETERMINE,&red->b);CHKERRQ(ierr); 204ca320bd4SBarry Smith ierr = VecDuplicate(red->b,&red->x);CHKERRQ(ierr); 205ca320bd4SBarry Smith ierr = MatGetVecs(pc->pmat,&tvec,PETSC_NULL);CHKERRQ(ierr); 206ca320bd4SBarry Smith ierr = VecScatterCreate(tvec,red->is,red->b,PETSC_NULL,&red->scatter);CHKERRQ(ierr); 207ca320bd4SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 208ca320bd4SBarry Smith ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_INITIAL_MATRIX,&tmat);CHKERRQ(ierr); 209ca320bd4SBarry Smith ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 210ca320bd4SBarry Smith ierr = MatDestroy(tmat);CHKERRQ(ierr); 2111d805cfdSBarry Smith } 212181dd334SBarry Smith 213181dd334SBarry Smith /* get diagonal portion of matrix */ 214181dd334SBarry Smith ierr = PetscMalloc(red->dcnt*sizeof(PetscScalar),&red->diag);CHKERRQ(ierr); 215181dd334SBarry Smith ierr = MatGetVecs(pc->pmat,&diag,PETSC_NULL);CHKERRQ(ierr); 216181dd334SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 2173649974fSBarry Smith ierr = VecGetArrayRead(diag,&d);CHKERRQ(ierr); 218181dd334SBarry Smith for (i=0; i<red->dcnt; i++) { 219181dd334SBarry Smith red->diag[i] = 1.0/d[red->drows[i]]; 220181dd334SBarry Smith } 2213649974fSBarry Smith ierr = VecRestoreArrayRead(diag,&d);CHKERRQ(ierr); 222181dd334SBarry Smith ierr = VecDestroy(diag);CHKERRQ(ierr); 223df826632SBarry Smith ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 224df826632SBarry Smith PetscFunctionReturn(0); 225df826632SBarry Smith } 226df826632SBarry Smith 227df826632SBarry Smith #undef __FUNCT__ 228df826632SBarry Smith #define __FUNCT__ "PCApply_Redistribute" 229df826632SBarry Smith static PetscErrorCode PCApply_Redistribute(PC pc,Vec b,Vec x) 230df826632SBarry Smith { 231df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 232df826632SBarry Smith PetscErrorCode ierr; 233181dd334SBarry Smith PetscInt dcnt = red->dcnt,i; 234181dd334SBarry Smith const PetscInt *drows = red->drows; 235181dd334SBarry Smith PetscScalar *xwork; 236181dd334SBarry Smith const PetscScalar *bwork,*diag = red->diag; 237df826632SBarry Smith 238df826632SBarry Smith PetscFunctionBegin; 239181dd334SBarry Smith if (!red->work) { 240181dd334SBarry Smith ierr = VecDuplicate(b,&red->work);CHKERRQ(ierr); 241181dd334SBarry Smith } 242181dd334SBarry Smith /* compute the rows of solution that have diagonal entries only */ 243181dd334SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); /* x = diag(A)^{-1} b */ 244181dd334SBarry Smith ierr = VecGetArray(x,&xwork);CHKERRQ(ierr); 2453649974fSBarry Smith ierr = VecGetArrayRead(b,&bwork);CHKERRQ(ierr); 246181dd334SBarry Smith for (i=0; i<dcnt; i++) { 247181dd334SBarry Smith xwork[drows[i]] = diag[i]*bwork[drows[i]]; 248181dd334SBarry Smith } 249181dd334SBarry Smith ierr = PetscLogFlops(dcnt);CHKERRQ(ierr); 250181dd334SBarry Smith ierr = VecRestoreArray(red->work,&xwork);CHKERRQ(ierr); 2513649974fSBarry Smith ierr = VecRestoreArrayRead(b,&bwork);CHKERRQ(ierr); 252181dd334SBarry Smith /* update the right hand side for the reduced system with diagonal rows (and corresponding columns) removed */ 253181dd334SBarry Smith ierr = MatMult(pc->pmat,x,red->work);CHKERRQ(ierr); 254181dd334SBarry Smith ierr = VecAYPX(red->work,-1.0,b);CHKERRQ(ierr); /* red->work = b - A x */ 255181dd334SBarry Smith 256181dd334SBarry Smith ierr = VecScatterBegin(red->scatter,red->work,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 257181dd334SBarry Smith ierr = VecScatterEnd(red->scatter,red->work,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 258df826632SBarry Smith ierr = KSPSolve(red->ksp,red->b,red->x);CHKERRQ(ierr); 259df826632SBarry Smith ierr = VecScatterBegin(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 260df826632SBarry Smith ierr = VecScatterEnd(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 261df826632SBarry Smith PetscFunctionReturn(0); 262df826632SBarry Smith } 263df826632SBarry Smith 264df826632SBarry Smith #undef __FUNCT__ 265df826632SBarry Smith #define __FUNCT__ "PCDestroy_Redistribute" 266df826632SBarry Smith static PetscErrorCode PCDestroy_Redistribute(PC pc) 267df826632SBarry Smith { 268df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 269df826632SBarry Smith PetscErrorCode ierr; 270df826632SBarry Smith 271df826632SBarry Smith PetscFunctionBegin; 272df826632SBarry Smith if (red->scatter) {ierr = VecScatterDestroy(red->scatter);CHKERRQ(ierr);} 273ca320bd4SBarry Smith if (red->is) {ierr = ISDestroy(red->is);CHKERRQ(ierr);} 274df826632SBarry Smith if (red->b) {ierr = VecDestroy(red->b);CHKERRQ(ierr);} 275df826632SBarry Smith if (red->x) {ierr = VecDestroy(red->x);CHKERRQ(ierr);} 276df826632SBarry Smith if (red->ksp) {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);} 277181dd334SBarry Smith if (red->work) {ierr = VecDestroy(red->work);CHKERRQ(ierr);} 278181dd334SBarry Smith ierr = PetscFree(red->drows); 279181dd334SBarry Smith ierr = PetscFree(red->diag); 280df826632SBarry Smith ierr = PetscFree(red);CHKERRQ(ierr); 281df826632SBarry Smith PetscFunctionReturn(0); 282df826632SBarry Smith } 283df826632SBarry Smith 284df826632SBarry Smith #undef __FUNCT__ 285df826632SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redistribute" 286df826632SBarry Smith static PetscErrorCode PCSetFromOptions_Redistribute(PC pc) 287df826632SBarry Smith { 288df826632SBarry Smith PetscErrorCode ierr; 289df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 290df826632SBarry Smith 291df826632SBarry Smith PetscFunctionBegin; 292df826632SBarry Smith ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 293df826632SBarry Smith PetscFunctionReturn(0); 294df826632SBarry Smith } 295df826632SBarry Smith 296*5e7ef714SBarry Smith #undef __FUNCT__ 297*5e7ef714SBarry Smith #define __FUNCT__ "PCRedistributeGetKSP" 298*5e7ef714SBarry Smith /*@ 299*5e7ef714SBarry Smith PCRedistributeGetKSP - Gets the KSP created by the PCREDISTRIBUTE 300*5e7ef714SBarry Smith 301*5e7ef714SBarry Smith Not Collective 302*5e7ef714SBarry Smith 303*5e7ef714SBarry Smith Input Parameter: 304*5e7ef714SBarry Smith . pc - the preconditioner context 305*5e7ef714SBarry Smith 306*5e7ef714SBarry Smith Output Parameter: 307*5e7ef714SBarry Smith . innerksp - the inner KSP 308*5e7ef714SBarry Smith 309*5e7ef714SBarry Smith Level: advanced 310*5e7ef714SBarry Smith 311*5e7ef714SBarry Smith .keywords: PC, redistribute solve 312*5e7ef714SBarry Smith @*/ 313*5e7ef714SBarry Smith PetscErrorCode PCRedistributeGetKSP(PC pc,KSP *innerksp) 314*5e7ef714SBarry Smith { 315*5e7ef714SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 316*5e7ef714SBarry Smith 317*5e7ef714SBarry Smith PetscFunctionBegin; 318*5e7ef714SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 319*5e7ef714SBarry Smith PetscValidPointer(innerksp,2); 320*5e7ef714SBarry Smith *innerksp = red->ksp; 321*5e7ef714SBarry Smith PetscFunctionReturn(0); 322*5e7ef714SBarry Smith } 323*5e7ef714SBarry Smith 324df826632SBarry Smith /* -------------------------------------------------------------------------------------*/ 325df826632SBarry Smith /*MC 326181dd334SBarry Smith PCREDISTRIBUTE - Redistributes a matrix for load balancing, removing the rows that only have a diagonal entry and then applys a KSP to that new matrix 327df826632SBarry Smith 328df826632SBarry Smith Options for the redistribute preconditioners can be set with -redistribute_ksp_xxx <values> and -redistribute_pc_xxx <values> 329df826632SBarry Smith 330181dd334SBarry Smith Notes: Usually run this with -ksp_type preonly 331181dd334SBarry Smith 332181dd334SBarry Smith If you have used MatZeroRows() to eliminate (for example, Dirichlet) boundary conditions for a symmetric problem then you can use, for example, -ksp_type preonly 333181dd334SBarry Smith -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc to take advantage of the symmetry. 334181dd334SBarry Smith 335df826632SBarry Smith Level: intermediate 336df826632SBarry Smith 337*5e7ef714SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedistributeGetKSP() 338df826632SBarry Smith M*/ 339df826632SBarry Smith 340df826632SBarry Smith EXTERN_C_BEGIN 341df826632SBarry Smith #undef __FUNCT__ 342df826632SBarry Smith #define __FUNCT__ "PCCreate_Redistribute" 3437087cfbeSBarry Smith PetscErrorCode PCCreate_Redistribute(PC pc) 344df826632SBarry Smith { 345df826632SBarry Smith PetscErrorCode ierr; 346df826632SBarry Smith PC_Redistribute *red; 347911f9fe8SBarry Smith const char *prefix; 348df826632SBarry Smith 349df826632SBarry Smith PetscFunctionBegin; 350df826632SBarry Smith ierr = PetscNewLog(pc,PC_Redistribute,&red);CHKERRQ(ierr); 351df826632SBarry Smith pc->data = (void*)red; 352df826632SBarry Smith 353df826632SBarry Smith pc->ops->apply = PCApply_Redistribute; 354df826632SBarry Smith pc->ops->applytranspose = 0; 355df826632SBarry Smith pc->ops->setup = PCSetUp_Redistribute; 356df826632SBarry Smith pc->ops->destroy = PCDestroy_Redistribute; 357df826632SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redistribute; 358df826632SBarry Smith pc->ops->view = PCView_Redistribute; 359911f9fe8SBarry Smith 360911f9fe8SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&red->ksp);CHKERRQ(ierr); 361911f9fe8SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 362911f9fe8SBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 363911f9fe8SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 364911f9fe8SBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 365911f9fe8SBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redistribute_");CHKERRQ(ierr); 366df826632SBarry Smith PetscFunctionReturn(0); 367df826632SBarry Smith } 368df826632SBarry Smith EXTERN_C_END 369