xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 65e19b508233067d6906a3416498cb9734b3efd5)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3eb9c0419SKris Buschelman /*
442c57489SHong Zhang   Defines projective product routines where A is a SeqAIJ matrix
5eb9c0419SKris Buschelman           C = P^T * A * P
6eb9c0419SKris Buschelman */
7eb9c0419SKris Buschelman 
87c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
97c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
1055a3bba9SHong Zhang #include "petscbt.h"
11eb9c0419SKris Buschelman 
12ff134f7aSHong Zhang #undef __FUNCT__
137ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ"
147ba1a0bfSKris Buschelman PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
157ba1a0bfSKris Buschelman {
167ba1a0bfSKris Buschelman   PetscErrorCode ierr;
177ba1a0bfSKris Buschelman 
187ba1a0bfSKris Buschelman   PetscFunctionBegin;
19*65e19b50SBarry 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);
207ba1a0bfSKris Buschelman   ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr);
217ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
227ba1a0bfSKris Buschelman }
237ba1a0bfSKris Buschelman 
247ba1a0bfSKris Buschelman #undef __FUNCT__
257ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_SeqAIJ"
267ba1a0bfSKris Buschelman PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C)
277ba1a0bfSKris Buschelman {
287ba1a0bfSKris Buschelman   PetscErrorCode ierr;
297ba1a0bfSKris Buschelman 
307ba1a0bfSKris Buschelman   PetscFunctionBegin;
31*65e19b50SBarry 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);
327ba1a0bfSKris Buschelman   ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr);
337ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
347ba1a0bfSKris Buschelman }
357ba1a0bfSKris Buschelman 
367ba1a0bfSKris Buschelman #undef __FUNCT__
379af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
3855a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
3955a3bba9SHong Zhang {
40dfbe8321SBarry Smith   PetscErrorCode     ierr;
41a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
42d20bfe6fSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
4355a3bba9SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
44f2b054eeSHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
45d0f46423SBarry Smith   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
46b8374ebeSBarry Smith   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
47d20bfe6fSHong Zhang   MatScalar          *ca;
4855a3bba9SHong Zhang   PetscBT            lnkbt;
49eb9c0419SKris Buschelman 
50eb9c0419SKris Buschelman   PetscFunctionBegin;
51d20bfe6fSHong Zhang   /* Get ij structure of P^T */
52eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
53d20bfe6fSHong Zhang   ptJ=ptj;
54eb9c0419SKris Buschelman 
55d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
56d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
5755a3bba9SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
58d20bfe6fSHong Zhang   ci[0] = 0;
59eb9c0419SKris Buschelman 
6055a3bba9SHong Zhang   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
6155a3bba9SHong Zhang   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
62d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
6355a3bba9SHong Zhang 
6455a3bba9SHong Zhang   /* create and initialize a linked list */
6555a3bba9SHong Zhang   nlnk = pn+1;
6655a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
67eb9c0419SKris Buschelman 
68f2b054eeSHong Zhang   /* Set initial free space to be fill*nnz(A). */
69d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
70f2b054eeSHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
71d20bfe6fSHong Zhang   current_space = free_space;
72d20bfe6fSHong Zhang 
73d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
74d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
75d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
76d20bfe6fSHong Zhang     ptanzi = 0;
77d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
78d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
79d20bfe6fSHong Zhang       arow = *ptJ++;
80d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
81d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
82d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
83d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
84d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
85d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
86d20bfe6fSHong Zhang         }
87d20bfe6fSHong Zhang       }
88d20bfe6fSHong Zhang     }
89d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
90d20bfe6fSHong Zhang     ptaj = ptasparserow;
91d20bfe6fSHong Zhang     cnzi   = 0;
92d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
93d20bfe6fSHong Zhang       prow = *ptaj++;
94d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
95d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
9655a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
9755a3bba9SHong Zhang       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
9855a3bba9SHong Zhang       cnzi += nlnk;
99d20bfe6fSHong Zhang     }
100d20bfe6fSHong Zhang 
101d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
102d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
103d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
1044238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
105f2b054eeSHong Zhang       nspacedouble++;
106d20bfe6fSHong Zhang     }
107d20bfe6fSHong Zhang 
108d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
10955a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
110d20bfe6fSHong Zhang     current_space->array           += cnzi;
111d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
112d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
113d20bfe6fSHong Zhang 
114d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
115d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
116d20bfe6fSHong Zhang     }
117d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
118d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
119d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
120d20bfe6fSHong Zhang   }
121d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
122d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
123d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
12455a3bba9SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
125a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
126d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
12755a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
128d20bfe6fSHong Zhang 
129d20bfe6fSHong Zhang   /* Allocate space for ca */
130d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
131d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
132d20bfe6fSHong Zhang 
133d20bfe6fSHong Zhang   /* put together the new matrix */
1347adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
135d20bfe6fSHong Zhang 
136d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
137d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
138d20bfe6fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
139e6b907acSBarry Smith   c->free_a  = PETSC_TRUE;
140e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
141d20bfe6fSHong Zhang   c->nonew   = 0;
142d20bfe6fSHong Zhang 
143d20bfe6fSHong Zhang   /* Clean up. */
144d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
145f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
146f2b054eeSHong Zhang   if (ci[pn] != 0) {
147f2b054eeSHong Zhang     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
148f2b054eeSHong Zhang     if (afill < 1.0) afill = 1.0;
149f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
150f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
151f2b054eeSHong Zhang   } else {
152f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
153f2b054eeSHong Zhang   }
154f2b054eeSHong Zhang #endif
155eb9c0419SKris Buschelman   PetscFunctionReturn(0);
156eb9c0419SKris Buschelman }
157eb9c0419SKris Buschelman 
158eb9c0419SKris Buschelman #undef __FUNCT__
1599af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
160dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
161dfbe8321SBarry Smith {
162dfbe8321SBarry Smith   PetscErrorCode ierr;
163dab953b1SSatish Balay   PetscLogDouble flops=0.0;
164d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
165d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
166d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
167b1d57f15SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
168b1d57f15SBarry Smith   PetscInt       *ci=c->i,*cj=c->j,*cjj;
169d0f46423SBarry Smith   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
170b1d57f15SBarry Smith   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
171d20bfe6fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
172eb9c0419SKris Buschelman 
173eb9c0419SKris Buschelman   PetscFunctionBegin;
174d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
175b1d57f15SBarry Smith   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
176b1d57f15SBarry Smith   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
177eb9c0419SKris Buschelman 
178b1d57f15SBarry Smith   apj      = (PetscInt *)(apa + cn);
179d20bfe6fSHong Zhang   apjdense = apj + cn;
180d20bfe6fSHong Zhang 
181d20bfe6fSHong Zhang   /* Clear old values in C */
182d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
183d20bfe6fSHong Zhang 
184d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
185d20bfe6fSHong Zhang     /* Form sparse row of A*P */
186d20bfe6fSHong Zhang     anzi  = ai[i+1] - ai[i];
187d20bfe6fSHong Zhang     apnzj = 0;
188d20bfe6fSHong Zhang     for (j=0;j<anzi;j++) {
189d20bfe6fSHong Zhang       prow = *aj++;
190d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
191d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
192d20bfe6fSHong Zhang       paj  = pa + pi[prow];
193d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
194d20bfe6fSHong Zhang         if (!apjdense[pjj[k]]) {
195d20bfe6fSHong Zhang           apjdense[pjj[k]] = -1;
196d20bfe6fSHong Zhang           apj[apnzj++]     = pjj[k];
197d20bfe6fSHong Zhang         }
198d20bfe6fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
199d20bfe6fSHong Zhang       }
200dab953b1SSatish Balay       flops += 2.0*pnzj;
201d20bfe6fSHong Zhang       aa++;
202d20bfe6fSHong Zhang     }
203d20bfe6fSHong Zhang 
204d20bfe6fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
205d6509036SKris Buschelman     /* Note: a array does not need sorting as it is in dense storage locations. */
206d20bfe6fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
207d20bfe6fSHong Zhang 
208d20bfe6fSHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
209d20bfe6fSHong Zhang     pnzi = pi[i+1] - pi[i];
210d20bfe6fSHong Zhang     for (j=0;j<pnzi;j++) {
211d20bfe6fSHong Zhang       nextap = 0;
212d20bfe6fSHong Zhang       crow   = *pJ++;
213d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
214d20bfe6fSHong Zhang       caj    = ca + ci[crow];
215d20bfe6fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
216d20bfe6fSHong Zhang       for (k=0;nextap<apnzj;k++) {
217f9f3b1caSBarry Smith #if defined(PETSC_USE_DEBUG)
218f9f3b1caSBarry Smith         if (k >= ci[crow+1] - ci[crow]) {
219e32f2f54SBarry Smith           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
220f9f3b1caSBarry Smith         }
221f9f3b1caSBarry Smith #endif
222d20bfe6fSHong Zhang         if (cjj[k]==apj[nextap]) {
223d20bfe6fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
224d20bfe6fSHong Zhang         }
225d20bfe6fSHong Zhang       }
226dab953b1SSatish Balay       flops += 2.0*apnzj;
227d20bfe6fSHong Zhang       pA++;
228d20bfe6fSHong Zhang     }
229d20bfe6fSHong Zhang 
230d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
231d20bfe6fSHong Zhang     for (j=0;j<apnzj;j++) {
232d20bfe6fSHong Zhang       apa[apj[j]]      = 0.;
233d20bfe6fSHong Zhang       apjdense[apj[j]] = 0;
234d20bfe6fSHong Zhang     }
235d20bfe6fSHong Zhang   }
236d20bfe6fSHong Zhang 
237d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
238d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
240d20bfe6fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
241d20bfe6fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
242d20bfe6fSHong Zhang 
243eb9c0419SKris Buschelman   PetscFunctionReturn(0);
244eb9c0419SKris Buschelman }
245