xref: /petsc/src/mat/impls/normal/normm.c (revision 9c30b7d2697335155d7490a7e085415ee7b4a02a)
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,*cols;
46   PetscScalar *diag,*work,*values;
47   PetscMap    cmap;
48 
49   PetscFunctionBegin;
50   ierr = PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr);
51   work = diag + A->N;
52   ierr = PetscMemzero(work,A->N*sizeof(PetscScalar));CHKERRQ(ierr);
53   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
54   for (i=rstart; i<rend; i++) {
55     ierr = MatGetRow(A,i,&nnz,&cols,&values);CHKERRQ(ierr);
56     for (j=0; j<nnz; j++) {
57       work[cols[j]] += values[j]*values[j];
58     }
59     ierr = MatRestoreRow(A,i,&nnz,&cols,&values);CHKERRQ(ierr);
60   }
61   ierr = MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr);
62   ierr = MatGetPetscMaps(A,PETSC_NULL,&cmap);CHKERRQ(ierr);
63   ierr = PetscMapGetLocalRange(cmap,&rstart,&rend);CHKERRQ(ierr);
64   ierr = VecGetArray(v,&values);CHKERRQ(ierr);
65   ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr);
66   ierr = VecRestoreArray(v,&values);CHKERRQ(ierr);
67   ierr = PetscFree(diag);CHKERRQ(ierr);
68   PetscFunctionReturn(0);
69 }
70 
71 #undef __FUNCT__
72 #define __FUNCT__ "MatCreateNormal"
73 /*@
74       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
75 
76    Collective on Mat
77 
78    Input Parameter:
79 .   A  - the (possibly rectangular) matrix
80 
81    Output Parameter:
82 .   N - the matrix that represents A'*A
83 
84    Notes: The product A'*A is NOT actually formed! Rather the new matrix
85           object performs the matrix-vector product by first multiplying by
86           A and then A'
87 @*/
88 int MatCreateNormal(Mat A,Mat *N)
89 {
90   int        ierr,m,n;
91   Mat_Normal *Na;
92 
93   PetscFunctionBegin;
94   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
95   ierr = MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr);
96   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
97 
98   ierr      = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr);
99   Na->A     = A;
100   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
101   (*N)->data = (void*) Na;
102 
103   ierr    = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
104   (*N)->ops->destroy     = MatDestroy_Normal;
105   (*N)->ops->mult        = MatMult_Normal;
106   (*N)->ops->getdiagonal = MatGetDiagonal_Normal;
107   (*N)->assembled        = PETSC_TRUE;
108   (*N)->N                = A->N;
109   (*N)->M                = A->N;
110   (*N)->n                = A->n;
111   (*N)->m                = A->n;
112   PetscFunctionReturn(0);
113 }
114 
115