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