xref: /petsc/src/mat/impls/normal/normm.c (revision c8a8475e04bcaa43590892a5c3e60c6f87bc31f7)
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 #undef __FUNCT__
39 #define __FUNCT__ "MatCreateNormal"
40 /*@
41       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
42 
43    Collective on Mat
44 
45    Input Parameter:
46 .   A  - the (possibly rectangular) matrix
47 
48    Output Parameter:
49 .   N - the matrix that represents A'*A
50 
51    Notes: The product A'*A is NOT actually formed! Rather the new matrix
52           object performs the matrix-vector product by first multiplying by
53           A and then A'
54 @*/
55 int MatCreateNormal(Mat A,Mat *N)
56 {
57   int        ierr,m,n;
58   Mat_Normal *Na;
59 
60   PetscFunctionBegin;
61   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
62   ierr = MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr);
63   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr);
64 
65   ierr      = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr);
66   Na->A     = A;
67   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
68   (*N)->data = (void*) Na;
69 
70   ierr    = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr);
71   (*N)->ops->destroy = MatDestroy_Normal;
72   (*N)->ops->mult    = MatMult_Normal;
73   (*N)->assembled    = PETSC_TRUE;
74   (*N)->N            = A->N;
75   (*N)->M            = A->N;
76   (*N)->n            = A->n;
77   (*N)->m            = A->n;
78   PetscFunctionReturn(0);
79 }
80 
81