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 Mat mat; 13df826632SBarry Smith VecScatter scatter; 14df826632SBarry Smith } PC_Redistribute; 15df826632SBarry Smith 16df826632SBarry Smith #undef __FUNCT__ 17df826632SBarry Smith #define __FUNCT__ "PCView_Redistribute" 18df826632SBarry Smith static PetscErrorCode PCView_Redistribute(PC pc,PetscViewer viewer) 19df826632SBarry Smith { 20df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 21df826632SBarry Smith PetscErrorCode ierr; 22df826632SBarry Smith PetscTruth iascii,isstring; 23df826632SBarry Smith 24df826632SBarry Smith PetscFunctionBegin; 25df826632SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 26df826632SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 27df826632SBarry Smith if (iascii) { 28df826632SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Redistribute preconditioner: \n");CHKERRQ(ierr); 29df826632SBarry Smith ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr); 30df826632SBarry Smith } else if (isstring) { 31df826632SBarry Smith ierr = PetscViewerStringSPrintf(viewer," Redistribute preconditioner");CHKERRQ(ierr); 32df826632SBarry Smith ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr); 33df826632SBarry Smith } else { 34df826632SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redistribute",((PetscObject)viewer)->type_name); 35df826632SBarry Smith } 36df826632SBarry Smith PetscFunctionReturn(0); 37df826632SBarry Smith } 38df826632SBarry Smith 39df826632SBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 40df826632SBarry Smith #undef __FUNCT__ 41df826632SBarry Smith #define __FUNCT__ "PCSetUp_Redistribute" 42df826632SBarry Smith static PetscErrorCode PCSetUp_Redistribute(PC pc) 43df826632SBarry Smith { 44df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 45df826632SBarry Smith PetscErrorCode ierr; 46df826632SBarry Smith const char *prefix; 47df826632SBarry Smith MPI_Comm comm; 48df826632SBarry Smith PetscInt rstart,rend,i,nz,cnt,*rows,ncnt; 49df826632SBarry Smith PetscMap *map,*nmap; 50df826632SBarry Smith PetscMPIInt size,rank,imdex,tag,n; 51e8dd6687SHong Zhang PetscInt *source = PETSC_NULL; 52df826632SBarry Smith PetscMPIInt *nprocs = PETSC_NULL,nrecvs; 53df826632SBarry Smith PetscInt j,nsends; 54df826632SBarry Smith PetscInt *owner = PETSC_NULL,*starts = PETSC_NULL,count,slen; 55e8dd6687SHong Zhang PetscInt *rvalues,*svalues,recvtotal; 56df826632SBarry Smith PetscMPIInt *onodes1,*olengths1; 57df826632SBarry Smith MPI_Request *send_waits = PETSC_NULL,*recv_waits = PETSC_NULL; 58df826632SBarry Smith MPI_Status recv_status,*send_status; 59df826632SBarry Smith 60df826632SBarry Smith PetscFunctionBegin; 61df826632SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 62df826632SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 63df826632SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 64*361c1e09SMatthew Knepley ierr = PetscObjectGetNewTag((PetscObject)pc,&tag);CHKERRQ(ierr); 65df826632SBarry Smith 66df826632SBarry Smith if (!pc->setupcalled) { 67df826632SBarry Smith if (!red->ksp) { 68df826632SBarry Smith ierr = KSPCreate(comm,&red->ksp);CHKERRQ(ierr); 69df826632SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 70df826632SBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 71df826632SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 72df826632SBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 73df826632SBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redistribute_");CHKERRQ(ierr); 74df826632SBarry Smith } 75df826632SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Cannot yet re-setup a redistribute KSP"); 76df826632SBarry Smith 77df826632SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&rstart,&rend);CHKERRQ(ierr); 78df826632SBarry Smith for (i=rstart; i<rend; i++) { 79df826632SBarry Smith ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 80df826632SBarry Smith if (nz > 1) cnt++; 81df826632SBarry Smith } 82df826632SBarry Smith ierr = PetscMalloc(cnt*sizeof(PetscInt),&rows);CHKERRQ(ierr); 83df826632SBarry Smith cnt = 0; 84df826632SBarry Smith for (i=rstart; i<rend; i++) { 85df826632SBarry Smith ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 86df826632SBarry Smith if (nz > 1) rows[cnt++] = i; 87df826632SBarry Smith } 88df826632SBarry Smith ierr = PetscMalloc(sizeof(PetscMap),&map);CHKERRQ(ierr); 89df826632SBarry Smith ierr = PetscMapInitialize(comm,map);CHKERRQ(ierr); 90df826632SBarry Smith ierr = PetscMapSetLocalSize(map,cnt);CHKERRQ(ierr); 91df826632SBarry Smith ierr = PetscMapSetBlockSize(map,1);CHKERRQ(ierr); 92df826632SBarry Smith ierr = PetscMapSetUp(map);CHKERRQ(ierr); 93df826632SBarry Smith rstart = map->rstart; 94df826632SBarry Smith rend = map->rend; 95df826632SBarry Smith 96df826632SBarry Smith ierr = PetscMalloc(sizeof(PetscMap),&nmap);CHKERRQ(ierr); 97df826632SBarry Smith ierr = PetscMapInitialize(comm,nmap);CHKERRQ(ierr); 98df826632SBarry Smith ierr = MPI_Allreduce(&cnt,&ncnt,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 99df826632SBarry Smith ierr = PetscMapSetSize(nmap,ncnt);CHKERRQ(ierr); 100df826632SBarry Smith ierr = PetscMapSetBlockSize(nmap,1);CHKERRQ(ierr); 101df826632SBarry Smith ierr = PetscMapSetUp(nmap);CHKERRQ(ierr); 102df826632SBarry Smith 103df826632SBarry Smith /* this code is taken from VecScatterCreate_PtoS() */ 104df826632SBarry Smith /* count number of contributors to each processor */ 105df826632SBarry Smith ierr = PetscMalloc2(size,PetscMPIInt,&nprocs,cnt,PetscInt,&owner);CHKERRQ(ierr); 106df826632SBarry Smith ierr = PetscMemzero(nprocs,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 107df826632SBarry Smith j = 0; 108df826632SBarry Smith nsends = 0; 109df826632SBarry Smith for (i=rstart; i<rend; i++) { 110df826632SBarry Smith if (i < nmap->range[j]) j = 0; 111df826632SBarry Smith for (; j<size; j++) { 112df826632SBarry Smith if (i < nmap->range[j+1]) { 113df826632SBarry Smith if (!nprocs[j]++) nsends++; 114df826632SBarry Smith owner[i] = j; 115df826632SBarry Smith break; 116df826632SBarry Smith } 117df826632SBarry Smith } 118df826632SBarry Smith } 119df826632SBarry Smith /* inform other processors of number of messages and max length*/ 120df826632SBarry Smith ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,nprocs,&nrecvs);CHKERRQ(ierr); 121df826632SBarry Smith ierr = PetscGatherMessageLengths(comm,nsends,nrecvs,nprocs,&onodes1,&olengths1);CHKERRQ(ierr); 122df826632SBarry Smith ierr = PetscSortMPIIntWithArray(nrecvs,onodes1,olengths1);CHKERRQ(ierr); 123df826632SBarry Smith recvtotal = 0; for (i=0; i<nrecvs; i++) recvtotal += olengths1[i]; 124df826632SBarry Smith 125df826632SBarry Smith /* post receives: rvalues - rows I will own; count - nu */ 126df826632SBarry Smith ierr = PetscMalloc3(recvtotal,PetscInt,&rvalues,nrecvs,PetscInt,&source,nrecvs,MPI_Request,&recv_waits);CHKERRQ(ierr); 127df826632SBarry Smith count = 0; 128df826632SBarry Smith for (i=0; i<nrecvs; i++) { 129df826632SBarry Smith ierr = MPI_Irecv((rvalues+count),olengths1[i],MPIU_INT,onodes1[i],tag,comm,recv_waits+i);CHKERRQ(ierr); 130df826632SBarry Smith count += olengths1[i]; 131df826632SBarry Smith } 132df826632SBarry Smith 133df826632SBarry Smith /* do sends: 134df826632SBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 135df826632SBarry Smith the ith processor 136df826632SBarry Smith */ 137df826632SBarry Smith ierr = PetscMalloc3(cnt,PetscInt,&svalues,nsends,MPI_Request,&send_waits,size+1,PetscInt,&starts);CHKERRQ(ierr); 138df826632SBarry Smith starts[0] = 0; 139df826632SBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 140df826632SBarry Smith for (i=0; i<cnt; i++) { 141df826632SBarry Smith if (owner[i] != rank) { 142df826632SBarry Smith svalues[starts[owner[i]]++] = rows[i]; 143df826632SBarry Smith } 144df826632SBarry Smith } 145df826632SBarry Smith 146df826632SBarry Smith starts[0] = 0; 147df826632SBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 148df826632SBarry Smith count = 0; 149df826632SBarry Smith for (i=0; i<size; i++) { 150df826632SBarry Smith if (nprocs[i]) { 151df826632SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 152df826632SBarry Smith } 153df826632SBarry Smith } 154df826632SBarry Smith 155df826632SBarry Smith /* wait on receives */ 156df826632SBarry Smith count = nrecvs; 157df826632SBarry Smith slen = 0; 158df826632SBarry Smith while (count) { 159df826632SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 160df826632SBarry Smith /* unpack receives into our local space */ 161df826632SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 162df826632SBarry Smith slen += n; 163df826632SBarry Smith count--; 164df826632SBarry Smith } 165df826632SBarry Smith if (slen != recvtotal) SETERRQ2(PETSC_ERR_PLIB,"Total message lengths %D not expected %D",slen,recvtotal); 166df826632SBarry Smith 167df826632SBarry Smith ierr = PetscFree(olengths1);CHKERRQ(ierr); 168df826632SBarry Smith ierr = PetscFree(onodes1);CHKERRQ(ierr); 169df826632SBarry Smith ierr = PetscFree3(rvalues,source,recv_waits);CHKERRQ(ierr); 170df826632SBarry Smith 171df826632SBarry Smith /* wait on sends */ 172df826632SBarry Smith if (nsends) { 173df826632SBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 174df826632SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 175df826632SBarry Smith ierr = PetscFree(send_status);CHKERRQ(ierr); 176df826632SBarry Smith } 177df826632SBarry Smith ierr = PetscFree3(svalues,send_waits,starts);CHKERRQ(ierr); 178df826632SBarry Smith 179df826632SBarry Smith 180df826632SBarry Smith ierr = PetscFree(map);CHKERRQ(ierr); 181df826632SBarry Smith ierr = PetscFree(nmap);CHKERRQ(ierr); 182df826632SBarry Smith 183df826632SBarry Smith ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 184df826632SBarry Smith PetscFunctionReturn(0); 185df826632SBarry Smith } 186df826632SBarry Smith 187df826632SBarry Smith #undef __FUNCT__ 188df826632SBarry Smith #define __FUNCT__ "PCApply_Redistribute" 189df826632SBarry Smith static PetscErrorCode PCApply_Redistribute(PC pc,Vec b,Vec x) 190df826632SBarry Smith { 191df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 192df826632SBarry Smith PetscErrorCode ierr; 193df826632SBarry Smith 194df826632SBarry Smith PetscFunctionBegin; 195df826632SBarry Smith ierr = VecScatterBegin(red->scatter,b,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 196df826632SBarry Smith ierr = VecScatterEnd(red->scatter,b,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 197df826632SBarry Smith ierr = KSPSolve(red->ksp,red->b,red->x);CHKERRQ(ierr); 198df826632SBarry Smith ierr = VecScatterBegin(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 199df826632SBarry Smith ierr = VecScatterEnd(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 200df826632SBarry Smith PetscFunctionReturn(0); 201df826632SBarry Smith } 202df826632SBarry Smith 203df826632SBarry Smith #undef __FUNCT__ 204df826632SBarry Smith #define __FUNCT__ "PCDestroy_Redistribute" 205df826632SBarry Smith static PetscErrorCode PCDestroy_Redistribute(PC pc) 206df826632SBarry Smith { 207df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 208df826632SBarry Smith PetscErrorCode ierr; 209df826632SBarry Smith 210df826632SBarry Smith PetscFunctionBegin; 211df826632SBarry Smith if (red->scatter) {ierr = VecScatterDestroy(red->scatter);CHKERRQ(ierr);} 212df826632SBarry Smith if (red->b) {ierr = VecDestroy(red->b);CHKERRQ(ierr);} 213df826632SBarry Smith if (red->x) {ierr = VecDestroy(red->x);CHKERRQ(ierr);} 214df826632SBarry Smith if (red->mat) {ierr = MatDestroy(red->mat);CHKERRQ(ierr);} 215df826632SBarry Smith if (red->ksp) {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);} 216df826632SBarry Smith ierr = PetscFree(red);CHKERRQ(ierr); 217df826632SBarry Smith PetscFunctionReturn(0); 218df826632SBarry Smith } 219df826632SBarry Smith 220df826632SBarry Smith #undef __FUNCT__ 221df826632SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redistribute" 222df826632SBarry Smith static PetscErrorCode PCSetFromOptions_Redistribute(PC pc) 223df826632SBarry Smith { 224df826632SBarry Smith PetscErrorCode ierr; 225df826632SBarry Smith PC_Redistribute *red = (PC_Redistribute*)pc->data; 226df826632SBarry Smith 227df826632SBarry Smith PetscFunctionBegin; 228df826632SBarry Smith ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 229df826632SBarry Smith PetscFunctionReturn(0); 230df826632SBarry Smith } 231df826632SBarry Smith 232df826632SBarry Smith /* -------------------------------------------------------------------------------------*/ 233df826632SBarry Smith /*MC 234df826632SBarry Smith PCREDISTRIBUTE - Redistributes a matrix for load balancing and then applys a KSP to that new matrix 235df826632SBarry Smith 236df826632SBarry Smith Options for the redistribute preconditioners can be set with -redistribute_ksp_xxx <values> and -redistribute_pc_xxx <values> 237df826632SBarry Smith 238df826632SBarry Smith Level: intermediate 239df826632SBarry Smith 240df826632SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types) 241df826632SBarry Smith M*/ 242df826632SBarry Smith 243df826632SBarry Smith EXTERN_C_BEGIN 244df826632SBarry Smith #undef __FUNCT__ 245df826632SBarry Smith #define __FUNCT__ "PCCreate_Redistribute" 246df826632SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redistribute(PC pc) 247df826632SBarry Smith { 248df826632SBarry Smith PetscErrorCode ierr; 249df826632SBarry Smith PC_Redistribute *red; 250df826632SBarry Smith 251df826632SBarry Smith PetscFunctionBegin; 252df826632SBarry Smith ierr = PetscNewLog(pc,PC_Redistribute,&red);CHKERRQ(ierr); 253df826632SBarry Smith pc->data = (void*)red; 254df826632SBarry Smith 255df826632SBarry Smith pc->ops->apply = PCApply_Redistribute; 256df826632SBarry Smith pc->ops->applytranspose = 0; 257df826632SBarry Smith pc->ops->setup = PCSetUp_Redistribute; 258df826632SBarry Smith pc->ops->destroy = PCDestroy_Redistribute; 259df826632SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redistribute; 260df826632SBarry Smith pc->ops->view = PCView_Redistribute; 261df826632SBarry Smith PetscFunctionReturn(0); 262df826632SBarry Smith } 263df826632SBarry Smith EXTERN_C_END 264