1 2 /* 3 Defines matrix-matrix product routines for pairs of MPIAIJ matrices 4 C = A * B 5 */ 6 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 7 #include <../src/mat/utils/freespace.h> 8 #include <../src/mat/impls/aij/mpi/mpiaij.h> 9 #include <petscbt.h> 10 #include <../src/mat/impls/dense/mpi/mpidense.h> 11 #include <petsc/private/vecimpl.h> 12 #include <petsc/private/vecscatterimpl.h> 13 14 #if defined(PETSC_HAVE_HYPRE) 15 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 16 #endif 17 18 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 19 { 20 PetscErrorCode ierr; 21 #if defined(PETSC_HAVE_HYPRE) 22 const char *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"}; 23 PetscInt nalg = 4; 24 #else 25 const char *algTypes[3] = {"scalable","nonscalable","seqmpi"}; 26 PetscInt nalg = 3; 27 #endif 28 PetscInt alg = 1; /* set nonscalable algorithm as default */ 29 MPI_Comm comm; 30 PetscBool flg; 31 32 PetscFunctionBegin; 33 if (scall == MAT_INITIAL_MATRIX) { 34 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 35 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend); 36 37 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatMatMult","Mat");CHKERRQ(ierr); 38 ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr); 39 ierr = PetscOptionsEnd();CHKERRQ(ierr); 40 41 if (!flg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */ 42 MatInfo Ainfo,Binfo; 43 PetscInt nz_local; 44 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 45 46 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 47 ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr); 48 nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated); 49 50 if (B->cmap->N > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 51 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 52 53 if (alg_scalable) { 54 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 55 ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,fill*nz_local);CHKERRQ(ierr); 56 } 57 } 58 59 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 60 switch (alg) { 61 case 1: 62 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr); 63 break; 64 case 2: 65 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);CHKERRQ(ierr); 66 break; 67 #if defined(PETSC_HAVE_HYPRE) 68 case 3: 69 ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr); 70 break; 71 #endif 72 default: 73 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 74 break; 75 } 76 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 77 78 if (alg == 0 || alg == 1) { 79 Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 80 Mat_APMPI *ap = c->ap; 81 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 82 ap->freestruct = PETSC_FALSE; 83 ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 84 ierr = PetscOptionsEnd();CHKERRQ(ierr); 85 } 86 } 87 88 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 89 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 90 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 95 { 96 PetscErrorCode ierr; 97 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 98 Mat_APMPI *ptap = a->ap; 99 100 PetscFunctionBegin; 101 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 102 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 103 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 104 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 105 ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 106 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 107 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 108 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 109 ierr = ptap->destroy(A);CHKERRQ(ierr); 110 ierr = PetscFree(ptap);CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C) 115 { 116 PetscErrorCode ierr; 117 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 118 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 119 Mat_SeqAIJ *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 120 PetscScalar *cda=cd->a,*coa=co->a; 121 Mat_SeqAIJ *p_loc,*p_oth; 122 PetscScalar *apa,*ca; 123 PetscInt cm =C->rmap->n; 124 Mat_APMPI *ptap=c->ap; 125 PetscInt *api,*apj,*apJ,i,k; 126 PetscInt cstart=C->cmap->rstart; 127 PetscInt cdnz,conz,k0,k1; 128 MPI_Comm comm; 129 PetscMPIInt size; 130 131 PetscFunctionBegin; 132 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 133 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 134 135 if (!ptap->P_oth && size>1) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 136 137 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 138 /*-----------------------------------------------------*/ 139 /* update numerical values of P_oth and P_loc */ 140 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 141 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 142 143 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 144 /*----------------------------------------------------------*/ 145 /* get data from symbolic products */ 146 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 147 p_oth = NULL; 148 if (size >1) { 149 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 150 } 151 152 /* get apa for storing dense row A[i,:]*P */ 153 apa = ptap->apa; 154 155 api = ptap->api; 156 apj = ptap->apj; 157 for (i=0; i<cm; i++) { 158 /* compute apa = A[i,:]*P */ 159 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 160 161 /* set values in C */ 162 apJ = apj + api[i]; 163 cdnz = cd->i[i+1] - cd->i[i]; 164 conz = co->i[i+1] - co->i[i]; 165 166 /* 1st off-diagoanl part of C */ 167 ca = coa + co->i[i]; 168 k = 0; 169 for (k0=0; k0<conz; k0++) { 170 if (apJ[k] >= cstart) break; 171 ca[k0] = apa[apJ[k]]; 172 apa[apJ[k++]] = 0.0; 173 } 174 175 /* diagonal part of C */ 176 ca = cda + cd->i[i]; 177 for (k1=0; k1<cdnz; k1++) { 178 ca[k1] = apa[apJ[k]]; 179 apa[apJ[k++]] = 0.0; 180 } 181 182 /* 2nd off-diagoanl part of C */ 183 ca = coa + co->i[i]; 184 for (; k0<conz; k0++) { 185 ca[k0] = apa[apJ[k]]; 186 apa[apJ[k++]] = 0.0; 187 } 188 } 189 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 190 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 191 192 if (ptap->freestruct) { 193 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 194 } 195 PetscFunctionReturn(0); 196 } 197 198 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C) 199 { 200 PetscErrorCode ierr; 201 MPI_Comm comm; 202 PetscMPIInt size; 203 Mat Cmpi; 204 Mat_APMPI *ptap; 205 PetscFreeSpaceList free_space=NULL,current_space=NULL; 206 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c; 207 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 208 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 209 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 210 PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 211 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 212 PetscBT lnkbt; 213 PetscScalar *apa; 214 PetscReal afill; 215 MatType mtype; 216 217 PetscFunctionBegin; 218 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 219 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 220 221 /* create struct Mat_APMPI and attached it to C later */ 222 ierr = PetscNew(&ptap);CHKERRQ(ierr); 223 224 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 225 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 226 227 /* get P_loc by taking all local rows of P */ 228 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 229 230 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 231 pi_loc = p_loc->i; pj_loc = p_loc->j; 232 if (size > 1) { 233 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 234 pi_oth = p_oth->i; pj_oth = p_oth->j; 235 } else { 236 p_oth = NULL; 237 pi_oth = NULL; pj_oth = NULL; 238 } 239 240 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 241 /*-------------------------------------------------------------------*/ 242 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 243 ptap->api = api; 244 api[0] = 0; 245 246 /* create and initialize a linked list */ 247 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 248 249 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 250 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 251 current_space = free_space; 252 253 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 254 for (i=0; i<am; i++) { 255 /* diagonal portion of A */ 256 nzi = adi[i+1] - adi[i]; 257 for (j=0; j<nzi; j++) { 258 row = *adj++; 259 pnz = pi_loc[row+1] - pi_loc[row]; 260 Jptr = pj_loc + pi_loc[row]; 261 /* add non-zero cols of P into the sorted linked list lnk */ 262 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 263 } 264 /* off-diagonal portion of A */ 265 nzi = aoi[i+1] - aoi[i]; 266 for (j=0; j<nzi; j++) { 267 row = *aoj++; 268 pnz = pi_oth[row+1] - pi_oth[row]; 269 Jptr = pj_oth + pi_oth[row]; 270 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 271 } 272 273 apnz = lnk[0]; 274 api[i+1] = api[i] + apnz; 275 276 /* if free space is not available, double the total space in the list */ 277 if (current_space->local_remaining<apnz) { 278 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 279 nspacedouble++; 280 } 281 282 /* Copy data into free space, then initialize lnk */ 283 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 284 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 285 286 current_space->array += apnz; 287 current_space->local_used += apnz; 288 current_space->local_remaining -= apnz; 289 } 290 291 /* Allocate space for apj, initialize apj, and */ 292 /* destroy list of free space and other temporary array(s) */ 293 ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 294 apj = ptap->apj; 295 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 296 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 297 298 /* malloc apa to store dense row A[i,:]*P */ 299 ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr); 300 301 ptap->apa = apa; 302 303 /* create and assemble symbolic parallel matrix Cmpi */ 304 /*----------------------------------------------------*/ 305 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 306 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 307 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 308 309 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 310 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 311 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 312 313 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 314 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 315 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 317 318 ptap->destroy = Cmpi->ops->destroy; 319 ptap->duplicate = Cmpi->ops->duplicate; 320 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 321 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 322 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 323 324 /* attach the supporting struct to Cmpi for reuse */ 325 c = (Mat_MPIAIJ*)Cmpi->data; 326 c->ap = ptap; 327 328 *C = Cmpi; 329 330 /* set MatInfo */ 331 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 332 if (afill < 1.0) afill = 1.0; 333 Cmpi->info.mallocs = nspacedouble; 334 Cmpi->info.fill_ratio_given = fill; 335 Cmpi->info.fill_ratio_needed = afill; 336 337 #if defined(PETSC_USE_INFO) 338 if (api[am]) { 339 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 340 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 341 } else { 342 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 343 } 344 #endif 345 PetscFunctionReturn(0); 346 } 347 348 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 349 { 350 PetscErrorCode ierr; 351 352 PetscFunctionBegin; 353 if (scall == MAT_INITIAL_MATRIX) { 354 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 355 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 356 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 357 } 358 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 359 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 360 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 typedef struct { 365 Mat workB; 366 PetscScalar *rvalues,*svalues; 367 MPI_Request *rwaits,*swaits; 368 } MPIAIJ_MPIDense; 369 370 PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx) 371 { 372 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 373 PetscErrorCode ierr; 374 375 PetscFunctionBegin; 376 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 377 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 378 ierr = PetscFree(contents);CHKERRQ(ierr); 379 PetscFunctionReturn(0); 380 } 381 382 /* 383 This is a "dummy function" that handles the case where matrix C was created as a dense matrix 384 directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option 385 386 It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C 387 */ 388 PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C) 389 { 390 PetscErrorCode ierr; 391 PetscBool flg; 392 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 393 PetscInt nz = aij->B->cmap->n,to_n,to_entries,from_n,from_entries; 394 PetscContainer container; 395 MPIAIJ_MPIDense *contents; 396 VecScatter ctx = aij->Mvctx; 397 398 PetscFunctionBegin; 399 ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr); 400 if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense"); 401 402 /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/ 403 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 404 if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ"); 405 406 C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 407 408 ierr = PetscNew(&contents);CHKERRQ(ierr); 409 /* Create work matrix used to store off processor rows of B needed for local product */ 410 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 411 /* Create work arrays needed */ 412 ierr = VecScatterGetRemoteCount_Private(ctx,PETSC_TRUE/*send*/,&to_n,&to_entries);CHKERRQ(ierr); 413 ierr = VecScatterGetRemoteCount_Private(ctx,PETSC_FALSE/*recv*/,&from_n,&from_entries);CHKERRQ(ierr); 414 ierr = PetscMalloc4(B->cmap->N*from_entries,&contents->rvalues,B->cmap->N*to_entries,&contents->svalues,from_n,&contents->rwaits,to_n,&contents->swaits);CHKERRQ(ierr); 415 416 ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 417 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 418 ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 419 ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr); 420 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 421 422 ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 427 { 428 PetscErrorCode ierr; 429 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 430 PetscInt nz = aij->B->cmap->n,to_n,to_entries,from_n,from_entries; 431 PetscContainer container; 432 MPIAIJ_MPIDense *contents; 433 VecScatter ctx = aij->Mvctx; 434 PetscInt m = A->rmap->n,n=B->cmap->n; 435 436 PetscFunctionBegin; 437 ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr); 438 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 439 ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 440 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 441 ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr); 442 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 443 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 444 445 (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 446 447 ierr = PetscNew(&contents);CHKERRQ(ierr); 448 /* Create work matrix used to store off processor rows of B needed for local product */ 449 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 450 /* Create work arrays needed */ 451 ierr = VecScatterGetRemoteCount_Private(ctx,PETSC_TRUE/*send*/,&to_n,&to_entries);CHKERRQ(ierr); 452 ierr = VecScatterGetRemoteCount_Private(ctx,PETSC_FALSE/*recv*/,&from_n,&from_entries);CHKERRQ(ierr); 453 ierr = PetscMalloc4(B->cmap->N*from_entries,&contents->rvalues,B->cmap->N*to_entries,&contents->svalues,from_n,&contents->rwaits,to_n,&contents->swaits);CHKERRQ(ierr); 454 455 ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 456 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 457 ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 458 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 459 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 460 PetscFunctionReturn(0); 461 } 462 463 /* 464 Performs an efficient scatter on the rows of B needed by this process; this is 465 a modification of the VecScatterBegin_() routines. 466 */ 467 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 468 { 469 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 470 PetscErrorCode ierr; 471 const PetscScalar *b; 472 PetscScalar *w,*svalues,*rvalues; 473 VecScatter ctx = aij->Mvctx; 474 PetscInt i,j,k; 475 const PetscInt *sindices,*sstarts,*rindices,*rstarts; 476 const PetscMPIInt *sprocs,*rprocs; 477 PetscInt nsends,nrecvs,nrecvs2; 478 MPI_Request *swaits,*rwaits; 479 MPI_Comm comm; 480 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n,nsends_mpi,nrecvs_mpi; 481 MPI_Status status; 482 MPIAIJ_MPIDense *contents; 483 PetscContainer container; 484 Mat workB; 485 486 PetscFunctionBegin; 487 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 488 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 489 if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 490 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 491 492 workB = *outworkB = contents->workB; 493 if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n); 494 ierr = VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);CHKERRQ(ierr); 495 ierr = VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,&rindices,&rprocs,NULL/*bs*/);CHKERRQ(ierr); 496 ierr = PetscMPIIntCast(nsends,&nsends_mpi);CHKERRQ(ierr); 497 ierr = PetscMPIIntCast(nrecvs,&nrecvs_mpi);CHKERRQ(ierr); 498 svalues = contents->svalues; 499 rvalues = contents->rvalues; 500 swaits = contents->swaits; 501 rwaits = contents->rwaits; 502 503 ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr); 504 ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr); 505 506 for (i=0; i<nrecvs; i++) { 507 ierr = MPI_Irecv(rvalues+ncols*(rstarts[i]-rstarts[0]),ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 508 } 509 510 for (i=0; i<nsends; i++) { 511 /* pack a message at a time */ 512 for (j=0; j<sstarts[i+1]-sstarts[i]; j++) { 513 for (k=0; k<ncols; k++) { 514 svalues[ncols*(sstarts[i]-sstarts[0]+j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 515 } 516 } 517 ierr = MPI_Isend(svalues+ncols*(sstarts[i]-sstarts[0]),ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 518 } 519 520 nrecvs2 = nrecvs; 521 while (nrecvs2) { 522 ierr = MPI_Waitany(nrecvs_mpi,rwaits,&imdex,&status);CHKERRQ(ierr); 523 nrecvs2--; 524 /* unpack a message at a time */ 525 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) { 526 for (k=0; k<ncols; k++) { 527 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex]-rstarts[0]+j) + k]; 528 } 529 } 530 } 531 if (nsends) {ierr = MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);CHKERRQ(ierr);} 532 533 ierr = VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);CHKERRQ(ierr); 534 ierr = VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,&rindices,&rprocs,NULL/*bs*/);CHKERRQ(ierr); 535 ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr); 536 ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr); 537 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 538 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 539 PetscFunctionReturn(0); 540 541 } 542 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 543 544 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 545 { 546 PetscErrorCode ierr; 547 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 548 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 549 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 550 Mat workB; 551 552 PetscFunctionBegin; 553 /* diagonal block of A times all local rows of B*/ 554 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 555 556 /* get off processor parts of B needed to complete the product */ 557 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 558 559 /* off-diagonal block of A times nonlocal rows of B */ 560 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 561 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 562 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 563 PetscFunctionReturn(0); 564 } 565 566 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 567 { 568 PetscErrorCode ierr; 569 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 570 Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 571 Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 572 PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj; 573 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 574 Mat_SeqAIJ *p_loc,*p_oth; 575 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 576 PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 577 PetscInt cm = C->rmap->n,anz,pnz; 578 Mat_APMPI *ptap = c->ap; 579 PetscScalar *apa_sparse; 580 PetscInt *api,*apj,*apJ,i,j,k,row; 581 PetscInt cstart = C->cmap->rstart; 582 PetscInt cdnz,conz,k0,k1,nextp; 583 MPI_Comm comm; 584 PetscMPIInt size; 585 586 PetscFunctionBegin; 587 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 588 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 589 590 if (!ptap->P_oth && size>1) { 591 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 592 } 593 apa_sparse = ptap->apa; 594 595 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 596 /*-----------------------------------------------------*/ 597 /* update numerical values of P_oth and P_loc */ 598 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 599 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 600 601 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 602 /*----------------------------------------------------------*/ 603 /* get data from symbolic products */ 604 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 605 pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a; 606 if (size >1) { 607 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 608 pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a; 609 } else { 610 p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL; 611 } 612 613 api = ptap->api; 614 apj = ptap->apj; 615 for (i=0; i<cm; i++) { 616 apJ = apj + api[i]; 617 618 /* diagonal portion of A */ 619 anz = adi[i+1] - adi[i]; 620 adj = ad->j + adi[i]; 621 ada = ad->a + adi[i]; 622 for (j=0; j<anz; j++) { 623 row = adj[j]; 624 pnz = pi_loc[row+1] - pi_loc[row]; 625 pj = pj_loc + pi_loc[row]; 626 pa = pa_loc + pi_loc[row]; 627 /* perform sparse axpy */ 628 valtmp = ada[j]; 629 nextp = 0; 630 for (k=0; nextp<pnz; k++) { 631 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 632 apa_sparse[k] += valtmp*pa[nextp++]; 633 } 634 } 635 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 636 } 637 638 /* off-diagonal portion of A */ 639 anz = aoi[i+1] - aoi[i]; 640 aoj = ao->j + aoi[i]; 641 aoa = ao->a + aoi[i]; 642 for (j=0; j<anz; j++) { 643 row = aoj[j]; 644 pnz = pi_oth[row+1] - pi_oth[row]; 645 pj = pj_oth + pi_oth[row]; 646 pa = pa_oth + pi_oth[row]; 647 /* perform sparse axpy */ 648 valtmp = aoa[j]; 649 nextp = 0; 650 for (k=0; nextp<pnz; k++) { 651 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 652 apa_sparse[k] += valtmp*pa[nextp++]; 653 } 654 } 655 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 656 } 657 658 /* set values in C */ 659 cdnz = cd->i[i+1] - cd->i[i]; 660 conz = co->i[i+1] - co->i[i]; 661 662 /* 1st off-diagoanl part of C */ 663 ca = coa + co->i[i]; 664 k = 0; 665 for (k0=0; k0<conz; k0++) { 666 if (apJ[k] >= cstart) break; 667 ca[k0] = apa_sparse[k]; 668 apa_sparse[k] = 0.0; 669 k++; 670 } 671 672 /* diagonal part of C */ 673 ca = cda + cd->i[i]; 674 for (k1=0; k1<cdnz; k1++) { 675 ca[k1] = apa_sparse[k]; 676 apa_sparse[k] = 0.0; 677 k++; 678 } 679 680 /* 2nd off-diagoanl part of C */ 681 ca = coa + co->i[i]; 682 for (; k0<conz; k0++) { 683 ca[k0] = apa_sparse[k]; 684 apa_sparse[k] = 0.0; 685 k++; 686 } 687 } 688 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 689 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 690 691 if (ptap->freestruct) { 692 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 693 } 694 PetscFunctionReturn(0); 695 } 696 697 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */ 698 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 699 { 700 PetscErrorCode ierr; 701 MPI_Comm comm; 702 PetscMPIInt size; 703 Mat Cmpi; 704 Mat_APMPI *ptap; 705 PetscFreeSpaceList free_space = NULL,current_space=NULL; 706 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c; 707 Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 708 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 709 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 710 PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 711 PetscInt am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20; 712 PetscReal afill; 713 PetscScalar *apa; 714 MatType mtype; 715 716 PetscFunctionBegin; 717 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 718 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 719 720 /* create struct Mat_APMPI and attached it to C later */ 721 ierr = PetscNew(&ptap);CHKERRQ(ierr); 722 723 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 724 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 725 726 /* get P_loc by taking all local rows of P */ 727 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 728 729 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 730 pi_loc = p_loc->i; pj_loc = p_loc->j; 731 if (size > 1) { 732 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 733 pi_oth = p_oth->i; pj_oth = p_oth->j; 734 } else { 735 p_oth = NULL; 736 pi_oth = NULL; pj_oth = NULL; 737 } 738 739 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 740 /*-------------------------------------------------------------------*/ 741 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 742 ptap->api = api; 743 api[0] = 0; 744 745 ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr); 746 747 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 748 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 749 current_space = free_space; 750 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 751 for (i=0; i<am; i++) { 752 /* diagonal portion of A */ 753 nzi = adi[i+1] - adi[i]; 754 for (j=0; j<nzi; j++) { 755 row = *adj++; 756 pnz = pi_loc[row+1] - pi_loc[row]; 757 Jptr = pj_loc + pi_loc[row]; 758 /* Expand list if it is not long enough */ 759 if (pnz+apnz_max > lsize) { 760 lsize = pnz+apnz_max; 761 ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 762 } 763 /* add non-zero cols of P into the sorted linked list lnk */ 764 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 765 apnz = *lnk; /* The first element in the list is the number of items in the list */ 766 api[i+1] = api[i] + apnz; 767 if (apnz > apnz_max) apnz_max = apnz; 768 } 769 /* off-diagonal portion of A */ 770 nzi = aoi[i+1] - aoi[i]; 771 for (j=0; j<nzi; j++) { 772 row = *aoj++; 773 pnz = pi_oth[row+1] - pi_oth[row]; 774 Jptr = pj_oth + pi_oth[row]; 775 /* Expand list if it is not long enough */ 776 if (pnz+apnz_max > lsize) { 777 lsize = pnz + apnz_max; 778 ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 779 } 780 /* add non-zero cols of P into the sorted linked list lnk */ 781 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 782 apnz = *lnk; /* The first element in the list is the number of items in the list */ 783 api[i+1] = api[i] + apnz; 784 if (apnz > apnz_max) apnz_max = apnz; 785 } 786 apnz = *lnk; 787 api[i+1] = api[i] + apnz; 788 if (apnz > apnz_max) apnz_max = apnz; 789 790 /* if free space is not available, double the total space in the list */ 791 if (current_space->local_remaining<apnz) { 792 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 793 nspacedouble++; 794 } 795 796 /* Copy data into free space, then initialize lnk */ 797 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 798 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 799 800 current_space->array += apnz; 801 current_space->local_used += apnz; 802 current_space->local_remaining -= apnz; 803 } 804 805 /* Allocate space for apj, initialize apj, and */ 806 /* destroy list of free space and other temporary array(s) */ 807 ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 808 apj = ptap->apj; 809 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 810 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 811 812 /* create and assemble symbolic parallel matrix Cmpi */ 813 /*----------------------------------------------------*/ 814 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 815 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 816 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 817 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 818 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 819 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 820 821 /* malloc apa for assembly Cmpi */ 822 ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr); 823 ptap->apa = apa; 824 825 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 826 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 827 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 828 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 829 830 ptap->destroy = Cmpi->ops->destroy; 831 ptap->duplicate = Cmpi->ops->duplicate; 832 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 833 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 834 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 835 836 /* attach the supporting struct to Cmpi for reuse */ 837 c = (Mat_MPIAIJ*)Cmpi->data; 838 c->ap = ptap; 839 *C = Cmpi; 840 841 /* set MatInfo */ 842 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 843 if (afill < 1.0) afill = 1.0; 844 Cmpi->info.mallocs = nspacedouble; 845 Cmpi->info.fill_ratio_given = fill; 846 Cmpi->info.fill_ratio_needed = afill; 847 848 #if defined(PETSC_USE_INFO) 849 if (api[am]) { 850 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 851 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 852 } else { 853 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 854 } 855 #endif 856 PetscFunctionReturn(0); 857 } 858 859 /* This function is needed for the seqMPI matrix-matrix multiplication. */ 860 /* Three input arrays are merged to one output array. The size of the */ 861 /* output array is also output. Duplicate entries only show up once. */ 862 static void Merge3SortedArrays(PetscInt size1, PetscInt *in1, 863 PetscInt size2, PetscInt *in2, 864 PetscInt size3, PetscInt *in3, 865 PetscInt *size4, PetscInt *out) 866 { 867 int i = 0, j = 0, k = 0, l = 0; 868 869 /* Traverse all three arrays */ 870 while (i<size1 && j<size2 && k<size3) { 871 if (in1[i] < in2[j] && in1[i] < in3[k]) { 872 out[l++] = in1[i++]; 873 } 874 else if(in2[j] < in1[i] && in2[j] < in3[k]) { 875 out[l++] = in2[j++]; 876 } 877 else if(in3[k] < in1[i] && in3[k] < in2[j]) { 878 out[l++] = in3[k++]; 879 } 880 else if(in1[i] == in2[j] && in1[i] < in3[k]) { 881 out[l++] = in1[i]; 882 i++, j++; 883 } 884 else if(in1[i] == in3[k] && in1[i] < in2[j]) { 885 out[l++] = in1[i]; 886 i++, k++; 887 } 888 else if(in3[k] == in2[j] && in2[j] < in1[i]) { 889 out[l++] = in2[j]; 890 k++, j++; 891 } 892 else if(in1[i] == in2[j] && in1[i] == in3[k]) { 893 out[l++] = in1[i]; 894 i++, j++, k++; 895 } 896 } 897 898 /* Traverse two remaining arrays */ 899 while (i<size1 && j<size2) { 900 if (in1[i] < in2[j]) { 901 out[l++] = in1[i++]; 902 } 903 else if(in1[i] > in2[j]) { 904 out[l++] = in2[j++]; 905 } 906 else { 907 out[l++] = in1[i]; 908 i++, j++; 909 } 910 } 911 912 while (i<size1 && k<size3) { 913 if (in1[i] < in3[k]) { 914 out[l++] = in1[i++]; 915 } 916 else if(in1[i] > in3[k]) { 917 out[l++] = in3[k++]; 918 } 919 else { 920 out[l++] = in1[i]; 921 i++, k++; 922 } 923 } 924 925 while (k<size3 && j<size2) { 926 if (in3[k] < in2[j]) { 927 out[l++] = in3[k++]; 928 } 929 else if(in3[k] > in2[j]) { 930 out[l++] = in2[j++]; 931 } 932 else { 933 out[l++] = in3[k]; 934 k++, j++; 935 } 936 } 937 938 /* Traverse one remaining array */ 939 while (i<size1) out[l++] = in1[i++]; 940 while (j<size2) out[l++] = in2[j++]; 941 while (k<size3) out[l++] = in3[k++]; 942 943 *size4 = l; 944 } 945 946 /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */ 947 /* adds up the products. Two of these three multiplications are performed with existing (sequential) */ 948 /* matrix-matrix multiplications. */ 949 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C) 950 { 951 PetscErrorCode ierr; 952 MPI_Comm comm; 953 PetscMPIInt size; 954 Mat Cmpi; 955 Mat_APMPI *ptap; 956 PetscFreeSpaceList free_space_diag=NULL, current_space=NULL; 957 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data; 958 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc; 959 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data; 960 Mat_MPIAIJ *c; 961 Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq; 962 PetscInt adponz, adpdnz; 963 PetscInt *pi_loc,*dnz,*onz; 964 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart; 965 PetscInt *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi, 966 *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj; 967 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend; 968 PetscBT lnkbt; 969 PetscScalar *apa; 970 PetscReal afill; 971 PetscMPIInt rank; 972 Mat adpd, aopoth; 973 MatType mtype; 974 const char *prefix; 975 976 PetscFunctionBegin; 977 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 978 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 979 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 980 ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend); CHKERRQ(ierr); 981 982 /* create struct Mat_APMPI and attached it to C later */ 983 ierr = PetscNew(&ptap);CHKERRQ(ierr); 984 985 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 986 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 987 988 /* get P_loc by taking all local rows of P */ 989 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 990 991 992 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 993 pi_loc = p_loc->i; 994 995 /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */ 996 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 997 ierr = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr); 998 999 adpoi[0] = 0; 1000 ptap->api = api; 1001 api[0] = 0; 1002 1003 /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */ 1004 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1005 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 1006 1007 /* Symbolic calc of A_loc_diag * P_loc_diag */ 1008 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1009 ierr = MatSetOptionsPrefix(a->A,prefix);CHKERRQ(ierr); 1010 ierr = MatAppendOptionsPrefix(a->A,"inner_diag_");CHKERRQ(ierr); 1011 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);CHKERRQ(ierr); 1012 adpd_seq = (Mat_SeqAIJ*)((adpd)->data); 1013 adpdi = adpd_seq->i; adpdj = adpd_seq->j; 1014 p_off = (Mat_SeqAIJ*)((p->B)->data); 1015 poff_i = p_off->i; poff_j = p_off->j; 1016 1017 /* j_temp stores indices of a result row before they are added to the linked list */ 1018 ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr); 1019 1020 1021 /* Symbolic calc of the A_diag * p_loc_off */ 1022 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 1023 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr); 1024 current_space = free_space_diag; 1025 1026 for (i=0; i<am; i++) { 1027 /* A_diag * P_loc_off */ 1028 nzi = adi[i+1] - adi[i]; 1029 for (j=0; j<nzi; j++) { 1030 row = *adj++; 1031 pnz = poff_i[row+1] - poff_i[row]; 1032 Jptr = poff_j + poff_i[row]; 1033 for(i1 = 0; i1 < pnz; i1++) { 1034 j_temp[i1] = p->garray[Jptr[i1]]; 1035 } 1036 /* add non-zero cols of P into the sorted linked list lnk */ 1037 ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr); 1038 } 1039 1040 adponz = lnk[0]; 1041 adpoi[i+1] = adpoi[i] + adponz; 1042 1043 /* if free space is not available, double the total space in the list */ 1044 if (current_space->local_remaining<adponz) { 1045 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1046 nspacedouble++; 1047 } 1048 1049 /* Copy data into free space, then initialize lnk */ 1050 ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1051 1052 current_space->array += adponz; 1053 current_space->local_used += adponz; 1054 current_space->local_remaining -= adponz; 1055 } 1056 1057 /* Symbolic calc of A_off * P_oth */ 1058 ierr = MatSetOptionsPrefix(a->B,prefix);CHKERRQ(ierr); 1059 ierr = MatAppendOptionsPrefix(a->B,"inner_offdiag_");CHKERRQ(ierr); 1060 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);CHKERRQ(ierr); 1061 aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data); 1062 aopothi = aopoth_seq->i; aopothj = aopoth_seq->j; 1063 1064 /* Allocate space for apj, adpj, aopj, ... */ 1065 /* destroy lists of free space and other temporary array(s) */ 1066 1067 ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr); 1068 ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr); 1069 1070 /* Copy from linked list to j-array */ 1071 ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr); 1072 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1073 1074 adpoJ = adpoj; 1075 adpdJ = adpdj; 1076 aopJ = aopothj; 1077 apj = ptap->apj; 1078 apJ = apj; /* still empty */ 1079 1080 /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */ 1081 /* A_diag * P_loc_diag to get A*P */ 1082 for (i = 0; i < am; i++) { 1083 aopnz = aopothi[i+1] - aopothi[i]; 1084 adponz = adpoi[i+1] - adpoi[i]; 1085 adpdnz = adpdi[i+1] - adpdi[i]; 1086 1087 /* Correct indices from A_diag*P_diag */ 1088 for(i1 = 0; i1 < adpdnz; i1++) { 1089 adpdJ[i1] += p_colstart; 1090 } 1091 /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */ 1092 Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ); 1093 ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz); CHKERRQ(ierr); 1094 1095 aopJ += aopnz; 1096 adpoJ += adponz; 1097 adpdJ += adpdnz; 1098 apJ += apnz; 1099 api[i+1] = api[i] + apnz; 1100 } 1101 1102 /* malloc apa to store dense row A[i,:]*P */ 1103 ierr = PetscCalloc1(pN+2,&apa);CHKERRQ(ierr); 1104 1105 ptap->apa = apa; 1106 /* create and assemble symbolic parallel matrix Cmpi */ 1107 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1108 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1109 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 1110 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1111 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1112 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1113 1114 1115 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 1116 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1117 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1118 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1119 1120 1121 ptap->destroy = Cmpi->ops->destroy; 1122 ptap->duplicate = Cmpi->ops->duplicate; 1123 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1124 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1125 1126 /* attach the supporting struct to Cmpi for reuse */ 1127 c = (Mat_MPIAIJ*)Cmpi->data; 1128 c->ap = ptap; 1129 *C = Cmpi; 1130 1131 /* set MatInfo */ 1132 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 1133 if (afill < 1.0) afill = 1.0; 1134 Cmpi->info.mallocs = nspacedouble; 1135 Cmpi->info.fill_ratio_given = fill; 1136 Cmpi->info.fill_ratio_needed = afill; 1137 1138 #if defined(PETSC_USE_INFO) 1139 if (api[am]) { 1140 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1141 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1142 } else { 1143 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1144 } 1145 #endif 1146 1147 ierr = MatDestroy(&aopoth);CHKERRQ(ierr); 1148 ierr = MatDestroy(&adpd);CHKERRQ(ierr); 1149 ierr = PetscFree(j_temp);CHKERRQ(ierr); 1150 ierr = PetscFree(adpoj);CHKERRQ(ierr); 1151 ierr = PetscFree(adpoi);CHKERRQ(ierr); 1152 PetscFunctionReturn(0); 1153 } 1154 1155 1156 /*-------------------------------------------------------------------------*/ 1157 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 1158 { 1159 PetscErrorCode ierr; 1160 const char *algTypes[3] = {"scalable","nonscalable","matmatmult"}; 1161 PetscInt aN=A->cmap->N,alg=1; /* set default algorithm */ 1162 PetscBool flg; 1163 1164 PetscFunctionBegin; 1165 if (scall == MAT_INITIAL_MATRIX) { 1166 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr); 1167 ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);CHKERRQ(ierr); 1168 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1169 1170 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 1171 switch (alg) { 1172 case 1: 1173 if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */ 1174 MatInfo Ainfo,Pinfo; 1175 PetscInt nz_local; 1176 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 1177 MPI_Comm comm; 1178 1179 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 1180 ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 1181 nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */ 1182 1183 if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 1184 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1185 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 1186 1187 if (alg_scalable) { 1188 alg = 0; /* scalable algorithm would slower than nonscalable algorithm */ 1189 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1190 break; 1191 } 1192 } 1193 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr); 1194 break; 1195 case 2: 1196 { 1197 Mat Pt; 1198 Mat_APMPI *ptap; 1199 Mat_MPIAIJ *c; 1200 ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 1201 ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 1202 c = (Mat_MPIAIJ*)(*C)->data; 1203 ptap = c->ap; 1204 if (ptap) { 1205 ptap->Pt = Pt; 1206 (*C)->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1207 } 1208 (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 1209 PetscFunctionReturn(0); 1210 } 1211 break; 1212 default: /* scalable algorithm */ 1213 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1214 break; 1215 } 1216 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 1217 1218 { 1219 Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 1220 Mat_APMPI *ap = c->ap; 1221 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 1222 ap->freestruct = PETSC_FALSE; 1223 ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 1224 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1225 } 1226 } 1227 1228 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 1229 ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 1230 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 1231 PetscFunctionReturn(0); 1232 } 1233 1234 /* This routine only works when scall=MAT_REUSE_MATRIX! */ 1235 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 1236 { 1237 PetscErrorCode ierr; 1238 Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data; 1239 Mat_APMPI *ptap= c->ap; 1240 Mat Pt; 1241 1242 PetscFunctionBegin; 1243 if (!ptap->Pt) { 1244 MPI_Comm comm; 1245 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1246 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1247 } 1248 1249 Pt=ptap->Pt; 1250 ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 1251 ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr); 1252 1253 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1254 if (ptap->freestruct) { 1255 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1256 } 1257 PetscFunctionReturn(0); 1258 } 1259 1260 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1261 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C) 1262 { 1263 PetscErrorCode ierr; 1264 Mat_APMPI *ptap; 1265 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1266 MPI_Comm comm; 1267 PetscMPIInt size,rank; 1268 Mat Cmpi; 1269 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1270 PetscInt pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n; 1271 PetscInt *lnk,i,k,nsend; 1272 PetscBT lnkbt; 1273 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 1274 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1275 PetscInt len,proc,*dnz,*onz,*owners,nzi; 1276 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1277 MPI_Request *swaits,*rwaits; 1278 MPI_Status *sstatus,rstatus; 1279 PetscLayout rowmap; 1280 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1281 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 1282 PetscInt *Jptr,*prmap=p->garray,con,j,Crmax; 1283 Mat_SeqAIJ *a_loc,*c_loc,*c_oth; 1284 PetscTable ta; 1285 MatType mtype; 1286 const char *prefix; 1287 1288 PetscFunctionBegin; 1289 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1290 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1291 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1292 1293 /* create symbolic parallel matrix Cmpi */ 1294 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1295 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1296 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1297 1298 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1299 1300 /* create struct Mat_APMPI and attached it to C later */ 1301 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1302 ptap->reuse = MAT_INITIAL_MATRIX; 1303 1304 /* (0) compute Rd = Pd^T, Ro = Po^T */ 1305 /* --------------------------------- */ 1306 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1307 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1308 1309 /* (1) compute symbolic A_loc */ 1310 /* ---------------------------*/ 1311 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1312 1313 /* (2-1) compute symbolic C_oth = Ro*A_loc */ 1314 /* ------------------------------------ */ 1315 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1316 ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 1317 ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 1318 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 1319 1320 /* (3) send coj of C_oth to other processors */ 1321 /* ------------------------------------------ */ 1322 /* determine row ownership */ 1323 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1324 rowmap->n = pn; 1325 rowmap->bs = 1; 1326 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1327 owners = rowmap->range; 1328 1329 /* determine the number of messages to send, their lengths */ 1330 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 1331 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1332 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1333 1334 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 1335 coi = c_oth->i; coj = c_oth->j; 1336 con = ptap->C_oth->rmap->n; 1337 proc = 0; 1338 for (i=0; i<con; i++) { 1339 while (prmap[i] >= owners[proc+1]) proc++; 1340 len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */ 1341 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1342 } 1343 1344 len = 0; /* max length of buf_si[], see (4) */ 1345 owners_co[0] = 0; 1346 nsend = 0; 1347 for (proc=0; proc<size; proc++) { 1348 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1349 if (len_s[proc]) { 1350 nsend++; 1351 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1352 len += len_si[proc]; 1353 } 1354 } 1355 1356 /* determine the number and length of messages to receive for coi and coj */ 1357 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1358 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 1359 1360 /* post the Irecv and Isend of coj */ 1361 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1362 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1363 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 1364 for (proc=0, k=0; proc<size; proc++) { 1365 if (!len_s[proc]) continue; 1366 i = owners_co[proc]; 1367 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1368 k++; 1369 } 1370 1371 /* (2-2) compute symbolic C_loc = Rd*A_loc */ 1372 /* ---------------------------------------- */ 1373 ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 1374 ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 1375 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 1376 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1377 1378 /* receives coj are complete */ 1379 for (i=0; i<nrecv; i++) { 1380 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1381 } 1382 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1383 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1384 1385 /* add received column indices into ta to update Crmax */ 1386 a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data; 1387 1388 /* create and initialize a linked list */ 1389 ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */ 1390 MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta); 1391 1392 for (k=0; k<nrecv; k++) {/* k-th received message */ 1393 Jptr = buf_rj[k]; 1394 for (j=0; j<len_r[k]; j++) { 1395 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1396 } 1397 } 1398 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 1399 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 1400 1401 /* (4) send and recv coi */ 1402 /*-----------------------*/ 1403 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1404 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1405 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1406 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1407 for (proc=0,k=0; proc<size; proc++) { 1408 if (!len_s[proc]) continue; 1409 /* form outgoing message for i-structure: 1410 buf_si[0]: nrows to be sent 1411 [1:nrows]: row index (global) 1412 [nrows+1:2*nrows+1]: i-structure index 1413 */ 1414 /*-------------------------------------------*/ 1415 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1416 buf_si_i = buf_si + nrows+1; 1417 buf_si[0] = nrows; 1418 buf_si_i[0] = 0; 1419 nrows = 0; 1420 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1421 nzi = coi[i+1] - coi[i]; 1422 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1423 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1424 nrows++; 1425 } 1426 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1427 k++; 1428 buf_si += len_si[proc]; 1429 } 1430 for (i=0; i<nrecv; i++) { 1431 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1432 } 1433 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1434 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1435 1436 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 1437 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1438 ierr = PetscFree(swaits);CHKERRQ(ierr); 1439 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1440 1441 /* (5) compute the local portion of Cmpi */ 1442 /* ------------------------------------------ */ 1443 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 1444 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 1445 current_space = free_space; 1446 1447 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1448 for (k=0; k<nrecv; k++) { 1449 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1450 nrows = *buf_ri_k[k]; 1451 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1452 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1453 } 1454 1455 ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr); 1456 ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1457 for (i=0; i<pn; i++) { 1458 /* add C_loc into Cmpi */ 1459 nzi = c_loc->i[i+1] - c_loc->i[i]; 1460 Jptr = c_loc->j + c_loc->i[i]; 1461 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1462 1463 /* add received col data into lnk */ 1464 for (k=0; k<nrecv; k++) { /* k-th received message */ 1465 if (i == *nextrow[k]) { /* i-th row */ 1466 nzi = *(nextci[k]+1) - *nextci[k]; 1467 Jptr = buf_rj[k] + *nextci[k]; 1468 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1469 nextrow[k]++; nextci[k]++; 1470 } 1471 } 1472 nzi = lnk[0]; 1473 1474 /* copy data into free space, then initialize lnk */ 1475 ierr = PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1476 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 1477 } 1478 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1479 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1480 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 1481 1482 /* local sizes and preallocation */ 1483 ierr = MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1484 if (P->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr);} 1485 if (A->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(Cmpi->cmap,A->cmap->bs);CHKERRQ(ierr);} 1486 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1487 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1488 1489 /* members in merge */ 1490 ierr = PetscFree(id_r);CHKERRQ(ierr); 1491 ierr = PetscFree(len_r);CHKERRQ(ierr); 1492 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1493 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1494 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1495 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1496 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 1497 1498 /* attach the supporting struct to Cmpi for reuse */ 1499 c = (Mat_MPIAIJ*)Cmpi->data; 1500 c->ap = ptap; 1501 ptap->destroy = Cmpi->ops->destroy; 1502 1503 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1504 Cmpi->assembled = PETSC_FALSE; 1505 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1506 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1507 1508 *C = Cmpi; 1509 PetscFunctionReturn(0); 1510 } 1511 1512 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 1513 { 1514 PetscErrorCode ierr; 1515 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1516 Mat_SeqAIJ *c_seq; 1517 Mat_APMPI *ptap = c->ap; 1518 Mat A_loc,C_loc,C_oth; 1519 PetscInt i,rstart,rend,cm,ncols,row; 1520 const PetscInt *cols; 1521 const PetscScalar *vals; 1522 1523 PetscFunctionBegin; 1524 if (!ptap->A_loc) { 1525 MPI_Comm comm; 1526 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1527 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1528 } 1529 1530 ierr = MatZeroEntries(C);CHKERRQ(ierr); 1531 1532 if (ptap->reuse == MAT_REUSE_MATRIX) { 1533 /* These matrices are obtained in MatTransposeMatMultSymbolic() */ 1534 /* 1) get R = Pd^T, Ro = Po^T */ 1535 /*----------------------------*/ 1536 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1537 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1538 1539 /* 2) compute numeric A_loc */ 1540 /*--------------------------*/ 1541 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1542 } 1543 1544 /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */ 1545 A_loc = ptap->A_loc; 1546 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr); 1547 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr); 1548 C_loc = ptap->C_loc; 1549 C_oth = ptap->C_oth; 1550 1551 /* add C_loc and Co to to C */ 1552 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1553 1554 /* C_loc -> C */ 1555 cm = C_loc->rmap->N; 1556 c_seq = (Mat_SeqAIJ*)C_loc->data; 1557 cols = c_seq->j; 1558 vals = c_seq->a; 1559 for (i=0; i<cm; i++) { 1560 ncols = c_seq->i[i+1] - c_seq->i[i]; 1561 row = rstart + i; 1562 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1563 cols += ncols; vals += ncols; 1564 } 1565 1566 /* Co -> C, off-processor part */ 1567 cm = C_oth->rmap->N; 1568 c_seq = (Mat_SeqAIJ*)C_oth->data; 1569 cols = c_seq->j; 1570 vals = c_seq->a; 1571 for (i=0; i<cm; i++) { 1572 ncols = c_seq->i[i+1] - c_seq->i[i]; 1573 row = p->garray[i]; 1574 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1575 cols += ncols; vals += ncols; 1576 } 1577 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1578 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1579 1580 ptap->reuse = MAT_REUSE_MATRIX; 1581 1582 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1583 if (ptap->freestruct) { 1584 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1585 } 1586 PetscFunctionReturn(0); 1587 } 1588 1589 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1590 { 1591 PetscErrorCode ierr; 1592 Mat_Merge_SeqsToMPI *merge; 1593 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1594 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1595 Mat_APMPI *ptap; 1596 PetscInt *adj; 1597 PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1598 MatScalar *ada,*ca,valtmp; 1599 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1600 MPI_Comm comm; 1601 PetscMPIInt size,rank,taga,*len_s; 1602 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1603 PetscInt **buf_ri,**buf_rj; 1604 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1605 MPI_Request *s_waits,*r_waits; 1606 MPI_Status *status; 1607 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1608 PetscInt *ai,*aj,*coi,*coj,*poJ,*pdJ; 1609 Mat A_loc; 1610 Mat_SeqAIJ *a_loc; 1611 1612 PetscFunctionBegin; 1613 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1614 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1615 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1616 1617 ptap = c->ap; 1618 if (!ptap->A_loc) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1619 merge = ptap->merge; 1620 1621 /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1622 /*------------------------------------------*/ 1623 /* get data from symbolic products */ 1624 coi = merge->coi; coj = merge->coj; 1625 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1626 bi = merge->bi; bj = merge->bj; 1627 owners = merge->rowmap->range; 1628 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1629 1630 /* get A_loc by taking all local rows of A */ 1631 A_loc = ptap->A_loc; 1632 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1633 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1634 ai = a_loc->i; 1635 aj = a_loc->j; 1636 1637 for (i=0; i<am; i++) { 1638 anz = ai[i+1] - ai[i]; 1639 adj = aj + ai[i]; 1640 ada = a_loc->a + ai[i]; 1641 1642 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1643 /*-------------------------------------------------------------*/ 1644 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1645 pnz = po->i[i+1] - po->i[i]; 1646 poJ = po->j + po->i[i]; 1647 pA = po->a + po->i[i]; 1648 for (j=0; j<pnz; j++) { 1649 row = poJ[j]; 1650 cj = coj + coi[row]; 1651 ca = coa + coi[row]; 1652 /* perform sparse axpy */ 1653 nexta = 0; 1654 valtmp = pA[j]; 1655 for (k=0; nexta<anz; k++) { 1656 if (cj[k] == adj[nexta]) { 1657 ca[k] += valtmp*ada[nexta]; 1658 nexta++; 1659 } 1660 } 1661 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1662 } 1663 1664 /* put the value into Cd (diagonal part) */ 1665 pnz = pd->i[i+1] - pd->i[i]; 1666 pdJ = pd->j + pd->i[i]; 1667 pA = pd->a + pd->i[i]; 1668 for (j=0; j<pnz; j++) { 1669 row = pdJ[j]; 1670 cj = bj + bi[row]; 1671 ca = ba + bi[row]; 1672 /* perform sparse axpy */ 1673 nexta = 0; 1674 valtmp = pA[j]; 1675 for (k=0; nexta<anz; k++) { 1676 if (cj[k] == adj[nexta]) { 1677 ca[k] += valtmp*ada[nexta]; 1678 nexta++; 1679 } 1680 } 1681 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1682 } 1683 } 1684 1685 /* 3) send and recv matrix values coa */ 1686 /*------------------------------------*/ 1687 buf_ri = merge->buf_ri; 1688 buf_rj = merge->buf_rj; 1689 len_s = merge->len_s; 1690 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1691 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1692 1693 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1694 for (proc=0,k=0; proc<size; proc++) { 1695 if (!len_s[proc]) continue; 1696 i = merge->owners_co[proc]; 1697 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1698 k++; 1699 } 1700 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1701 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1702 1703 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1704 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1705 ierr = PetscFree(coa);CHKERRQ(ierr); 1706 1707 /* 4) insert local Cseq and received values into Cmpi */ 1708 /*----------------------------------------------------*/ 1709 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1710 for (k=0; k<merge->nrecv; k++) { 1711 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1712 nrows = *(buf_ri_k[k]); 1713 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1714 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1715 } 1716 1717 for (i=0; i<cm; i++) { 1718 row = owners[rank] + i; /* global row index of C_seq */ 1719 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1720 ba_i = ba + bi[i]; 1721 bnz = bi[i+1] - bi[i]; 1722 /* add received vals into ba */ 1723 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1724 /* i-th row */ 1725 if (i == *nextrow[k]) { 1726 cnz = *(nextci[k]+1) - *nextci[k]; 1727 cj = buf_rj[k] + *(nextci[k]); 1728 ca = abuf_r[k] + *(nextci[k]); 1729 nextcj = 0; 1730 for (j=0; nextcj<cnz; j++) { 1731 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1732 ba_i[j] += ca[nextcj++]; 1733 } 1734 } 1735 nextrow[k]++; nextci[k]++; 1736 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1737 } 1738 } 1739 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1740 } 1741 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1742 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1743 1744 ierr = PetscFree(ba);CHKERRQ(ierr); 1745 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1746 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1747 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1748 1749 if (ptap->freestruct) { 1750 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1751 } 1752 PetscFunctionReturn(0); 1753 } 1754 1755 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1756 { 1757 PetscErrorCode ierr; 1758 Mat Cmpi,A_loc,POt,PDt; 1759 Mat_APMPI *ptap; 1760 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1761 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c; 1762 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1763 PetscInt nnz; 1764 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1765 PetscInt am =A->rmap->n,pn=P->cmap->n; 1766 MPI_Comm comm; 1767 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1768 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1769 PetscInt len,proc,*dnz,*onz,*owners; 1770 PetscInt nzi,*bi,*bj; 1771 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1772 MPI_Request *swaits,*rwaits; 1773 MPI_Status *sstatus,rstatus; 1774 Mat_Merge_SeqsToMPI *merge; 1775 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1776 PetscReal afill =1.0,afill_tmp; 1777 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax; 1778 PetscScalar *vals; 1779 Mat_SeqAIJ *a_loc,*pdt,*pot; 1780 PetscTable ta; 1781 MatType mtype; 1782 1783 PetscFunctionBegin; 1784 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1785 /* check if matrix local sizes are compatible */ 1786 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 1787 1788 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1789 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1790 1791 /* create struct Mat_APMPI and attached it to C later */ 1792 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1793 1794 /* get A_loc by taking all local rows of A */ 1795 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1796 1797 ptap->A_loc = A_loc; 1798 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1799 ai = a_loc->i; 1800 aj = a_loc->j; 1801 1802 /* determine symbolic Co=(p->B)^T*A - send to others */ 1803 /*----------------------------------------------------*/ 1804 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1805 pdt = (Mat_SeqAIJ*)PDt->data; 1806 pdti = pdt->i; pdtj = pdt->j; 1807 1808 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1809 pot = (Mat_SeqAIJ*)POt->data; 1810 poti = pot->i; potj = pot->j; 1811 1812 /* then, compute symbolic Co = (p->B)^T*A */ 1813 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1814 >= (num of nonzero rows of C_seq) - pn */ 1815 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1816 coi[0] = 0; 1817 1818 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1819 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am])); 1820 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1821 current_space = free_space; 1822 1823 /* create and initialize a linked list */ 1824 ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr); 1825 MatRowMergeMax_SeqAIJ(a_loc,am,ta); 1826 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1827 1828 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 1829 1830 for (i=0; i<pon; i++) { 1831 pnz = poti[i+1] - poti[i]; 1832 ptJ = potj + poti[i]; 1833 for (j=0; j<pnz; j++) { 1834 row = ptJ[j]; /* row of A_loc == col of Pot */ 1835 anz = ai[row+1] - ai[row]; 1836 Jptr = aj + ai[row]; 1837 /* add non-zero cols of AP into the sorted linked list lnk */ 1838 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1839 } 1840 nnz = lnk[0]; 1841 1842 /* If free space is not available, double the total space in the list */ 1843 if (current_space->local_remaining<nnz) { 1844 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1845 nspacedouble++; 1846 } 1847 1848 /* Copy data into free space, and zero out denserows */ 1849 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1850 1851 current_space->array += nnz; 1852 current_space->local_used += nnz; 1853 current_space->local_remaining -= nnz; 1854 1855 coi[i+1] = coi[i] + nnz; 1856 } 1857 1858 ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1859 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1860 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */ 1861 1862 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1863 if (afill_tmp > afill) afill = afill_tmp; 1864 1865 /* send j-array (coj) of Co to other processors */ 1866 /*----------------------------------------------*/ 1867 /* determine row ownership */ 1868 ierr = PetscNew(&merge);CHKERRQ(ierr); 1869 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1870 1871 merge->rowmap->n = pn; 1872 merge->rowmap->bs = 1; 1873 1874 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1875 owners = merge->rowmap->range; 1876 1877 /* determine the number of messages to send, their lengths */ 1878 ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1879 ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 1880 1881 len_s = merge->len_s; 1882 merge->nsend = 0; 1883 1884 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1885 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1886 1887 proc = 0; 1888 for (i=0; i<pon; i++) { 1889 while (prmap[i] >= owners[proc+1]) proc++; 1890 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1891 len_s[proc] += coi[i+1] - coi[i]; 1892 } 1893 1894 len = 0; /* max length of buf_si[] */ 1895 owners_co[0] = 0; 1896 for (proc=0; proc<size; proc++) { 1897 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1898 if (len_si[proc]) { 1899 merge->nsend++; 1900 len_si[proc] = 2*(len_si[proc] + 1); 1901 len += len_si[proc]; 1902 } 1903 } 1904 1905 /* determine the number and length of messages to receive for coi and coj */ 1906 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1907 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1908 1909 /* post the Irecv and Isend of coj */ 1910 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1911 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1912 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1913 for (proc=0, k=0; proc<size; proc++) { 1914 if (!len_s[proc]) continue; 1915 i = owners_co[proc]; 1916 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1917 k++; 1918 } 1919 1920 /* receives and sends of coj are complete */ 1921 ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1922 for (i=0; i<merge->nrecv; i++) { 1923 PetscMPIInt icompleted; 1924 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1925 } 1926 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1927 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1928 1929 /* add received column indices into table to update Armax */ 1930 /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */ 1931 for (k=0; k<merge->nrecv; k++) {/* k-th received message */ 1932 Jptr = buf_rj[k]; 1933 for (j=0; j<merge->len_r[k]; j++) { 1934 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1935 } 1936 } 1937 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1938 /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */ 1939 1940 /* send and recv coi */ 1941 /*-------------------*/ 1942 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1943 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1944 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1945 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1946 for (proc=0,k=0; proc<size; proc++) { 1947 if (!len_s[proc]) continue; 1948 /* form outgoing message for i-structure: 1949 buf_si[0]: nrows to be sent 1950 [1:nrows]: row index (global) 1951 [nrows+1:2*nrows+1]: i-structure index 1952 */ 1953 /*-------------------------------------------*/ 1954 nrows = len_si[proc]/2 - 1; 1955 buf_si_i = buf_si + nrows+1; 1956 buf_si[0] = nrows; 1957 buf_si_i[0] = 0; 1958 nrows = 0; 1959 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1960 nzi = coi[i+1] - coi[i]; 1961 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1962 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1963 nrows++; 1964 } 1965 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1966 k++; 1967 buf_si += len_si[proc]; 1968 } 1969 i = merge->nrecv; 1970 while (i--) { 1971 PetscMPIInt icompleted; 1972 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1973 } 1974 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1975 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1976 ierr = PetscFree(len_si);CHKERRQ(ierr); 1977 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1978 ierr = PetscFree(swaits);CHKERRQ(ierr); 1979 ierr = PetscFree(sstatus);CHKERRQ(ierr); 1980 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1981 1982 /* compute the local portion of C (mpi mat) */ 1983 /*------------------------------------------*/ 1984 /* allocate bi array and free space for accumulating nonzero column info */ 1985 ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1986 bi[0] = 0; 1987 1988 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1989 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am]))); 1990 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1991 current_space = free_space; 1992 1993 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1994 for (k=0; k<merge->nrecv; k++) { 1995 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1996 nrows = *buf_ri_k[k]; 1997 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1998 nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */ 1999 } 2000 2001 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 2002 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2003 rmax = 0; 2004 for (i=0; i<pn; i++) { 2005 /* add pdt[i,:]*AP into lnk */ 2006 pnz = pdti[i+1] - pdti[i]; 2007 ptJ = pdtj + pdti[i]; 2008 for (j=0; j<pnz; j++) { 2009 row = ptJ[j]; /* row of AP == col of Pt */ 2010 anz = ai[row+1] - ai[row]; 2011 Jptr = aj + ai[row]; 2012 /* add non-zero cols of AP into the sorted linked list lnk */ 2013 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2014 } 2015 2016 /* add received col data into lnk */ 2017 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 2018 if (i == *nextrow[k]) { /* i-th row */ 2019 nzi = *(nextci[k]+1) - *nextci[k]; 2020 Jptr = buf_rj[k] + *nextci[k]; 2021 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2022 nextrow[k]++; nextci[k]++; 2023 } 2024 } 2025 nnz = lnk[0]; 2026 2027 /* if free space is not available, make more free space */ 2028 if (current_space->local_remaining<nnz) { 2029 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 2030 nspacedouble++; 2031 } 2032 /* copy data into free space, then initialize lnk */ 2033 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2034 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2035 2036 current_space->array += nnz; 2037 current_space->local_used += nnz; 2038 current_space->local_remaining -= nnz; 2039 2040 bi[i+1] = bi[i] + nnz; 2041 if (nnz > rmax) rmax = nnz; 2042 } 2043 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2044 2045 ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 2046 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2047 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 2048 if (afill_tmp > afill) afill = afill_tmp; 2049 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 2050 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 2051 2052 ierr = MatDestroy(&POt);CHKERRQ(ierr); 2053 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 2054 2055 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 2056 /*----------------------------------------------------------------------------------*/ 2057 ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 2058 2059 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 2060 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 2061 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 2062 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 2063 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 2064 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 2065 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2066 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 2067 for (i=0; i<pn; i++) { 2068 row = i + rstart; 2069 nnz = bi[i+1] - bi[i]; 2070 Jptr = bj + bi[i]; 2071 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 2072 } 2073 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2074 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2075 ierr = PetscFree(vals);CHKERRQ(ierr); 2076 2077 merge->bi = bi; 2078 merge->bj = bj; 2079 merge->coi = coi; 2080 merge->coj = coj; 2081 merge->buf_ri = buf_ri; 2082 merge->buf_rj = buf_rj; 2083 merge->owners_co = owners_co; 2084 2085 /* attach the supporting struct to Cmpi for reuse */ 2086 c = (Mat_MPIAIJ*)Cmpi->data; 2087 2088 c->ap = ptap; 2089 ptap->api = NULL; 2090 ptap->apj = NULL; 2091 ptap->merge = merge; 2092 ptap->apa = NULL; 2093 ptap->destroy = Cmpi->ops->destroy; 2094 ptap->duplicate = Cmpi->ops->duplicate; 2095 2096 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 2097 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 2098 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 2099 2100 *C = Cmpi; 2101 #if defined(PETSC_USE_INFO) 2102 if (bi[pn] != 0) { 2103 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 2104 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 2105 } else { 2106 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 2107 } 2108 #endif 2109 PetscFunctionReturn(0); 2110 } 2111