xref: /petsc/src/mat/impls/normal/normm.c (revision 7cfd0b8c5eefbe860503f3f6a36a2b6430f202f5)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2c8a8475eSBarry Smith 
37c4f633dSBarry Smith #include "private/matimpl.h"          /*I "petscmat.h" I*/
4c8a8475eSBarry Smith 
5c8a8475eSBarry Smith typedef struct {
6*7cfd0b8cSBarry Smith   Mat         A;
7*7cfd0b8cSBarry Smith   Vec         w,left,right,leftwork,rightwork;
8b20f02adSBarry Smith   PetscScalar scale;
9c8a8475eSBarry Smith } Mat_Normal;
10c8a8475eSBarry Smith 
11c8a8475eSBarry Smith #undef __FUNCT__
12*7cfd0b8cSBarry Smith #define __FUNCT__ "MatDiagonalScale_Normal"
13*7cfd0b8cSBarry Smith PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right)
14*7cfd0b8cSBarry Smith {
15*7cfd0b8cSBarry Smith   Mat_Normal     *a = (Mat_Normal*)inA->data;
16*7cfd0b8cSBarry Smith   PetscErrorCode ierr;
17*7cfd0b8cSBarry Smith 
18*7cfd0b8cSBarry Smith   PetscFunctionBegin;
19*7cfd0b8cSBarry Smith   if (left) {
20*7cfd0b8cSBarry Smith     if (!a->left) {
21*7cfd0b8cSBarry Smith       ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr);
22*7cfd0b8cSBarry Smith       ierr = VecCopy(left,a->left);CHKERRQ(ierr);
23*7cfd0b8cSBarry Smith     } else {
24*7cfd0b8cSBarry Smith       ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr);
25*7cfd0b8cSBarry Smith     }
26*7cfd0b8cSBarry Smith   }
27*7cfd0b8cSBarry Smith   if (right) {
28*7cfd0b8cSBarry Smith     if (!a->right) {
29*7cfd0b8cSBarry Smith       ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr);
30*7cfd0b8cSBarry Smith       ierr = VecCopy(right,a->right);CHKERRQ(ierr);
31*7cfd0b8cSBarry Smith     } else {
32*7cfd0b8cSBarry Smith       ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr);
33*7cfd0b8cSBarry Smith     }
34*7cfd0b8cSBarry Smith   }
35*7cfd0b8cSBarry Smith   PetscFunctionReturn(0);
36*7cfd0b8cSBarry Smith }
37*7cfd0b8cSBarry Smith 
38*7cfd0b8cSBarry Smith #undef __FUNCT__
39c8a8475eSBarry Smith #define __FUNCT__ "MatMult_Normal"
40dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
41c8a8475eSBarry Smith {
42c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
43dfbe8321SBarry Smith   PetscErrorCode ierr;
44*7cfd0b8cSBarry Smith   Vec            in;
45c8a8475eSBarry Smith 
46c8a8475eSBarry Smith   PetscFunctionBegin;
47*7cfd0b8cSBarry Smith   in = x;
48*7cfd0b8cSBarry Smith   if (Na->right) {
49*7cfd0b8cSBarry Smith     if (!Na->rightwork) {
50*7cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr);
51*7cfd0b8cSBarry Smith     }
52*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr);
53*7cfd0b8cSBarry Smith     in   = Na->rightwork;
54*7cfd0b8cSBarry Smith   }
55*7cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
56c8a8475eSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
57*7cfd0b8cSBarry Smith   if (Na->left) {
58*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(y,Na->left,y);CHKERRQ(ierr);
59*7cfd0b8cSBarry Smith   }
60b20f02adSBarry Smith   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
61c8a8475eSBarry Smith   PetscFunctionReturn(0);
62c8a8475eSBarry Smith }
63c8a8475eSBarry Smith 
64c8a8475eSBarry Smith #undef __FUNCT__
65db090513SMatthew Knepley #define __FUNCT__ "MatMultAdd_Normal"
66db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
67db090513SMatthew Knepley {
68db090513SMatthew Knepley   Mat_Normal     *Na = (Mat_Normal*)N->data;
69db090513SMatthew Knepley   PetscErrorCode ierr;
70*7cfd0b8cSBarry Smith   Vec            in;
71db090513SMatthew Knepley 
72db090513SMatthew Knepley   PetscFunctionBegin;
73*7cfd0b8cSBarry Smith   in = v1;
74*7cfd0b8cSBarry Smith   if (Na->right) {
75*7cfd0b8cSBarry Smith     if (!Na->rightwork) {
76*7cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr);
77*7cfd0b8cSBarry Smith     }
78*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr);
79*7cfd0b8cSBarry Smith     in   = Na->rightwork;
80*7cfd0b8cSBarry Smith   }
81*7cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
82b20f02adSBarry Smith   ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr);
83*7cfd0b8cSBarry Smith   if (Na->left) {
84*7cfd0b8cSBarry Smith     ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr);
85*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(v3,Na->left,v3);CHKERRQ(ierr);
86*7cfd0b8cSBarry Smith     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
87*7cfd0b8cSBarry Smith   } else {
88b20f02adSBarry Smith     ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr);
89*7cfd0b8cSBarry Smith   }
90b20f02adSBarry Smith   PetscFunctionReturn(0);
91b20f02adSBarry Smith }
92b20f02adSBarry Smith 
93b20f02adSBarry Smith #undef __FUNCT__
94b20f02adSBarry Smith #define __FUNCT__ "MatMultTranspose_Normal"
95b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y)
96b20f02adSBarry Smith {
97b20f02adSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
98b20f02adSBarry Smith   PetscErrorCode ierr;
99*7cfd0b8cSBarry Smith   Vec            in;
100b20f02adSBarry Smith 
101b20f02adSBarry Smith   PetscFunctionBegin;
102*7cfd0b8cSBarry Smith   in = x;
103*7cfd0b8cSBarry Smith   if (Na->left) {
104*7cfd0b8cSBarry Smith     if (!Na->leftwork) {
105*7cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr);
106*7cfd0b8cSBarry Smith     }
107*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr);
108*7cfd0b8cSBarry Smith     in   = Na->leftwork;
109*7cfd0b8cSBarry Smith   }
110*7cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
111b20f02adSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
112*7cfd0b8cSBarry Smith   if (Na->right) {
113*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(y,Na->right,y);CHKERRQ(ierr);
114*7cfd0b8cSBarry Smith   }
115b20f02adSBarry Smith   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
116b20f02adSBarry Smith   PetscFunctionReturn(0);
117b20f02adSBarry Smith }
118b20f02adSBarry Smith 
119b20f02adSBarry Smith #undef __FUNCT__
120b20f02adSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Normal"
121b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
122b20f02adSBarry Smith {
123b20f02adSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
124b20f02adSBarry Smith   PetscErrorCode ierr;
125*7cfd0b8cSBarry Smith   Vec            in;
126b20f02adSBarry Smith 
127b20f02adSBarry Smith   PetscFunctionBegin;
128*7cfd0b8cSBarry Smith   in = v1;
129*7cfd0b8cSBarry Smith   if (Na->left) {
130*7cfd0b8cSBarry Smith     if (!Na->leftwork) {
131*7cfd0b8cSBarry Smith       ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr);
132*7cfd0b8cSBarry Smith     }
133*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr);
134*7cfd0b8cSBarry Smith     in   = Na->leftwork;
135*7cfd0b8cSBarry Smith   }
136*7cfd0b8cSBarry Smith   ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr);
137b20f02adSBarry Smith   ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr);
138*7cfd0b8cSBarry Smith   if (Na->right) {
139*7cfd0b8cSBarry Smith     ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr);
140*7cfd0b8cSBarry Smith     ierr = VecPointwiseMult(v3,Na->right,v3);CHKERRQ(ierr);
141*7cfd0b8cSBarry Smith     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
142*7cfd0b8cSBarry Smith   } else {
143db090513SMatthew Knepley     ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr);
144*7cfd0b8cSBarry Smith   }
145db090513SMatthew Knepley   PetscFunctionReturn(0);
146db090513SMatthew Knepley }
147db090513SMatthew Knepley 
148db090513SMatthew Knepley #undef __FUNCT__
149c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal"
150dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N)
151c8a8475eSBarry Smith {
152c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
153dfbe8321SBarry Smith   PetscErrorCode ierr;
154c8a8475eSBarry Smith 
155c8a8475eSBarry Smith   PetscFunctionBegin;
156c3122656SLisandro Dalcin   if (Na->A) { ierr = MatDestroy(Na->A);CHKERRQ(ierr); }
157c3122656SLisandro Dalcin   if (Na->w) { ierr = VecDestroy(Na->w);CHKERRQ(ierr); }
158*7cfd0b8cSBarry Smith   if (Na->left) { ierr = VecDestroy(Na->left);CHKERRQ(ierr); }
159*7cfd0b8cSBarry Smith   if (Na->right) { ierr = VecDestroy(Na->right);CHKERRQ(ierr); }
160*7cfd0b8cSBarry Smith   if (Na->leftwork) { ierr = VecDestroy(Na->leftwork);CHKERRQ(ierr); }
161*7cfd0b8cSBarry Smith   if (Na->rightwork) { ierr = VecDestroy(Na->rightwork);CHKERRQ(ierr); }
162c8a8475eSBarry Smith   ierr = PetscFree(Na);CHKERRQ(ierr);
163c8a8475eSBarry Smith   PetscFunctionReturn(0);
164c8a8475eSBarry Smith }
165c8a8475eSBarry Smith 
16617a6fd94SBarry Smith /*
16717a6fd94SBarry Smith       Slow, nonscalable version
16817a6fd94SBarry Smith */
16917a6fd94SBarry Smith #undef __FUNCT__
17017a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal"
171dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
17217a6fd94SBarry Smith {
17317a6fd94SBarry Smith   Mat_Normal        *Na = (Mat_Normal*)N->data;
17417a6fd94SBarry Smith   Mat               A = Na->A;
175dfbe8321SBarry Smith   PetscErrorCode    ierr;
176b24ad042SBarry Smith   PetscInt          i,j,rstart,rend,nnz;
177b24ad042SBarry Smith   const PetscInt    *cols;
17817a6fd94SBarry Smith   PetscScalar       *diag,*work,*values;
179b3cc6726SBarry Smith   const PetscScalar *mvalues;
18017a6fd94SBarry Smith 
18117a6fd94SBarry Smith   PetscFunctionBegin;
182d0f46423SBarry Smith   ierr = PetscMalloc(2*A->cmap->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr);
183d0f46423SBarry Smith   work = diag + A->cmap->N;
184d0f46423SBarry Smith   ierr = PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
18517a6fd94SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
18617a6fd94SBarry Smith   for (i=rstart; i<rend; i++) {
187b3cc6726SBarry Smith     ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
18817a6fd94SBarry Smith     for (j=0; j<nnz; j++) {
189b3cc6726SBarry Smith       work[cols[j]] += mvalues[j]*mvalues[j];
19017a6fd94SBarry Smith     }
191b3cc6726SBarry Smith     ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
19217a6fd94SBarry Smith   }
193d0f46423SBarry Smith   ierr = MPI_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPI_SUM,((PetscObject)N)->comm);CHKERRQ(ierr);
194d0f46423SBarry Smith   rstart = N->cmap->rstart;
195d0f46423SBarry Smith   rend   = N->cmap->rend;
19617a6fd94SBarry Smith   ierr = VecGetArray(v,&values);CHKERRQ(ierr);
19717a6fd94SBarry Smith   ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr);
19817a6fd94SBarry Smith   ierr = VecRestoreArray(v,&values);CHKERRQ(ierr);
19917a6fd94SBarry Smith   ierr = PetscFree(diag);CHKERRQ(ierr);
200b20f02adSBarry Smith   ierr = VecScale(v,Na->scale);CHKERRQ(ierr);
20117a6fd94SBarry Smith   PetscFunctionReturn(0);
20217a6fd94SBarry Smith }
203c8a8475eSBarry Smith 
204c8a8475eSBarry Smith #undef __FUNCT__
205c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal"
206c8a8475eSBarry Smith /*@
207c8a8475eSBarry Smith       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
208c8a8475eSBarry Smith 
209c8a8475eSBarry Smith    Collective on Mat
210c8a8475eSBarry Smith 
211c8a8475eSBarry Smith    Input Parameter:
212c8a8475eSBarry Smith .   A  - the (possibly rectangular) matrix
213c8a8475eSBarry Smith 
214c8a8475eSBarry Smith    Output Parameter:
215c8a8475eSBarry Smith .   N - the matrix that represents A'*A
216c8a8475eSBarry Smith 
217c30e7cdbSSatish Balay    Level: intermediate
218c30e7cdbSSatish Balay 
219c8a8475eSBarry Smith    Notes: The product A'*A is NOT actually formed! Rather the new matrix
220c8a8475eSBarry Smith           object performs the matrix-vector product by first multiplying by
221c8a8475eSBarry Smith           A and then A'
222c8a8475eSBarry Smith @*/
223be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat A,Mat *N)
224c8a8475eSBarry Smith {
225dfbe8321SBarry Smith   PetscErrorCode ierr;
226b24ad042SBarry Smith   PetscInt       m,n;
227c8a8475eSBarry Smith   Mat_Normal     *Na;
228c8a8475eSBarry Smith 
229c8a8475eSBarry Smith   PetscFunctionBegin;
230c8a8475eSBarry Smith   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
2317adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,N);CHKERRQ(ierr);
232f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
233c8a8475eSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
234c8a8475eSBarry Smith 
23538f2d2fdSLisandro Dalcin   ierr      = PetscNewLog(*N,Mat_Normal,&Na);CHKERRQ(ierr);
23638f2d2fdSLisandro Dalcin   (*N)->data = (void*) Na;
237c8a8475eSBarry Smith   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
238c3122656SLisandro Dalcin   Na->A     = A;
239b20f02adSBarry Smith   Na->scale = 1.0;
240c8a8475eSBarry Smith 
2417adad957SLisandro Dalcin   ierr    = VecCreateMPI(((PetscObject)A)->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
242c8a8475eSBarry Smith   (*N)->ops->destroy          = MatDestroy_Normal;
243c8a8475eSBarry Smith   (*N)->ops->mult             = MatMult_Normal;
244b20f02adSBarry Smith   (*N)->ops->multtranspose    = MatMultTranspose_Normal;
245b20f02adSBarry Smith   (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal;
246db090513SMatthew Knepley   (*N)->ops->multadd          = MatMultAdd_Normal;
24717a6fd94SBarry Smith   (*N)->ops->getdiagonal      = MatGetDiagonal_Normal;
248c8a8475eSBarry Smith   (*N)->assembled             = PETSC_TRUE;
249d0f46423SBarry Smith   (*N)->cmap->N               = A->cmap->N;
250d0f46423SBarry Smith   (*N)->rmap->N               = A->cmap->N;
251d0f46423SBarry Smith   (*N)->cmap->n               = A->cmap->n;
252d0f46423SBarry Smith   (*N)->rmap->n               = A->cmap->n;
253c8a8475eSBarry Smith   PetscFunctionReturn(0);
254c8a8475eSBarry Smith }
255c8a8475eSBarry Smith 
256