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