xref: /petsc/src/mat/impls/normal/normm.c (revision 958c9bcce270cdc364c8832dcd56bbb1c4a38e24)
1 
2 #include "src/mat/matimpl.h"          /*I "petscmat.h" I*/
3 
4 typedef struct {
5   Mat A;
6   Vec w;
7 } Mat_Normal;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatMult_Normal"
11 int MatMult_Normal(Mat N,Vec x,Vec y)
12 {
13   Mat_Normal *Na = (Mat_Normal*)N->data;
14   int        ierr;
15 
16   PetscFunctionBegin;
17   ierr = MatMult(Na->A,x,Na->w);CHKERRQ(ierr);
18   ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr);
19   PetscFunctionReturn(0);
20 }
21 
22 #undef __FUNCT__
23 #define __FUNCT__ "MatDestroy_Normal"
24 int MatDestroy_Normal(Mat N)
25 {
26   Mat_Normal *Na = (Mat_Normal*)N->data;
27   int        ierr;
28 
29   PetscFunctionBegin;
30   ierr = PetscObjectDereference((PetscObject)Na->A);CHKERRQ(ierr);
31   ierr = VecDestroy(Na->w);CHKERRQ(ierr);
32   ierr = PetscFree(Na);CHKERRQ(ierr);
33   PetscFunctionReturn(0);
34 }
35 
36 /*
37       Slow, nonscalable version
38 */
39 #undef __FUNCT__
40 #define __FUNCT__ "MatGetDiagonal_Normal"
41 int MatGetDiagonal_Normal(Mat N,Vec v)
42 {
43   Mat_Normal        *Na = (Mat_Normal*)N->data;
44   Mat               A = Na->A;
45   int               ierr,i,j,rstart,rend,nnz;
46   const int         *cols;
47   PetscScalar       *diag,*work,*values;
48   const PetscScalar *mvalues;
49   PetscMap          cmap;
50 
51   PetscFunctionBegin;
52   ierr = PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr);
53   work = diag + A->N;
54   ierr = PetscMemzero(work,A->N*sizeof(PetscScalar));CHKERRQ(ierr);
55   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
56   for (i=rstart; i<rend; i++) {
57     ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
58     for (j=0; j<nnz; j++) {
59       work[cols[j]] += mvalues[j]*mvalues[j];
60     }
61     ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr);
62   }
63   ierr = MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr);
64   ierr = MatGetPetscMaps(A,PETSC_NULL,&cmap);CHKERRQ(ierr);
65   ierr = PetscMapGetLocalRange(cmap,&rstart,&rend);CHKERRQ(ierr);
66   ierr = VecGetArray(v,&values);CHKERRQ(ierr);
67   ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr);
68   ierr = VecRestoreArray(v,&values);CHKERRQ(ierr);
69   ierr = PetscFree(diag);CHKERRQ(ierr);
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "MatCreateNormal"
75 /*@
76       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
77 
78    Collective on Mat
79 
80    Input Parameter:
81 .   A  - the (possibly rectangular) matrix
82 
83    Output Parameter:
84 .   N - the matrix that represents A'*A
85 
86    Notes: The product A'*A is NOT actually formed! Rather the new matrix
87           object performs the matrix-vector product by first multiplying by
88           A and then A'
89 @*/
90 int MatCreateNormal(Mat A,Mat *N)
91 {
92   int        ierr,m,n;
93   Mat_Normal *Na;
94 
95   PetscFunctionBegin;
96   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
97   ierr = MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr);
98   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
99 
100   ierr      = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr);
101   Na->A     = A;
102   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
103   (*N)->data = (void*) Na;
104 
105   ierr    = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
106   (*N)->ops->destroy     = MatDestroy_Normal;
107   (*N)->ops->mult        = MatMult_Normal;
108   (*N)->ops->getdiagonal = MatGetDiagonal_Normal;
109   (*N)->assembled        = PETSC_TRUE;
110   (*N)->N                = A->N;
111   (*N)->M                = A->N;
112   (*N)->n                = A->n;
113   (*N)->m                = A->n;
114   PetscFunctionReturn(0);
115 }
116 
117