1 2 /* 3 Defines projective product routines where A is a SeqAIJ matrix 4 C = P^T * A * P 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/utils/freespace.h> 9 #include <petscbt.h> 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ" 13 PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 14 { 15 PetscErrorCode ierr; 16 17 PetscFunctionBegin; 18 if (!P->ops->ptapsymbolic_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name); 19 ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr); 20 PetscFunctionReturn(0); 21 } 22 23 #undef __FUNCT__ 24 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ" 25 PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C) 26 { 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 if (!P->ops->ptapnumeric_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name); 31 ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 37 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 38 { 39 PetscErrorCode ierr; 40 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 41 Mat_PtAP *ptap = a->ptap; 42 43 PetscFunctionBegin; 44 /* free ptap, then A */ 45 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 46 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 47 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 48 ierr = (ptap->destroy)(A);CHKERRQ(ierr); 49 ierr = PetscFree(ptap);CHKERRQ(ierr); 50 PetscFunctionReturn(0); 51 } 52 53 #undef __FUNCT__ 54 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2" 55 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,PetscReal fill,Mat *C) 56 { 57 PetscErrorCode ierr; 58 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 59 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 60 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 61 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 62 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N; 63 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 64 MatScalar *ca; 65 PetscBT lnkbt; 66 67 PetscFunctionBegin; 68 /* Get ij structure of P^T */ 69 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 70 ptJ=ptj; 71 72 /* Allocate ci array, arrays for fill computation and */ 73 /* free space for accumulating nonzero column info */ 74 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 75 ci[0] = 0; 76 77 ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 78 ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 79 ptasparserow = ptadenserow + an; 80 81 /* create and initialize a linked list */ 82 nlnk = pn+1; 83 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 84 85 /* Set initial free space to be fill*nnz(A). */ 86 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 87 ierr = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space); 88 current_space = free_space; 89 90 /* Determine symbolic info for each row of C: */ 91 for (i=0;i<pn;i++) { 92 ptnzi = pti[i+1] - pti[i]; 93 ptanzi = 0; 94 /* Determine symbolic row of PtA: */ 95 for (j=0;j<ptnzi;j++) { 96 arow = *ptJ++; 97 anzj = ai[arow+1] - ai[arow]; 98 ajj = aj + ai[arow]; 99 for (k=0;k<anzj;k++) { 100 if (!ptadenserow[ajj[k]]) { 101 ptadenserow[ajj[k]] = -1; 102 ptasparserow[ptanzi++] = ajj[k]; 103 } 104 } 105 } 106 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 107 ptaj = ptasparserow; 108 cnzi = 0; 109 for (j=0;j<ptanzi;j++) { 110 prow = *ptaj++; 111 pnzj = pi[prow+1] - pi[prow]; 112 pjj = pj + pi[prow]; 113 /* add non-zero cols of P into the sorted linked list lnk */ 114 ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 115 cnzi += nlnk; 116 } 117 118 /* If free space is not available, make more free space */ 119 /* Double the amount of total space in the list */ 120 if (current_space->local_remaining<cnzi) { 121 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 122 nspacedouble++; 123 } 124 125 /* Copy data into free space, and zero out denserows */ 126 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 127 current_space->array += cnzi; 128 current_space->local_used += cnzi; 129 current_space->local_remaining -= cnzi; 130 131 for (j=0;j<ptanzi;j++) { 132 ptadenserow[ptasparserow[j]] = 0; 133 } 134 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 135 /* For now, we will recompute what is needed. */ 136 ci[i+1] = ci[i] + cnzi; 137 } 138 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 139 /* Allocate space for cj, initialize cj, and */ 140 /* destroy list of free space and other temporary array(s) */ 141 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 142 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 143 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 144 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 145 146 /* Allocate space for ca */ 147 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 148 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 149 150 /* put together the new matrix */ 151 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 152 153 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 154 /* Since these are PETSc arrays, change flags to free them as necessary. */ 155 c = (Mat_SeqAIJ *)((*C)->data); 156 c->free_a = PETSC_TRUE; 157 c->free_ij = PETSC_TRUE; 158 c->nonew = 0; 159 A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2; /* should use *C->ops until PtAP insterface is updated to double dispatch as MatMatMult() */ 160 161 /* Clean up. */ 162 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 163 #if defined(PETSC_USE_INFO) 164 if (ci[pn] != 0) { 165 PetscReal afill = ((PetscReal)ci[pn])/ai[am]; 166 if (afill < 1.0) afill = 1.0; 167 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 168 ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 169 } else { 170 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 171 } 172 #endif 173 PetscFunctionReturn(0); 174 } 175 176 #undef __FUNCT__ 177 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2" 178 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,Mat C) 179 { 180 PetscErrorCode ierr; 181 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 182 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 183 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 184 PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 185 PetscInt *ci=c->i,*cj=c->j,*cjj; 186 PetscInt am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 187 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 188 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 189 190 PetscFunctionBegin; 191 /* Allocate temporary array for storage of one row of A*P */ 192 ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr); 193 194 apjdense = (PetscInt *)(apa + cn); 195 apj = apjdense + cn; 196 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 197 198 /* Clear old values in C */ 199 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 200 201 for (i=0; i<am; i++) { 202 /* Form sparse row of A*P */ 203 anzi = ai[i+1] - ai[i]; 204 apnzj = 0; 205 for (j=0; j<anzi; j++) { 206 prow = *aj++; 207 pnzj = pi[prow+1] - pi[prow]; 208 pjj = pj + pi[prow]; 209 paj = pa + pi[prow]; 210 for (k=0;k<pnzj;k++) { 211 if (!apjdense[pjj[k]]) { 212 apjdense[pjj[k]] = -1; 213 apj[apnzj++] = pjj[k]; 214 } 215 apa[pjj[k]] += (*aa)*paj[k]; 216 } 217 ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr); 218 aa++; 219 } 220 221 /* Sort the j index array for quick sparse axpy. */ 222 /* Note: a array does not need sorting as it is in dense storage locations. */ 223 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 224 225 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 226 pnzi = pi[i+1] - pi[i]; 227 for (j=0; j<pnzi; j++) { 228 nextap = 0; 229 crow = *pJ++; 230 cjj = cj + ci[crow]; 231 caj = ca + ci[crow]; 232 /* Perform sparse axpy operation. Note cjj includes apj. */ 233 for (k=0;nextap<apnzj;k++) { 234 #if defined(PETSC_USE_DEBUG) 235 if (k >= ci[crow+1] - ci[crow]) { 236 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 237 } 238 #endif 239 if (cjj[k]==apj[nextap]) { 240 caj[k] += (*pA)*apa[apj[nextap++]]; 241 } 242 } 243 ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 244 pA++; 245 } 246 247 /* Zero the current row info for A*P */ 248 for (j=0; j<apnzj; j++) { 249 apa[apj[j]] = 0.; 250 apjdense[apj[j]] = 0; 251 } 252 } 253 254 /* Assemble the final matrix and clean up */ 255 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 256 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 257 ierr = PetscFree(apa);CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 #undef __FUNCT__ 262 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 263 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 264 { 265 PetscErrorCode ierr; 266 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c,*ap; 267 PetscInt *ai=a->i,*pi=p->i,*api,*apj; 268 PetscInt *ci,*cj; 269 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 270 MatScalar *ca; 271 Mat_PtAP *ptap; 272 PetscInt sparse_axpy=0; 273 Mat Pt,AP; 274 275 PetscFunctionBegin; 276 /* flag 'sparse_axpy' determines which implementations to be used: 277 0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; (default) 278 1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops 279 (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz; 280 2: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */ 281 ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr); 282 if (sparse_axpy == 2){ 283 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(A,P,fill,C);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 /* Get symbolic Pt = P^T */ 288 ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr); 289 290 /* Get symbolic AP = A*P */ 291 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 292 ap = (Mat_SeqAIJ*)AP->data; 293 api = ap->i; 294 apj = ap->j; 295 ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */ 296 297 /* Get C = Pt*AP */ 298 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 299 c = (Mat_SeqAIJ*)(*C)->data; 300 ci = c->i; 301 cj = c->j; 302 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 303 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 304 c->a = ca; 305 c->free_a = PETSC_TRUE; 306 307 /* Create a supporting struct for reuse by MatPtAPNumeric() */ 308 ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr); 309 c->ptap = ptap; 310 ptap->destroy = (*C)->ops->destroy; 311 (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 312 313 /* Allocate temporary array for storage of one row of A*P */ 314 ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr); 315 ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr); 316 if (sparse_axpy == 1){ 317 A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 318 } else { 319 A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 320 } 321 ptap->api = api; 322 ptap->apj = apj; 323 324 /* Clean up. */ 325 ierr = MatDestroy(&Pt);CHKERRQ(ierr); 326 ierr = MatDestroy(&AP);CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 330 /* #define PROFILE_MatPtAPNumeric */ 331 #undef __FUNCT__ 332 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 333 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 334 { 335 PetscErrorCode ierr; 336 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 337 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 338 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 339 const PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 340 const PetscScalar *aa=a->a,*pa=p->a,*pval; 341 const PetscInt *apj,*pcol,*cjj; 342 const PetscInt am=A->rmap->N,cm=C->rmap->N; 343 PetscInt i,j,k,anz,apnz,pnz,prow,crow,cnz; 344 PetscScalar *apa,*ca=c->a,*caj,pvalj; 345 Mat_PtAP *ptap = c->ptap; 346 #if defined(PROFILE_MatPtAPNumeric) 347 PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0; 348 PetscInt flops0=0,flops1=0; 349 #endif 350 351 PetscFunctionBegin; 352 /* Get temporary array for storage of one row of A*P */ 353 apa = ptap->apa; 354 355 /* Clear old values in C */ 356 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 357 358 for (i=0;i<am;i++) { 359 /* Form sparse row of AP[i,:] = A[i,:]*P */ 360 #if defined(PROFILE_MatPtAPNumeric) 361 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 362 #endif 363 anz = ai[i+1] - ai[i]; 364 apnz = 0; 365 for (j=0; j<anz; j++) { 366 prow = aj[j]; 367 pnz = pi[prow+1] - pi[prow]; 368 pcol = pj + pi[prow]; 369 pval = pa + pi[prow]; 370 for (k=0; k<pnz; k++) { 371 apa[pcol[k]] += aa[j]*pval[k]; 372 } 373 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 374 #if defined(PROFILE_MatPtAPNumeric) 375 flops0 += 2.0*pnz; 376 #endif 377 } 378 aj += anz; aa += anz; 379 #if defined(PROFILE_MatPtAPNumeric) 380 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 381 time_Cseq0 += tf - t0; 382 #endif 383 384 /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 385 #if defined(PROFILE_MatPtAPNumeric) 386 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 387 #endif 388 apj = ptap->apj + ptap->api[i]; 389 apnz = ptap->api[i+1] - ptap->api[i]; 390 pnz = pi[i+1] - pi[i]; 391 pcol = pj + pi[i]; 392 pval = pa + pi[i]; 393 394 /* Perform dense axpy */ 395 for (j=0; j<pnz; j++) { 396 crow = pcol[j]; 397 cjj = cj + ci[crow]; 398 caj = ca + ci[crow]; 399 pvalj = pval[j]; 400 cnz = ci[crow+1] - ci[crow]; 401 for (k=0; k<cnz; k++){ 402 caj[k] += pvalj*apa[cjj[k]]; 403 } 404 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 405 #if defined(PROFILE_MatPtAPNumeric) 406 flops1 += 2.0*cnz; 407 #endif 408 } 409 #if defined(PROFILE_MatPtAPNumeric) 410 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 411 time_Cseq1 += tf - t0; 412 #endif 413 414 /* Zero the current row info for A*P */ 415 for (j=0; j<apnz; j++) apa[apj[j]] = 0.0; 416 } 417 418 /* Assemble the final matrix and clean up */ 419 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 421 #if defined(PROFILE_MatPtAPNumeric) 422 printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1); 423 #endif 424 PetscFunctionReturn(0); 425 } 426 427 #undef __FUNCT__ 428 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 429 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C) 430 { 431 PetscErrorCode ierr; 432 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 433 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 434 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 435 const PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 436 const PetscScalar *aa=a->a,*pa=p->a,*pval; 437 const PetscInt *apj,*pcol,*cjj; 438 const PetscInt am=A->rmap->N,cm=C->rmap->N; 439 PetscInt i,j,k,anz,apnz,pnz,prow,crow,apcol,nextap; 440 PetscScalar *apa,*ca=c->a,*caj,pvalj; 441 Mat_PtAP *ptap = c->ptap; 442 #if defined(PROFILE_MatPtAPNumeric) 443 PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0; 444 PetscInt flops0=0,flops1=0; 445 #endif 446 447 PetscFunctionBegin; 448 /* Get temporary array for storage of one row of A*P */ 449 apa = ptap->apa; 450 451 /* Clear old values in C */ 452 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 453 454 for (i=0;i<am;i++) { 455 /* Form sparse row of AP[i,:] = A[i,:]*P */ 456 #if defined(PROFILE_MatPtAPNumeric) 457 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 458 #endif 459 anz = ai[i+1] - ai[i]; 460 apnz = 0; 461 for (j=0; j<anz; j++) { 462 prow = aj[j]; 463 pnz = pi[prow+1] - pi[prow]; 464 pcol = pj + pi[prow]; 465 pval = pa + pi[prow]; 466 for (k=0; k<pnz; k++) { 467 apa[pcol[k]] += aa[j]*pval[k]; 468 } 469 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 470 #if defined(PROFILE_MatPtAPNumeric) 471 flops0 += 2.0*pnz; 472 #endif 473 } 474 aj += anz; aa += anz; 475 #if defined(PROFILE_MatPtAPNumeric) 476 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 477 time_Cseq0 += tf - t0; 478 #endif 479 480 /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 481 #if defined(PROFILE_MatPtAPNumeric) 482 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 483 #endif 484 apj = ptap->apj + ptap->api[i]; 485 apnz = ptap->api[i+1] - ptap->api[i]; 486 pnz = pi[i+1] - pi[i]; 487 pcol = pj + pi[i]; 488 pval = pa + pi[i]; 489 490 /* Perform sparse axpy */ 491 for (j=0; j<pnz; j++) { 492 crow = pcol[j]; 493 cjj = cj + ci[crow]; 494 caj = ca + ci[crow]; 495 pvalj = pval[j]; 496 nextap = 1; 497 apcol = apj[nextap]; 498 for (k=0; nextap<apnz; k++) { 499 if (cjj[k] == apcol) { 500 caj[k] += pvalj*apa[apcol]; 501 apcol = apj[nextap++]; 502 } 503 } 504 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 505 #if defined(PROFILE_MatPtAPNumeric) 506 flops1 += 2.0*apnz; 507 #endif 508 } 509 #if defined(PROFILE_MatPtAPNumeric) 510 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 511 time_Cseq1 += tf - t0; 512 #endif 513 514 /* Zero the current row info for A*P */ 515 for (j=0; j<apnz; j++) { 516 apcol = apj[j]; 517 apa[apcol] = 0.; 518 } 519 } 520 521 /* Assemble the final matrix and clean up */ 522 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 523 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 524 #if defined(PROFILE_MatPtAPNumeric) 525 printf("MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1); 526 #endif 527 PetscFunctionReturn(0); 528 } 529