1 /* 2 Defines projective product routines where A is a AIJ 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 "src/mat/impls/aij/mpi/mpiaij.h" 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "MatPtAP" 12 /*@ 13 MatPtAP - Creates the matrix projection C = P^T * A * P 14 15 Collective on Mat 16 17 Input Parameters: 18 + A - the matrix 19 . P - the projection matrix 20 . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 21 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)) 22 23 Output Parameters: 24 . C - the product matrix 25 26 Notes: 27 C will be created and must be destroyed by the user with MatDestroy(). 28 29 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 30 which inherit from SeqAIJ. C will be of type MATSEQAIJ. 31 32 Level: intermediate 33 34 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 35 @*/ 36 PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) { 37 PetscErrorCode ierr; 38 PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *); 39 PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *); 40 41 PetscFunctionBegin; 42 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 43 PetscValidType(A,1); 44 MatPreallocated(A); 45 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 46 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 47 PetscValidHeaderSpecific(P,MAT_COOKIE,2); 48 PetscValidType(P,2); 49 MatPreallocated(P); 50 if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 51 if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 52 PetscValidPointer(C,3); 53 54 if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 55 56 if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 57 58 /* For now, we do not dispatch based on the type of A and P */ 59 /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 60 fA = A->ops->ptap; 61 if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name); 62 fP = P->ops->ptap; 63 if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name); 64 if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 65 66 ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 67 ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr); 68 ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 74 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 75 { 76 PetscErrorCode ierr; 77 Mat C_mpi,AP_seq,P_seq,P_subseq,*psubseq; 78 Mat_MPIAIJ *p = (Mat_MPIAIJ*)P->data; 79 Mat_MatMatMultMPI *mult; 80 int i,prow,prstart,prend,m=P->m,pncols; 81 const int *pcols; 82 const PetscScalar *pvals; 83 int rank; 84 85 PetscFunctionBegin; 86 ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 87 88 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,P,fill,&C_mpi);CHKERRQ(ierr); 89 mult = (Mat_MatMatMultMPI*)C_mpi->spptr; 90 P_seq = mult->bseq[0]; 91 AP_seq = mult->C_seq; 92 prstart = mult->brstart; 93 prend = prstart + m; 94 ierr = PetscPrintf(PETSC_COMM_SELF," [%d] prstart: %d, prend: %d, dim of P_seq: %d, %d\n",rank,prstart,prend,P_seq->m,P_seq->n); 95 96 /* get seq matrix P_subseq by taking local rows of P */ 97 IS isrow,iscol; 98 ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr); 99 ierr = ISCreateStride(PETSC_COMM_SELF,P_seq->n,0,1,&iscol);CHKERRQ(ierr); 100 ierr = MatGetSubMatrices(P_seq,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psubseq);CHKERRQ(ierr); 101 P_subseq = psubseq[0]; 102 ierr = PetscPrintf(PETSC_COMM_SELF," [%d] dim of P_subseq: %d, %d\n",rank,P_subseq->m,P_subseq->n); 103 104 /* Compute P_subseq^T*C_seq using outer product (P_loc^T)[:,i]*C_seq[i,:]. */ 105 for (i=0;i<m;i++) { 106 prow = prstart + i; 107 ierr = MatGetRow(P_seq,prow,&pncols,&pcols,&pvals);CHKERRQ(ierr); 108 ierr = MatRestoreRow(P_seq,prow,&pncols,&pcols,&pvals);CHKERRQ(ierr); 109 } 110 111 *C = C_mpi; /* to be removed! */ 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 117 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 118 { 119 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 120 PetscErrorCode ierr; 121 122 PetscFunctionBegin; 123 PetscFunctionReturn(0); 124 } 125 126 #undef __FUNCT__ 127 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 128 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 129 { 130 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 131 PetscErrorCode ierr; 132 133 PetscFunctionBegin; 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 139 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 140 { 141 PetscErrorCode ierr; 142 PetscFunctionBegin; 143 if (scall == MAT_INITIAL_MATRIX){ 144 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 145 } 146 ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 147 PetscFunctionReturn(0); 148 } 149 150 #undef __FUNCT__ 151 #define __FUNCT__ "MatPtAPSymbolic" 152 /* 153 MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 154 155 Collective on Mat 156 157 Input Parameters: 158 + A - the matrix 159 - P - the projection matrix 160 161 Output Parameters: 162 . C - the (i,j) structure of the product matrix 163 164 Notes: 165 C will be created and must be destroyed by the user with MatDestroy(). 166 167 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 168 which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 169 this (i,j) structure by calling MatPtAPNumeric(). 170 171 Level: intermediate 172 173 .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 174 */ 175 PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) { 176 PetscErrorCode ierr; 177 PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); 178 PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*); 179 180 PetscFunctionBegin; 181 182 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 183 PetscValidType(A,1); 184 MatPreallocated(A); 185 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 186 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 187 188 PetscValidHeaderSpecific(P,MAT_COOKIE,2); 189 PetscValidType(P,2); 190 MatPreallocated(P); 191 if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 192 if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 193 194 PetscValidPointer(C,3); 195 196 if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 197 if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 198 199 /* For now, we do not dispatch based on the type of A and P */ 200 /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 201 fA = A->ops->ptapsymbolic; 202 if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); 203 fP = P->ops->ptapsymbolic; 204 if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); 205 if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 206 207 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 208 ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); 209 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 210 211 PetscFunctionReturn(0); 212 } 213 214 typedef struct { 215 Mat symAP; 216 } Mat_PtAPstruct; 217 218 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 219 220 #undef __FUNCT__ 221 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 222 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 223 { 224 PetscErrorCode ierr; 225 Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; 226 227 PetscFunctionBegin; 228 ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); 229 ierr = PetscFree(ptap);CHKERRQ(ierr); 230 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 236 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) { 237 PetscErrorCode ierr; 238 int *pti,*ptj; 239 Mat Pt,AP; 240 Mat_PtAPstruct *ptap; 241 242 PetscFunctionBegin; 243 /* create symbolic Pt */ 244 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 245 ierr = MatCreateSeqAIJWithArrays(P->comm,P->N,P->M,pti,ptj,PETSC_NULL,&Pt);CHKERRQ(ierr); 246 247 /* get symbolic AP=A*P and C=Pt*AP */ 248 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 249 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 250 251 /* clean up */ 252 ierr = MatRestoreSymbolicTranspose_SeqAIJ(Pt,&pti,&ptj);CHKERRQ(ierr); 253 ierr = MatDestroy(Pt);CHKERRQ(ierr); 254 255 /* save symbolic AP - to be used by MatPtAPNumeric_SeqAIJ_SeqAIJ() */ 256 ierr = PetscNew(Mat_PtAPstruct,&ptap);CHKERRQ(ierr); 257 ptap->symAP = AP; 258 (*C)->spptr = (void*)ptap; 259 (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 260 261 PetscFunctionReturn(0); 262 } 263 264 #include "src/mat/impls/maij/maij.h" 265 EXTERN_C_BEGIN 266 #undef __FUNCT__ 267 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 268 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) { 269 /* This routine requires testing -- I don't think it works. */ 270 PetscErrorCode ierr; 271 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 272 Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 273 Mat P=pp->AIJ; 274 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 275 int *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 276 int *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 277 int an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 278 int i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 279 MatScalar *ca; 280 281 PetscFunctionBegin; 282 /* Start timer */ 283 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 284 285 /* Get ij structure of P^T */ 286 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 287 288 /* Allocate ci array, arrays for fill computation and */ 289 /* free space for accumulating nonzero column info */ 290 ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr); 291 ci[0] = 0; 292 293 ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr); 294 ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr); 295 ptasparserow = ptadenserow + an; 296 denserow = ptasparserow + an; 297 sparserow = denserow + pn; 298 299 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 300 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 301 ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 302 current_space = free_space; 303 304 /* Determine symbolic info for each row of C: */ 305 for (i=0;i<pn/ppdof;i++) { 306 ptnzi = pti[i+1] - pti[i]; 307 ptanzi = 0; 308 ptJ = ptj + pti[i]; 309 for (dof=0;dof<ppdof;dof++) { 310 /* Determine symbolic row of PtA: */ 311 for (j=0;j<ptnzi;j++) { 312 arow = ptJ[j] + dof; 313 anzj = ai[arow+1] - ai[arow]; 314 ajj = aj + ai[arow]; 315 for (k=0;k<anzj;k++) { 316 if (!ptadenserow[ajj[k]]) { 317 ptadenserow[ajj[k]] = -1; 318 ptasparserow[ptanzi++] = ajj[k]; 319 } 320 } 321 } 322 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 323 ptaj = ptasparserow; 324 cnzi = 0; 325 for (j=0;j<ptanzi;j++) { 326 pdof = *ptaj%dof; 327 prow = (*ptaj++)/dof; 328 pnzj = pi[prow+1] - pi[prow]; 329 pjj = pj + pi[prow]; 330 for (k=0;k<pnzj;k++) { 331 if (!denserow[pjj[k]+pdof]) { 332 denserow[pjj[k]+pdof] = -1; 333 sparserow[cnzi++] = pjj[k]+pdof; 334 } 335 } 336 } 337 338 /* sort sparserow */ 339 ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 340 341 /* If free space is not available, make more free space */ 342 /* Double the amount of total space in the list */ 343 if (current_space->local_remaining<cnzi) { 344 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 345 } 346 347 /* Copy data into free space, and zero out denserows */ 348 ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr); 349 current_space->array += cnzi; 350 current_space->local_used += cnzi; 351 current_space->local_remaining -= cnzi; 352 353 for (j=0;j<ptanzi;j++) { 354 ptadenserow[ptasparserow[j]] = 0; 355 } 356 for (j=0;j<cnzi;j++) { 357 denserow[sparserow[j]] = 0; 358 } 359 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 360 /* For now, we will recompute what is needed. */ 361 ci[i+1+dof] = ci[i+dof] + cnzi; 362 } 363 } 364 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 365 /* Allocate space for cj, initialize cj, and */ 366 /* destroy list of free space and other temporary array(s) */ 367 ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr); 368 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 369 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 370 371 /* Allocate space for ca */ 372 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 373 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 374 375 /* put together the new matrix */ 376 ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 377 378 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 379 /* Since these are PETSc arrays, change flags to free them as necessary. */ 380 c = (Mat_SeqAIJ *)((*C)->data); 381 c->freedata = PETSC_TRUE; 382 c->nonew = 0; 383 384 /* Clean up. */ 385 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 386 387 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 388 PetscFunctionReturn(0); 389 } 390 EXTERN_C_END 391 392 #undef __FUNCT__ 393 #define __FUNCT__ "MatPtAPNumeric" 394 /* 395 MatPtAPNumeric - Computes the matrix projection C = P^T * A * P 396 397 Collective on Mat 398 399 Input Parameters: 400 + A - the matrix 401 - P - the projection matrix 402 403 Output Parameters: 404 . C - the product matrix 405 406 Notes: 407 C must have been created by calling MatPtAPSymbolic and must be destroyed by 408 the user using MatDeatroy(). 409 410 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 411 which inherit from SeqAIJ. C will be of type MATSEQAIJ. 412 413 Level: intermediate 414 415 .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() 416 */ 417 PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) { 418 PetscErrorCode ierr; 419 PetscErrorCode (*fA)(Mat,Mat,Mat); 420 PetscErrorCode (*fP)(Mat,Mat,Mat); 421 422 PetscFunctionBegin; 423 424 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 425 PetscValidType(A,1); 426 MatPreallocated(A); 427 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 428 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 429 430 PetscValidHeaderSpecific(P,MAT_COOKIE,2); 431 PetscValidType(P,2); 432 MatPreallocated(P); 433 if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 434 if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 435 436 PetscValidHeaderSpecific(C,MAT_COOKIE,3); 437 PetscValidType(C,3); 438 MatPreallocated(C); 439 if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 440 if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 441 442 if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); 443 if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 444 if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 445 if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); 446 447 /* For now, we do not dispatch based on the type of A and P */ 448 /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 449 fA = A->ops->ptapnumeric; 450 if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); 451 fP = P->ops->ptapnumeric; 452 if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); 453 if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 454 455 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 456 ierr = (*fA)(A,P,C);CHKERRQ(ierr); 457 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 458 459 PetscFunctionReturn(0); 460 } 461 462 #undef __FUNCT__ 463 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 464 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 465 { 466 PetscErrorCode ierr; 467 int flops=0; 468 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 469 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 470 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 471 int *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 472 int *ci=c->i,*cj=c->j,*cjj; 473 int am=A->M,cn=C->N,cm=C->M; 474 int i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 475 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 476 Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)C->spptr; 477 Mat_SeqAIJ *ap = (Mat_SeqAIJ *)(ptap->symAP)->data; 478 int *api=ap->i,*apj=ap->j,apj_nextap; 479 480 PetscFunctionBegin; 481 /* Allocate temporary array for storage of one row of A*P */ 482 ierr = PetscMalloc(cn*sizeof(MatScalar),&apa);CHKERRQ(ierr); 483 ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr); 484 485 /* Clear old values in C */ 486 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 487 488 for (i=0;i<am;i++) { 489 /* Get sparse values of A*P[i,:] */ 490 anzi = ai[i+1] - ai[i]; 491 apnzj = 0; 492 for (j=0;j<anzi;j++) { 493 prow = *aj++; 494 pnzj = pi[prow+1] - pi[prow]; 495 pjj = pj + pi[prow]; 496 paj = pa + pi[prow]; 497 for (k=0;k<pnzj;k++) { 498 apa[pjj[k]] += (*aa)*paj[k]; 499 } 500 flops += 2*pnzj; 501 aa++; 502 } 503 504 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 505 apj = ap->j + api[i]; 506 apnzj = api[i+1] - api[i]; 507 pnzi = pi[i+1] - pi[i]; 508 for (j=0;j<pnzi;j++) { 509 nextap = 0; 510 crow = *pJ++; 511 cjj = cj + ci[crow]; 512 caj = ca + ci[crow]; 513 /* Perform sparse axpy operation. Note cjj includes apj. */ 514 for (k=0; nextap<apnzj; k++) { 515 apj_nextap = *(apj+nextap); 516 if (cjj[k]==apj_nextap) { 517 caj[k] += (*pA)*apa[apj_nextap]; 518 nextap++; 519 } 520 } 521 flops += 2*apnzj; 522 pA++; 523 } 524 525 /* Zero the current row values for A*P */ 526 for (j=0;j<apnzj;j++) apa[apj[j]] = 0.0; 527 } 528 529 /* Assemble the final matrix and clean up */ 530 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 531 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 532 ierr = PetscFree(apa);CHKERRQ(ierr); 533 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 534 535 PetscFunctionReturn(0); 536 } 537