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