xref: /petsc/src/mat/impls/normal/normm.c (revision b24ad0428a60bb3d8d9be05f5b432634d28bad96)
1c8a8475eSBarry Smith 
2c8a8475eSBarry Smith #include "src/mat/matimpl.h"          /*I "petscmat.h" I*/
3c8a8475eSBarry Smith 
4c8a8475eSBarry Smith typedef struct {
5c8a8475eSBarry Smith   Mat A;
6c8a8475eSBarry Smith   Vec w;
7c8a8475eSBarry Smith } Mat_Normal;
8c8a8475eSBarry Smith 
9c8a8475eSBarry Smith #undef __FUNCT__
10c8a8475eSBarry Smith #define __FUNCT__ "MatMult_Normal"
11dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
12c8a8475eSBarry Smith {
13c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
14dfbe8321SBarry Smith   PetscErrorCode ierr;
15c8a8475eSBarry Smith 
16c8a8475eSBarry Smith   PetscFunctionBegin;
17c8a8475eSBarry Smith   ierr = MatMult(Na->A,x,Na->w);CHKERRQ(ierr);
18c8a8475eSBarry Smith   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
19c8a8475eSBarry Smith   PetscFunctionReturn(0);
20c8a8475eSBarry Smith }
21c8a8475eSBarry Smith 
22c8a8475eSBarry Smith #undef __FUNCT__
23c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal"
24dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N)
25c8a8475eSBarry Smith {
26c8a8475eSBarry Smith   Mat_Normal     *Na = (Mat_Normal*)N->data;
27dfbe8321SBarry Smith   PetscErrorCode ierr;
28c8a8475eSBarry Smith 
29c8a8475eSBarry Smith   PetscFunctionBegin;
30c8a8475eSBarry Smith   ierr = PetscObjectDereference((PetscObject)Na->A);CHKERRQ(ierr);
31c8a8475eSBarry Smith   ierr = VecDestroy(Na->w);CHKERRQ(ierr);
32c8a8475eSBarry Smith   ierr = PetscFree(Na);CHKERRQ(ierr);
33c8a8475eSBarry Smith   PetscFunctionReturn(0);
34c8a8475eSBarry Smith }
35c8a8475eSBarry Smith 
3617a6fd94SBarry Smith /*
3717a6fd94SBarry Smith       Slow, nonscalable version
3817a6fd94SBarry Smith */
3917a6fd94SBarry Smith #undef __FUNCT__
4017a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal"
41dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
4217a6fd94SBarry Smith {
4317a6fd94SBarry Smith   Mat_Normal        *Na = (Mat_Normal*)N->data;
4417a6fd94SBarry Smith   Mat               A = Na->A;
45dfbe8321SBarry Smith   PetscErrorCode    ierr;
46*b24ad042SBarry Smith   PetscInt          i,j,rstart,rend,nnz;
47*b24ad042SBarry Smith   const PetscInt    *cols;
4817a6fd94SBarry Smith   PetscScalar       *diag,*work,*values;
49b3cc6726SBarry Smith   const PetscScalar *mvalues;
5017a6fd94SBarry Smith   PetscMap          cmap;
5117a6fd94SBarry Smith 
5217a6fd94SBarry Smith   PetscFunctionBegin;
5317a6fd94SBarry Smith   ierr = PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr);
5417a6fd94SBarry Smith   work = diag + A->N;
5517a6fd94SBarry Smith   ierr = PetscMemzero(work,A->N*sizeof(PetscScalar));CHKERRQ(ierr);
5617a6fd94SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
5717a6fd94SBarry Smith   for (i=rstart; i<rend; i++) {
58b3cc6726SBarry Smith     ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
5917a6fd94SBarry Smith     for (j=0; j<nnz; j++) {
60b3cc6726SBarry Smith       work[cols[j]] += mvalues[j]*mvalues[j];
6117a6fd94SBarry Smith     }
62b3cc6726SBarry Smith     ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
6317a6fd94SBarry Smith   }
6417a6fd94SBarry Smith   ierr = MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr);
6517a6fd94SBarry Smith   ierr = MatGetPetscMaps(A,PETSC_NULL,&cmap);CHKERRQ(ierr);
6617a6fd94SBarry Smith   ierr = PetscMapGetLocalRange(cmap,&rstart,&rend);CHKERRQ(ierr);
6717a6fd94SBarry Smith   ierr = VecGetArray(v,&values);CHKERRQ(ierr);
6817a6fd94SBarry Smith   ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr);
6917a6fd94SBarry Smith   ierr = VecRestoreArray(v,&values);CHKERRQ(ierr);
7017a6fd94SBarry Smith   ierr = PetscFree(diag);CHKERRQ(ierr);
7117a6fd94SBarry Smith   PetscFunctionReturn(0);
7217a6fd94SBarry Smith }
73c8a8475eSBarry Smith 
74c8a8475eSBarry Smith #undef __FUNCT__
75c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal"
76c8a8475eSBarry Smith /*@
77c8a8475eSBarry Smith       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
78c8a8475eSBarry Smith 
79c8a8475eSBarry Smith    Collective on Mat
80c8a8475eSBarry Smith 
81c8a8475eSBarry Smith    Input Parameter:
82c8a8475eSBarry Smith .   A  - the (possibly rectangular) matrix
83c8a8475eSBarry Smith 
84c8a8475eSBarry Smith    Output Parameter:
85c8a8475eSBarry Smith .   N - the matrix that represents A'*A
86c8a8475eSBarry Smith 
87c8a8475eSBarry Smith    Notes: The product A'*A is NOT actually formed! Rather the new matrix
88c8a8475eSBarry Smith           object performs the matrix-vector product by first multiplying by
89c8a8475eSBarry Smith           A and then A'
90c8a8475eSBarry Smith @*/
91dfbe8321SBarry Smith PetscErrorCode MatCreateNormal(Mat A,Mat *N)
92c8a8475eSBarry Smith {
93dfbe8321SBarry Smith   PetscErrorCode ierr;
94*b24ad042SBarry Smith   PetscInt       m,n;
95c8a8475eSBarry Smith   Mat_Normal     *Na;
96c8a8475eSBarry Smith 
97c8a8475eSBarry Smith   PetscFunctionBegin;
98c8a8475eSBarry Smith   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
99c8a8475eSBarry Smith   ierr = MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr);
100c8a8475eSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
101c8a8475eSBarry Smith 
102c8a8475eSBarry Smith   ierr      = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr);
103c8a8475eSBarry Smith   Na->A     = A;
104c8a8475eSBarry Smith   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
105c8a8475eSBarry Smith   (*N)->data = (void*) Na;
106c8a8475eSBarry Smith 
107c8a8475eSBarry Smith   ierr    = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
108c8a8475eSBarry Smith   (*N)->ops->destroy     = MatDestroy_Normal;
109c8a8475eSBarry Smith   (*N)->ops->mult        = MatMult_Normal;
11017a6fd94SBarry Smith   (*N)->ops->getdiagonal = MatGetDiagonal_Normal;
111c8a8475eSBarry Smith   (*N)->assembled        = PETSC_TRUE;
112c8a8475eSBarry Smith   (*N)->N                = A->N;
113c8a8475eSBarry Smith   (*N)->M                = A->N;
114c8a8475eSBarry Smith   (*N)->n                = A->n;
115c8a8475eSBarry Smith   (*N)->m                = A->n;
116c8a8475eSBarry Smith   PetscFunctionReturn(0);
117c8a8475eSBarry Smith }
118c8a8475eSBarry Smith 
119