1 /* 2 Defines projective product routines where A is a SeqAIJ matrix 3 C = P^T * A * P 4 */ 5 6 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 7 #include "src/mat/utils/freespace.h" 8 #include "petscbt.h" 9 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 = MatPtAPNumeric_SeqAIJ_SeqAIJ(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__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 31 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 32 { 33 PetscErrorCode ierr; 34 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 35 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 36 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 37 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; 38 PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 39 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 40 MatScalar *ca; 41 PetscBT lnkbt; 42 43 PetscFunctionBegin; 44 /* Get ij structure of P^T */ 45 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 46 ptJ=ptj; 47 48 /* Allocate ci array, arrays for fill computation and */ 49 /* free space for accumulating nonzero column info */ 50 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 51 ci[0] = 0; 52 53 ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 54 ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 55 ptasparserow = ptadenserow + an; 56 57 /* create and initialize a linked list */ 58 nlnk = pn+1; 59 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 60 61 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 62 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 63 ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 64 current_space = free_space; 65 66 /* Determine symbolic info for each row of C: */ 67 for (i=0;i<pn;i++) { 68 ptnzi = pti[i+1] - pti[i]; 69 ptanzi = 0; 70 /* Determine symbolic row of PtA: */ 71 for (j=0;j<ptnzi;j++) { 72 arow = *ptJ++; 73 anzj = ai[arow+1] - ai[arow]; 74 ajj = aj + ai[arow]; 75 for (k=0;k<anzj;k++) { 76 if (!ptadenserow[ajj[k]]) { 77 ptadenserow[ajj[k]] = -1; 78 ptasparserow[ptanzi++] = ajj[k]; 79 } 80 } 81 } 82 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 83 ptaj = ptasparserow; 84 cnzi = 0; 85 for (j=0;j<ptanzi;j++) { 86 prow = *ptaj++; 87 pnzj = pi[prow+1] - pi[prow]; 88 pjj = pj + pi[prow]; 89 /* add non-zero cols of P into the sorted linked list lnk */ 90 ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 91 cnzi += nlnk; 92 } 93 94 /* If free space is not available, make more free space */ 95 /* Double the amount of total space in the list */ 96 if (current_space->local_remaining<cnzi) { 97 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 98 } 99 100 /* Copy data into free space, and zero out denserows */ 101 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 102 current_space->array += cnzi; 103 current_space->local_used += cnzi; 104 current_space->local_remaining -= cnzi; 105 106 for (j=0;j<ptanzi;j++) { 107 ptadenserow[ptasparserow[j]] = 0; 108 } 109 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 110 /* For now, we will recompute what is needed. */ 111 ci[i+1] = ci[i] + cnzi; 112 } 113 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 114 /* Allocate space for cj, initialize cj, and */ 115 /* destroy list of free space and other temporary array(s) */ 116 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 117 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 118 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 119 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 120 121 /* Allocate space for ca */ 122 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 123 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 124 125 /* put together the new matrix */ 126 ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 127 128 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 129 /* Since these are PETSc arrays, change flags to free them as necessary. */ 130 c = (Mat_SeqAIJ *)((*C)->data); 131 c->freedata = PETSC_TRUE; 132 c->nonew = 0; 133 134 /* Clean up. */ 135 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 136 137 PetscFunctionReturn(0); 138 } 139 140 #include "src/mat/impls/maij/maij.h" 141 #undef __FUNCT__ 142 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 143 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) 144 { 145 /* This routine requires testing -- but it's getting better. */ 146 PetscErrorCode ierr; 147 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 148 Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 149 Mat P=pp->AIJ; 150 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 151 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 152 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*denserow,*sparserow,*ptaj; 153 PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof,cn; 154 PetscInt i,j,k,dof,pshift,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 155 MatScalar *ca; 156 157 PetscFunctionBegin; 158 /* Start timer */ 159 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 160 161 /* Get ij structure of P^T */ 162 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 163 164 cn = pn*ppdof; 165 /* Allocate ci array, arrays for fill computation and */ 166 /* free space for accumulating nonzero column info */ 167 ierr = PetscMalloc((cn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 168 ci[0] = 0; 169 170 /* Work arrays for rows of P^T*A */ 171 ierr = PetscMalloc((2*cn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 172 ierr = PetscMemzero(ptadenserow,(2*cn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 173 ptasparserow = ptadenserow + an; 174 denserow = ptasparserow + an; 175 sparserow = denserow + cn; 176 177 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 178 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 179 /* Note, aspect ratio of P is the same as the aspect ratio of SeqAIJ inside P */ 180 ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 181 current_space = free_space; 182 183 /* Determine symbolic info for each row of C: */ 184 for (i=0;i<pn;i++) { 185 ptnzi = pti[i+1] - pti[i]; 186 ptanzi = 0; 187 ptJ = ptj + pti[i]; 188 for (dof=0;dof<ppdof;dof++) { 189 /* Determine symbolic row of PtA: */ 190 for (j=0;j<ptnzi;j++) { 191 /* Expand ptJ[j] by block size and shift by dof to get the right row of A */ 192 arow = ptJ[j]*ppdof + dof; 193 /* Nonzeros of P^T*A will be in same locations as any element of A in that row */ 194 anzj = ai[arow+1] - ai[arow]; 195 ajj = aj + ai[arow]; 196 for (k=0;k<anzj;k++) { 197 if (!ptadenserow[ajj[k]]) { 198 ptadenserow[ajj[k]] = -1; 199 ptasparserow[ptanzi++] = ajj[k]; 200 } 201 } 202 } 203 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 204 ptaj = ptasparserow; 205 cnzi = 0; 206 for (j=0;j<ptanzi;j++) { 207 /* Get offset within block of P */ 208 pshift = *ptaj%ppdof; 209 /* Get block row of P */ 210 prow = (*ptaj++)/ppdof; 211 /* P has same number of nonzeros per row as the compressed form */ 212 pnzj = pi[prow+1] - pi[prow]; 213 pjj = pj + pi[prow]; 214 for (k=0;k<pnzj;k++) { 215 /* Locations in C are shifted by the offset within the block */ 216 if (!denserow[pjj[k]+pshift]) { 217 denserow[pjj[k]+pshift] = -1; 218 sparserow[cnzi++] = pjj[k]+pshift; 219 } 220 } 221 } 222 223 /* sort sparserow */ 224 ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 225 226 /* If free space is not available, make more free space */ 227 /* Double the amount of total space in the list */ 228 if (current_space->local_remaining<cnzi) { 229 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 230 } 231 232 /* Copy data into free space, and zero out denserows */ 233 ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr); 234 current_space->array += cnzi; 235 current_space->local_used += cnzi; 236 current_space->local_remaining -= cnzi; 237 238 for (j=0;j<ptanzi;j++) { 239 ptadenserow[ptasparserow[j]] = 0; 240 } 241 for (j=0;j<cnzi;j++) { 242 denserow[sparserow[j]] = 0; 243 } 244 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 245 /* For now, we will recompute what is needed. */ 246 ci[i*ppdof+1+dof] = ci[i*ppdof+dof] + cnzi; 247 } 248 } 249 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 250 /* Allocate space for cj, initialize cj, and */ 251 /* destroy list of free space and other temporary array(s) */ 252 ierr = PetscMalloc((ci[cn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 253 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 254 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 255 256 /* Allocate space for ca */ 257 ierr = PetscMalloc((ci[cn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 258 ierr = PetscMemzero(ca,(ci[cn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 259 260 /* put together the new matrix */ 261 ierr = MatCreateSeqAIJWithArrays(A->comm,cn,cn,ci,cj,ca,C);CHKERRQ(ierr); 262 263 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 264 /* Since these are PETSc arrays, change flags to free them as necessary. */ 265 c = (Mat_SeqAIJ *)((*C)->data); 266 c->freedata = PETSC_TRUE; 267 c->nonew = 0; 268 269 /* Clean up. */ 270 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 271 272 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 278 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 279 { 280 PetscErrorCode ierr; 281 PetscInt flops=0; 282 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 283 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 284 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 285 PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 286 PetscInt *ci=c->i,*cj=c->j,*cjj; 287 PetscInt am=A->M,cn=C->N,cm=C->M; 288 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 289 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 290 291 PetscFunctionBegin; 292 /* Allocate temporary array for storage of one row of A*P */ 293 ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 294 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 295 296 apj = (PetscInt *)(apa + cn); 297 apjdense = apj + cn; 298 299 /* Clear old values in C */ 300 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 301 302 for (i=0;i<am;i++) { 303 /* Form sparse row of A*P */ 304 anzi = ai[i+1] - ai[i]; 305 apnzj = 0; 306 for (j=0;j<anzi;j++) { 307 prow = *aj++; 308 pnzj = pi[prow+1] - pi[prow]; 309 pjj = pj + pi[prow]; 310 paj = pa + pi[prow]; 311 for (k=0;k<pnzj;k++) { 312 if (!apjdense[pjj[k]]) { 313 apjdense[pjj[k]] = -1; 314 apj[apnzj++] = pjj[k]; 315 } 316 apa[pjj[k]] += (*aa)*paj[k]; 317 } 318 flops += 2*pnzj; 319 aa++; 320 } 321 322 /* Sort the j index array for quick sparse axpy. */ 323 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 324 325 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 326 pnzi = pi[i+1] - pi[i]; 327 for (j=0;j<pnzi;j++) { 328 nextap = 0; 329 crow = *pJ++; 330 cjj = cj + ci[crow]; 331 caj = ca + ci[crow]; 332 /* Perform sparse axpy operation. Note cjj includes apj. */ 333 for (k=0;nextap<apnzj;k++) { 334 if (cjj[k]==apj[nextap]) { 335 caj[k] += (*pA)*apa[apj[nextap++]]; 336 } 337 } 338 flops += 2*apnzj; 339 pA++; 340 } 341 342 /* Zero the current row info for A*P */ 343 for (j=0;j<apnzj;j++) { 344 apa[apj[j]] = 0.; 345 apjdense[apj[j]] = 0; 346 } 347 } 348 349 /* Assemble the final matrix and clean up */ 350 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352 ierr = PetscFree(apa);CHKERRQ(ierr); 353 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 354 355 PetscFunctionReturn(0); 356 } 357