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 = PetscLLAdd(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)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 193 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 194 195 apj = (PetscInt *)(apa + cn); 196 apjdense = apj + cn; 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; 267 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*api,*apj; 268 PetscInt *ci,*cj,ndouble_ap,ndouble_ptap; 269 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N; 270 MatScalar *ca; 271 Mat_PtAP *ptap; 272 PetscInt sparse_axpy=0; 273 274 PetscFunctionBegin; 275 /* flag 'sparse_axpy' determines which implementations to be used: 276 0: do dense axpy in MatPtAPNumeric() - fastest, but requires more storage; 277 1: do one sparse axpy in MatPtAPNumeric() 278 2: do two sparse axpy in MatPtAPNumeric() - slowest, but does not store structure of A*P */ 279 ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr); 280 if (sparse_axpy == 2){ 281 /* Implementation used in petsc-3.2: 282 doese not store structure of A*P, requires two sparse axpy in each step 283 of loop (i=0;i<am;i++). Slow, but use less memory */ 284 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(A,P,fill,C);CHKERRQ(ierr); 285 PetscFunctionReturn(0); 286 } 287 288 /* Get ij structure of Pt = P^T */ 289 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 290 ptJ=ptj; 291 292 /* Get structure of AP = A*P */ 293 ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,an,pn,pi,pj,fill,&api,&apj,&ndouble_ap);CHKERRQ(ierr); 294 295 /* Get structure of C = Pt*AP */ 296 ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(pn,pti,ptj,am,pn,api,apj,fill,&ci,&cj,&ndouble_ptap);CHKERRQ(ierr); 297 298 /* Allocate space for ca */ 299 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 300 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 301 302 /* put together the new matrix */ 303 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 304 305 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 306 /* Since these are PETSc arrays, change flags to free them as necessary. */ 307 c = (Mat_SeqAIJ *)(*C)->data; 308 c->free_a = PETSC_TRUE; 309 c->free_ij = PETSC_TRUE; 310 c->nonew = 0; 311 312 /* create a supporting struct for reuse by MatPtAPNumeric() */ 313 ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr); 314 c->ptap = ptap; 315 ptap->destroy = (*C)->ops->destroy; 316 (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 317 318 /* Allocate temporary array for storage of one row of A*P */ 319 if (sparse_axpy == 1){ 320 ierr = PetscMalloc((c->rmax+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr); 321 ierr = PetscMemzero(ptap->apa,(c->rmax+1)*sizeof(MatScalar));CHKERRQ(ierr); 322 A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 323 } else { 324 ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr); 325 ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(MatScalar));CHKERRQ(ierr); 326 A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 327 } 328 ptap->api = api; 329 ptap->apj = apj; 330 331 /* Clean up. */ 332 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 333 #if defined(PETSC_USE_INFO) 334 if (ci[pn] != 0) { 335 PetscReal apfill,ptapfill; 336 apfill = ((PetscReal)api[am])/(ai[am]+pi[an]); 337 if (apfill < 1.0) apfill = 1.0; 338 ierr = PetscInfo3((*C),"A*P: Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble_ap,fill,apfill);CHKERRQ(ierr); 339 ptapfill = ((PetscReal)ci[pn])/(pi[an]+api[am]); 340 if (ptapfill < 1.0) ptapfill = 1.0; 341 ierr = PetscInfo3((*C),"Pt*AP: Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble_ptap,fill,ptapfill);CHKERRQ(ierr); 342 343 ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",PetscMax(apfill,ptapfill));CHKERRQ(ierr); 344 ierr = PetscInfo4((*C),"nonzeros: A %D, P %D, A*P %D, C=PtAP %D\n",ai[am],pi[an],api[am],ci[pn]);CHKERRQ(ierr); 345 } else { 346 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 347 } 348 #endif 349 PetscFunctionReturn(0); 350 } 351 352 #undef __FUNCT__ 353 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 354 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 355 { 356 PetscErrorCode ierr; 357 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 358 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 359 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 360 PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 361 PetscScalar *aa=a->a,*pa=p->a; 362 PetscInt *apj,*pcol,*cjj,cnz; 363 PetscInt am=A->rmap->N,cm=C->rmap->N; 364 PetscInt i,j,k,anz,apnz,pnz,prow,crow,apcol; 365 PetscScalar *apa,*pval,*ca=c->a,*caj; 366 Mat_PtAP *ptap = c->ptap; 367 368 PetscFunctionBegin; 369 /* Get temporary array for storage of one row of A*P */ 370 apa = ptap->apa; 371 372 /* Clear old values in C */ 373 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 374 375 for (i=0;i<am;i++) { 376 /* Form sparse row of AP[i,:] = A[i,:]*P */ 377 anz = ai[i+1] - ai[i]; 378 apnz = 0; 379 for (j=0; j<anz; j++) { 380 prow = aj[j]; 381 pnz = pi[prow+1] - pi[prow]; 382 pcol = pj + pi[prow]; 383 pval = pa + pi[prow]; 384 for (k=0; k<pnz; k++) { 385 apa[pcol[k]] += aa[j]*pval[k]; 386 } 387 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 388 } 389 aj += anz; aa += anz; 390 391 /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 392 apj = ptap->apj + ptap->api[i]; 393 apnz = ptap->api[i+1] - ptap->api[i]; 394 pnz = pi[i+1] - pi[i]; 395 pcol = pj + pi[i]; 396 pval = pa + pi[i]; 397 398 /* Perform dense axpy */ 399 for (j=0; j<pnz; j++) { 400 crow = pcol[j]; 401 cjj = cj + ci[crow]; 402 caj = ca + ci[crow]; 403 cnz = ci[crow+1] - ci[crow]; 404 for (k=0; k<cnz; k++){ 405 caj[k] += pval[j]*apa[cjj[k]]; 406 } 407 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 408 } 409 410 /* Zero the current row info for A*P */ 411 for (j=0; j<apnz; j++) { 412 apcol = apj[j]; 413 apa[apcol] = 0.; 414 } 415 } 416 417 /* Assemble the final matrix and clean up */ 418 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 419 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420 PetscFunctionReturn(0); 421 } 422 423 #undef __FUNCT__ 424 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 425 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C) 426 { 427 PetscErrorCode ierr; 428 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 429 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 430 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 431 PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 432 PetscScalar *aa=a->a,*pa=p->a; 433 PetscInt *apj,*pcol,*cjj; 434 PetscInt am=A->rmap->N,cm=C->rmap->N; 435 PetscInt i,j,k,anz,apnz,pnz,prow,crow,apcol,nextap; 436 PetscScalar *apa,*pval,*ca=c->a,*caj; 437 Mat_PtAP *ptap = c->ptap; 438 439 PetscFunctionBegin; 440 /* Get temporary array for storage of one row of A*P */ 441 apa = ptap->apa; 442 443 /* Clear old values in C */ 444 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 445 446 for (i=0;i<am;i++) { 447 /* Form sparse row of AP[i,:] = A[i,:]*P */ 448 anz = ai[i+1] - ai[i]; 449 apnz = 0; 450 for (j=0; j<anz; j++) { 451 prow = aj[j]; 452 pnz = pi[prow+1] - pi[prow]; 453 pcol = pj + pi[prow]; 454 pval = pa + pi[prow]; 455 for (k=0; k<pnz; k++) { 456 apa[pcol[k]] += aa[j]*pval[k]; 457 } 458 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 459 } 460 aj += anz; aa += anz; 461 462 /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 463 apj = ptap->apj + ptap->api[i]; 464 apnz = ptap->api[i+1] - ptap->api[i]; 465 pnz = pi[i+1] - pi[i]; 466 pcol = pj + pi[i]; 467 pval = pa + pi[i]; 468 469 /* Perform sparse axpy */ 470 for (j=0; j<pnz; j++) { 471 crow = pcol[j]; 472 cjj = cj + ci[crow]; 473 caj = ca + ci[crow]; 474 nextap = 0; 475 apcol = apj[nextap]; 476 for (k=0; nextap<apnz; k++) { 477 if (cjj[k] == apcol) { 478 caj[k] += pval[j]*apa[apcol]; 479 apcol = apj[++nextap]; 480 } 481 } 482 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 483 } 484 485 /* Zero the current row info for A*P */ 486 for (j=0; j<apnz; j++) { 487 apcol = apj[j]; 488 apa[apcol] = 0.; 489 } 490 } 491 492 /* Assemble the final matrix and clean up */ 493 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 494 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 495 PetscFunctionReturn(0); 496 } 497