xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 78844af44ce4cb7a2f055fa18fb2484121bec20a)
1be1d678aSKris Buschelman 
2eb9c0419SKris Buschelman /*
342c57489SHong Zhang   Defines projective product routines where A is a SeqAIJ matrix
4eb9c0419SKris Buschelman           C = P^T * A * P
5eb9c0419SKris Buschelman */
6eb9c0419SKris Buschelman 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <petscbt.h>
10eb9c0419SKris Buschelman 
11ff134f7aSHong Zhang #undef __FUNCT__
127ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ"
137ba1a0bfSKris Buschelman PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
147ba1a0bfSKris Buschelman {
157ba1a0bfSKris Buschelman   PetscErrorCode ierr;
167ba1a0bfSKris Buschelman 
177ba1a0bfSKris Buschelman   PetscFunctionBegin;
1865e19b50SBarry Smith   if (!P->ops->ptapsymbolic_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
197ba1a0bfSKris Buschelman   ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr);
207ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
217ba1a0bfSKris Buschelman }
227ba1a0bfSKris Buschelman 
237ba1a0bfSKris Buschelman #undef __FUNCT__
247ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_SeqAIJ"
257ba1a0bfSKris Buschelman PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C)
267ba1a0bfSKris Buschelman {
277ba1a0bfSKris Buschelman   PetscErrorCode ierr;
287ba1a0bfSKris Buschelman 
297ba1a0bfSKris Buschelman   PetscFunctionBegin;
3065e19b50SBarry Smith   if (!P->ops->ptapnumeric_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
317ba1a0bfSKris Buschelman   ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr);
327ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
337ba1a0bfSKris Buschelman }
347ba1a0bfSKris Buschelman 
357ba1a0bfSKris Buschelman #undef __FUNCT__
369af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
3755a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
3855a3bba9SHong Zhang {
39dfbe8321SBarry Smith   PetscErrorCode     ierr;
40a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
41d20bfe6fSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
4255a3bba9SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
43f2b054eeSHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
44d0f46423SBarry Smith   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
45b8374ebeSBarry Smith   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
46d20bfe6fSHong Zhang   MatScalar          *ca;
4755a3bba9SHong Zhang   PetscBT            lnkbt;
48eb9c0419SKris Buschelman 
49eb9c0419SKris Buschelman   PetscFunctionBegin;
50d20bfe6fSHong Zhang   /* Get ij structure of P^T */
51eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
52d20bfe6fSHong Zhang   ptJ=ptj;
53eb9c0419SKris Buschelman 
54d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
55d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
5655a3bba9SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
57d20bfe6fSHong Zhang   ci[0] = 0;
58eb9c0419SKris Buschelman 
5955a3bba9SHong Zhang   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
6055a3bba9SHong Zhang   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
61d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
6255a3bba9SHong Zhang 
6355a3bba9SHong Zhang   /* create and initialize a linked list */
6455a3bba9SHong Zhang   nlnk = pn+1;
6555a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
66eb9c0419SKris Buschelman 
67f2b054eeSHong Zhang   /* Set initial free space to be fill*nnz(A). */
68d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
69f2b054eeSHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
70d20bfe6fSHong Zhang   current_space = free_space;
71d20bfe6fSHong Zhang 
72d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
73d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
74d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
75d20bfe6fSHong Zhang     ptanzi = 0;
76d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
77d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
78d20bfe6fSHong Zhang       arow = *ptJ++;
79d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
80d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
81d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
82d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
83d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
84d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
85d20bfe6fSHong Zhang         }
86d20bfe6fSHong Zhang       }
87d20bfe6fSHong Zhang     }
88d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
89d20bfe6fSHong Zhang     ptaj = ptasparserow;
90d20bfe6fSHong Zhang     cnzi   = 0;
91d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
92d20bfe6fSHong Zhang       prow = *ptaj++;
93d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
94d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
9555a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
9655a3bba9SHong Zhang       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
9755a3bba9SHong Zhang       cnzi += nlnk;
98d20bfe6fSHong Zhang     }
99d20bfe6fSHong Zhang 
100d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
101d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
102d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
1034238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
104f2b054eeSHong Zhang       nspacedouble++;
105d20bfe6fSHong Zhang     }
106d20bfe6fSHong Zhang 
107d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
10855a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
109d20bfe6fSHong Zhang     current_space->array           += cnzi;
110d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
111d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
112d20bfe6fSHong Zhang 
113d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
114d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
115d20bfe6fSHong Zhang     }
116d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
117d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
118d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
119d20bfe6fSHong Zhang   }
120d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
121d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
122d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
12355a3bba9SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
124a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
125d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
12655a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
127d20bfe6fSHong Zhang 
128d20bfe6fSHong Zhang   /* Allocate space for ca */
129d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
130d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
131d20bfe6fSHong Zhang 
132d20bfe6fSHong Zhang   /* put together the new matrix */
1337adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
134d20bfe6fSHong Zhang 
135d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
136d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
137d20bfe6fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
138e6b907acSBarry Smith   c->free_a  = PETSC_TRUE;
139e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
140d20bfe6fSHong Zhang   c->nonew   = 0;
141d20bfe6fSHong Zhang 
142d20bfe6fSHong Zhang   /* Clean up. */
143d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
144f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
145f2b054eeSHong Zhang   if (ci[pn] != 0) {
146f2b054eeSHong Zhang     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
147f2b054eeSHong Zhang     if (afill < 1.0) afill = 1.0;
148f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
149f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
150f2b054eeSHong Zhang   } else {
151f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
152f2b054eeSHong Zhang   }
153f2b054eeSHong Zhang #endif
154eb9c0419SKris Buschelman   PetscFunctionReturn(0);
155eb9c0419SKris Buschelman }
156eb9c0419SKris Buschelman 
157eb9c0419SKris Buschelman #undef __FUNCT__
1589af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
159dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
160dfbe8321SBarry Smith {
161dfbe8321SBarry Smith   PetscErrorCode ierr;
162d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
163d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
164d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
16578844af4SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pcol;
16678844af4SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj,cnz;
1677dca3203SBarry Smith   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
16878844af4SHong Zhang   PetscInt       i,j,k,anz,apnz,pnz,prow,crow,apcol,nextap;
16978844af4SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pval,*ca=c->a,*caj;
17078844af4SHong Zhang   PetscBool      sparse_axpy=PETSC_FALSE;
171eb9c0419SKris Buschelman 
172eb9c0419SKris Buschelman   PetscFunctionBegin;
173d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
174*b5b29cf7SHong Zhang   ierr = PetscMalloc3(cn,MatScalar,&apa,cn,PetscInt,&apj,cn,PetscInt,&apjdense);CHKERRQ(ierr);
175*b5b29cf7SHong Zhang   ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr);
176*b5b29cf7SHong Zhang   ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr);
177d20bfe6fSHong Zhang 
178d20bfe6fSHong Zhang   /* Clear old values in C */
179d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
180d20bfe6fSHong Zhang 
181d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
18278844af4SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
18378844af4SHong Zhang     anz  = ai[i+1] - ai[i];
18478844af4SHong Zhang     apnz = 0;
18578844af4SHong Zhang     for (j=0; j<anz; j++) {
18678844af4SHong Zhang       prow = aj[j];
18778844af4SHong Zhang       pnz  = pi[prow+1] - pi[prow];
18878844af4SHong Zhang       pcol = pj + pi[prow];
18978844af4SHong Zhang       pval = pa + pi[prow];
19078844af4SHong Zhang       for (k=0; k<pnz; k++) {
19178844af4SHong Zhang         if (!apjdense[pcol[k]]) {
19278844af4SHong Zhang           apjdense[pcol[k]] = -1;
19378844af4SHong Zhang           apj[apnz++]       = pcol[k];
194d20bfe6fSHong Zhang         }
19578844af4SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
196d20bfe6fSHong Zhang       }
19778844af4SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
19878844af4SHong Zhang     }
19978844af4SHong Zhang     aj += anz; aa += anz;
20078844af4SHong Zhang 
20178844af4SHong Zhang     if (sparse_axpy){
20278844af4SHong Zhang       ierr = PetscSortInt(apnz,apj);CHKERRQ(ierr);
203d20bfe6fSHong Zhang     }
204d20bfe6fSHong Zhang 
20578844af4SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
20678844af4SHong Zhang     pnz  = pi[i+1] - pi[i];
20778844af4SHong Zhang     pcol = pj + pi[i];
20878844af4SHong Zhang     pval = pa + pi[i];
20978844af4SHong Zhang     for (j=0; j<pnz; j++) {
21078844af4SHong Zhang       crow = pcol[j];
211d20bfe6fSHong Zhang       cjj  = cj + ci[crow];
212d20bfe6fSHong Zhang       caj  = ca + ci[crow];
21378844af4SHong Zhang 
21478844af4SHong Zhang       if (sparse_axpy){  /* Perform sparse axpy */
21578844af4SHong Zhang         nextap = 0;
216*b5b29cf7SHong Zhang         apcol = apj[nextap];
21778844af4SHong Zhang         for (k=0; nextap<apnz; k++) {
218*b5b29cf7SHong Zhang           if (cjj[k] == apcol) {
21978844af4SHong Zhang             caj[k] += pval[j]*apa[apcol];
220*b5b29cf7SHong Zhang             apcol   = apj[++nextap];
221d20bfe6fSHong Zhang           }
222d20bfe6fSHong Zhang         }
22378844af4SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
22478844af4SHong Zhang       } else { /* Perform dense axpy */
22578844af4SHong Zhang         cnz  = ci[crow+1] - ci[crow];
22678844af4SHong Zhang         for (k=0; k<cnz; k++){
22778844af4SHong Zhang           caj[k] += pval[j]*apa[cjj[k]];
22878844af4SHong Zhang         }
22978844af4SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
23078844af4SHong Zhang       }
231d20bfe6fSHong Zhang     }
232d20bfe6fSHong Zhang 
233d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
23478844af4SHong Zhang     for (j=0; j<apnz; j++) {
235*b5b29cf7SHong Zhang       apcol           = apj[j];
236*b5b29cf7SHong Zhang       apa[apcol]      = 0.;
237*b5b29cf7SHong Zhang       apjdense[apcol] = 0;
238d20bfe6fSHong Zhang     }
239d20bfe6fSHong Zhang   }
240d20bfe6fSHong Zhang 
241d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
242d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
244*b5b29cf7SHong Zhang   ierr = PetscFree3(apa,apj,apjdense);CHKERRQ(ierr);
245eb9c0419SKris Buschelman   PetscFunctionReturn(0);
246eb9c0419SKris Buschelman }
247