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