xref: /petsc/src/mat/impls/lrc/lrc.c (revision 0575051bb73012f2ba419badb39151cb6298121b)
1 
2 #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
3 
4 typedef struct {
5   Mat         A;           /* sparse matrix */
6   Mat         U,V;         /* dense tall-skinny matrices */
7   Vec         c;           /* sequential vector containing the diagonal of C */
8   Vec         work1,work2; /* sequential vectors that hold partial products */
9   PetscMPIInt nwork;       /* length of work vectors */
10   Vec         xl,yl;       /* auxiliary sequential vectors for matmult operation */
11 } Mat_LRC;
12 
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "MatMult_LRC"
16 PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y)
17 {
18   Mat_LRC           *Na = (Mat_LRC*)N->data;
19   PetscErrorCode    ierr;
20   PetscScalar       *w1,*w2;
21   const PetscScalar *a;
22 
23   PetscFunctionBegin;
24   ierr = VecGetArrayRead(x,&a);CHKERRQ(ierr);
25   ierr = VecPlaceArray(Na->xl,a);CHKERRQ(ierr);
26   ierr = VecGetLocalVector(y,Na->yl);CHKERRQ(ierr);
27 
28   /* multiply the local part of V with the local part of x */
29 #if defined(PETSC_USE_COMPLEX)
30   ierr = MatMultHermitianTranspose(Na->V,Na->xl,Na->work1);CHKERRQ(ierr);
31 #else
32   ierr = MatMultTranspose(Na->V,Na->xl,Na->work1);CHKERRQ(ierr);
33 #endif
34 
35   /* form the sum of all the local multiplies: this is work2 = V'*x =
36      sum_{all processors} work1 */
37   ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr);
38   ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr);
39   ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr);
40   ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr);
41   ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr);
42 
43   if (Na->c) {  /* work2 = C*work2 */
44     ierr = VecPointwiseMult(Na->work2,Na->c,Na->work2);CHKERRQ(ierr);
45   }
46 
47   if (Na->A) {
48     /* form y = A*x */
49     ierr = MatMult(Na->A,x,y);CHKERRQ(ierr);
50     /* multiply-add y = y + U*work2 */
51     ierr = MatMultAdd(Na->U,Na->work2,Na->yl,Na->yl);CHKERRQ(ierr);
52   } else {
53     /* multiply y = U*work2 */
54     ierr = MatMult(Na->U,Na->work2,Na->yl);CHKERRQ(ierr);
55   }
56 
57   ierr = VecRestoreArrayRead(x,&a);CHKERRQ(ierr);
58   ierr = VecResetArray(Na->xl);CHKERRQ(ierr);
59   ierr = VecRestoreLocalVector(y,Na->yl);CHKERRQ(ierr);
60   PetscFunctionReturn(0);
61 }
62 
63 #undef __FUNCT__
64 #define __FUNCT__ "MatDestroy_LRC"
65 PetscErrorCode MatDestroy_LRC(Mat N)
66 {
67   Mat_LRC        *Na = (Mat_LRC*)N->data;
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
72   ierr = MatDestroy(&Na->U);CHKERRQ(ierr);
73   ierr = MatDestroy(&Na->V);CHKERRQ(ierr);
74   ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
75   ierr = VecDestroy(&Na->work1);CHKERRQ(ierr);
76   ierr = VecDestroy(&Na->work2);CHKERRQ(ierr);
77   ierr = VecDestroy(&Na->xl);CHKERRQ(ierr);
78   ierr = VecDestroy(&Na->yl);CHKERRQ(ierr);
79   ierr = PetscFree(N->data);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "MatCreateLRC"
85 /*@
86    MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V'
87 
88    Collective on Mat
89 
90    Input Parameters:
91 +  A  - the (sparse) matrix (can be NULL)
92 .  U, V - two dense rectangular (tall and skinny) matrices
93 -  c  - a sequential vector containing the diagonal of C (can be NULL)
94 
95    Output Parameter:
96 .  N - the matrix that represents A + U*C*V'
97 
98    Notes:
99    The matrix A + U*C*V' is not formed! Rather the new matrix
100    object performs the matrix-vector product by first multiplying by
101    A and then adding the other term.
102 
103    C is a diagonal matrix (represented as a vector) of order k,
104    where k is the number of columns of both U and V.
105 
106    If A is NULL then the new object behaves like a low-rank matrix U*C*V'.
107 
108    Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'.
109 
110    If c is NULL then the low-rank correction is just U*V'.
111 
112    Level: intermediate
113 @*/
114 PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N)
115 {
116   PetscErrorCode ierr;
117   PetscBool      match;
118   PetscInt       m,n,k,m1,n1,k1;
119   Mat_LRC        *Na;
120 
121   PetscFunctionBegin;
122   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1);
123   PetscValidHeaderSpecific(U,MAT_CLASSID,2);
124   if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3);
125   if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,4);
126   else V=U;
127   if (A) PetscCheckSameComm(A,1,U,2);
128   PetscCheckSameComm(U,2,V,4);
129 
130   ierr = PetscObjectTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
131   if (!match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense");
132   if (V) {
133     ierr = PetscObjectTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
134     if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Matrix V must be of type dense");
135   }
136 
137   ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr);
138   ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr);
139   if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns");
140   ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr);
141   ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr);
142   if (A) {
143     ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr);
144     if (m!=m1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U and A do not match");
145     if (n!=n1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V and A do not match");
146   }
147   if (c) {
148     ierr = VecGetSize(c,&k1);CHKERRQ(ierr);
149     if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"The length of c does not match the number of columns of U and V");
150     ierr = VecGetLocalSize(c,&k1);CHKERRQ(ierr);
151     if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"c must be a sequential vector");
152   }
153 
154   ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr);
155   ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
156   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr);
157 
158   ierr       = PetscNewLog(*N,&Na);CHKERRQ(ierr);
159   (*N)->data = (void*)Na;
160   Na->A      = A;
161   Na->c      = c;
162 
163   ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr);
164   ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr);
165   if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); }
166   ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr);
167   ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr);
168   if (c) { ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); }
169 
170   ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr);
171   ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr);
172   ierr = PetscMPIIntCast(U->cmap->N,&Na->nwork);CHKERRQ(ierr);
173 
174   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,V->rmap->n,NULL,&Na->xl);CHKERRQ(ierr);
175   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,U->rmap->n,NULL,&Na->yl);CHKERRQ(ierr);
176 
177   (*N)->ops->destroy = MatDestroy_LRC;
178   (*N)->ops->mult    = MatMult_LRC;
179   (*N)->assembled    = PETSC_TRUE;
180   (*N)->preallocated = PETSC_TRUE;
181   (*N)->cmap->N      = V->rmap->N;
182   (*N)->rmap->N      = U->rmap->N;
183   (*N)->cmap->n      = V->rmap->n;
184   (*N)->rmap->n      = U->rmap->n;
185   PetscFunctionReturn(0);
186 }
187 
188