xref: /petsc/src/mat/impls/normal/normm.c (revision 74ed9c26ab880fdcd1d3cdbb5b1e39d1b833147d)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2c8a8475eSBarry Smith 
37c4f633dSBarry Smith #include "private/matimpl.h"          /*I "petscmat.h" I*/
4c8a8475eSBarry Smith 
5c8a8475eSBarry Smith typedef struct {
67cfd0b8cSBarry Smith   Mat         A;
77cfd0b8cSBarry Smith   Vec         w,left,right,leftwork,rightwork;
8b20f02adSBarry Smith   PetscScalar scale;
9c8a8475eSBarry Smith } Mat_Normal;
10c8a8475eSBarry Smith 
11c8a8475eSBarry Smith #undef __FUNCT__
1269ef1043SBarry Smith #define __FUNCT__ "MatScale_Normal"
1369ef1043SBarry Smith PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale)
1469ef1043SBarry Smith {
1569ef1043SBarry Smith   Mat_Normal     *a = (Mat_Normal*)inA->data;
1669ef1043SBarry Smith 
1769ef1043SBarry Smith   PetscFunctionBegin;
1869ef1043SBarry Smith   a->scale *= scale;
1969ef1043SBarry Smith   PetscFunctionReturn(0);
2069ef1043SBarry Smith }
2169ef1043SBarry Smith 
2269ef1043SBarry Smith #undef __FUNCT__
237cfd0b8cSBarry Smith #define __FUNCT__ "MatDiagonalScale_Normal"
247cfd0b8cSBarry Smith PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right)
257cfd0b8cSBarry Smith {
267cfd0b8cSBarry Smith   Mat_Normal     *a = (Mat_Normal*)inA->data;
277cfd0b8cSBarry Smith   PetscErrorCode ierr;
287cfd0b8cSBarry Smith 
297cfd0b8cSBarry Smith   PetscFunctionBegin;
307cfd0b8cSBarry Smith   if (left) {
317cfd0b8cSBarry Smith     if (!a->left) {
327cfd0b8cSBarry Smith       ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr);
337cfd0b8cSBarry Smith       ierr = VecCopy(left,a->left);CHKERRQ(ierr);
347cfd0b8cSBarry Smith     } else {
357cfd0b8cSBarry Smith       ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr);
367cfd0b8cSBarry Smith     }
377cfd0b8cSBarry Smith   }
387cfd0b8cSBarry Smith   if (right) {
397cfd0b8cSBarry Smith     if (!a->right) {
407cfd0b8cSBarry Smith       ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr);
417cfd0b8cSBarry Smith       ierr = VecCopy(right,a->right);CHKERRQ(ierr);
427cfd0b8cSBarry Smith     } else {
437cfd0b8cSBarry Smith       ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr);
447cfd0b8cSBarry Smith     }
457cfd0b8cSBarry Smith   }
467cfd0b8cSBarry Smith   PetscFunctionReturn(0);
477cfd0b8cSBarry Smith }
487cfd0b8cSBarry Smith 
497cfd0b8cSBarry Smith #undef __FUNCT__
50c8a8475eSBarry Smith #define __FUNCT__ "MatMult_Normal"
51dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
52c8a8475eSBarry Smith {
53c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
557cfd0b8cSBarry Smith   Vec            in;
56c8a8475eSBarry Smith 
57c8a8475eSBarry Smith   PetscFunctionBegin;
587cfd0b8cSBarry Smith   in = x;
597cfd0b8cSBarry Smith   if (Na->right) {
607cfd0b8cSBarry Smith     if (!Na->rightwork) {
617cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr);
627cfd0b8cSBarry Smith     }
637cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr);
647cfd0b8cSBarry Smith     in   = Na->rightwork;
657cfd0b8cSBarry Smith   }
667cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
67c8a8475eSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
687cfd0b8cSBarry Smith   if (Na->left) {
697cfd0b8cSBarry Smith     ierr = VecPointwiseMult(y,Na->left,y);CHKERRQ(ierr);
707cfd0b8cSBarry Smith   }
71b20f02adSBarry Smith   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
72c8a8475eSBarry Smith   PetscFunctionReturn(0);
73c8a8475eSBarry Smith }
74c8a8475eSBarry Smith 
75c8a8475eSBarry Smith #undef __FUNCT__
76db090513SMatthew Knepley #define __FUNCT__ "MatMultAdd_Normal"
77db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
78db090513SMatthew Knepley {
79db090513SMatthew Knepley   Mat_Normal     *Na = (Mat_Normal*)N->data;
80db090513SMatthew Knepley   PetscErrorCode ierr;
817cfd0b8cSBarry Smith   Vec            in;
82db090513SMatthew Knepley 
83db090513SMatthew Knepley   PetscFunctionBegin;
847cfd0b8cSBarry Smith   in = v1;
857cfd0b8cSBarry Smith   if (Na->right) {
867cfd0b8cSBarry Smith     if (!Na->rightwork) {
877cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr);
887cfd0b8cSBarry Smith     }
897cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr);
907cfd0b8cSBarry Smith     in   = Na->rightwork;
917cfd0b8cSBarry Smith   }
927cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
93b20f02adSBarry Smith   ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr);
947cfd0b8cSBarry Smith   if (Na->left) {
957cfd0b8cSBarry Smith     ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr);
967cfd0b8cSBarry Smith     ierr = VecPointwiseMult(v3,Na->left,v3);CHKERRQ(ierr);
977cfd0b8cSBarry Smith     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
987cfd0b8cSBarry Smith   } else {
99b20f02adSBarry Smith     ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr);
1007cfd0b8cSBarry Smith   }
101b20f02adSBarry Smith   PetscFunctionReturn(0);
102b20f02adSBarry Smith }
103b20f02adSBarry Smith 
104b20f02adSBarry Smith #undef __FUNCT__
105b20f02adSBarry Smith #define __FUNCT__ "MatMultTranspose_Normal"
106b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y)
107b20f02adSBarry Smith {
108b20f02adSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
109b20f02adSBarry Smith   PetscErrorCode ierr;
1107cfd0b8cSBarry Smith   Vec            in;
111b20f02adSBarry Smith 
112b20f02adSBarry Smith   PetscFunctionBegin;
1137cfd0b8cSBarry Smith   in = x;
1147cfd0b8cSBarry Smith   if (Na->left) {
1157cfd0b8cSBarry Smith     if (!Na->leftwork) {
1167cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr);
1177cfd0b8cSBarry Smith     }
1187cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr);
1197cfd0b8cSBarry Smith     in   = Na->leftwork;
1207cfd0b8cSBarry Smith   }
1217cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
122b20f02adSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
1237cfd0b8cSBarry Smith   if (Na->right) {
1247cfd0b8cSBarry Smith     ierr = VecPointwiseMult(y,Na->right,y);CHKERRQ(ierr);
1257cfd0b8cSBarry Smith   }
126b20f02adSBarry Smith   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
127b20f02adSBarry Smith   PetscFunctionReturn(0);
128b20f02adSBarry Smith }
129b20f02adSBarry Smith 
130b20f02adSBarry Smith #undef __FUNCT__
131b20f02adSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Normal"
132b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
133b20f02adSBarry Smith {
134b20f02adSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
135b20f02adSBarry Smith   PetscErrorCode ierr;
1367cfd0b8cSBarry Smith   Vec            in;
137b20f02adSBarry Smith 
138b20f02adSBarry Smith   PetscFunctionBegin;
1397cfd0b8cSBarry Smith   in = v1;
1407cfd0b8cSBarry Smith   if (Na->left) {
1417cfd0b8cSBarry Smith     if (!Na->leftwork) {
1427cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr);
1437cfd0b8cSBarry Smith     }
1447cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr);
1457cfd0b8cSBarry Smith     in   = Na->leftwork;
1467cfd0b8cSBarry Smith   }
1477cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
148b20f02adSBarry Smith   ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr);
1497cfd0b8cSBarry Smith   if (Na->right) {
1507cfd0b8cSBarry Smith     ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr);
1517cfd0b8cSBarry Smith     ierr = VecPointwiseMult(v3,Na->right,v3);CHKERRQ(ierr);
1527cfd0b8cSBarry Smith     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
1537cfd0b8cSBarry Smith   } else {
154db090513SMatthew Knepley     ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr);
1557cfd0b8cSBarry Smith   }
156db090513SMatthew Knepley   PetscFunctionReturn(0);
157db090513SMatthew Knepley }
158db090513SMatthew Knepley 
159db090513SMatthew Knepley #undef __FUNCT__
160c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal"
161dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N)
162c8a8475eSBarry Smith {
163c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
164dfbe8321SBarry Smith   PetscErrorCode ierr;
165c8a8475eSBarry Smith 
166c8a8475eSBarry Smith   PetscFunctionBegin;
167c3122656SLisandro Dalcin   if (Na->A) { ierr = MatDestroy(Na->A);CHKERRQ(ierr); }
168c3122656SLisandro Dalcin   if (Na->w) { ierr = VecDestroy(Na->w);CHKERRQ(ierr); }
1697cfd0b8cSBarry Smith   if (Na->left) { ierr = VecDestroy(Na->left);CHKERRQ(ierr); }
1707cfd0b8cSBarry Smith   if (Na->right) { ierr = VecDestroy(Na->right);CHKERRQ(ierr); }
1717cfd0b8cSBarry Smith   if (Na->leftwork) { ierr = VecDestroy(Na->leftwork);CHKERRQ(ierr); }
1727cfd0b8cSBarry Smith   if (Na->rightwork) { ierr = VecDestroy(Na->rightwork);CHKERRQ(ierr); }
173c8a8475eSBarry Smith   ierr = PetscFree(Na);CHKERRQ(ierr);
174c8a8475eSBarry Smith   PetscFunctionReturn(0);
175c8a8475eSBarry Smith }
176c8a8475eSBarry Smith 
17717a6fd94SBarry Smith /*
17817a6fd94SBarry Smith       Slow, nonscalable version
17917a6fd94SBarry Smith */
18017a6fd94SBarry Smith #undef __FUNCT__
18117a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal"
182dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
18317a6fd94SBarry Smith {
18417a6fd94SBarry Smith   Mat_Normal        *Na = (Mat_Normal*)N->data;
18517a6fd94SBarry Smith   Mat               A = Na->A;
186dfbe8321SBarry Smith   PetscErrorCode    ierr;
187b24ad042SBarry Smith   PetscInt          i,j,rstart,rend,nnz;
188b24ad042SBarry Smith   const PetscInt    *cols;
18917a6fd94SBarry Smith   PetscScalar       *diag,*work,*values;
190b3cc6726SBarry Smith   const PetscScalar *mvalues;
19117a6fd94SBarry Smith 
19217a6fd94SBarry Smith   PetscFunctionBegin;
193*74ed9c26SBarry Smith   ierr = PetscMalloc2(A->cmap->N,PetscScalar,&diag,A->cmap->N,PetscScalar,&work);CHKERRQ(ierr);
194d0f46423SBarry Smith   ierr = PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
19517a6fd94SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
19617a6fd94SBarry Smith   for (i=rstart; i<rend; i++) {
197b3cc6726SBarry Smith     ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
19817a6fd94SBarry Smith     for (j=0; j<nnz; j++) {
199b3cc6726SBarry Smith       work[cols[j]] += mvalues[j]*mvalues[j];
20017a6fd94SBarry Smith     }
201b3cc6726SBarry Smith     ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
20217a6fd94SBarry Smith   }
203d0f46423SBarry Smith   ierr = MPI_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPI_SUM,((PetscObject)N)->comm);CHKERRQ(ierr);
204d0f46423SBarry Smith   rstart = N->cmap->rstart;
205d0f46423SBarry Smith   rend   = N->cmap->rend;
20617a6fd94SBarry Smith   ierr = VecGetArray(v,&values);CHKERRQ(ierr);
20717a6fd94SBarry Smith   ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr);
20817a6fd94SBarry Smith   ierr = VecRestoreArray(v,&values);CHKERRQ(ierr);
209*74ed9c26SBarry Smith   ierr = PetscFree2(diag,work);CHKERRQ(ierr);
210b20f02adSBarry Smith   ierr = VecScale(v,Na->scale);CHKERRQ(ierr);
21117a6fd94SBarry Smith   PetscFunctionReturn(0);
21217a6fd94SBarry Smith }
213c8a8475eSBarry Smith 
214c8a8475eSBarry Smith #undef __FUNCT__
215c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal"
216c8a8475eSBarry Smith /*@
217c8a8475eSBarry Smith       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
218c8a8475eSBarry Smith 
219c8a8475eSBarry Smith    Collective on Mat
220c8a8475eSBarry Smith 
221c8a8475eSBarry Smith    Input Parameter:
222c8a8475eSBarry Smith .   A  - the (possibly rectangular) matrix
223c8a8475eSBarry Smith 
224c8a8475eSBarry Smith    Output Parameter:
225c8a8475eSBarry Smith .   N - the matrix that represents A'*A
226c8a8475eSBarry Smith 
227c30e7cdbSSatish Balay    Level: intermediate
228c30e7cdbSSatish Balay 
229c8a8475eSBarry Smith    Notes: The product A'*A is NOT actually formed! Rather the new matrix
230c8a8475eSBarry Smith           object performs the matrix-vector product by first multiplying by
231c8a8475eSBarry Smith           A and then A'
232c8a8475eSBarry Smith @*/
233be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat A,Mat *N)
234c8a8475eSBarry Smith {
235dfbe8321SBarry Smith   PetscErrorCode ierr;
236b24ad042SBarry Smith   PetscInt       m,n;
237c8a8475eSBarry Smith   Mat_Normal     *Na;
238c8a8475eSBarry Smith 
239c8a8475eSBarry Smith   PetscFunctionBegin;
240c8a8475eSBarry Smith   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
2417adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,N);CHKERRQ(ierr);
242f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
243c8a8475eSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
244c8a8475eSBarry Smith 
24538f2d2fdSLisandro Dalcin   ierr      = PetscNewLog(*N,Mat_Normal,&Na);CHKERRQ(ierr);
24638f2d2fdSLisandro Dalcin   (*N)->data = (void*) Na;
247c8a8475eSBarry Smith   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
248c3122656SLisandro Dalcin   Na->A     = A;
249b20f02adSBarry Smith   Na->scale = 1.0;
250c8a8475eSBarry Smith 
2517adad957SLisandro Dalcin   ierr    = VecCreateMPI(((PetscObject)A)->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
252c8a8475eSBarry Smith   (*N)->ops->destroy          = MatDestroy_Normal;
253c8a8475eSBarry Smith   (*N)->ops->mult             = MatMult_Normal;
254b20f02adSBarry Smith   (*N)->ops->multtranspose    = MatMultTranspose_Normal;
255b20f02adSBarry Smith   (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal;
256db090513SMatthew Knepley   (*N)->ops->multadd          = MatMultAdd_Normal;
25717a6fd94SBarry Smith   (*N)->ops->getdiagonal      = MatGetDiagonal_Normal;
25869ef1043SBarry Smith   (*N)->ops->scale            = MatScale_Normal;
25969ef1043SBarry Smith   (*N)->ops->diagonalscale    = MatDiagonalScale_Normal;
260c8a8475eSBarry Smith   (*N)->assembled             = PETSC_TRUE;
261d0f46423SBarry Smith   (*N)->cmap->N               = A->cmap->N;
262d0f46423SBarry Smith   (*N)->rmap->N               = A->cmap->N;
263d0f46423SBarry Smith   (*N)->cmap->n               = A->cmap->n;
264d0f46423SBarry Smith   (*N)->rmap->n               = A->cmap->n;
265c8a8475eSBarry Smith   PetscFunctionReturn(0);
266c8a8475eSBarry Smith }
267c8a8475eSBarry Smith 
268