xref: /petsc/src/mat/impls/normal/normm.c (revision b20f02adef4958dda258787f341a9f0049e7d5c9)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2c8a8475eSBarry Smith 
37c4f633dSBarry Smith #include "private/matimpl.h"          /*I "petscmat.h" I*/
4c8a8475eSBarry Smith 
5c8a8475eSBarry Smith typedef struct {
6*b20f02adSBarry Smith   Mat         A,left,right;   /* left and right scaling not yet implemented */
7c8a8475eSBarry Smith   Vec         w;
8*b20f02adSBarry Smith   PetscScalar scale;
9c8a8475eSBarry Smith } Mat_Normal;
10c8a8475eSBarry Smith 
11c8a8475eSBarry Smith #undef __FUNCT__
12c8a8475eSBarry Smith #define __FUNCT__ "MatMult_Normal"
13dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
14c8a8475eSBarry Smith {
15c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
16dfbe8321SBarry Smith   PetscErrorCode ierr;
17c8a8475eSBarry Smith 
18c8a8475eSBarry Smith   PetscFunctionBegin;
19c8a8475eSBarry Smith   ierr = MatMult(Na->A,x,Na->w);CHKERRQ(ierr);
20c8a8475eSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
21*b20f02adSBarry Smith   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
22c8a8475eSBarry Smith   PetscFunctionReturn(0);
23c8a8475eSBarry Smith }
24c8a8475eSBarry Smith 
25c8a8475eSBarry Smith #undef __FUNCT__
26db090513SMatthew Knepley #define __FUNCT__ "MatMultAdd_Normal"
27db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
28db090513SMatthew Knepley {
29db090513SMatthew Knepley   Mat_Normal     *Na = (Mat_Normal*)N->data;
30db090513SMatthew Knepley   PetscErrorCode ierr;
31db090513SMatthew Knepley 
32db090513SMatthew Knepley   PetscFunctionBegin;
33db090513SMatthew Knepley   ierr = MatMult(Na->A,v1,Na->w);CHKERRQ(ierr);
34*b20f02adSBarry Smith   ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr);
35*b20f02adSBarry Smith   ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr);
36*b20f02adSBarry Smith   PetscFunctionReturn(0);
37*b20f02adSBarry Smith }
38*b20f02adSBarry Smith 
39*b20f02adSBarry Smith #undef __FUNCT__
40*b20f02adSBarry Smith #define __FUNCT__ "MatMultTranspose_Normal"
41*b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y)
42*b20f02adSBarry Smith {
43*b20f02adSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
44*b20f02adSBarry Smith   PetscErrorCode ierr;
45*b20f02adSBarry Smith 
46*b20f02adSBarry Smith   PetscFunctionBegin;
47*b20f02adSBarry Smith   ierr = MatMult(Na->A,x,Na->w);CHKERRQ(ierr);
48*b20f02adSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
49*b20f02adSBarry Smith   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
50*b20f02adSBarry Smith   PetscFunctionReturn(0);
51*b20f02adSBarry Smith }
52*b20f02adSBarry Smith 
53*b20f02adSBarry Smith #undef __FUNCT__
54*b20f02adSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Normal"
55*b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
56*b20f02adSBarry Smith {
57*b20f02adSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
58*b20f02adSBarry Smith   PetscErrorCode ierr;
59*b20f02adSBarry Smith 
60*b20f02adSBarry Smith   PetscFunctionBegin;
61*b20f02adSBarry Smith   ierr = MatMult(Na->A,v1,Na->w);CHKERRQ(ierr);
62*b20f02adSBarry Smith   ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr);
63db090513SMatthew Knepley   ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr);
64db090513SMatthew Knepley   PetscFunctionReturn(0);
65db090513SMatthew Knepley }
66db090513SMatthew Knepley 
67db090513SMatthew Knepley #undef __FUNCT__
68c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal"
69dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N)
70c8a8475eSBarry Smith {
71c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
72dfbe8321SBarry Smith   PetscErrorCode ierr;
73c8a8475eSBarry Smith 
74c8a8475eSBarry Smith   PetscFunctionBegin;
75c3122656SLisandro Dalcin   if (Na->A) { ierr = MatDestroy(Na->A);CHKERRQ(ierr); }
76c3122656SLisandro Dalcin   if (Na->w) { ierr = VecDestroy(Na->w);CHKERRQ(ierr); }
77c8a8475eSBarry Smith   ierr = PetscFree(Na);CHKERRQ(ierr);
78c8a8475eSBarry Smith   PetscFunctionReturn(0);
79c8a8475eSBarry Smith }
80c8a8475eSBarry Smith 
8117a6fd94SBarry Smith /*
8217a6fd94SBarry Smith       Slow, nonscalable version
8317a6fd94SBarry Smith */
8417a6fd94SBarry Smith #undef __FUNCT__
8517a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal"
86dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
8717a6fd94SBarry Smith {
8817a6fd94SBarry Smith   Mat_Normal        *Na = (Mat_Normal*)N->data;
8917a6fd94SBarry Smith   Mat               A = Na->A;
90dfbe8321SBarry Smith   PetscErrorCode    ierr;
91b24ad042SBarry Smith   PetscInt          i,j,rstart,rend,nnz;
92b24ad042SBarry Smith   const PetscInt    *cols;
9317a6fd94SBarry Smith   PetscScalar       *diag,*work,*values;
94b3cc6726SBarry Smith   const PetscScalar *mvalues;
9517a6fd94SBarry Smith 
9617a6fd94SBarry Smith   PetscFunctionBegin;
97d0f46423SBarry Smith   ierr = PetscMalloc(2*A->cmap->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr);
98d0f46423SBarry Smith   work = diag + A->cmap->N;
99d0f46423SBarry Smith   ierr = PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
10017a6fd94SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
10117a6fd94SBarry Smith   for (i=rstart; i<rend; i++) {
102b3cc6726SBarry Smith     ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
10317a6fd94SBarry Smith     for (j=0; j<nnz; j++) {
104b3cc6726SBarry Smith       work[cols[j]] += mvalues[j]*mvalues[j];
10517a6fd94SBarry Smith     }
106b3cc6726SBarry Smith     ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
10717a6fd94SBarry Smith   }
108d0f46423SBarry Smith   ierr = MPI_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPI_SUM,((PetscObject)N)->comm);CHKERRQ(ierr);
109d0f46423SBarry Smith   rstart = N->cmap->rstart;
110d0f46423SBarry Smith   rend   = N->cmap->rend;
11117a6fd94SBarry Smith   ierr = VecGetArray(v,&values);CHKERRQ(ierr);
11217a6fd94SBarry Smith   ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr);
11317a6fd94SBarry Smith   ierr = VecRestoreArray(v,&values);CHKERRQ(ierr);
11417a6fd94SBarry Smith   ierr = PetscFree(diag);CHKERRQ(ierr);
115*b20f02adSBarry Smith   ierr = VecScale(v,Na->scale);CHKERRQ(ierr);
11617a6fd94SBarry Smith   PetscFunctionReturn(0);
11717a6fd94SBarry Smith }
118c8a8475eSBarry Smith 
119c8a8475eSBarry Smith #undef __FUNCT__
120c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal"
121c8a8475eSBarry Smith /*@
122c8a8475eSBarry Smith       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
123c8a8475eSBarry Smith 
124c8a8475eSBarry Smith    Collective on Mat
125c8a8475eSBarry Smith 
126c8a8475eSBarry Smith    Input Parameter:
127c8a8475eSBarry Smith .   A  - the (possibly rectangular) matrix
128c8a8475eSBarry Smith 
129c8a8475eSBarry Smith    Output Parameter:
130c8a8475eSBarry Smith .   N - the matrix that represents A'*A
131c8a8475eSBarry Smith 
132c30e7cdbSSatish Balay    Level: intermediate
133c30e7cdbSSatish Balay 
134c8a8475eSBarry Smith    Notes: The product A'*A is NOT actually formed! Rather the new matrix
135c8a8475eSBarry Smith           object performs the matrix-vector product by first multiplying by
136c8a8475eSBarry Smith           A and then A'
137c8a8475eSBarry Smith @*/
138be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat A,Mat *N)
139c8a8475eSBarry Smith {
140dfbe8321SBarry Smith   PetscErrorCode ierr;
141b24ad042SBarry Smith   PetscInt       m,n;
142c8a8475eSBarry Smith   Mat_Normal     *Na;
143c8a8475eSBarry Smith 
144c8a8475eSBarry Smith   PetscFunctionBegin;
145c8a8475eSBarry Smith   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1467adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,N);CHKERRQ(ierr);
147f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
148c8a8475eSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
149c8a8475eSBarry Smith 
15038f2d2fdSLisandro Dalcin   ierr      = PetscNewLog(*N,Mat_Normal,&Na);CHKERRQ(ierr);
15138f2d2fdSLisandro Dalcin   (*N)->data = (void*) Na;
152c8a8475eSBarry Smith   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
153c3122656SLisandro Dalcin   Na->A     = A;
154*b20f02adSBarry Smith   Na->scale = 1.0;
155c8a8475eSBarry Smith 
1567adad957SLisandro Dalcin   ierr    = VecCreateMPI(((PetscObject)A)->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
157c8a8475eSBarry Smith   (*N)->ops->destroy          = MatDestroy_Normal;
158c8a8475eSBarry Smith   (*N)->ops->mult             = MatMult_Normal;
159*b20f02adSBarry Smith   (*N)->ops->multtranspose    = MatMultTranspose_Normal;
160*b20f02adSBarry Smith   (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal;
161db090513SMatthew Knepley   (*N)->ops->multadd          = MatMultAdd_Normal;
16217a6fd94SBarry Smith   (*N)->ops->getdiagonal      = MatGetDiagonal_Normal;
163c8a8475eSBarry Smith   (*N)->assembled             = PETSC_TRUE;
164d0f46423SBarry Smith   (*N)->cmap->N               = A->cmap->N;
165d0f46423SBarry Smith   (*N)->rmap->N               = A->cmap->N;
166d0f46423SBarry Smith   (*N)->cmap->n               = A->cmap->n;
167d0f46423SBarry Smith   (*N)->rmap->n               = A->cmap->n;
168c8a8475eSBarry Smith   PetscFunctionReturn(0);
169c8a8475eSBarry Smith }
170c8a8475eSBarry Smith 
171