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