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 #undef __FUNCT__ 11 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ" 12 PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 13 { 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 if (!P->ops->ptapsymbolic_seqaij) { 18 SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 19 } 20 ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr); 21 PetscFunctionReturn(0); 22 } 23 24 #undef __FUNCT__ 25 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ" 26 PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C) 27 { 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 if (!P->ops->ptapnumeric_seqaij) { 32 SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 33 } 34 ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr); 35 PetscFunctionReturn(0); 36 } 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 40 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 41 { 42 PetscErrorCode ierr; 43 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 44 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 45 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 46 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; 47 PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 48 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 49 MatScalar *ca; 50 PetscBT lnkbt; 51 52 PetscFunctionBegin; 53 /* Get ij structure of P^T */ 54 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 55 ptJ=ptj; 56 57 /* Allocate ci array, arrays for fill computation and */ 58 /* free space for accumulating nonzero column info */ 59 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 60 ci[0] = 0; 61 62 ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 63 ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 64 ptasparserow = ptadenserow + an; 65 66 /* create and initialize a linked list */ 67 nlnk = pn+1; 68 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 69 70 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 71 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 72 ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 73 current_space = free_space; 74 75 /* Determine symbolic info for each row of C: */ 76 for (i=0;i<pn;i++) { 77 ptnzi = pti[i+1] - pti[i]; 78 ptanzi = 0; 79 /* Determine symbolic row of PtA: */ 80 for (j=0;j<ptnzi;j++) { 81 arow = *ptJ++; 82 anzj = ai[arow+1] - ai[arow]; 83 ajj = aj + ai[arow]; 84 for (k=0;k<anzj;k++) { 85 if (!ptadenserow[ajj[k]]) { 86 ptadenserow[ajj[k]] = -1; 87 ptasparserow[ptanzi++] = ajj[k]; 88 } 89 } 90 } 91 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 92 ptaj = ptasparserow; 93 cnzi = 0; 94 for (j=0;j<ptanzi;j++) { 95 prow = *ptaj++; 96 pnzj = pi[prow+1] - pi[prow]; 97 pjj = pj + pi[prow]; 98 /* add non-zero cols of P into the sorted linked list lnk */ 99 ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 100 cnzi += nlnk; 101 } 102 103 /* If free space is not available, make more free space */ 104 /* Double the amount of total space in the list */ 105 if (current_space->local_remaining<cnzi) { 106 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 107 } 108 109 /* Copy data into free space, and zero out denserows */ 110 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 111 current_space->array += cnzi; 112 current_space->local_used += cnzi; 113 current_space->local_remaining -= cnzi; 114 115 for (j=0;j<ptanzi;j++) { 116 ptadenserow[ptasparserow[j]] = 0; 117 } 118 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 119 /* For now, we will recompute what is needed. */ 120 ci[i+1] = ci[i] + cnzi; 121 } 122 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 123 /* Allocate space for cj, initialize cj, and */ 124 /* destroy list of free space and other temporary array(s) */ 125 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 126 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 127 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 128 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 129 130 /* Allocate space for ca */ 131 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 132 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 133 134 /* put together the new matrix */ 135 ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 136 137 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 138 /* Since these are PETSc arrays, change flags to free them as necessary. */ 139 c = (Mat_SeqAIJ *)((*C)->data); 140 c->freedata = PETSC_TRUE; 141 c->nonew = 0; 142 143 /* Clean up. */ 144 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 145 146 PetscFunctionReturn(0); 147 } 148 149 #undef __FUNCT__ 150 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 151 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 152 { 153 PetscErrorCode ierr; 154 PetscInt flops=0; 155 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 156 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 157 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 158 PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 159 PetscInt *ci=c->i,*cj=c->j,*cjj; 160 PetscInt am=A->M,cn=C->N,cm=C->M; 161 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 162 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 163 164 PetscFunctionBegin; 165 /* Allocate temporary array for storage of one row of A*P */ 166 ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 167 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 168 169 apj = (PetscInt *)(apa + cn); 170 apjdense = apj + cn; 171 172 /* Clear old values in C */ 173 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 174 175 for (i=0;i<am;i++) { 176 /* Form sparse row of A*P */ 177 anzi = ai[i+1] - ai[i]; 178 apnzj = 0; 179 for (j=0;j<anzi;j++) { 180 prow = *aj++; 181 pnzj = pi[prow+1] - pi[prow]; 182 pjj = pj + pi[prow]; 183 paj = pa + pi[prow]; 184 for (k=0;k<pnzj;k++) { 185 if (!apjdense[pjj[k]]) { 186 apjdense[pjj[k]] = -1; 187 apj[apnzj++] = pjj[k]; 188 } 189 apa[pjj[k]] += (*aa)*paj[k]; 190 } 191 flops += 2*pnzj; 192 aa++; 193 } 194 195 /* Sort the j index array for quick sparse axpy. */ 196 /* Note: a array does not need sorting as it is in dense storage locations. */ 197 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 198 199 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 200 pnzi = pi[i+1] - pi[i]; 201 for (j=0;j<pnzi;j++) { 202 nextap = 0; 203 crow = *pJ++; 204 cjj = cj + ci[crow]; 205 caj = ca + ci[crow]; 206 /* Perform sparse axpy operation. Note cjj includes apj. */ 207 for (k=0;nextap<apnzj;k++) { 208 if (cjj[k]==apj[nextap]) { 209 caj[k] += (*pA)*apa[apj[nextap++]]; 210 } 211 } 212 flops += 2*apnzj; 213 pA++; 214 } 215 216 /* Zero the current row info for A*P */ 217 for (j=0;j<apnzj;j++) { 218 apa[apj[j]] = 0.; 219 apjdense[apj[j]] = 0; 220 } 221 } 222 223 /* Assemble the final matrix and clean up */ 224 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 226 ierr = PetscFree(apa);CHKERRQ(ierr); 227 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 228 229 PetscFunctionReturn(0); 230 } 231