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