1eb9c0419SKris Buschelman /* 29af31e4aSHong Zhang Defines projective product routines where A is a AIJ matrix 3eb9c0419SKris Buschelman C = P^T * A * P 4eb9c0419SKris Buschelman */ 5eb9c0419SKris Buschelman 6231952e2SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 7eb9c0419SKris Buschelman #include "src/mat/utils/freespace.h" 89af31e4aSHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 9eb9c0419SKris Buschelman 107e44a18eSKris Buschelman EXTERN int RegisterMatMatMultRoutines_Private(Mat); 11eb9c0419SKris Buschelman 129af31e4aSHong Zhang static int MAT_PtAPSymbolic = 0; 139af31e4aSHong Zhang static int MAT_PtAPNumeric = 0; 14eb9c0419SKris Buschelman 15eb9c0419SKris Buschelman #undef __FUNCT__ 169af31e4aSHong Zhang #define __FUNCT__ "MatPtAP" 174d3841fdSKris Buschelman /*@ 189af31e4aSHong Zhang MatPtAP - Creates the matrix projection C = P^T * A * P 194d3841fdSKris Buschelman 204d3841fdSKris Buschelman Collective on Mat 214d3841fdSKris Buschelman 224d3841fdSKris Buschelman Input Parameters: 234d3841fdSKris Buschelman + A - the matrix 24*f747e1aeSHong Zhang . P - the projection matrix 25*f747e1aeSHong Zhang . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 26*f747e1aeSHong Zhang - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)) 274d3841fdSKris Buschelman 284d3841fdSKris Buschelman Output Parameters: 294d3841fdSKris Buschelman . C - the product matrix 304d3841fdSKris Buschelman 314d3841fdSKris Buschelman Notes: 324d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 334d3841fdSKris Buschelman 344d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 354d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. 364d3841fdSKris Buschelman 374d3841fdSKris Buschelman Level: intermediate 384d3841fdSKris Buschelman 399af31e4aSHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 404d3841fdSKris Buschelman @*/ 419af31e4aSHong Zhang int MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) { 42eb9c0419SKris Buschelman int ierr; 43eb9c0419SKris Buschelman 44eb9c0419SKris Buschelman PetscFunctionBegin; 459af31e4aSHong Zhang PetscValidHeaderSpecific(A,MAT_COOKIE,1); 469af31e4aSHong Zhang PetscValidType(A,1); 479af31e4aSHong Zhang MatPreallocated(A); 489af31e4aSHong Zhang if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 499af31e4aSHong Zhang if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 509af31e4aSHong Zhang PetscValidHeaderSpecific(P,MAT_COOKIE,2); 519af31e4aSHong Zhang PetscValidType(P,2); 529af31e4aSHong Zhang MatPreallocated(P); 539af31e4aSHong Zhang if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 549af31e4aSHong Zhang if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 559af31e4aSHong Zhang PetscValidPointer(C,3); 569af31e4aSHong Zhang if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 57eb9c0419SKris Buschelman 589af31e4aSHong Zhang if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 59eb9c0419SKris Buschelman 609af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 619af31e4aSHong Zhang ierr = (*A->ops->matptap)(A,P,scall,fill,C);CHKERRQ(ierr); 629af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 63eb9c0419SKris Buschelman PetscFunctionReturn(0); 64eb9c0419SKris Buschelman } 65eb9c0419SKris Buschelman 66eb9c0419SKris Buschelman #undef __FUNCT__ 679af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 689af31e4aSHong Zhang int MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 699af31e4aSHong Zhang { 709af31e4aSHong Zhang int ierr; 719af31e4aSHong Zhang PetscFunctionBegin; 729af31e4aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 739af31e4aSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 749af31e4aSHong Zhang } 759af31e4aSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 769af31e4aSHong Zhang PetscFunctionReturn(0); 779af31e4aSHong Zhang } 789af31e4aSHong Zhang 799af31e4aSHong Zhang #undef __FUNCT__ 809af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic" 814d3841fdSKris Buschelman /*@ 829af31e4aSHong Zhang MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 834d3841fdSKris Buschelman 844d3841fdSKris Buschelman Collective on Mat 854d3841fdSKris Buschelman 864d3841fdSKris Buschelman Input Parameters: 874d3841fdSKris Buschelman + A - the matrix 884d3841fdSKris Buschelman - P - the projection matrix 894d3841fdSKris Buschelman 904d3841fdSKris Buschelman Output Parameters: 914d3841fdSKris Buschelman . C - the (i,j) structure of the product matrix 924d3841fdSKris Buschelman 934d3841fdSKris Buschelman Notes: 944d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 954d3841fdSKris Buschelman 964d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 974d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 989af31e4aSHong Zhang this (i,j) structure by calling MatPtAPNumeric(). 994d3841fdSKris Buschelman 1004d3841fdSKris Buschelman Level: intermediate 1014d3841fdSKris Buschelman 1029af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 1034d3841fdSKris Buschelman @*/ 1049af31e4aSHong Zhang int MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) { 105eb9c0419SKris Buschelman int ierr; 106eb9c0419SKris Buschelman char funct[80]; 107bf812060SBarry Smith int (*f)(Mat,Mat,Mat*); 108eb9c0419SKris Buschelman 109eb9c0419SKris Buschelman PetscFunctionBegin; 110eb9c0419SKris Buschelman 1114482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 112c9780b6fSBarry Smith PetscValidType(A,1); 113eb9c0419SKris Buschelman MatPreallocated(A); 114eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 115eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 116eb9c0419SKris Buschelman 1174482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 118c9780b6fSBarry Smith PetscValidType(P,2); 119eb9c0419SKris Buschelman MatPreallocated(P); 120eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 121eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 122eb9c0419SKris Buschelman 1234482741eSBarry Smith PetscValidPointer(C,3); 1244482741eSBarry Smith 125eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 126eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 127eb9c0419SKris Buschelman 1284d3841fdSKris Buschelman /* Currently only _seqaij_seqaij is implemented, so just query for it. */ 1294d3841fdSKris Buschelman /* When other implementations exist, attack the multiple dispatch problem. */ 1309af31e4aSHong Zhang ierr = PetscStrcpy(funct,"MatPtAPSymbolic_seqaij_seqaij");CHKERRQ(ierr); 1314d3841fdSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)P,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); 1329af31e4aSHong Zhang if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic is not supported for P of type %s",P->type_name); 1334d3841fdSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); 1349af31e4aSHong Zhang if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic is not supported for A of type %s",A->type_name); 1354d3841fdSKris Buschelman 136cd34a33cSKris Buschelman ierr = (*f)(A,P,C);CHKERRQ(ierr); 137eb9c0419SKris Buschelman 138eb9c0419SKris Buschelman PetscFunctionReturn(0); 139eb9c0419SKris Buschelman } 140eb9c0419SKris Buschelman 141*f747e1aeSHong Zhang typedef struct { 142*f747e1aeSHong Zhang Mat symAP; 143*f747e1aeSHong Zhang } Mat_PtAPstruct; 144*f747e1aeSHong Zhang 145*f747e1aeSHong Zhang #undef __FUNCT__ 146*f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 147*f747e1aeSHong Zhang int MatDestroy_SeqAIJ_PtAP(Mat A) 148*f747e1aeSHong Zhang { 149*f747e1aeSHong Zhang int ierr; 150*f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; 151*f747e1aeSHong Zhang 152*f747e1aeSHong Zhang PetscFunctionBegin; 153*f747e1aeSHong Zhang ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); 154*f747e1aeSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 155*f747e1aeSHong Zhang 156*f747e1aeSHong Zhang ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 157*f747e1aeSHong Zhang PetscFunctionReturn(0); 158*f747e1aeSHong Zhang } 159*f747e1aeSHong Zhang 160eb9c0419SKris Buschelman #undef __FUNCT__ 1619af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 1629af31e4aSHong Zhang int MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) { 163*f747e1aeSHong Zhang int ierr,*pti,*ptj; 164*f747e1aeSHong Zhang Mat Pt,AP; 165*f747e1aeSHong Zhang Mat_PtAPstruct *ptap; 166eb9c0419SKris Buschelman 167eb9c0419SKris Buschelman PetscFunctionBegin; 168*f747e1aeSHong Zhang /* create symbolic Pt */ 169eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 170*f747e1aeSHong Zhang ierr = MatCreateSeqAIJWithArrays(P->comm,P->N,P->M,pti,ptj,PETSC_NULL,&Pt);CHKERRQ(ierr); 171eb9c0419SKris Buschelman 172*f747e1aeSHong Zhang /* get symbolic AP=A*P and C=Pt*AP */ 173*f747e1aeSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 174*f747e1aeSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 175eb9c0419SKris Buschelman 176*f747e1aeSHong Zhang /* clean up */ 177*f747e1aeSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(Pt,&pti,&ptj);CHKERRQ(ierr); 178*f747e1aeSHong Zhang ierr = MatDestroy(Pt);CHKERRQ(ierr); 179eb9c0419SKris Buschelman 180*f747e1aeSHong Zhang /* save symbolic AP - to be used by MatPtAPNumeric_SeqAIJ_SeqAIJ() */ 181*f747e1aeSHong Zhang ierr = PetscNew(Mat_PtAPstruct,&ptap);CHKERRQ(ierr); 182*f747e1aeSHong Zhang ptap->symAP = AP; 183*f747e1aeSHong Zhang (*C)->spptr = (void*)ptap; 184*f747e1aeSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 185eb9c0419SKris Buschelman 186eb9c0419SKris Buschelman PetscFunctionReturn(0); 187eb9c0419SKris Buschelman } 188eb9c0419SKris Buschelman 1893985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h" 1903985e5eaSKris Buschelman EXTERN_C_BEGIN 1913985e5eaSKris Buschelman #undef __FUNCT__ 1929af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 1939af31e4aSHong Zhang int MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) { 1945c66b693SKris Buschelman /* This routine requires testing -- I don't think it works. */ 1953985e5eaSKris Buschelman int ierr; 1963985e5eaSKris Buschelman FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1973985e5eaSKris Buschelman Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 1983985e5eaSKris Buschelman Mat P=pp->AIJ; 1993985e5eaSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 2003985e5eaSKris Buschelman int *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 2013985e5eaSKris Buschelman int *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 2023985e5eaSKris Buschelman int an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 203fe05a634SKris Buschelman int i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 2043985e5eaSKris Buschelman MatScalar *ca; 2053985e5eaSKris Buschelman 2063985e5eaSKris Buschelman PetscFunctionBegin; 2073985e5eaSKris Buschelman /* Start timer */ 2089af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 2093985e5eaSKris Buschelman 2103985e5eaSKris Buschelman /* Get ij structure of P^T */ 2113985e5eaSKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 2123985e5eaSKris Buschelman 2133985e5eaSKris Buschelman /* Allocate ci array, arrays for fill computation and */ 2143985e5eaSKris Buschelman /* free space for accumulating nonzero column info */ 2153985e5eaSKris Buschelman ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr); 2163985e5eaSKris Buschelman ci[0] = 0; 2173985e5eaSKris Buschelman 2183985e5eaSKris Buschelman ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr); 2193985e5eaSKris Buschelman ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr); 2203985e5eaSKris Buschelman ptasparserow = ptadenserow + an; 2213985e5eaSKris Buschelman denserow = ptasparserow + an; 2223985e5eaSKris Buschelman sparserow = denserow + pn; 2233985e5eaSKris Buschelman 2243985e5eaSKris Buschelman /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 2253985e5eaSKris Buschelman /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 2263985e5eaSKris Buschelman ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 2273985e5eaSKris Buschelman current_space = free_space; 2283985e5eaSKris Buschelman 2293985e5eaSKris Buschelman /* Determine symbolic info for each row of C: */ 2303985e5eaSKris Buschelman for (i=0;i<pn/ppdof;i++) { 2313985e5eaSKris Buschelman ptnzi = pti[i+1] - pti[i]; 2323985e5eaSKris Buschelman ptanzi = 0; 2333985e5eaSKris Buschelman ptJ = ptj + pti[i]; 2343985e5eaSKris Buschelman for (dof=0;dof<ppdof;dof++) { 2353985e5eaSKris Buschelman /* Determine symbolic row of PtA: */ 2363985e5eaSKris Buschelman for (j=0;j<ptnzi;j++) { 2373985e5eaSKris Buschelman arow = ptJ[j] + dof; 2383985e5eaSKris Buschelman anzj = ai[arow+1] - ai[arow]; 2393985e5eaSKris Buschelman ajj = aj + ai[arow]; 2403985e5eaSKris Buschelman for (k=0;k<anzj;k++) { 2413985e5eaSKris Buschelman if (!ptadenserow[ajj[k]]) { 2423985e5eaSKris Buschelman ptadenserow[ajj[k]] = -1; 2433985e5eaSKris Buschelman ptasparserow[ptanzi++] = ajj[k]; 2443985e5eaSKris Buschelman } 2453985e5eaSKris Buschelman } 2463985e5eaSKris Buschelman } 2473985e5eaSKris Buschelman /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 2483985e5eaSKris Buschelman ptaj = ptasparserow; 2493985e5eaSKris Buschelman cnzi = 0; 2503985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 251fe05a634SKris Buschelman pdof = *ptaj%dof; 2523985e5eaSKris Buschelman prow = (*ptaj++)/dof; 2533985e5eaSKris Buschelman pnzj = pi[prow+1] - pi[prow]; 2543985e5eaSKris Buschelman pjj = pj + pi[prow]; 2553985e5eaSKris Buschelman for (k=0;k<pnzj;k++) { 256fe05a634SKris Buschelman if (!denserow[pjj[k]+pdof]) { 257fe05a634SKris Buschelman denserow[pjj[k]+pdof] = -1; 258fe05a634SKris Buschelman sparserow[cnzi++] = pjj[k]+pdof; 2593985e5eaSKris Buschelman } 2603985e5eaSKris Buschelman } 2613985e5eaSKris Buschelman } 2623985e5eaSKris Buschelman 2633985e5eaSKris Buschelman /* sort sparserow */ 2643985e5eaSKris Buschelman ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 2653985e5eaSKris Buschelman 2663985e5eaSKris Buschelman /* If free space is not available, make more free space */ 2673985e5eaSKris Buschelman /* Double the amount of total space in the list */ 2683985e5eaSKris Buschelman if (current_space->local_remaining<cnzi) { 2693985e5eaSKris Buschelman ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2703985e5eaSKris Buschelman } 2713985e5eaSKris Buschelman 2723985e5eaSKris Buschelman /* Copy data into free space, and zero out denserows */ 2733985e5eaSKris Buschelman ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr); 2743985e5eaSKris Buschelman current_space->array += cnzi; 2753985e5eaSKris Buschelman current_space->local_used += cnzi; 2763985e5eaSKris Buschelman current_space->local_remaining -= cnzi; 2773985e5eaSKris Buschelman 2783985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 2793985e5eaSKris Buschelman ptadenserow[ptasparserow[j]] = 0; 2803985e5eaSKris Buschelman } 2813985e5eaSKris Buschelman for (j=0;j<cnzi;j++) { 2823985e5eaSKris Buschelman denserow[sparserow[j]] = 0; 2833985e5eaSKris Buschelman } 2843985e5eaSKris Buschelman /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 2853985e5eaSKris Buschelman /* For now, we will recompute what is needed. */ 2863985e5eaSKris Buschelman ci[i+1+dof] = ci[i+dof] + cnzi; 2873985e5eaSKris Buschelman } 2883985e5eaSKris Buschelman } 2893985e5eaSKris Buschelman /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 2903985e5eaSKris Buschelman /* Allocate space for cj, initialize cj, and */ 2913985e5eaSKris Buschelman /* destroy list of free space and other temporary array(s) */ 2923985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr); 2933985e5eaSKris Buschelman ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 2943985e5eaSKris Buschelman ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 2953985e5eaSKris Buschelman 2963985e5eaSKris Buschelman /* Allocate space for ca */ 2973985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 2983985e5eaSKris Buschelman ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 2993985e5eaSKris Buschelman 3003985e5eaSKris Buschelman /* put together the new matrix */ 3013985e5eaSKris Buschelman ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 3023985e5eaSKris Buschelman 3033985e5eaSKris Buschelman /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 3043985e5eaSKris Buschelman /* Since these are PETSc arrays, change flags to free them as necessary. */ 3053985e5eaSKris Buschelman c = (Mat_SeqAIJ *)((*C)->data); 3063985e5eaSKris Buschelman c->freedata = PETSC_TRUE; 3073985e5eaSKris Buschelman c->nonew = 0; 3083985e5eaSKris Buschelman 3093985e5eaSKris Buschelman /* Clean up. */ 3103985e5eaSKris Buschelman ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 3113985e5eaSKris Buschelman 3129af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 3133985e5eaSKris Buschelman PetscFunctionReturn(0); 3143985e5eaSKris Buschelman } 3153985e5eaSKris Buschelman EXTERN_C_END 3163985e5eaSKris Buschelman 317eb9c0419SKris Buschelman #undef __FUNCT__ 3189af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric" 3194d3841fdSKris Buschelman /*@ 3209af31e4aSHong Zhang MatPtAPNumeric - Computes the matrix projection C = P^T * A * P 3214d3841fdSKris Buschelman 3224d3841fdSKris Buschelman Collective on Mat 3234d3841fdSKris Buschelman 3244d3841fdSKris Buschelman Input Parameters: 3254d3841fdSKris Buschelman + A - the matrix 3264d3841fdSKris Buschelman - P - the projection matrix 3274d3841fdSKris Buschelman 3284d3841fdSKris Buschelman Output Parameters: 3294d3841fdSKris Buschelman . C - the product matrix 3304d3841fdSKris Buschelman 3314d3841fdSKris Buschelman Notes: 3329af31e4aSHong Zhang C must have been created by calling MatPtAPSymbolic and must be destroyed by 3334d3841fdSKris Buschelman the user using MatDeatroy(). 3344d3841fdSKris Buschelman 3354d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 3364d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. 3374d3841fdSKris Buschelman 3384d3841fdSKris Buschelman Level: intermediate 3394d3841fdSKris Buschelman 3409af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() 3414d3841fdSKris Buschelman @*/ 3429af31e4aSHong Zhang int MatPtAPNumeric(Mat A,Mat P,Mat C) { 343eb9c0419SKris Buschelman int ierr; 344eb9c0419SKris Buschelman char funct[80]; 3454d3841fdSKris Buschelman int (*f)(Mat,Mat,Mat); 346eb9c0419SKris Buschelman 347eb9c0419SKris Buschelman PetscFunctionBegin; 348eb9c0419SKris Buschelman 3494482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 350c9780b6fSBarry Smith PetscValidType(A,1); 351eb9c0419SKris Buschelman MatPreallocated(A); 352eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 353eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 354eb9c0419SKris Buschelman 3554482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 356c9780b6fSBarry Smith PetscValidType(P,2); 357eb9c0419SKris Buschelman MatPreallocated(P); 358eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 359eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 360eb9c0419SKris Buschelman 3614482741eSBarry Smith PetscValidHeaderSpecific(C,MAT_COOKIE,3); 362c9780b6fSBarry Smith PetscValidType(C,3); 363eb9c0419SKris Buschelman MatPreallocated(C); 364eb9c0419SKris Buschelman if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 365eb9c0419SKris Buschelman if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 366eb9c0419SKris Buschelman 367eb9c0419SKris Buschelman if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); 368eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 369eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 370eb9c0419SKris Buschelman if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); 371eb9c0419SKris Buschelman 3724d3841fdSKris Buschelman /* Currently only _seqaij_seqaij is implemented, so just query for it. */ 3734d3841fdSKris Buschelman /* When other implementations exist, attack the multiple dispatch problem. */ 3749af31e4aSHong Zhang ierr = PetscStrcpy(funct,"MatPtAPNumeric_seqaij_seqaij");CHKERRQ(ierr); 3754d3841fdSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)P,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); 3769af31e4aSHong Zhang if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric is not supported for P of type %s",P->type_name); 3774d3841fdSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); 3789af31e4aSHong Zhang if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric is not supported for A of type %s",A->type_name); 3794d3841fdSKris Buschelman 3804d3841fdSKris Buschelman ierr = (*f)(A,P,C);CHKERRQ(ierr); 381eb9c0419SKris Buschelman 382eb9c0419SKris Buschelman PetscFunctionReturn(0); 383eb9c0419SKris Buschelman } 384eb9c0419SKris Buschelman 385eb9c0419SKris Buschelman #undef __FUNCT__ 3869af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 3879af31e4aSHong Zhang int MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) { 388eb9c0419SKris Buschelman int ierr,flops=0; 389eb9c0419SKris Buschelman Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 390eb9c0419SKris Buschelman Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 391eb9c0419SKris Buschelman Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 392*f747e1aeSHong Zhang int *ai=a->i,*aj=a->j,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 393eb9c0419SKris Buschelman int *ci=c->i,*cj=c->j,*cjj; 394eb9c0419SKris Buschelman int am=A->M,cn=C->N,cm=C->M; 395eb9c0419SKris Buschelman int i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 396eb9c0419SKris Buschelman MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 397*f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)C->spptr; 398*f747e1aeSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ *)(ptap->symAP)->data; 399*f747e1aeSHong Zhang int *api=ap->i,*apj=ap->j,apj_nextap; 400eb9c0419SKris Buschelman 401eb9c0419SKris Buschelman PetscFunctionBegin; 402eb9c0419SKris Buschelman /* Allocate temporary array for storage of one row of A*P */ 403*f747e1aeSHong Zhang ierr = PetscMalloc(cn*sizeof(MatScalar),&apa);CHKERRQ(ierr); 404*f747e1aeSHong Zhang ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr); 405eb9c0419SKris Buschelman 406eb9c0419SKris Buschelman /* Clear old values in C */ 407eb9c0419SKris Buschelman ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 408eb9c0419SKris Buschelman 409eb9c0419SKris Buschelman for (i=0;i<am;i++) { 410*f747e1aeSHong Zhang /* Get sparse values of A*P[i,:] */ 411eb9c0419SKris Buschelman anzi = ai[i+1] - ai[i]; 412eb9c0419SKris Buschelman apnzj = 0; 413eb9c0419SKris Buschelman for (j=0;j<anzi;j++) { 414eb9c0419SKris Buschelman prow = *aj++; 415eb9c0419SKris Buschelman pnzj = pi[prow+1] - pi[prow]; 416eb9c0419SKris Buschelman pjj = pj + pi[prow]; 417eb9c0419SKris Buschelman paj = pa + pi[prow]; 418eb9c0419SKris Buschelman for (k=0;k<pnzj;k++) { 419eb9c0419SKris Buschelman apa[pjj[k]] += (*aa)*paj[k]; 420eb9c0419SKris Buschelman } 421eb9c0419SKris Buschelman flops += 2*pnzj; 422eb9c0419SKris Buschelman aa++; 423eb9c0419SKris Buschelman } 424eb9c0419SKris Buschelman 425eb9c0419SKris Buschelman /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 426*f747e1aeSHong Zhang apj = ap->j + api[i]; 427*f747e1aeSHong Zhang apnzj = api[i+1] - api[i]; 428eb9c0419SKris Buschelman pnzi = pi[i+1] - pi[i]; 429eb9c0419SKris Buschelman for (j=0;j<pnzi;j++) { 430eb9c0419SKris Buschelman nextap = 0; 431eb9c0419SKris Buschelman crow = *pJ++; 432eb9c0419SKris Buschelman cjj = cj + ci[crow]; 433eb9c0419SKris Buschelman caj = ca + ci[crow]; 434eb9c0419SKris Buschelman /* Perform sparse axpy operation. Note cjj includes apj. */ 435eb9c0419SKris Buschelman for (k=0; nextap<apnzj; k++) { 436*f747e1aeSHong Zhang apj_nextap = *(apj+nextap); 437*f747e1aeSHong Zhang if (cjj[k]==apj_nextap) { 438*f747e1aeSHong Zhang caj[k] += (*pA)*apa[apj_nextap]; 439*f747e1aeSHong Zhang nextap++; 440eb9c0419SKris Buschelman } 441eb9c0419SKris Buschelman } 442eb9c0419SKris Buschelman flops += 2*apnzj; 443eb9c0419SKris Buschelman pA++; 444eb9c0419SKris Buschelman } 445eb9c0419SKris Buschelman 446*f747e1aeSHong Zhang /* Zero the current row values for A*P */ 447*f747e1aeSHong Zhang for (j=0;j<apnzj;j++) apa[apj[j]] = 0.0; 448eb9c0419SKris Buschelman } 449eb9c0419SKris Buschelman 450eb9c0419SKris Buschelman /* Assemble the final matrix and clean up */ 451eb9c0419SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 452eb9c0419SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 453eb9c0419SKris Buschelman ierr = PetscFree(apa);CHKERRQ(ierr); 454eb9c0419SKris Buschelman ierr = PetscLogFlops(flops);CHKERRQ(ierr); 455eb9c0419SKris Buschelman 456eb9c0419SKris Buschelman PetscFunctionReturn(0); 457eb9c0419SKris Buschelman } 458eb9c0419SKris Buschelman 459eb9c0419SKris Buschelman #undef __FUNCT__ 4609af31e4aSHong Zhang #define __FUNCT__ "RegisterPtAPRoutines_Private" 4619af31e4aSHong Zhang int RegisterPtAPRoutines_Private(Mat A) { 462eb9c0419SKris Buschelman int ierr; 463eb9c0419SKris Buschelman 464eb9c0419SKris Buschelman PetscFunctionBegin; 4659af31e4aSHong Zhang if (!MAT_PtAPSymbolic) { 4669af31e4aSHong Zhang ierr = PetscLogEventRegister(&MAT_PtAPSymbolic,"MatPtAPSymbolic",MAT_COOKIE);CHKERRQ(ierr); 467eb9c0419SKris Buschelman } 4689af31e4aSHong Zhang if (!MAT_PtAPNumeric) { 4699af31e4aSHong Zhang ierr = PetscLogEventRegister(&MAT_PtAPNumeric,"MatPtAPNumeric",MAT_COOKIE);CHKERRQ(ierr); 470eb9c0419SKris Buschelman } 4715c66b693SKris Buschelman ierr = RegisterMatMatMultRoutines_Private(A);CHKERRQ(ierr); 4729af31e4aSHong Zhang 473eb9c0419SKris Buschelman PetscFunctionReturn(0); 474eb9c0419SKris Buschelman } 475