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 18 PetscFunctionBegin; 19 if (scall == MAT_INITIAL_MATRIX) { 20 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 21 } 22 ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 23 PetscFunctionReturn(0); 24 } 25 26 #undef __FUNCT__ 27 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 28 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 29 { 30 PetscErrorCode ierr; 31 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 32 Mat_PtAP *ptap = a->ptap; 33 34 PetscFunctionBegin; 35 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 36 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 37 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 38 ierr = (ptap->destroy)(A);CHKERRQ(ierr); 39 ierr = PetscFree(ptap);CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy" 45 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C) 46 { 47 PetscErrorCode ierr; 48 PetscFreeSpaceList free_space=NULL,current_space=NULL; 49 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 50 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 51 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 52 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 53 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 54 MatScalar *ca; 55 PetscBT lnkbt; 56 PetscReal afill; 57 58 PetscFunctionBegin; 59 /* Get ij structure of P^T */ 60 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 61 ptJ = ptj; 62 63 /* Allocate ci array, arrays for fill computation and */ 64 /* free space for accumulating nonzero column info */ 65 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 66 ci[0] = 0; 67 68 ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 69 ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 70 ptasparserow = ptadenserow + an; 71 72 /* create and initialize a linked list */ 73 nlnk = pn+1; 74 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 75 76 /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 77 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+pi[pm])),&free_space);CHKERRQ(ierr); 78 current_space = free_space; 79 80 /* Determine symbolic info for each row of C: */ 81 for (i=0; i<pn; i++) { 82 ptnzi = pti[i+1] - pti[i]; 83 ptanzi = 0; 84 /* Determine symbolic row of PtA: */ 85 for (j=0; j<ptnzi; j++) { 86 arow = *ptJ++; 87 anzj = ai[arow+1] - ai[arow]; 88 ajj = aj + ai[arow]; 89 for (k=0; k<anzj; k++) { 90 if (!ptadenserow[ajj[k]]) { 91 ptadenserow[ajj[k]] = -1; 92 ptasparserow[ptanzi++] = ajj[k]; 93 } 94 } 95 } 96 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 97 ptaj = ptasparserow; 98 cnzi = 0; 99 for (j=0; j<ptanzi; j++) { 100 prow = *ptaj++; 101 pnzj = pi[prow+1] - pi[prow]; 102 pjj = pj + pi[prow]; 103 /* add non-zero cols of P into the sorted linked list lnk */ 104 ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 105 cnzi += nlnk; 106 } 107 108 /* If free space is not available, make more free space */ 109 /* Double the amount of total space in the list */ 110 if (current_space->local_remaining<cnzi) { 111 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 112 nspacedouble++; 113 } 114 115 /* Copy data into free space, and zero out denserows */ 116 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 117 118 current_space->array += cnzi; 119 current_space->local_used += cnzi; 120 current_space->local_remaining -= cnzi; 121 122 for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 123 124 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 125 /* For now, we will recompute what is needed. */ 126 ci[i+1] = ci[i] + cnzi; 127 } 128 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 129 /* Allocate space for cj, initialize cj, and */ 130 /* destroy list of free space and other temporary array(s) */ 131 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 132 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 133 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 134 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 135 136 /* Allocate space for ca */ 137 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 138 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 139 140 /* put together the new matrix */ 141 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 142 143 (*C)->rmap->bs = P->cmap->bs; 144 (*C)->cmap->bs = P->cmap->bs; 145 146 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 147 /* Since these are PETSc arrays, change flags to free them as necessary. */ 148 c = (Mat_SeqAIJ*)((*C)->data); 149 c->free_a = PETSC_TRUE; 150 c->free_ij = PETSC_TRUE; 151 c->nonew = 0; 152 (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 153 154 /* set MatInfo */ 155 afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5); 156 if (afill < 1.0) afill = 1.0; 157 c->maxnz = ci[pn]; 158 c->nz = ci[pn]; 159 (*C)->info.mallocs = nspacedouble; 160 (*C)->info.fill_ratio_given = fill; 161 (*C)->info.fill_ratio_needed = afill; 162 163 /* Clean up. */ 164 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 165 #if defined(PETSC_USE_INFO) 166 if (ci[pn] != 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_SparseAxpy" 178 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(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 (cn: non-scalable) */ 192 ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr); 193 194 apjdense = (PetscInt*)(apa + cn); 195 apj = apjdense + cn; 196 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 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]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 236 #endif 237 if (cjj[k]==apj[nextap]) { 238 caj[k] += (*pA)*apa[apj[nextap++]]; 239 } 240 } 241 ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 242 pA++; 243 } 244 245 /* Zero the current row info for A*P */ 246 for (j=0; j<apnzj; j++) { 247 apa[apj[j]] = 0.; 248 apjdense[apj[j]] = 0; 249 } 250 } 251 252 /* Assemble the final matrix and clean up */ 253 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 254 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255 256 ierr = PetscFree(apa);CHKERRQ(ierr); 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 262 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 263 { 264 PetscErrorCode ierr; 265 Mat_SeqAIJ *ap,*c; 266 PetscInt *api,*apj,*ci,pn=P->cmap->N; 267 MatScalar *ca; 268 Mat_PtAP *ptap; 269 Mat Pt,AP; 270 PetscBool sparse_axpy=PETSC_TRUE; 271 272 PetscFunctionBegin; 273 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 274 /* flag 'sparse_axpy' determines which implementations to be used: 275 0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; 276 1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */ 277 ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,NULL);CHKERRQ(ierr); 278 ierr = PetscOptionsEnd();CHKERRQ(ierr); 279 if (sparse_axpy) { 280 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 284 /* Get symbolic Pt = P^T */ 285 ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr); 286 287 /* Get symbolic AP = A*P */ 288 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 289 290 ap = (Mat_SeqAIJ*)AP->data; 291 api = ap->i; 292 apj = ap->j; 293 ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */ 294 295 /* Get C = Pt*AP */ 296 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 297 298 c = (Mat_SeqAIJ*)(*C)->data; 299 ci = c->i; 300 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 301 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 302 c->a = ca; 303 c->free_a = PETSC_TRUE; 304 305 /* Create a supporting struct for reuse by MatPtAPNumeric() */ 306 ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr); 307 308 c->ptap = ptap; 309 ptap->destroy = (*C)->ops->destroy; 310 (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 311 312 /* Allocate temporary array for storage of one row of A*P */ 313 ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr); 314 ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));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 = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);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