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 PetscCheckFalse(n <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Column block size %" PetscInt_FMT " must be positive",n); 615 616 for (i=0; i<BN; i+=n) { 617 ierr = MatDenseGetSubMatrix(B,i,PetscMin(i+n,BN),&Bb);CHKERRQ(ierr); 618 ierr = MatDenseGetSubMatrix(C,i,PetscMin(i+n,BN),&Cb);CHKERRQ(ierr); 619 620 /* get off processor parts of B needed to complete C=A*B */ 621 ierr = MatMPIDenseScatter(A,Bb,i+n>BN,C,&workB);CHKERRQ(ierr); 622 623 /* off-diagonal block of A times nonlocal rows of B */ 624 cdense = (Mat_MPIDense*)Cb->data; 625 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);CHKERRQ(ierr); 626 627 ierr = MatDenseRestoreSubMatrix(B,&Bb);CHKERRQ(ierr); 628 ierr = MatDenseRestoreSubMatrix(C,&Cb);CHKERRQ(ierr); 629 } 630 } 631 PetscFunctionReturn(0); 632 } 633 634 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 635 { 636 PetscErrorCode ierr; 637 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 638 Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 639 Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 640 PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj; 641 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 642 Mat_SeqAIJ *p_loc,*p_oth; 643 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 644 PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 645 PetscInt cm = C->rmap->n,anz,pnz; 646 Mat_APMPI *ptap; 647 PetscScalar *apa_sparse; 648 const PetscScalar *dummy; 649 PetscInt *api,*apj,*apJ,i,j,k,row; 650 PetscInt cstart = C->cmap->rstart; 651 PetscInt cdnz,conz,k0,k1,nextp; 652 MPI_Comm comm; 653 PetscMPIInt size; 654 655 PetscFunctionBegin; 656 MatCheckProduct(C,3); 657 ptap = (Mat_APMPI*)C->product->data; 658 PetscCheckFalse(!ptap,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 659 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 660 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 661 PetscCheckFalse(!ptap->P_oth && size>1,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()"); 662 663 /* flag CPU mask for C */ 664 #if defined(PETSC_HAVE_DEVICE) 665 if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU; 666 if (c->A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->A->offloadmask = PETSC_OFFLOAD_CPU; 667 if (c->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->B->offloadmask = PETSC_OFFLOAD_CPU; 668 #endif 669 apa_sparse = ptap->apa; 670 671 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 672 /*-----------------------------------------------------*/ 673 /* update numerical values of P_oth and P_loc */ 674 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 675 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 676 677 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 678 /*----------------------------------------------------------*/ 679 /* get data from symbolic products */ 680 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 681 pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a; 682 if (size >1) { 683 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 684 pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a; 685 } else { 686 p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL; 687 } 688 689 /* trigger copy to CPU */ 690 ierr = MatSeqAIJGetArrayRead(a->A,&dummy);CHKERRQ(ierr); 691 ierr = MatSeqAIJRestoreArrayRead(a->A,&dummy);CHKERRQ(ierr); 692 ierr = MatSeqAIJGetArrayRead(a->B,&dummy);CHKERRQ(ierr); 693 ierr = MatSeqAIJRestoreArrayRead(a->B,&dummy);CHKERRQ(ierr); 694 api = ptap->api; 695 apj = ptap->apj; 696 for (i=0; i<cm; i++) { 697 apJ = apj + api[i]; 698 699 /* diagonal portion of A */ 700 anz = adi[i+1] - adi[i]; 701 adj = ad->j + adi[i]; 702 ada = ad->a + adi[i]; 703 for (j=0; j<anz; j++) { 704 row = adj[j]; 705 pnz = pi_loc[row+1] - pi_loc[row]; 706 pj = pj_loc + pi_loc[row]; 707 pa = pa_loc + pi_loc[row]; 708 /* perform sparse axpy */ 709 valtmp = ada[j]; 710 nextp = 0; 711 for (k=0; nextp<pnz; k++) { 712 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 713 apa_sparse[k] += valtmp*pa[nextp++]; 714 } 715 } 716 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 717 } 718 719 /* off-diagonal portion of A */ 720 anz = aoi[i+1] - aoi[i]; 721 aoj = ao->j + aoi[i]; 722 aoa = ao->a + aoi[i]; 723 for (j=0; j<anz; j++) { 724 row = aoj[j]; 725 pnz = pi_oth[row+1] - pi_oth[row]; 726 pj = pj_oth + pi_oth[row]; 727 pa = pa_oth + pi_oth[row]; 728 /* perform sparse axpy */ 729 valtmp = aoa[j]; 730 nextp = 0; 731 for (k=0; nextp<pnz; k++) { 732 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 733 apa_sparse[k] += valtmp*pa[nextp++]; 734 } 735 } 736 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 737 } 738 739 /* set values in C */ 740 cdnz = cd->i[i+1] - cd->i[i]; 741 conz = co->i[i+1] - co->i[i]; 742 743 /* 1st off-diagonal part of C */ 744 ca = coa + co->i[i]; 745 k = 0; 746 for (k0=0; k0<conz; k0++) { 747 if (apJ[k] >= cstart) break; 748 ca[k0] = apa_sparse[k]; 749 apa_sparse[k] = 0.0; 750 k++; 751 } 752 753 /* diagonal part of C */ 754 ca = cda + cd->i[i]; 755 for (k1=0; k1<cdnz; k1++) { 756 ca[k1] = apa_sparse[k]; 757 apa_sparse[k] = 0.0; 758 k++; 759 } 760 761 /* 2nd off-diagonal part of C */ 762 ca = coa + co->i[i]; 763 for (; k0<conz; k0++) { 764 ca[k0] = apa_sparse[k]; 765 apa_sparse[k] = 0.0; 766 k++; 767 } 768 } 769 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 770 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 771 PetscFunctionReturn(0); 772 } 773 774 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */ 775 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat C) 776 { 777 PetscErrorCode ierr; 778 MPI_Comm comm; 779 PetscMPIInt size; 780 Mat_APMPI *ptap; 781 PetscFreeSpaceList free_space = NULL,current_space=NULL; 782 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 783 Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 784 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 785 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 786 PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=1; 787 PetscInt am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20; 788 PetscReal afill; 789 MatType mtype; 790 791 PetscFunctionBegin; 792 MatCheckProduct(C,4); 793 PetscCheckFalse(C->product->data,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 794 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 795 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 796 797 /* create struct Mat_APMPI and attached it to C later */ 798 ierr = PetscNew(&ptap);CHKERRQ(ierr); 799 800 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 801 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 802 803 /* get P_loc by taking all local rows of P */ 804 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 805 806 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 807 pi_loc = p_loc->i; pj_loc = p_loc->j; 808 if (size > 1) { 809 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 810 pi_oth = p_oth->i; pj_oth = p_oth->j; 811 } else { 812 p_oth = NULL; 813 pi_oth = NULL; pj_oth = NULL; 814 } 815 816 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 817 /*-------------------------------------------------------------------*/ 818 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 819 ptap->api = api; 820 api[0] = 0; 821 822 ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr); 823 824 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 825 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 826 current_space = free_space; 827 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 828 for (i=0; i<am; i++) { 829 /* diagonal portion of A */ 830 nzi = adi[i+1] - adi[i]; 831 for (j=0; j<nzi; j++) { 832 row = *adj++; 833 pnz = pi_loc[row+1] - pi_loc[row]; 834 Jptr = pj_loc + pi_loc[row]; 835 /* Expand list if it is not long enough */ 836 if (pnz+apnz_max > lsize) { 837 lsize = pnz+apnz_max; 838 ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 839 } 840 /* add non-zero cols of P into the sorted linked list lnk */ 841 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 842 apnz = *lnk; /* The first element in the list is the number of items in the list */ 843 api[i+1] = api[i] + apnz; 844 if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */ 845 } 846 /* off-diagonal portion of A */ 847 nzi = aoi[i+1] - aoi[i]; 848 for (j=0; j<nzi; j++) { 849 row = *aoj++; 850 pnz = pi_oth[row+1] - pi_oth[row]; 851 Jptr = pj_oth + pi_oth[row]; 852 /* Expand list if it is not long enough */ 853 if (pnz+apnz_max > lsize) { 854 lsize = pnz + apnz_max; 855 ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 856 } 857 /* add non-zero cols of P into the sorted linked list lnk */ 858 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 859 apnz = *lnk; /* The first element in the list is the number of items in the list */ 860 api[i+1] = api[i] + apnz; 861 if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */ 862 } 863 864 /* add missing diagonal entry */ 865 if (C->force_diagonals) { 866 j = i + rstart; /* column index */ 867 ierr = PetscLLCondensedAddSorted_Scalable(1,&j,lnk);CHKERRQ(ierr); 868 } 869 870 apnz = *lnk; 871 api[i+1] = api[i] + apnz; 872 if (apnz > apnz_max) apnz_max = apnz; 873 874 /* if free space is not available, double the total space in the list */ 875 if (current_space->local_remaining<apnz) { 876 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 877 nspacedouble++; 878 } 879 880 /* Copy data into free space, then initialize lnk */ 881 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 882 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 883 884 current_space->array += apnz; 885 current_space->local_used += apnz; 886 current_space->local_remaining -= apnz; 887 } 888 889 /* Allocate space for apj, initialize apj, and */ 890 /* destroy list of free space and other temporary array(s) */ 891 ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 892 apj = ptap->apj; 893 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 894 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 895 896 /* create and assemble symbolic parallel matrix C */ 897 /*----------------------------------------------------*/ 898 ierr = MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 899 ierr = MatSetBlockSizesFromMats(C,A,P);CHKERRQ(ierr); 900 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 901 ierr = MatSetType(C,mtype);CHKERRQ(ierr); 902 ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr); 903 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 904 905 /* malloc apa for assembly C */ 906 ierr = PetscCalloc1(apnz_max,&ptap->apa);CHKERRQ(ierr); 907 908 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);CHKERRQ(ierr); 909 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 910 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 911 ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 912 913 C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 914 C->ops->productnumeric = MatProductNumeric_AB; 915 916 /* attach the supporting struct to C for reuse */ 917 C->product->data = ptap; 918 C->product->destroy = MatDestroy_MPIAIJ_MatMatMult; 919 920 /* set MatInfo */ 921 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 922 if (afill < 1.0) afill = 1.0; 923 C->info.mallocs = nspacedouble; 924 C->info.fill_ratio_given = fill; 925 C->info.fill_ratio_needed = afill; 926 927 #if defined(PETSC_USE_INFO) 928 if (api[am]) { 929 ierr = PetscInfo(C,"Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 930 ierr = PetscInfo(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 931 } else { 932 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 933 } 934 #endif 935 PetscFunctionReturn(0); 936 } 937 938 /* This function is needed for the seqMPI matrix-matrix multiplication. */ 939 /* Three input arrays are merged to one output array. The size of the */ 940 /* output array is also output. Duplicate entries only show up once. */ 941 static void Merge3SortedArrays(PetscInt size1, PetscInt *in1, 942 PetscInt size2, PetscInt *in2, 943 PetscInt size3, PetscInt *in3, 944 PetscInt *size4, PetscInt *out) 945 { 946 int i = 0, j = 0, k = 0, l = 0; 947 948 /* Traverse all three arrays */ 949 while (i<size1 && j<size2 && k<size3) { 950 if (in1[i] < in2[j] && in1[i] < in3[k]) { 951 out[l++] = in1[i++]; 952 } 953 else if (in2[j] < in1[i] && in2[j] < in3[k]) { 954 out[l++] = in2[j++]; 955 } 956 else if (in3[k] < in1[i] && in3[k] < in2[j]) { 957 out[l++] = in3[k++]; 958 } 959 else if (in1[i] == in2[j] && in1[i] < in3[k]) { 960 out[l++] = in1[i]; 961 i++, j++; 962 } 963 else if (in1[i] == in3[k] && in1[i] < in2[j]) { 964 out[l++] = in1[i]; 965 i++, k++; 966 } 967 else if (in3[k] == in2[j] && in2[j] < in1[i]) { 968 out[l++] = in2[j]; 969 k++, j++; 970 } 971 else if (in1[i] == in2[j] && in1[i] == in3[k]) { 972 out[l++] = in1[i]; 973 i++, j++, k++; 974 } 975 } 976 977 /* Traverse two remaining arrays */ 978 while (i<size1 && j<size2) { 979 if (in1[i] < in2[j]) { 980 out[l++] = in1[i++]; 981 } 982 else if (in1[i] > in2[j]) { 983 out[l++] = in2[j++]; 984 } 985 else { 986 out[l++] = in1[i]; 987 i++, j++; 988 } 989 } 990 991 while (i<size1 && k<size3) { 992 if (in1[i] < in3[k]) { 993 out[l++] = in1[i++]; 994 } 995 else if (in1[i] > in3[k]) { 996 out[l++] = in3[k++]; 997 } 998 else { 999 out[l++] = in1[i]; 1000 i++, k++; 1001 } 1002 } 1003 1004 while (k<size3 && j<size2) { 1005 if (in3[k] < in2[j]) { 1006 out[l++] = in3[k++]; 1007 } 1008 else if (in3[k] > in2[j]) { 1009 out[l++] = in2[j++]; 1010 } 1011 else { 1012 out[l++] = in3[k]; 1013 k++, j++; 1014 } 1015 } 1016 1017 /* Traverse one remaining array */ 1018 while (i<size1) out[l++] = in1[i++]; 1019 while (j<size2) out[l++] = in2[j++]; 1020 while (k<size3) out[l++] = in3[k++]; 1021 1022 *size4 = l; 1023 } 1024 1025 /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */ 1026 /* adds up the products. Two of these three multiplications are performed with existing (sequential) */ 1027 /* matrix-matrix multiplications. */ 1028 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat C) 1029 { 1030 PetscErrorCode ierr; 1031 MPI_Comm comm; 1032 PetscMPIInt size; 1033 Mat_APMPI *ptap; 1034 PetscFreeSpaceList free_space_diag=NULL, current_space=NULL; 1035 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data; 1036 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc; 1037 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data; 1038 Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq; 1039 PetscInt adponz, adpdnz; 1040 PetscInt *pi_loc,*dnz,*onz; 1041 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart; 1042 PetscInt *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi, 1043 *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj; 1044 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend; 1045 PetscBT lnkbt; 1046 PetscReal afill; 1047 PetscMPIInt rank; 1048 Mat adpd, aopoth; 1049 MatType mtype; 1050 const char *prefix; 1051 1052 PetscFunctionBegin; 1053 MatCheckProduct(C,4); 1054 PetscCheckFalse(C->product->data,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 1055 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1056 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 1057 ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 1058 ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);CHKERRQ(ierr); 1059 1060 /* create struct Mat_APMPI and attached it to C later */ 1061 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1062 1063 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1064 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1065 1066 /* get P_loc by taking all local rows of P */ 1067 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1068 1069 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1070 pi_loc = p_loc->i; 1071 1072 /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */ 1073 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 1074 ierr = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr); 1075 1076 adpoi[0] = 0; 1077 ptap->api = api; 1078 api[0] = 0; 1079 1080 /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */ 1081 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1082 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 1083 1084 /* Symbolic calc of A_loc_diag * P_loc_diag */ 1085 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1086 ierr = MatProductCreate(a->A,p->A,NULL,&adpd);CHKERRQ(ierr); 1087 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1088 ierr = MatSetOptionsPrefix(adpd,prefix);CHKERRQ(ierr); 1089 ierr = MatAppendOptionsPrefix(adpd,"inner_diag_");CHKERRQ(ierr); 1090 1091 ierr = MatProductSetType(adpd,MATPRODUCT_AB);CHKERRQ(ierr); 1092 ierr = MatProductSetAlgorithm(adpd,"sorted");CHKERRQ(ierr); 1093 ierr = MatProductSetFill(adpd,fill);CHKERRQ(ierr); 1094 ierr = MatProductSetFromOptions(adpd);CHKERRQ(ierr); 1095 1096 adpd->force_diagonals = C->force_diagonals; 1097 ierr = MatProductSymbolic(adpd);CHKERRQ(ierr); 1098 1099 adpd_seq = (Mat_SeqAIJ*)((adpd)->data); 1100 adpdi = adpd_seq->i; adpdj = adpd_seq->j; 1101 p_off = (Mat_SeqAIJ*)((p->B)->data); 1102 poff_i = p_off->i; poff_j = p_off->j; 1103 1104 /* j_temp stores indices of a result row before they are added to the linked list */ 1105 ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr); 1106 1107 /* Symbolic calc of the A_diag * p_loc_off */ 1108 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 1109 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr); 1110 current_space = free_space_diag; 1111 1112 for (i=0; i<am; i++) { 1113 /* A_diag * P_loc_off */ 1114 nzi = adi[i+1] - adi[i]; 1115 for (j=0; j<nzi; j++) { 1116 row = *adj++; 1117 pnz = poff_i[row+1] - poff_i[row]; 1118 Jptr = poff_j + poff_i[row]; 1119 for (i1 = 0; i1 < pnz; i1++) { 1120 j_temp[i1] = p->garray[Jptr[i1]]; 1121 } 1122 /* add non-zero cols of P into the sorted linked list lnk */ 1123 ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr); 1124 } 1125 1126 adponz = lnk[0]; 1127 adpoi[i+1] = adpoi[i] + adponz; 1128 1129 /* if free space is not available, double the total space in the list */ 1130 if (current_space->local_remaining<adponz) { 1131 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1132 nspacedouble++; 1133 } 1134 1135 /* Copy data into free space, then initialize lnk */ 1136 ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1137 1138 current_space->array += adponz; 1139 current_space->local_used += adponz; 1140 current_space->local_remaining -= adponz; 1141 } 1142 1143 /* Symbolic calc of A_off * P_oth */ 1144 ierr = MatSetOptionsPrefix(a->B,prefix);CHKERRQ(ierr); 1145 ierr = MatAppendOptionsPrefix(a->B,"inner_offdiag_");CHKERRQ(ierr); 1146 ierr = MatCreate(PETSC_COMM_SELF,&aopoth);CHKERRQ(ierr); 1147 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, aopoth);CHKERRQ(ierr); 1148 aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data); 1149 aopothi = aopoth_seq->i; aopothj = aopoth_seq->j; 1150 1151 /* Allocate space for apj, adpj, aopj, ... */ 1152 /* destroy lists of free space and other temporary array(s) */ 1153 1154 ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr); 1155 ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr); 1156 1157 /* Copy from linked list to j-array */ 1158 ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr); 1159 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1160 1161 adpoJ = adpoj; 1162 adpdJ = adpdj; 1163 aopJ = aopothj; 1164 apj = ptap->apj; 1165 apJ = apj; /* still empty */ 1166 1167 /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */ 1168 /* A_diag * P_loc_diag to get A*P */ 1169 for (i = 0; i < am; i++) { 1170 aopnz = aopothi[i+1] - aopothi[i]; 1171 adponz = adpoi[i+1] - adpoi[i]; 1172 adpdnz = adpdi[i+1] - adpdi[i]; 1173 1174 /* Correct indices from A_diag*P_diag */ 1175 for (i1 = 0; i1 < adpdnz; i1++) { 1176 adpdJ[i1] += p_colstart; 1177 } 1178 /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */ 1179 Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ); 1180 ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);CHKERRQ(ierr); 1181 1182 aopJ += aopnz; 1183 adpoJ += adponz; 1184 adpdJ += adpdnz; 1185 apJ += apnz; 1186 api[i+1] = api[i] + apnz; 1187 } 1188 1189 /* malloc apa to store dense row A[i,:]*P */ 1190 ierr = PetscCalloc1(pN+2,&ptap->apa);CHKERRQ(ierr); 1191 1192 /* create and assemble symbolic parallel matrix C */ 1193 ierr = MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1194 ierr = MatSetBlockSizesFromMats(C,A,P);CHKERRQ(ierr); 1195 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1196 ierr = MatSetType(C,mtype);CHKERRQ(ierr); 1197 ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr); 1198 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1199 1200 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);CHKERRQ(ierr); 1201 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1202 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1203 ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1204 1205 C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1206 C->ops->productnumeric = MatProductNumeric_AB; 1207 1208 /* attach the supporting struct to C for reuse */ 1209 C->product->data = ptap; 1210 C->product->destroy = MatDestroy_MPIAIJ_MatMatMult; 1211 1212 /* set MatInfo */ 1213 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 1214 if (afill < 1.0) afill = 1.0; 1215 C->info.mallocs = nspacedouble; 1216 C->info.fill_ratio_given = fill; 1217 C->info.fill_ratio_needed = afill; 1218 1219 #if defined(PETSC_USE_INFO) 1220 if (api[am]) { 1221 ierr = PetscInfo(C,"Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1222 ierr = PetscInfo(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1223 } else { 1224 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 1225 } 1226 #endif 1227 1228 ierr = MatDestroy(&aopoth);CHKERRQ(ierr); 1229 ierr = MatDestroy(&adpd);CHKERRQ(ierr); 1230 ierr = PetscFree(j_temp);CHKERRQ(ierr); 1231 ierr = PetscFree(adpoj);CHKERRQ(ierr); 1232 ierr = PetscFree(adpoi);CHKERRQ(ierr); 1233 PetscFunctionReturn(0); 1234 } 1235 1236 /*-------------------------------------------------------------------------*/ 1237 /* This routine only works when scall=MAT_REUSE_MATRIX! */ 1238 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 1239 { 1240 PetscErrorCode ierr; 1241 Mat_APMPI *ptap; 1242 Mat Pt; 1243 1244 PetscFunctionBegin; 1245 MatCheckProduct(C,3); 1246 ptap = (Mat_APMPI*)C->product->data; 1247 PetscCheckFalse(!ptap,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 1248 PetscCheckFalse(!ptap->Pt,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()"); 1249 1250 Pt = ptap->Pt; 1251 ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 1252 ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(Pt,A,C);CHKERRQ(ierr); 1253 PetscFunctionReturn(0); 1254 } 1255 1256 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1257 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat C) 1258 { 1259 PetscErrorCode ierr; 1260 Mat_APMPI *ptap; 1261 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data; 1262 MPI_Comm comm; 1263 PetscMPIInt size,rank; 1264 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1265 PetscInt pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n; 1266 PetscInt *lnk,i,k,nsend,rstart; 1267 PetscBT lnkbt; 1268 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,nrecv; 1269 PETSC_UNUSED PetscMPIInt icompleted=0; 1270 PetscInt **buf_rj,**buf_ri,**buf_ri_k,row,ncols,*cols; 1271 PetscInt len,proc,*dnz,*onz,*owners,nzi; 1272 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1273 MPI_Request *swaits,*rwaits; 1274 MPI_Status *sstatus,rstatus; 1275 PetscLayout rowmap; 1276 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1277 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 1278 PetscInt *Jptr,*prmap=p->garray,con,j,Crmax; 1279 Mat_SeqAIJ *a_loc,*c_loc,*c_oth; 1280 PetscTable ta; 1281 MatType mtype; 1282 const char *prefix; 1283 1284 PetscFunctionBegin; 1285 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1286 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 1287 ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 1288 1289 /* create symbolic parallel matrix C */ 1290 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1291 ierr = MatSetType(C,mtype);CHKERRQ(ierr); 1292 1293 C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1294 1295 /* create struct Mat_APMPI and attached it to C later */ 1296 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1297 ptap->reuse = MAT_INITIAL_MATRIX; 1298 1299 /* (0) compute Rd = Pd^T, Ro = Po^T */ 1300 /* --------------------------------- */ 1301 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1302 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1303 1304 /* (1) compute symbolic A_loc */ 1305 /* ---------------------------*/ 1306 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1307 1308 /* (2-1) compute symbolic C_oth = Ro*A_loc */ 1309 /* ------------------------------------ */ 1310 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1311 ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 1312 ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 1313 ierr = MatCreate(PETSC_COMM_SELF,&ptap->C_oth);CHKERRQ(ierr); 1314 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,ptap->C_oth);CHKERRQ(ierr); 1315 1316 /* (3) send coj of C_oth to other processors */ 1317 /* ------------------------------------------ */ 1318 /* determine row ownership */ 1319 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1320 rowmap->n = pn; 1321 rowmap->bs = 1; 1322 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1323 owners = rowmap->range; 1324 1325 /* determine the number of messages to send, their lengths */ 1326 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 1327 ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr); 1328 ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr); 1329 1330 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 1331 coi = c_oth->i; coj = c_oth->j; 1332 con = ptap->C_oth->rmap->n; 1333 proc = 0; 1334 for (i=0; i<con; i++) { 1335 while (prmap[i] >= owners[proc+1]) proc++; 1336 len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */ 1337 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1338 } 1339 1340 len = 0; /* max length of buf_si[], see (4) */ 1341 owners_co[0] = 0; 1342 nsend = 0; 1343 for (proc=0; proc<size; proc++) { 1344 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1345 if (len_s[proc]) { 1346 nsend++; 1347 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1348 len += len_si[proc]; 1349 } 1350 } 1351 1352 /* determine the number and length of messages to receive for coi and coj */ 1353 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1354 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 1355 1356 /* post the Irecv and Isend of coj */ 1357 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1358 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1359 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 1360 for (proc=0, k=0; proc<size; proc++) { 1361 if (!len_s[proc]) continue; 1362 i = owners_co[proc]; 1363 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRMPI(ierr); 1364 k++; 1365 } 1366 1367 /* (2-2) compute symbolic C_loc = Rd*A_loc */ 1368 /* ---------------------------------------- */ 1369 ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 1370 ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 1371 ierr = MatCreate(PETSC_COMM_SELF,&ptap->C_loc);CHKERRQ(ierr); 1372 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,ptap->C_loc);CHKERRQ(ierr); 1373 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1374 1375 /* receives coj are complete */ 1376 for (i=0; i<nrecv; i++) { 1377 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 1378 } 1379 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1380 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);} 1381 1382 /* add received column indices into ta to update Crmax */ 1383 a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data; 1384 1385 /* create and initialize a linked list */ 1386 ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */ 1387 MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta); 1388 1389 for (k=0; k<nrecv; k++) {/* k-th received message */ 1390 Jptr = buf_rj[k]; 1391 for (j=0; j<len_r[k]; j++) { 1392 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1393 } 1394 } 1395 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 1396 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 1397 1398 /* (4) send and recv coi */ 1399 /*-----------------------*/ 1400 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1401 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1402 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1403 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1404 for (proc=0,k=0; proc<size; proc++) { 1405 if (!len_s[proc]) continue; 1406 /* form outgoing message for i-structure: 1407 buf_si[0]: nrows to be sent 1408 [1:nrows]: row index (global) 1409 [nrows+1:2*nrows+1]: i-structure index 1410 */ 1411 /*-------------------------------------------*/ 1412 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1413 buf_si_i = buf_si + nrows+1; 1414 buf_si[0] = nrows; 1415 buf_si_i[0] = 0; 1416 nrows = 0; 1417 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1418 nzi = coi[i+1] - coi[i]; 1419 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1420 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1421 nrows++; 1422 } 1423 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRMPI(ierr); 1424 k++; 1425 buf_si += len_si[proc]; 1426 } 1427 for (i=0; i<nrecv; i++) { 1428 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 1429 } 1430 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1431 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);} 1432 1433 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 1434 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1435 ierr = PetscFree(swaits);CHKERRQ(ierr); 1436 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1437 1438 /* (5) compute the local portion of C */ 1439 /* ------------------------------------------ */ 1440 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of C */ 1441 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 1442 current_space = free_space; 1443 1444 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1445 for (k=0; k<nrecv; k++) { 1446 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1447 nrows = *buf_ri_k[k]; 1448 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1449 nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recved i-structure */ 1450 } 1451 1452 ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr); 1453 ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1454 for (i=0; i<pn; i++) { /* for each local row of C */ 1455 /* add C_loc into C */ 1456 nzi = c_loc->i[i+1] - c_loc->i[i]; 1457 Jptr = c_loc->j + c_loc->i[i]; 1458 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1459 1460 /* add received col data into lnk */ 1461 for (k=0; k<nrecv; k++) { /* k-th received message */ 1462 if (i == *nextrow[k]) { /* i-th row */ 1463 nzi = *(nextci[k]+1) - *nextci[k]; 1464 Jptr = buf_rj[k] + *nextci[k]; 1465 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1466 nextrow[k]++; nextci[k]++; 1467 } 1468 } 1469 1470 /* add missing diagonal entry */ 1471 if (C->force_diagonals) { 1472 k = i + owners[rank]; /* column index */ 1473 ierr = PetscLLCondensedAddSorted(1,&k,lnk,lnkbt);CHKERRQ(ierr); 1474 } 1475 1476 nzi = lnk[0]; 1477 1478 /* copy data into free space, then initialize lnk */ 1479 ierr = PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1480 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 1481 } 1482 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1483 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1484 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 1485 1486 /* local sizes and preallocation */ 1487 ierr = MatSetSizes(C,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1488 if (P->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(C->rmap,P->cmap->bs);CHKERRQ(ierr);} 1489 if (A->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(C->cmap,A->cmap->bs);CHKERRQ(ierr);} 1490 ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr); 1491 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1492 1493 /* add C_loc and C_oth to C */ 1494 ierr = MatGetOwnershipRange(C,&rstart,NULL);CHKERRQ(ierr); 1495 for (i=0; i<pn; i++) { 1496 ncols = c_loc->i[i+1] - c_loc->i[i]; 1497 cols = c_loc->j + c_loc->i[i]; 1498 row = rstart + i; 1499 ierr = MatSetValues(C,1,(const PetscInt*)&row,ncols,(const PetscInt*)cols,NULL,INSERT_VALUES);CHKERRQ(ierr); 1500 1501 if (C->force_diagonals) { 1502 ierr = MatSetValues(C,1,(const PetscInt*)&row,1,(const PetscInt*)&row,NULL,INSERT_VALUES);CHKERRQ(ierr); 1503 } 1504 } 1505 for (i=0; i<con; i++) { 1506 ncols = c_oth->i[i+1] - c_oth->i[i]; 1507 cols = c_oth->j + c_oth->i[i]; 1508 row = prmap[i]; 1509 ierr = MatSetValues(C,1,(const PetscInt*)&row,ncols,(const PetscInt*)cols,NULL,INSERT_VALUES);CHKERRQ(ierr); 1510 } 1511 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1512 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1513 ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1514 1515 /* members in merge */ 1516 ierr = PetscFree(id_r);CHKERRQ(ierr); 1517 ierr = PetscFree(len_r);CHKERRQ(ierr); 1518 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1519 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1520 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1521 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1522 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 1523 1524 /* attach the supporting struct to C for reuse */ 1525 C->product->data = ptap; 1526 C->product->destroy = MatDestroy_MPIAIJ_PtAP; 1527 PetscFunctionReturn(0); 1528 } 1529 1530 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 1531 { 1532 PetscErrorCode ierr; 1533 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data; 1534 Mat_SeqAIJ *c_seq; 1535 Mat_APMPI *ptap; 1536 Mat A_loc,C_loc,C_oth; 1537 PetscInt i,rstart,rend,cm,ncols,row; 1538 const PetscInt *cols; 1539 const PetscScalar *vals; 1540 1541 PetscFunctionBegin; 1542 MatCheckProduct(C,3); 1543 ptap = (Mat_APMPI*)C->product->data; 1544 PetscCheckFalse(!ptap,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 1545 PetscCheckFalse(!ptap->A_loc,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()"); 1546 ierr = MatZeroEntries(C);CHKERRQ(ierr); 1547 1548 if (ptap->reuse == MAT_REUSE_MATRIX) { 1549 /* These matrices are obtained in MatTransposeMatMultSymbolic() */ 1550 /* 1) get R = Pd^T, Ro = Po^T */ 1551 /*----------------------------*/ 1552 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1553 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1554 1555 /* 2) compute numeric A_loc */ 1556 /*--------------------------*/ 1557 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1558 } 1559 1560 /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */ 1561 A_loc = ptap->A_loc; 1562 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr); 1563 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr); 1564 C_loc = ptap->C_loc; 1565 C_oth = ptap->C_oth; 1566 1567 /* add C_loc and C_oth to C */ 1568 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1569 1570 /* C_loc -> C */ 1571 cm = C_loc->rmap->N; 1572 c_seq = (Mat_SeqAIJ*)C_loc->data; 1573 cols = c_seq->j; 1574 vals = c_seq->a; 1575 for (i=0; i<cm; i++) { 1576 ncols = c_seq->i[i+1] - c_seq->i[i]; 1577 row = rstart + i; 1578 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1579 cols += ncols; vals += ncols; 1580 } 1581 1582 /* Co -> C, off-processor part */ 1583 cm = C_oth->rmap->N; 1584 c_seq = (Mat_SeqAIJ*)C_oth->data; 1585 cols = c_seq->j; 1586 vals = c_seq->a; 1587 for (i=0; i<cm; i++) { 1588 ncols = c_seq->i[i+1] - c_seq->i[i]; 1589 row = p->garray[i]; 1590 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1591 cols += ncols; vals += ncols; 1592 } 1593 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1594 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1595 ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1596 1597 ptap->reuse = MAT_REUSE_MATRIX; 1598 PetscFunctionReturn(0); 1599 } 1600 1601 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1602 { 1603 PetscErrorCode ierr; 1604 Mat_Merge_SeqsToMPI *merge; 1605 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data; 1606 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1607 Mat_APMPI *ptap; 1608 PetscInt *adj; 1609 PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1610 MatScalar *ada,*ca,valtmp; 1611 PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1612 MPI_Comm comm; 1613 PetscMPIInt size,rank,taga,*len_s; 1614 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1615 PetscInt **buf_ri,**buf_rj; 1616 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1617 MPI_Request *s_waits,*r_waits; 1618 MPI_Status *status; 1619 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1620 const PetscScalar *dummy; 1621 PetscInt *ai,*aj,*coi,*coj,*poJ,*pdJ; 1622 Mat A_loc; 1623 Mat_SeqAIJ *a_loc; 1624 1625 PetscFunctionBegin; 1626 MatCheckProduct(C,3); 1627 ptap = (Mat_APMPI*)C->product->data; 1628 PetscCheckFalse(!ptap,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 1629 PetscCheckFalse(!ptap->A_loc,PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()"); 1630 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1631 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 1632 ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 1633 1634 merge = ptap->merge; 1635 1636 /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1637 /*------------------------------------------*/ 1638 /* get data from symbolic products */ 1639 coi = merge->coi; coj = merge->coj; 1640 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1641 bi = merge->bi; bj = merge->bj; 1642 owners = merge->rowmap->range; 1643 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1644 1645 /* get A_loc by taking all local rows of A */ 1646 A_loc = ptap->A_loc; 1647 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1648 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1649 ai = a_loc->i; 1650 aj = a_loc->j; 1651 1652 /* trigger copy to CPU */ 1653 ierr = MatSeqAIJGetArrayRead(p->A,&dummy);CHKERRQ(ierr); 1654 ierr = MatSeqAIJRestoreArrayRead(p->A,&dummy);CHKERRQ(ierr); 1655 ierr = MatSeqAIJGetArrayRead(p->B,&dummy);CHKERRQ(ierr); 1656 ierr = MatSeqAIJRestoreArrayRead(p->B,&dummy);CHKERRQ(ierr); 1657 for (i=0; i<am; i++) { 1658 anz = ai[i+1] - ai[i]; 1659 adj = aj + ai[i]; 1660 ada = a_loc->a + ai[i]; 1661 1662 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1663 /*-------------------------------------------------------------*/ 1664 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1665 pnz = po->i[i+1] - po->i[i]; 1666 poJ = po->j + po->i[i]; 1667 pA = po->a + po->i[i]; 1668 for (j=0; j<pnz; j++) { 1669 row = poJ[j]; 1670 cj = coj + coi[row]; 1671 ca = coa + coi[row]; 1672 /* perform sparse axpy */ 1673 nexta = 0; 1674 valtmp = pA[j]; 1675 for (k=0; nexta<anz; k++) { 1676 if (cj[k] == adj[nexta]) { 1677 ca[k] += valtmp*ada[nexta]; 1678 nexta++; 1679 } 1680 } 1681 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1682 } 1683 1684 /* put the value into Cd (diagonal part) */ 1685 pnz = pd->i[i+1] - pd->i[i]; 1686 pdJ = pd->j + pd->i[i]; 1687 pA = pd->a + pd->i[i]; 1688 for (j=0; j<pnz; j++) { 1689 row = pdJ[j]; 1690 cj = bj + bi[row]; 1691 ca = ba + bi[row]; 1692 /* perform sparse axpy */ 1693 nexta = 0; 1694 valtmp = pA[j]; 1695 for (k=0; nexta<anz; k++) { 1696 if (cj[k] == adj[nexta]) { 1697 ca[k] += valtmp*ada[nexta]; 1698 nexta++; 1699 } 1700 } 1701 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1702 } 1703 } 1704 1705 /* 3) send and recv matrix values coa */ 1706 /*------------------------------------*/ 1707 buf_ri = merge->buf_ri; 1708 buf_rj = merge->buf_rj; 1709 len_s = merge->len_s; 1710 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1711 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1712 1713 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1714 for (proc=0,k=0; proc<size; proc++) { 1715 if (!len_s[proc]) continue; 1716 i = merge->owners_co[proc]; 1717 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRMPI(ierr); 1718 k++; 1719 } 1720 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRMPI(ierr);} 1721 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRMPI(ierr);} 1722 1723 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1724 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1725 ierr = PetscFree(coa);CHKERRQ(ierr); 1726 1727 /* 4) insert local Cseq and received values into Cmpi */ 1728 /*----------------------------------------------------*/ 1729 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1730 for (k=0; k<merge->nrecv; k++) { 1731 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1732 nrows = *(buf_ri_k[k]); 1733 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1734 nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recved i-structure */ 1735 } 1736 1737 for (i=0; i<cm; i++) { 1738 row = owners[rank] + i; /* global row index of C_seq */ 1739 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1740 ba_i = ba + bi[i]; 1741 bnz = bi[i+1] - bi[i]; 1742 /* add received vals into ba */ 1743 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1744 /* i-th row */ 1745 if (i == *nextrow[k]) { 1746 cnz = *(nextci[k]+1) - *nextci[k]; 1747 cj = buf_rj[k] + *(nextci[k]); 1748 ca = abuf_r[k] + *(nextci[k]); 1749 nextcj = 0; 1750 for (j=0; nextcj<cnz; j++) { 1751 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1752 ba_i[j] += ca[nextcj++]; 1753 } 1754 } 1755 nextrow[k]++; nextci[k]++; 1756 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1757 } 1758 } 1759 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1760 } 1761 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1762 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1763 1764 ierr = PetscFree(ba);CHKERRQ(ierr); 1765 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1766 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1767 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1768 PetscFunctionReturn(0); 1769 } 1770 1771 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat C) 1772 { 1773 PetscErrorCode ierr; 1774 Mat A_loc; 1775 Mat_APMPI *ptap; 1776 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1777 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data; 1778 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1779 PetscInt nnz; 1780 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1781 PetscInt am =A->rmap->n,pn=P->cmap->n; 1782 MPI_Comm comm; 1783 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1784 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1785 PetscInt len,proc,*dnz,*onz,*owners; 1786 PetscInt nzi,*bi,*bj; 1787 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1788 MPI_Request *swaits,*rwaits; 1789 MPI_Status *sstatus,rstatus; 1790 Mat_Merge_SeqsToMPI *merge; 1791 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1792 PetscReal afill =1.0,afill_tmp; 1793 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax; 1794 Mat_SeqAIJ *a_loc; 1795 PetscTable ta; 1796 MatType mtype; 1797 1798 PetscFunctionBegin; 1799 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1800 /* check if matrix local sizes are compatible */ 1801 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); 1802 1803 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 1804 ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 1805 1806 /* create struct Mat_APMPI and attached it to C later */ 1807 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1808 1809 /* get A_loc by taking all local rows of A */ 1810 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1811 1812 ptap->A_loc = A_loc; 1813 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1814 ai = a_loc->i; 1815 aj = a_loc->j; 1816 1817 /* determine symbolic Co=(p->B)^T*A - send to others */ 1818 /*----------------------------------------------------*/ 1819 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 1820 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 1821 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1822 >= (num of nonzero rows of C_seq) - pn */ 1823 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1824 coi[0] = 0; 1825 1826 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1827 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am])); 1828 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1829 current_space = free_space; 1830 1831 /* create and initialize a linked list */ 1832 ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr); 1833 MatRowMergeMax_SeqAIJ(a_loc,am,ta); 1834 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1835 1836 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 1837 1838 for (i=0; i<pon; i++) { 1839 pnz = poti[i+1] - poti[i]; 1840 ptJ = potj + poti[i]; 1841 for (j=0; j<pnz; j++) { 1842 row = ptJ[j]; /* row of A_loc == col of Pot */ 1843 anz = ai[row+1] - ai[row]; 1844 Jptr = aj + ai[row]; 1845 /* add non-zero cols of AP into the sorted linked list lnk */ 1846 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1847 } 1848 nnz = lnk[0]; 1849 1850 /* If free space is not available, double the total space in the list */ 1851 if (current_space->local_remaining<nnz) { 1852 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1853 nspacedouble++; 1854 } 1855 1856 /* Copy data into free space, and zero out denserows */ 1857 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1858 1859 current_space->array += nnz; 1860 current_space->local_used += nnz; 1861 current_space->local_remaining -= nnz; 1862 1863 coi[i+1] = coi[i] + nnz; 1864 } 1865 1866 ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1867 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1868 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */ 1869 1870 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1871 if (afill_tmp > afill) afill = afill_tmp; 1872 1873 /* send j-array (coj) of Co to other processors */ 1874 /*----------------------------------------------*/ 1875 /* determine row ownership */ 1876 ierr = PetscNew(&merge);CHKERRQ(ierr); 1877 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1878 1879 merge->rowmap->n = pn; 1880 merge->rowmap->bs = 1; 1881 1882 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1883 owners = merge->rowmap->range; 1884 1885 /* determine the number of messages to send, their lengths */ 1886 ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1887 ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 1888 1889 len_s = merge->len_s; 1890 merge->nsend = 0; 1891 1892 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1893 1894 proc = 0; 1895 for (i=0; i<pon; i++) { 1896 while (prmap[i] >= owners[proc+1]) proc++; 1897 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1898 len_s[proc] += coi[i+1] - coi[i]; 1899 } 1900 1901 len = 0; /* max length of buf_si[] */ 1902 owners_co[0] = 0; 1903 for (proc=0; proc<size; proc++) { 1904 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1905 if (len_si[proc]) { 1906 merge->nsend++; 1907 len_si[proc] = 2*(len_si[proc] + 1); 1908 len += len_si[proc]; 1909 } 1910 } 1911 1912 /* determine the number and length of messages to receive for coi and coj */ 1913 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1914 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1915 1916 /* post the Irecv and Isend of coj */ 1917 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1918 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1919 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1920 for (proc=0, k=0; proc<size; proc++) { 1921 if (!len_s[proc]) continue; 1922 i = owners_co[proc]; 1923 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRMPI(ierr); 1924 k++; 1925 } 1926 1927 /* receives and sends of coj are complete */ 1928 ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1929 for (i=0; i<merge->nrecv; i++) { 1930 PETSC_UNUSED PetscMPIInt icompleted; 1931 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 1932 } 1933 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1934 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRMPI(ierr);} 1935 1936 /* add received column indices into table to update Armax */ 1937 /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/tutorials/ex56.c! */ 1938 for (k=0; k<merge->nrecv; k++) {/* k-th received message */ 1939 Jptr = buf_rj[k]; 1940 for (j=0; j<merge->len_r[k]; j++) { 1941 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1942 } 1943 } 1944 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1945 1946 /* send and recv coi */ 1947 /*-------------------*/ 1948 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1949 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1950 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1951 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1952 for (proc=0,k=0; proc<size; proc++) { 1953 if (!len_s[proc]) continue; 1954 /* form outgoing message for i-structure: 1955 buf_si[0]: nrows to be sent 1956 [1:nrows]: row index (global) 1957 [nrows+1:2*nrows+1]: i-structure index 1958 */ 1959 /*-------------------------------------------*/ 1960 nrows = len_si[proc]/2 - 1; 1961 buf_si_i = buf_si + nrows+1; 1962 buf_si[0] = nrows; 1963 buf_si_i[0] = 0; 1964 nrows = 0; 1965 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1966 nzi = coi[i+1] - coi[i]; 1967 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1968 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1969 nrows++; 1970 } 1971 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRMPI(ierr); 1972 k++; 1973 buf_si += len_si[proc]; 1974 } 1975 i = merge->nrecv; 1976 while (i--) { 1977 PETSC_UNUSED PetscMPIInt icompleted; 1978 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 1979 } 1980 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1981 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRMPI(ierr);} 1982 ierr = PetscFree(len_si);CHKERRQ(ierr); 1983 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1984 ierr = PetscFree(swaits);CHKERRQ(ierr); 1985 ierr = PetscFree(sstatus);CHKERRQ(ierr); 1986 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1987 1988 /* compute the local portion of C (mpi mat) */ 1989 /*------------------------------------------*/ 1990 /* allocate bi array and free space for accumulating nonzero column info */ 1991 ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1992 bi[0] = 0; 1993 1994 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1995 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am]))); 1996 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1997 current_space = free_space; 1998 1999 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 2000 for (k=0; k<merge->nrecv; k++) { 2001 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 2002 nrows = *buf_ri_k[k]; 2003 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 2004 nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th received i-structure */ 2005 } 2006 2007 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 2008 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2009 rmax = 0; 2010 for (i=0; i<pn; i++) { 2011 /* add pdt[i,:]*AP into lnk */ 2012 pnz = pdti[i+1] - pdti[i]; 2013 ptJ = pdtj + pdti[i]; 2014 for (j=0; j<pnz; j++) { 2015 row = ptJ[j]; /* row of AP == col of Pt */ 2016 anz = ai[row+1] - ai[row]; 2017 Jptr = aj + ai[row]; 2018 /* add non-zero cols of AP into the sorted linked list lnk */ 2019 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2020 } 2021 2022 /* add received col data into lnk */ 2023 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 2024 if (i == *nextrow[k]) { /* i-th row */ 2025 nzi = *(nextci[k]+1) - *nextci[k]; 2026 Jptr = buf_rj[k] + *nextci[k]; 2027 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2028 nextrow[k]++; nextci[k]++; 2029 } 2030 } 2031 2032 /* add missing diagonal entry */ 2033 if (C->force_diagonals) { 2034 k = i + owners[rank]; /* column index */ 2035 ierr = PetscLLCondensedAddSorted_Scalable(1,&k,lnk);CHKERRQ(ierr); 2036 } 2037 2038 nnz = lnk[0]; 2039 2040 /* if free space is not available, make more free space */ 2041 if (current_space->local_remaining<nnz) { 2042 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 2043 nspacedouble++; 2044 } 2045 /* copy data into free space, then initialize lnk */ 2046 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2047 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2048 2049 current_space->array += nnz; 2050 current_space->local_used += nnz; 2051 current_space->local_remaining -= nnz; 2052 2053 bi[i+1] = bi[i] + nnz; 2054 if (nnz > rmax) rmax = nnz; 2055 } 2056 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2057 2058 ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 2059 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2060 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 2061 if (afill_tmp > afill) afill = afill_tmp; 2062 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 2063 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 2064 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 2065 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 2066 2067 /* create symbolic parallel matrix C - why cannot be assembled in Numeric part */ 2068 /*-------------------------------------------------------------------------------*/ 2069 ierr = MatSetSizes(C,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 2070 ierr = MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 2071 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 2072 ierr = MatSetType(C,mtype);CHKERRQ(ierr); 2073 ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr); 2074 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2075 ierr = MatSetBlockSize(C,1);CHKERRQ(ierr); 2076 ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 2077 for (i=0; i<pn; i++) { 2078 row = i + rstart; 2079 nnz = bi[i+1] - bi[i]; 2080 Jptr = bj + bi[i]; 2081 ierr = MatSetValues(C,1,&row,nnz,Jptr,NULL,INSERT_VALUES);CHKERRQ(ierr); 2082 } 2083 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2084 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2085 ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 2086 merge->bi = bi; 2087 merge->bj = bj; 2088 merge->coi = coi; 2089 merge->coj = coj; 2090 merge->buf_ri = buf_ri; 2091 merge->buf_rj = buf_rj; 2092 merge->owners_co = owners_co; 2093 2094 /* attach the supporting struct to C for reuse */ 2095 C->product->data = ptap; 2096 C->product->destroy = MatDestroy_MPIAIJ_PtAP; 2097 ptap->merge = merge; 2098 2099 C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 2100 2101 #if defined(PETSC_USE_INFO) 2102 if (bi[pn] != 0) { 2103 ierr = PetscInfo(C,"Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 2104 ierr = PetscInfo(C,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 2105 } else { 2106 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 2107 } 2108 #endif 2109 PetscFunctionReturn(0); 2110 } 2111 2112 /* ---------------------------------------------------------------- */ 2113 static PetscErrorCode MatProductSymbolic_AtB_MPIAIJ_MPIAIJ(Mat C) 2114 { 2115 PetscErrorCode ierr; 2116 Mat_Product *product = C->product; 2117 Mat A=product->A,B=product->B; 2118 PetscReal fill=product->fill; 2119 PetscBool flg; 2120 2121 PetscFunctionBegin; 2122 /* scalable */ 2123 ierr = PetscStrcmp(product->alg,"scalable",&flg);CHKERRQ(ierr); 2124 if (flg) { 2125 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 2126 goto next; 2127 } 2128 2129 /* nonscalable */ 2130 ierr = PetscStrcmp(product->alg,"nonscalable",&flg);CHKERRQ(ierr); 2131 if (flg) { 2132 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr); 2133 goto next; 2134 } 2135 2136 /* matmatmult */ 2137 ierr = PetscStrcmp(product->alg,"at*b",&flg);CHKERRQ(ierr); 2138 if (flg) { 2139 Mat At; 2140 Mat_APMPI *ptap; 2141 2142 ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr); 2143 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(At,B,fill,C);CHKERRQ(ierr); 2144 ptap = (Mat_APMPI*)C->product->data; 2145 if (ptap) { 2146 ptap->Pt = At; 2147 C->product->destroy = MatDestroy_MPIAIJ_PtAP; 2148 } 2149 C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 2150 goto next; 2151 } 2152 2153 /* backend general code */ 2154 ierr = PetscStrcmp(product->alg,"backend",&flg);CHKERRQ(ierr); 2155 if (flg) { 2156 ierr = MatProductSymbolic_MPIAIJBACKEND(C);CHKERRQ(ierr); 2157 PetscFunctionReturn(0); 2158 } 2159 2160 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProduct type is not supported"); 2161 2162 next: 2163 C->ops->productnumeric = MatProductNumeric_AtB; 2164 PetscFunctionReturn(0); 2165 } 2166 2167 /* ---------------------------------------------------------------- */ 2168 /* Set options for MatMatMultxxx_MPIAIJ_MPIAIJ */ 2169 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AB(Mat C) 2170 { 2171 PetscErrorCode ierr; 2172 Mat_Product *product = C->product; 2173 Mat A=product->A,B=product->B; 2174 #if defined(PETSC_HAVE_HYPRE) 2175 const char *algTypes[5] = {"scalable","nonscalable","seqmpi","backend","hypre"}; 2176 PetscInt nalg = 5; 2177 #else 2178 const char *algTypes[4] = {"scalable","nonscalable","seqmpi","backend",}; 2179 PetscInt nalg = 4; 2180 #endif 2181 PetscInt alg = 1; /* set nonscalable algorithm as default */ 2182 PetscBool flg; 2183 MPI_Comm comm; 2184 2185 PetscFunctionBegin; 2186 /* Check matrix local sizes */ 2187 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 2188 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); 2189 2190 /* Set "nonscalable" as default algorithm */ 2191 ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr); 2192 if (flg) { 2193 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2194 2195 /* Set "scalable" as default if BN and local nonzeros of A and B are large */ 2196 if (B->cmap->N > 100000) { /* may switch to scalable algorithm as default */ 2197 MatInfo Ainfo,Binfo; 2198 PetscInt nz_local; 2199 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 2200 2201 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 2202 ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr); 2203 nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated); 2204 2205 if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE; 2206 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRMPI(ierr); 2207 2208 if (alg_scalable) { 2209 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 2210 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2211 ierr = PetscInfo(B,"Use scalable algorithm, BN %" PetscInt_FMT ", fill*nz_allocated %g\n",B->cmap->N,(double)(product->fill*nz_local));CHKERRQ(ierr); 2212 } 2213 } 2214 } 2215 2216 /* Get runtime option */ 2217 if (product->api_user) { 2218 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");CHKERRQ(ierr); 2219 ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2220 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2221 } else { 2222 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");CHKERRQ(ierr); 2223 ierr = PetscOptionsEList("-mat_product_algorithm","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2224 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2225 } 2226 if (flg) { 2227 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2228 } 2229 2230 C->ops->productsymbolic = MatProductSymbolic_AB_MPIAIJ_MPIAIJ; 2231 PetscFunctionReturn(0); 2232 } 2233 2234 /* Set options for MatTransposeMatMultXXX_MPIAIJ_MPIAIJ */ 2235 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AtB(Mat C) 2236 { 2237 PetscErrorCode ierr; 2238 Mat_Product *product = C->product; 2239 Mat A=product->A,B=product->B; 2240 const char *algTypes[4] = {"scalable","nonscalable","at*b","backend"}; 2241 PetscInt nalg = 4; 2242 PetscInt alg = 1; /* set default algorithm */ 2243 PetscBool flg; 2244 MPI_Comm comm; 2245 2246 PetscFunctionBegin; 2247 /* Check matrix local sizes */ 2248 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 2249 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); 2250 2251 /* Set default algorithm */ 2252 ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr); 2253 if (flg) { 2254 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2255 } 2256 2257 /* Set "scalable" as default if BN and local nonzeros of A and B are large */ 2258 if (alg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */ 2259 MatInfo Ainfo,Binfo; 2260 PetscInt nz_local; 2261 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 2262 2263 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 2264 ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr); 2265 nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated); 2266 2267 if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE; 2268 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRMPI(ierr); 2269 2270 if (alg_scalable) { 2271 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 2272 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2273 ierr = PetscInfo(B,"Use scalable algorithm, BN %" PetscInt_FMT ", fill*nz_allocated %g\n",B->cmap->N,(double)(product->fill*nz_local));CHKERRQ(ierr); 2274 } 2275 } 2276 2277 /* Get runtime option */ 2278 if (product->api_user) { 2279 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr); 2280 ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2281 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2282 } else { 2283 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");CHKERRQ(ierr); 2284 ierr = PetscOptionsEList("-mat_product_algorithm","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2285 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2286 } 2287 if (flg) { 2288 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2289 } 2290 2291 C->ops->productsymbolic = MatProductSymbolic_AtB_MPIAIJ_MPIAIJ; 2292 PetscFunctionReturn(0); 2293 } 2294 2295 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_PtAP(Mat C) 2296 { 2297 PetscErrorCode ierr; 2298 Mat_Product *product = C->product; 2299 Mat A=product->A,P=product->B; 2300 MPI_Comm comm; 2301 PetscBool flg; 2302 PetscInt alg=1; /* set default algorithm */ 2303 #if !defined(PETSC_HAVE_HYPRE) 2304 const char *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","backend"}; 2305 PetscInt nalg=5; 2306 #else 2307 const char *algTypes[6] = {"scalable","nonscalable","allatonce","allatonce_merged","backend","hypre"}; 2308 PetscInt nalg=6; 2309 #endif 2310 PetscInt pN=P->cmap->N; 2311 2312 PetscFunctionBegin; 2313 /* Check matrix local sizes */ 2314 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 2315 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); 2316 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); 2317 2318 /* Set "nonscalable" as default algorithm */ 2319 ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr); 2320 if (flg) { 2321 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2322 2323 /* Set "scalable" as default if BN and local nonzeros of A and B are large */ 2324 if (pN > 100000) { 2325 MatInfo Ainfo,Pinfo; 2326 PetscInt nz_local; 2327 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 2328 2329 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 2330 ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 2331 nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 2332 2333 if (pN > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE; 2334 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRMPI(ierr); 2335 2336 if (alg_scalable) { 2337 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 2338 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2339 } 2340 } 2341 } 2342 2343 /* Get runtime option */ 2344 if (product->api_user) { 2345 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 2346 ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2347 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2348 } else { 2349 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");CHKERRQ(ierr); 2350 ierr = PetscOptionsEList("-mat_product_algorithm","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2351 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2352 } 2353 if (flg) { 2354 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2355 } 2356 2357 C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ; 2358 PetscFunctionReturn(0); 2359 } 2360 2361 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_RARt(Mat C) 2362 { 2363 Mat_Product *product = C->product; 2364 Mat A = product->A,R=product->B; 2365 2366 PetscFunctionBegin; 2367 /* Check matrix local sizes */ 2368 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); 2369 2370 C->ops->productsymbolic = MatProductSymbolic_RARt_MPIAIJ_MPIAIJ; 2371 PetscFunctionReturn(0); 2372 } 2373 2374 /* 2375 Set options for ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm 2376 */ 2377 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABC(Mat C) 2378 { 2379 PetscErrorCode ierr; 2380 Mat_Product *product = C->product; 2381 PetscBool flg = PETSC_FALSE; 2382 PetscInt alg = 1; /* default algorithm */ 2383 const char *algTypes[3] = {"scalable","nonscalable","seqmpi"}; 2384 PetscInt nalg = 3; 2385 2386 PetscFunctionBegin; 2387 /* Set default algorithm */ 2388 ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr); 2389 if (flg) { 2390 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2391 } 2392 2393 /* Get runtime option */ 2394 if (product->api_user) { 2395 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");CHKERRQ(ierr); 2396 ierr = PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2397 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2398 } else { 2399 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");CHKERRQ(ierr); 2400 ierr = PetscOptionsEList("-mat_product_algorithm","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2401 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2402 } 2403 if (flg) { 2404 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2405 } 2406 2407 C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ; 2408 C->ops->productsymbolic = MatProductSymbolic_ABC; 2409 PetscFunctionReturn(0); 2410 } 2411 2412 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat C) 2413 { 2414 PetscErrorCode ierr; 2415 Mat_Product *product = C->product; 2416 2417 PetscFunctionBegin; 2418 switch (product->type) { 2419 case MATPRODUCT_AB: 2420 ierr = MatProductSetFromOptions_MPIAIJ_AB(C);CHKERRQ(ierr); 2421 break; 2422 case MATPRODUCT_AtB: 2423 ierr = MatProductSetFromOptions_MPIAIJ_AtB(C);CHKERRQ(ierr); 2424 break; 2425 case MATPRODUCT_PtAP: 2426 ierr = MatProductSetFromOptions_MPIAIJ_PtAP(C);CHKERRQ(ierr); 2427 break; 2428 case MATPRODUCT_RARt: 2429 ierr = MatProductSetFromOptions_MPIAIJ_RARt(C);CHKERRQ(ierr); 2430 break; 2431 case MATPRODUCT_ABC: 2432 ierr = MatProductSetFromOptions_MPIAIJ_ABC(C);CHKERRQ(ierr); 2433 break; 2434 default: 2435 break; 2436 } 2437 PetscFunctionReturn(0); 2438 } 2439