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