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 #if defined(PETSC_HAVE_HYPRE) 13 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat); 14 #endif 15 16 PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C) 17 { 18 PetscErrorCode ierr; 19 Mat_Product *product = C->product; 20 Mat A=product->A,P=product->B; 21 MatProductAlgorithm alg=product->alg; 22 PetscReal fill=product->fill; 23 PetscBool flg; 24 Mat Pt; 25 26 PetscFunctionBegin; 27 /* "scalable" */ 28 ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr); 29 if (flg) { 30 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 31 C->ops->productnumeric = MatProductNumeric_PtAP; 32 PetscFunctionReturn(0); 33 } 34 35 /* "rap" */ 36 ierr = PetscStrcmp(alg,"rap",&flg);CHKERRQ(ierr); 37 if (flg) { 38 Mat_MatTransMatMult *atb; 39 40 ierr = PetscNew(&atb);CHKERRQ(ierr); 41 ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 42 ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt,A,P,fill,C);CHKERRQ(ierr); 43 44 atb->At = Pt; 45 atb->data = C->product->data; 46 atb->destroy = C->product->destroy; 47 C->product->data = atb; 48 C->product->destroy = MatDestroy_SeqAIJ_MatTransMatMult; 49 C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 50 C->ops->productnumeric = MatProductNumeric_PtAP; 51 PetscFunctionReturn(0); 52 } 53 54 /* hypre */ 55 #if defined(PETSC_HAVE_HYPRE) 56 ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr); 57 if (flg) { 58 ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 59 PetscFunctionReturn(0); 60 } 61 #endif 62 63 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType is not supported"); 64 PetscFunctionReturn(0); 65 } 66 67 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat C) 68 { 69 PetscErrorCode ierr; 70 PetscFreeSpaceList free_space=NULL,current_space=NULL; 71 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 72 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 73 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 74 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 75 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 76 MatScalar *ca; 77 PetscBT lnkbt; 78 PetscReal afill; 79 80 PetscFunctionBegin; 81 /* Get ij structure of P^T */ 82 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 83 ptJ = ptj; 84 85 /* Allocate ci array, arrays for fill computation and */ 86 /* free space for accumulating nonzero column info */ 87 ierr = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr); 88 ci[0] = 0; 89 90 ierr = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr); 91 ptasparserow = ptadenserow + an; 92 93 /* create and initialize a linked list */ 94 nlnk = pn+1; 95 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 96 97 /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 98 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr); 99 current_space = free_space; 100 101 /* Determine symbolic info for each row of C: */ 102 for (i=0; i<pn; i++) { 103 ptnzi = pti[i+1] - pti[i]; 104 ptanzi = 0; 105 /* Determine symbolic row of PtA: */ 106 for (j=0; j<ptnzi; j++) { 107 arow = *ptJ++; 108 anzj = ai[arow+1] - ai[arow]; 109 ajj = aj + ai[arow]; 110 for (k=0; k<anzj; k++) { 111 if (!ptadenserow[ajj[k]]) { 112 ptadenserow[ajj[k]] = -1; 113 ptasparserow[ptanzi++] = ajj[k]; 114 } 115 } 116 } 117 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 118 ptaj = ptasparserow; 119 cnzi = 0; 120 for (j=0; j<ptanzi; j++) { 121 prow = *ptaj++; 122 pnzj = pi[prow+1] - pi[prow]; 123 pjj = pj + pi[prow]; 124 /* add non-zero cols of P into the sorted linked list lnk */ 125 ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 126 cnzi += nlnk; 127 } 128 129 /* If free space is not available, make more free space */ 130 /* Double the amount of total space in the list */ 131 if (current_space->local_remaining<cnzi) { 132 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 133 nspacedouble++; 134 } 135 136 /* Copy data into free space, and zero out denserows */ 137 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 138 139 current_space->array += cnzi; 140 current_space->local_used += cnzi; 141 current_space->local_remaining -= cnzi; 142 143 for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 144 145 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 146 /* For now, we will recompute what is needed. */ 147 ci[i+1] = ci[i] + cnzi; 148 } 149 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 150 /* Allocate space for cj, initialize cj, and */ 151 /* destroy list of free space and other temporary array(s) */ 152 ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr); 153 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 154 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 155 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 156 157 ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr); 158 159 /* put together the new matrix */ 160 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 161 ierr = MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 162 163 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 164 /* Since these are PETSc arrays, change flags to free them as necessary. */ 165 c = (Mat_SeqAIJ*)((C)->data); 166 c->free_a = PETSC_TRUE; 167 c->free_ij = PETSC_TRUE; 168 c->nonew = 0; 169 170 C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 171 172 /* set MatInfo */ 173 afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5); 174 if (afill < 1.0) afill = 1.0; 175 c->maxnz = ci[pn]; 176 c->nz = ci[pn]; 177 C->info.mallocs = nspacedouble; 178 C->info.fill_ratio_given = fill; 179 C->info.fill_ratio_needed = afill; 180 181 /* Clean up. */ 182 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 183 #if defined(PETSC_USE_INFO) 184 if (ci[pn] != 0) { 185 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 186 ierr = PetscInfo1(C,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 187 } else { 188 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 189 } 190 #endif 191 PetscFunctionReturn(0); 192 } 193 194 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C) 195 { 196 PetscErrorCode ierr; 197 Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 198 Mat_SeqAIJ *p = (Mat_SeqAIJ*) P->data; 199 Mat_SeqAIJ *c = (Mat_SeqAIJ*) C->data; 200 PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 201 PetscInt *ci=c->i,*cj=c->j,*cjj; 202 PetscInt am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 203 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 204 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 205 206 PetscFunctionBegin; 207 /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 208 ierr = PetscCalloc2(cn,&apa,cn,&apjdense);CHKERRQ(ierr); 209 ierr = PetscMalloc1(cn,&apj);CHKERRQ(ierr); 210 211 /* Clear old values in C */ 212 ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr); 213 214 for (i=0; i<am; i++) { 215 /* Form sparse row of A*P */ 216 anzi = ai[i+1] - ai[i]; 217 apnzj = 0; 218 for (j=0; j<anzi; j++) { 219 prow = *aj++; 220 pnzj = pi[prow+1] - pi[prow]; 221 pjj = pj + pi[prow]; 222 paj = pa + pi[prow]; 223 for (k=0; k<pnzj; k++) { 224 if (!apjdense[pjj[k]]) { 225 apjdense[pjj[k]] = -1; 226 apj[apnzj++] = pjj[k]; 227 } 228 apa[pjj[k]] += (*aa)*paj[k]; 229 } 230 ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr); 231 aa++; 232 } 233 234 /* Sort the j index array for quick sparse axpy. */ 235 /* Note: a array does not need sorting as it is in dense storage locations. */ 236 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 237 238 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 239 pnzi = pi[i+1] - pi[i]; 240 for (j=0; j<pnzi; j++) { 241 nextap = 0; 242 crow = *pJ++; 243 cjj = cj + ci[crow]; 244 caj = ca + ci[crow]; 245 /* Perform sparse axpy operation. Note cjj includes apj. */ 246 for (k=0; nextap<apnzj; k++) { 247 if (PetscUnlikelyDebug(k >= ci[crow+1] - ci[crow])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 248 if (cjj[k]==apj[nextap]) { 249 caj[k] += (*pA)*apa[apj[nextap++]]; 250 } 251 } 252 ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 253 pA++; 254 } 255 256 /* Zero the current row info for A*P */ 257 for (j=0; j<apnzj; j++) { 258 apa[apj[j]] = 0.; 259 apjdense[apj[j]] = 0; 260 } 261 } 262 263 /* Assemble the final matrix and clean up */ 264 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 265 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 266 267 ierr = PetscFree2(apa,apjdense);CHKERRQ(ierr); 268 ierr = PetscFree(apj);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 273 { 274 PetscErrorCode ierr; 275 Mat_MatTransMatMult *atb; 276 277 PetscFunctionBegin; 278 MatCheckProduct(C,3); 279 atb = (Mat_MatTransMatMult*)C->product->data; 280 if (!atb) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing data structure"); 281 ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&atb->At);CHKERRQ(ierr); 282 if (!C->ops->matmultnumeric) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing numeric operation"); 283 /* when using rap, MatMatMatMultSymbolic used a different data */ 284 if (atb->data) C->product->data = atb->data; 285 ierr = (*C->ops->matmatmultnumeric)(atb->At,A,P,C);CHKERRQ(ierr); 286 C->product->data = atb; 287 PetscFunctionReturn(0); 288 } 289