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 #include "petscbt.h" 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatPtAP" 13 /*@ 14 MatPtAP - Creates the matrix projection C = P^T * A * P 15 16 Collective on Mat 17 18 Input Parameters: 19 + A - the matrix 20 . P - the projection matrix 21 . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 22 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)) 23 24 Output Parameters: 25 . C - the product matrix 26 27 Notes: 28 C will be created and must be destroyed by the user with MatDestroy(). 29 30 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 31 which inherit from SeqAIJ. C will be of type MATSEQAIJ. 32 33 Level: intermediate 34 35 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 36 @*/ 37 PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 38 { 39 PetscErrorCode ierr; 40 PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *); 41 PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *); 42 43 PetscFunctionBegin; 44 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 45 PetscValidType(A,1); 46 MatPreallocated(A); 47 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 48 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 49 PetscValidHeaderSpecific(P,MAT_COOKIE,2); 50 PetscValidType(P,2); 51 MatPreallocated(P); 52 if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 53 if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 54 PetscValidPointer(C,3); 55 56 if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 57 58 if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 59 60 /* For now, we do not dispatch based on the type of A and P */ 61 /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 62 fA = A->ops->ptap; 63 if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name); 64 fP = P->ops->ptap; 65 if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name); 66 if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 67 68 ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 69 ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr); 70 ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 EXTERN PetscErrorCode MatPtAPReducedPt_MPIAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,PetscInt,PetscInt,Mat*); 75 EXTERN PetscErrorCode MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(Mat,Mat,PetscReal,PetscInt,PetscInt,Mat*); 76 EXTERN PetscErrorCode MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(Mat,Mat,PetscInt,PetscInt,Mat); 77 78 #undef __FUNCT__ 79 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 80 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 81 { 82 PetscErrorCode ierr; 83 84 PetscFunctionBegin; 85 if (scall == MAT_INITIAL_MATRIX){ 86 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 87 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ 88 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 89 } else if (scall == MAT_REUSE_MATRIX){ 90 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 91 ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 92 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 93 } else { 94 SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall); 95 } 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 101 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 102 { 103 PetscErrorCode ierr; 104 Mat P_seq,C_seq; 105 PetscInt prstart,prend,m=P->m; 106 107 PetscFunctionBegin; 108 /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */ 109 ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_seq);CHKERRQ(ierr); 110 111 /* compute C_seq = P_loc^T * A_loc * P_seq */ 112 prend = prstart + m; 113 ierr = MatPtAPReducedPt_MPIAIJ_SeqAIJ(A,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr); 114 115 ierr = MatDestroy(P_seq);CHKERRQ(ierr); 116 117 /* add C_seq into mpi C */ 118 ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 119 120 PetscFunctionReturn(0); 121 } 122 123 #undef __FUNCT__ 124 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 125 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 126 { 127 PetscErrorCode ierr; 128 Mat P_seq,C_seq; 129 PetscInt prstart,prend,m=P->m; 130 Mat_Merge_SeqsToMPI *merge; 131 PetscObjectContainer container; 132 133 PetscFunctionBegin; 134 ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr); 135 if (container) { 136 ierr = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr); 137 } else { 138 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container"); 139 } 140 141 /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */ 142 ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_seq);CHKERRQ(ierr); 143 144 /* compute C_seq = P_loc^T * A_loc * P_seq */ 145 prend = prstart + m; 146 C_seq = merge->C_seq; 147 ierr = MatPtAPReducedPt_MPIAIJ_SeqAIJ(A,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr); 148 ierr = MatDestroy(P_seq);CHKERRQ(ierr); 149 150 /* add C_seq into mpi C */ 151 ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 152 153 PetscFunctionReturn(0); 154 } 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 158 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 159 { 160 PetscErrorCode ierr; 161 162 PetscFunctionBegin; 163 if (scall == MAT_INITIAL_MATRIX){ 164 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 165 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 166 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 167 } 168 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 169 ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 170 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "MatPtAPSymbolic" 176 /* 177 MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 178 179 Collective on Mat 180 181 Input Parameters: 182 + A - the matrix 183 - P - the projection matrix 184 185 Output Parameters: 186 . C - the (i,j) structure of the product matrix 187 188 Notes: 189 C will be created and must be destroyed by the user with MatDestroy(). 190 191 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 192 which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 193 this (i,j) structure by calling MatPtAPNumeric(). 194 195 Level: intermediate 196 197 .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 198 */ 199 PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) 200 { 201 PetscErrorCode ierr; 202 PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); 203 PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*); 204 205 PetscFunctionBegin; 206 207 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 208 PetscValidType(A,1); 209 MatPreallocated(A); 210 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 211 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 212 213 PetscValidHeaderSpecific(P,MAT_COOKIE,2); 214 PetscValidType(P,2); 215 MatPreallocated(P); 216 if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 217 if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 218 219 PetscValidPointer(C,3); 220 221 if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 222 if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 223 224 /* For now, we do not dispatch based on the type of A and P */ 225 /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 226 fA = A->ops->ptapsymbolic; 227 if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); 228 fP = P->ops->ptapsymbolic; 229 if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); 230 if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 231 232 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 233 ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); 234 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 235 236 PetscFunctionReturn(0); 237 } 238 239 typedef struct { 240 Mat symAP; 241 } Mat_PtAPstruct; 242 243 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 244 245 #undef __FUNCT__ 246 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 247 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 248 { 249 PetscErrorCode ierr; 250 Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; 251 252 PetscFunctionBegin; 253 ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); 254 ierr = PetscFree(ptap);CHKERRQ(ierr); 255 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 256 PetscFunctionReturn(0); 257 } 258 259 #undef __FUNCT__ 260 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 261 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 262 { 263 PetscErrorCode ierr; 264 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 265 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 266 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 267 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; 268 PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 269 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 270 MatScalar *ca; 271 PetscBT lnkbt; 272 273 PetscFunctionBegin; 274 /* Get ij structure of P^T */ 275 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 276 ptJ=ptj; 277 278 /* Allocate ci array, arrays for fill computation and */ 279 /* free space for accumulating nonzero column info */ 280 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 281 ci[0] = 0; 282 283 ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 284 ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 285 ptasparserow = ptadenserow + an; 286 287 /* create and initialize a linked list */ 288 nlnk = pn+1; 289 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 290 291 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 292 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 293 ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 294 current_space = free_space; 295 296 /* Determine symbolic info for each row of C: */ 297 for (i=0;i<pn;i++) { 298 ptnzi = pti[i+1] - pti[i]; 299 ptanzi = 0; 300 /* Determine symbolic row of PtA: */ 301 for (j=0;j<ptnzi;j++) { 302 arow = *ptJ++; 303 anzj = ai[arow+1] - ai[arow]; 304 ajj = aj + ai[arow]; 305 for (k=0;k<anzj;k++) { 306 if (!ptadenserow[ajj[k]]) { 307 ptadenserow[ajj[k]] = -1; 308 ptasparserow[ptanzi++] = ajj[k]; 309 } 310 } 311 } 312 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 313 ptaj = ptasparserow; 314 cnzi = 0; 315 for (j=0;j<ptanzi;j++) { 316 prow = *ptaj++; 317 pnzj = pi[prow+1] - pi[prow]; 318 pjj = pj + pi[prow]; 319 /* add non-zero cols of P into the sorted linked list lnk */ 320 ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 321 cnzi += nlnk; 322 } 323 324 /* If free space is not available, make more free space */ 325 /* Double the amount of total space in the list */ 326 if (current_space->local_remaining<cnzi) { 327 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 328 } 329 330 /* Copy data into free space, and zero out denserows */ 331 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 332 current_space->array += cnzi; 333 current_space->local_used += cnzi; 334 current_space->local_remaining -= cnzi; 335 336 for (j=0;j<ptanzi;j++) { 337 ptadenserow[ptasparserow[j]] = 0; 338 } 339 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 340 /* For now, we will recompute what is needed. */ 341 ci[i+1] = ci[i] + cnzi; 342 } 343 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 344 /* Allocate space for cj, initialize cj, and */ 345 /* destroy list of free space and other temporary array(s) */ 346 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 347 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 348 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 349 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 350 351 /* Allocate space for ca */ 352 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 353 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 354 355 /* put together the new matrix */ 356 ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 357 358 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 359 /* Since these are PETSc arrays, change flags to free them as necessary. */ 360 c = (Mat_SeqAIJ *)((*C)->data); 361 c->freedata = PETSC_TRUE; 362 c->nonew = 0; 363 364 /* Clean up. */ 365 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 366 367 PetscFunctionReturn(0); 368 } 369 370 #include "src/mat/impls/maij/maij.h" 371 EXTERN_C_BEGIN 372 #undef __FUNCT__ 373 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 374 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) 375 { 376 /* This routine requires testing -- I don't think it works. */ 377 PetscErrorCode ierr; 378 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 379 Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 380 Mat P=pp->AIJ; 381 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 382 PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 383 PetscInt *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 384 PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 385 PetscInt i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 386 MatScalar *ca; 387 388 PetscFunctionBegin; 389 /* Start timer */ 390 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 391 392 /* Get ij structure of P^T */ 393 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 394 395 /* Allocate ci array, arrays for fill computation and */ 396 /* free space for accumulating nonzero column info */ 397 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 398 ci[0] = 0; 399 400 ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 401 ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 402 ptasparserow = ptadenserow + an; 403 denserow = ptasparserow + an; 404 sparserow = denserow + pn; 405 406 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 407 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 408 ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 409 current_space = free_space; 410 411 /* Determine symbolic info for each row of C: */ 412 for (i=0;i<pn/ppdof;i++) { 413 ptnzi = pti[i+1] - pti[i]; 414 ptanzi = 0; 415 ptJ = ptj + pti[i]; 416 for (dof=0;dof<ppdof;dof++) { 417 /* Determine symbolic row of PtA: */ 418 for (j=0;j<ptnzi;j++) { 419 arow = ptJ[j] + dof; 420 anzj = ai[arow+1] - ai[arow]; 421 ajj = aj + ai[arow]; 422 for (k=0;k<anzj;k++) { 423 if (!ptadenserow[ajj[k]]) { 424 ptadenserow[ajj[k]] = -1; 425 ptasparserow[ptanzi++] = ajj[k]; 426 } 427 } 428 } 429 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 430 ptaj = ptasparserow; 431 cnzi = 0; 432 for (j=0;j<ptanzi;j++) { 433 pdof = *ptaj%dof; 434 prow = (*ptaj++)/dof; 435 pnzj = pi[prow+1] - pi[prow]; 436 pjj = pj + pi[prow]; 437 for (k=0;k<pnzj;k++) { 438 if (!denserow[pjj[k]+pdof]) { 439 denserow[pjj[k]+pdof] = -1; 440 sparserow[cnzi++] = pjj[k]+pdof; 441 } 442 } 443 } 444 445 /* sort sparserow */ 446 ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 447 448 /* If free space is not available, make more free space */ 449 /* Double the amount of total space in the list */ 450 if (current_space->local_remaining<cnzi) { 451 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 452 } 453 454 /* Copy data into free space, and zero out denserows */ 455 ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr); 456 current_space->array += cnzi; 457 current_space->local_used += cnzi; 458 current_space->local_remaining -= cnzi; 459 460 for (j=0;j<ptanzi;j++) { 461 ptadenserow[ptasparserow[j]] = 0; 462 } 463 for (j=0;j<cnzi;j++) { 464 denserow[sparserow[j]] = 0; 465 } 466 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 467 /* For now, we will recompute what is needed. */ 468 ci[i+1+dof] = ci[i+dof] + cnzi; 469 } 470 } 471 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 472 /* Allocate space for cj, initialize cj, and */ 473 /* destroy list of free space and other temporary array(s) */ 474 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 475 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 476 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 477 478 /* Allocate space for ca */ 479 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 480 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 481 482 /* put together the new matrix */ 483 ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 484 485 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 486 /* Since these are PETSc arrays, change flags to free them as necessary. */ 487 c = (Mat_SeqAIJ *)((*C)->data); 488 c->freedata = PETSC_TRUE; 489 c->nonew = 0; 490 491 /* Clean up. */ 492 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 493 494 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 495 PetscFunctionReturn(0); 496 } 497 EXTERN_C_END 498 499 #undef __FUNCT__ 500 #define __FUNCT__ "MatPtAPNumeric" 501 /* 502 MatPtAPNumeric - Computes the matrix projection C = P^T * A * P 503 504 Collective on Mat 505 506 Input Parameters: 507 + A - the matrix 508 - P - the projection matrix 509 510 Output Parameters: 511 . C - the product matrix 512 513 Notes: 514 C must have been created by calling MatPtAPSymbolic and must be destroyed by 515 the user using MatDeatroy(). 516 517 This routine is currently only implemented for pairs of AIJ matrices and classes 518 which inherit from AIJ. C will be of type MATAIJ. 519 520 Level: intermediate 521 522 .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() 523 */ 524 PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) 525 { 526 PetscErrorCode ierr; 527 PetscErrorCode (*fA)(Mat,Mat,Mat); 528 PetscErrorCode (*fP)(Mat,Mat,Mat); 529 530 PetscFunctionBegin; 531 532 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 533 PetscValidType(A,1); 534 MatPreallocated(A); 535 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 536 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 537 538 PetscValidHeaderSpecific(P,MAT_COOKIE,2); 539 PetscValidType(P,2); 540 MatPreallocated(P); 541 if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 542 if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 543 544 PetscValidHeaderSpecific(C,MAT_COOKIE,3); 545 PetscValidType(C,3); 546 MatPreallocated(C); 547 if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 548 if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 549 550 if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); 551 if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 552 if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 553 if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); 554 555 /* For now, we do not dispatch based on the type of A and P */ 556 /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 557 fA = A->ops->ptapnumeric; 558 if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); 559 fP = P->ops->ptapnumeric; 560 if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); 561 if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 562 563 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 564 ierr = (*fA)(A,P,C);CHKERRQ(ierr); 565 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 566 567 PetscFunctionReturn(0); 568 } 569 570 #undef __FUNCT__ 571 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 572 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 573 { 574 PetscErrorCode ierr; 575 PetscInt flops=0; 576 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 577 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 578 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 579 PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 580 PetscInt *ci=c->i,*cj=c->j,*cjj; 581 PetscInt am=A->M,cn=C->N,cm=C->M; 582 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 583 MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 584 585 PetscFunctionBegin; 586 /* Allocate temporary array for storage of one row of A*P */ 587 ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 588 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 589 590 apj = (PetscInt *)(apa + cn); 591 apjdense = apj + cn; 592 593 /* Clear old values in C */ 594 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 595 596 for (i=0;i<am;i++) { 597 /* Form sparse row of A*P */ 598 anzi = ai[i+1] - ai[i]; 599 apnzj = 0; 600 for (j=0;j<anzi;j++) { 601 prow = *aj++; 602 pnzj = pi[prow+1] - pi[prow]; 603 pjj = pj + pi[prow]; 604 paj = pa + pi[prow]; 605 for (k=0;k<pnzj;k++) { 606 if (!apjdense[pjj[k]]) { 607 apjdense[pjj[k]] = -1; 608 apj[apnzj++] = pjj[k]; 609 } 610 apa[pjj[k]] += (*aa)*paj[k]; 611 } 612 flops += 2*pnzj; 613 aa++; 614 } 615 616 /* Sort the j index array for quick sparse axpy. */ 617 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 618 619 /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 620 pnzi = pi[i+1] - pi[i]; 621 for (j=0;j<pnzi;j++) { 622 nextap = 0; 623 crow = *pJ++; 624 cjj = cj + ci[crow]; 625 caj = ca + ci[crow]; 626 /* Perform sparse axpy operation. Note cjj includes apj. */ 627 for (k=0;nextap<apnzj;k++) { 628 if (cjj[k]==apj[nextap]) { 629 caj[k] += (*pA)*apa[apj[nextap++]]; 630 } 631 } 632 flops += 2*apnzj; 633 pA++; 634 } 635 636 /* Zero the current row info for A*P */ 637 for (j=0;j<apnzj;j++) { 638 apa[apj[j]] = 0.; 639 apjdense[apj[j]] = 0; 640 } 641 } 642 643 /* Assemble the final matrix and clean up */ 644 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 645 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 646 ierr = PetscFree(apa);CHKERRQ(ierr); 647 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 648 649 PetscFunctionReturn(0); 650 } 651 652 /* Compute local C = P[rstart:rend,:]^T * A * P - used by MatPtAP_MPIAIJ_MPIAIJ() */ 653 654 #undef __FUNCT__ 655 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt" 656 /*@C 657 MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices, 658 used by MatPtAP_MPIAIJ_MPIAIJ() 659 660 Collective on Mat 661 662 Input Parameters: 663 + A - the matrix in mpiaij format 664 . P - the projection matrix in seqaij format 665 . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 666 . fill - expected fill (not being used when scall==MAT_REUSE_MATRIX) 667 . prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose 668 669 Output Parameters: 670 . C - the product matrix in seqaij format 671 672 Level: developer 673 674 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 675 @*/ 676 677 #undef __FUNCT__ 678 #define __FUNCT__ "MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ" 679 PetscErrorCode MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(Mat A,Mat P,PetscInt prstart,PetscInt prend,Mat C) 680 { 681 PetscErrorCode ierr; 682 PetscInt flops=0; 683 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 684 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 685 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 686 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 687 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,*apj,*apjdense,cstart=a->cstart,cend=a->cend,col; 688 PetscInt *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj; 689 PetscInt *ci=c->i,*cj=c->j,*cjj; 690 PetscInt am=A->m,cn=C->N,cm=C->M; 691 PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 692 MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 693 694 PetscFunctionBegin; 695 pA=p->a+pi[prstart]; 696 /* Allocate temporary array for storage of one row of A*P */ 697 ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 698 ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 699 700 apj = (PetscInt *)(apa + cn); 701 apjdense = apj + cn; 702 703 /* Clear old values in C */ 704 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 705 706 for (i=0;i<am;i++) { 707 /* Form sparse row of A*P */ 708 /* diagonal portion of A */ 709 anzi = adi[i+1] - adi[i]; 710 apnzj = 0; 711 for (j=0;j<anzi;j++) { 712 prow = *adj + prstart; 713 adj++; 714 pnzj = pi[prow+1] - pi[prow]; 715 pjj = pj + pi[prow]; 716 paj = pa + pi[prow]; 717 for (k=0;k<pnzj;k++) { 718 if (!apjdense[pjj[k]]) { 719 apjdense[pjj[k]] = -1; 720 apj[apnzj++] = pjj[k]; 721 } 722 apa[pjj[k]] += (*ada)*paj[k]; 723 } 724 flops += 2*pnzj; 725 ada++; 726 } 727 /* off-diagonal portion of A */ 728 anzi = aoi[i+1] - aoi[i]; 729 for (j=0;j<anzi;j++) { 730 col = a->garray[*aoj]; 731 if (col < cstart){ 732 prow = *aoj; 733 } else if (col >= cend){ 734 prow = *aoj + prend-prstart; 735 } else { 736 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map"); 737 } 738 aoj++; 739 pnzj = pi[prow+1] - pi[prow]; 740 pjj = pj + pi[prow]; 741 paj = pa + pi[prow]; 742 for (k=0;k<pnzj;k++) { 743 if (!apjdense[pjj[k]]) { 744 apjdense[pjj[k]] = -1; 745 apj[apnzj++] = pjj[k]; 746 } 747 apa[pjj[k]] += (*aoa)*paj[k]; 748 } 749 flops += 2*pnzj; 750 aoa++; 751 } 752 753 /* Sort the j index array for quick sparse axpy. */ 754 ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 755 756 /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */ 757 pnzi = pi[i+1+prstart] - pi[i+prstart]; 758 for (j=0;j<pnzi;j++) { 759 nextap = 0; 760 crow = *pJ++; 761 cjj = cj + ci[crow]; 762 caj = ca + ci[crow]; 763 /* Perform sparse axpy operation. Note cjj includes apj. */ 764 for (k=0;nextap<apnzj;k++) { 765 if (cjj[k]==apj[nextap]) { 766 caj[k] += (*pA)*apa[apj[nextap++]]; 767 } 768 } 769 flops += 2*apnzj; 770 pA++; 771 } 772 773 /* Zero the current row info for A*P */ 774 for (j=0;j<apnzj;j++) { 775 apa[apj[j]] = 0.; 776 apjdense[apj[j]] = 0; 777 } 778 } 779 780 /* Assemble the final matrix and clean up */ 781 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 782 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 783 ierr = PetscFree(apa);CHKERRQ(ierr); 784 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 785 786 PetscFunctionReturn(0); 787 } 788 789 #undef __FUNCT__ 790 #define __FUNCT__ "MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ" 791 PetscErrorCode MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C) 792 { 793 PetscErrorCode ierr; 794 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 795 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 796 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 797 Mat_SeqAIJ *p=(Mat_SeqAIJ*)P->data,*c; 798 PetscInt *pti,*ptj,*ptJ,*ajj,*pi=p->i,*pj=p->j,*pjj; 799 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz,cstart=a->cstart,cend=a->cend,col; 800 PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; 801 PetscInt an=A->N,am=A->m,pn=P->N,pm=P->M; 802 PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 803 PetscInt m=prend-prstart,nlnk,*lnk; 804 MatScalar *ca; 805 PetscBT lnkbt; 806 807 PetscFunctionBegin; 808 809 /* Get ij structure of P[rstart:rend,:]^T */ 810 ierr = MatGetSymbolicTransposeReduced_SeqAIJ(P,prstart,prend,&pti,&ptj);CHKERRQ(ierr); 811 ptJ=ptj; 812 813 /* Allocate ci array, arrays for fill computation and */ 814 /* free space for accumulating nonzero column info */ 815 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 816 ci[0] = 0; 817 818 ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 819 ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 820 ptasparserow = ptadenserow + an; 821 822 /* create and initialize a linked list */ 823 nlnk = pn+1; 824 ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 825 826 /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 827 /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 828 nnz = adi[am] + aoi[am]; 829 ierr = GetMoreSpace((PetscInt)(fill*nnz/pm)*pn,&free_space); 830 current_space = free_space; 831 832 /* Determine symbolic info for each row of C: */ 833 for (i=0;i<pn;i++) { 834 ptnzi = pti[i+1] - pti[i]; 835 ptanzi = 0; 836 /* Determine symbolic row of PtA_reduced: */ 837 for (j=0;j<ptnzi;j++) { 838 arow = *ptJ++; 839 /* diagonal portion of A */ 840 anzj = adi[arow+1] - adi[arow]; 841 ajj = adj + adi[arow]; 842 for (k=0;k<anzj;k++) { 843 col = ajj[k]+prstart; 844 if (!ptadenserow[col]) { 845 ptadenserow[col] = -1; 846 ptasparserow[ptanzi++] = col; 847 } 848 } 849 /* off-diagonal portion of A */ 850 anzj = aoi[arow+1] - aoi[arow]; 851 ajj = aoj + aoi[arow]; 852 for (k=0;k<anzj;k++) { 853 col = a->garray[ajj[k]]; /* global col */ 854 if (col < cstart){ 855 col = ajj[k]; 856 } else if (col >= cend){ 857 col = ajj[k] + m; 858 } else { 859 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map"); 860 } 861 if (!ptadenserow[col]) { 862 ptadenserow[col] = -1; 863 ptasparserow[ptanzi++] = col; 864 } 865 } 866 } 867 /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 868 ptaj = ptasparserow; 869 cnzi = 0; 870 for (j=0;j<ptanzi;j++) { 871 prow = *ptaj++; 872 pnzj = pi[prow+1] - pi[prow]; 873 pjj = pj + pi[prow]; 874 /* if (rank==1) printf(" prow: %d, pnzj: %d\n",prow,pnzj); */ 875 /* add non-zero cols of P into the sorted linked list lnk */ 876 ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 877 cnzi += nlnk; 878 } 879 880 /* If free space is not available, make more free space */ 881 /* Double the amount of total space in the list */ 882 if (current_space->local_remaining<cnzi) { 883 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 884 } 885 886 /* Copy data into free space, and zero out denserows */ 887 ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 888 current_space->array += cnzi; 889 current_space->local_used += cnzi; 890 current_space->local_remaining -= cnzi; 891 892 for (j=0;j<ptanzi;j++) { 893 ptadenserow[ptasparserow[j]] = 0; 894 } 895 896 /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 897 /* For now, we will recompute what is needed. */ 898 ci[i+1] = ci[i] + cnzi; 899 } 900 /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 901 /* Allocate space for cj, initialize cj, and */ 902 /* destroy list of free space and other temporary array(s) */ 903 ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 904 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 905 ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 906 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 907 908 /* Allocate space for ca */ 909 ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 910 ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 911 912 /* put together the new matrix */ 913 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 914 915 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 916 /* Since these are PETSc arrays, change flags to free them as necessary. */ 917 c = (Mat_SeqAIJ *)((*C)->data); 918 c->freedata = PETSC_TRUE; 919 c->nonew = 0; 920 921 /* Clean up. */ 922 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 923 924 PetscFunctionReturn(0); 925 } 926 927 static PetscEvent logkey_matptapreducedpt = 0; 928 #undef __FUNCT__ 929 #define __FUNCT__ "MatPtAPReducedPt_MPIAIJ_SeqAIJ" 930 PetscErrorCode MatPtAPReducedPt_MPIAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C) 931 { 932 PetscErrorCode ierr; 933 934 PetscFunctionBegin; 935 if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart); 936 if (prend-prstart > P->m) SETERRQ2(PETSC_ERR_ARG_SIZ," prend-prstart %d cannot be larger than P->m %d",prend-prstart,P->m); 937 if (!logkey_matptapreducedpt) { 938 ierr = PetscLogEventRegister(&logkey_matptapreducedpt,"MatPtAP_ReducedPt",MAT_COOKIE); 939 } 940 PetscLogEventBegin(logkey_matptapreducedpt,A,P,0,0);CHKERRQ(ierr); 941 942 if (scall == MAT_INITIAL_MATRIX){ 943 ierr = MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(A,P,fill,prstart,prend,C);CHKERRQ(ierr); 944 } 945 ierr = MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(A,P,prstart,prend,*C);CHKERRQ(ierr); 946 ierr = PetscLogEventEnd(logkey_matptapreducedpt,A,P,0,0);CHKERRQ(ierr); 947 948 PetscFunctionReturn(0); 949 } 950 951