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 = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr); 87 ci[0] = 0; 88 89 ierr = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr); 90 ptasparserow = ptadenserow + an; 91 92 /* create and initialize a linked list */ 93 nlnk = pn+1; 94 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 95 96 /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 97 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+pi[pm])),&free_space);CHKERRQ(ierr); 98 current_space = free_space; 99 100 /* Determine symbolic info for each row of C: */ 101 for (i=0; i<pn; i++) { 102 ptnzi = pti[i+1] - pti[i]; 103 ptanzi = 0; 104 /* Determine symbolic row of PtA: */ 105 for (j=0; j<ptnzi; j++) { 106 arow = *ptJ++; 107 anzj = ai[arow+1] - ai[arow]; 108 ajj = aj + ai[arow]; 109 for (k=0; k<anzj; k++) { 110 if (!ptadenserow[ajj[k]]) { 111 ptadenserow[ajj[k]] = -1; 112 ptasparserow[ptanzi++] = ajj[k]; 113 } 114 } 115 } 116 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 117 ptaj = ptasparserow; 118 cnzi = 0; 119 for (j=0; j<ptanzi; j++) { 120 prow = *ptaj++; 121 pnzj = pi[prow+1] - pi[prow]; 122 pjj = pj + pi[prow]; 123 /* add non-zero cols of P into the sorted linked list lnk */ 124 ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 125 cnzi += nlnk; 126 } 127 128 /* If free space is not available, make more free space */ 129 /* Double the amount of total space in the list */ 130 if (current_space->local_remaining<cnzi) { 131 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 132 nspacedouble++; 133 } 134 135 /* Copy data into free space, and zero out denserows */ 136 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 137 138 current_space->array += cnzi; 139 current_space->local_used += cnzi; 140 current_space->local_remaining -= cnzi; 141 142 for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 143 144 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 145 /* For now, we will recompute what is needed. */ 146 ci[i+1] = ci[i] + cnzi; 147 } 148 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 149 /* Allocate space for cj, initialize cj, and */ 150 /* destroy list of free space and other temporary array(s) */ 151 ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr); 152 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 153 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 154 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 155 156 ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr); 157 158 /* put together the new matrix */ 159 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 160 ierr = MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 161 162 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 163 /* Since these are PETSc arrays, change flags to free them as necessary. */ 164 c = (Mat_SeqAIJ*)((*C)->data); 165 c->free_a = PETSC_TRUE; 166 c->free_ij = PETSC_TRUE; 167 c->nonew = 0; 168 (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 169 170 /* set MatInfo */ 171 afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5); 172 if (afill < 1.0) afill = 1.0; 173 c->maxnz = ci[pn]; 174 c->nz = ci[pn]; 175 (*C)->info.mallocs = nspacedouble; 176 (*C)->info.fill_ratio_given = fill; 177 (*C)->info.fill_ratio_needed = afill; 178 179 /* Clean up. */ 180 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 181 #if defined(PETSC_USE_INFO) 182 if (ci[pn] != 0) { 183 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 184 ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 185 } else { 186 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 187 } 188 #endif 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 194 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C) 195 { 196 PetscErrorCode ierr; 197 Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 198 Mat_SeqAIJ *p = (Mat_SeqAIJ*) P->data; 199 Mat_SeqAIJ *c = (Mat_SeqAIJ*) C->data; 200 PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 201 PetscInt *ci=c->i,*cj=c->j,*cjj; 202 PetscInt am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 203 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 204 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 205 206 PetscFunctionBegin; 207 /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 208 ierr = PetscMalloc3(cn,&apa,cn,&apjdense,c->rmax,&apj);CHKERRQ(ierr); 209 ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr); 210 ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr); 211 212 /* Clear old values in C */ 213 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 214 215 for (i=0; i<am; i++) { 216 /* Form sparse row of A*P */ 217 anzi = ai[i+1] - ai[i]; 218 apnzj = 0; 219 for (j=0; j<anzi; j++) { 220 prow = *aj++; 221 pnzj = pi[prow+1] - pi[prow]; 222 pjj = pj + pi[prow]; 223 paj = pa + pi[prow]; 224 for (k=0; k<pnzj; k++) { 225 if (!apjdense[pjj[k]]) { 226 apjdense[pjj[k]] = -1; 227 apj[apnzj++] = pjj[k]; 228 } 229 apa[pjj[k]] += (*aa)*paj[k]; 230 } 231 ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr); 232 aa++; 233 } 234 235 /* Sort the j index array for quick sparse axpy. */ 236 /* Note: a array does not need sorting as it is in dense storage locations. */ 237 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 238 239 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 240 pnzi = pi[i+1] - pi[i]; 241 for (j=0; j<pnzi; j++) { 242 nextap = 0; 243 crow = *pJ++; 244 cjj = cj + ci[crow]; 245 caj = ca + ci[crow]; 246 /* Perform sparse axpy operation. Note cjj includes apj. */ 247 for (k=0; nextap<apnzj; k++) { 248 #if defined(PETSC_USE_DEBUG) 249 if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 250 #endif 251 if (cjj[k]==apj[nextap]) { 252 caj[k] += (*pA)*apa[apj[nextap++]]; 253 } 254 } 255 ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 256 pA++; 257 } 258 259 /* Zero the current row info for A*P */ 260 for (j=0; j<apnzj; j++) { 261 apa[apj[j]] = 0.; 262 apjdense[apj[j]] = 0; 263 } 264 } 265 266 /* Assemble the final matrix and clean up */ 267 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 268 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 269 270 ierr = PetscFree3(apa,apjdense,apj);CHKERRQ(ierr); 271 PetscFunctionReturn(0); 272 } 273 274 #undef __FUNCT__ 275 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy" 276 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat A,Mat P,PetscReal fill,Mat *C) 277 { 278 PetscErrorCode ierr; 279 Mat_SeqAIJ *ap,*c; 280 PetscInt *api,*apj,*ci,pn=P->cmap->N; 281 MatScalar *ca; 282 Mat_PtAP *ptap; 283 Mat Pt,AP; 284 285 PetscFunctionBegin; 286 /* Get symbolic Pt = P^T */ 287 ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr); 288 289 /* Get symbolic AP = A*P */ 290 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 291 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 300 c = (Mat_SeqAIJ*)(*C)->data; 301 ci = c->i; 302 ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr); 303 c->a = ca; 304 c->free_a = PETSC_TRUE; 305 306 /* Create a supporting struct for reuse by MatPtAPNumeric() */ 307 ierr = PetscNew(&ptap);CHKERRQ(ierr); 308 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 = PetscCalloc1(pn+1,&ptap->apa);CHKERRQ(ierr); 315 316 (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 317 318 ptap->api = api; 319 ptap->apj = apj; 320 321 /* Clean up. */ 322 ierr = MatDestroy(&Pt);CHKERRQ(ierr); 323 ierr = MatDestroy(&AP);CHKERRQ(ierr); 324 #if defined(PETSC_USE_INFO) 325 ierr = PetscInfo1((*C),"given fill %g\n",(double)fill);CHKERRQ(ierr); 326 #endif 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 = PetscTime(&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 = PetscTime(&tf);CHKERRQ(ierr); 381 382 time_Cseq0 += tf - t0; 383 #endif 384 385 /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 386 #if defined(PROFILE_MatPtAPNumeric) 387 ierr = PetscTime(&t0);CHKERRQ(ierr); 388 #endif 389 apj = ptap->apj + ptap->api[i]; 390 apnz = ptap->api[i+1] - ptap->api[i]; 391 pnz = pi[i+1] - pi[i]; 392 pcol = pj + pi[i]; 393 pval = pa + pi[i]; 394 395 /* Perform dense axpy */ 396 for (j=0; j<pnz; j++) { 397 crow = pcol[j]; 398 cjj = cj + ci[crow]; 399 caj = ca + ci[crow]; 400 pvalj = pval[j]; 401 cnz = ci[crow+1] - ci[crow]; 402 for (k=0; k<cnz; k++) caj[k] += pvalj*apa[cjj[k]]; 403 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 404 #if defined(PROFILE_MatPtAPNumeric) 405 flops1 += 2.0*cnz; 406 #endif 407 } 408 #if defined(PROFILE_MatPtAPNumeric) 409 ierr = PetscTime(&tf);CHKERRQ(ierr); 410 time_Cseq1 += tf - t0; 411 #endif 412 413 /* Zero the current row info for A*P */ 414 for (j=0; j<apnz; j++) apa[apj[j]] = 0.0; 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 #if defined(PROFILE_MatPtAPNumeric) 421 printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1); 422 #endif 423 PetscFunctionReturn(0); 424 } 425