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