xref: /petsc/src/mat/impls/normal/normmh.c (revision a48507d3cf03db7c0f24e6f3ea654929f15734d0)
1ebda224cSfranklin5 
2ebda224cSfranklin5 #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
3ebda224cSfranklin5 
4ebda224cSfranklin5 typedef struct {
5ebda224cSfranklin5   Mat         A;
6*a48507d3SJose E. Roman   Mat         D; /* local submatrix for diagonal part */
7ebda224cSfranklin5   Vec         w,left,right,leftwork,rightwork;
8ebda224cSfranklin5   PetscScalar scale;
9ebda224cSfranklin5 } Mat_Normal;
10ebda224cSfranklin5 
11ebda224cSfranklin5 PetscErrorCode MatScaleHermitian_Normal(Mat inA,PetscScalar scale)
12ebda224cSfranklin5 {
13ebda224cSfranklin5   Mat_Normal *a = (Mat_Normal*)inA->data;
14ebda224cSfranklin5 
15ebda224cSfranklin5   PetscFunctionBegin;
16ebda224cSfranklin5   a->scale *= scale;
17ebda224cSfranklin5   PetscFunctionReturn(0);
18ebda224cSfranklin5 }
19ebda224cSfranklin5 
20ebda224cSfranklin5 PetscErrorCode MatDiagonalScaleHermitian_Normal(Mat inA,Vec left,Vec right)
21ebda224cSfranklin5 {
22ebda224cSfranklin5   Mat_Normal     *a = (Mat_Normal*)inA->data;
23ebda224cSfranklin5 
24ebda224cSfranklin5   PetscFunctionBegin;
25ebda224cSfranklin5   if (left) {
26ebda224cSfranklin5     if (!a->left) {
279566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left,&a->left));
289566063dSJacob Faibussowitsch       PetscCall(VecCopy(left,a->left));
29ebda224cSfranklin5     } else {
309566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left,left,a->left));
31ebda224cSfranklin5     }
32ebda224cSfranklin5   }
33ebda224cSfranklin5   if (right) {
34ebda224cSfranklin5     if (!a->right) {
359566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right,&a->right));
369566063dSJacob Faibussowitsch       PetscCall(VecCopy(right,a->right));
37ebda224cSfranklin5     } else {
389566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right,right,a->right));
39ebda224cSfranklin5     }
40ebda224cSfranklin5   }
41ebda224cSfranklin5   PetscFunctionReturn(0);
42ebda224cSfranklin5 }
43ebda224cSfranklin5 
44ebda224cSfranklin5 PetscErrorCode MatMultHermitian_Normal(Mat N,Vec x,Vec y)
45ebda224cSfranklin5 {
46ebda224cSfranklin5   Mat_Normal     *Na = (Mat_Normal*)N->data;
47ebda224cSfranklin5   Vec            in;
48ebda224cSfranklin5 
49ebda224cSfranklin5   PetscFunctionBegin;
50ebda224cSfranklin5   in = x;
51ebda224cSfranklin5   if (Na->right) {
52ebda224cSfranklin5     if (!Na->rightwork) {
539566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(Na->right,&Na->rightwork));
54ebda224cSfranklin5     }
559566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in));
56ebda224cSfranklin5     in   = Na->rightwork;
57ebda224cSfranklin5   }
589566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A,in,Na->w));
599566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y));
60ebda224cSfranklin5   if (Na->left) {
619566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(y,Na->left,y));
62ebda224cSfranklin5   }
639566063dSJacob Faibussowitsch   PetscCall(VecScale(y,Na->scale));
64ebda224cSfranklin5   PetscFunctionReturn(0);
65ebda224cSfranklin5 }
66ebda224cSfranklin5 
67ebda224cSfranklin5 PetscErrorCode MatMultHermitianAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
68ebda224cSfranklin5 {
69ebda224cSfranklin5   Mat_Normal     *Na = (Mat_Normal*)N->data;
70ebda224cSfranklin5   Vec            in;
71ebda224cSfranklin5 
72ebda224cSfranklin5   PetscFunctionBegin;
73ebda224cSfranklin5   in = v1;
74ebda224cSfranklin5   if (Na->right) {
75ebda224cSfranklin5     if (!Na->rightwork) {
769566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(Na->right,&Na->rightwork));
77ebda224cSfranklin5     }
789566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in));
79ebda224cSfranklin5     in   = Na->rightwork;
80ebda224cSfranklin5   }
819566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A,in,Na->w));
829566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w,Na->scale));
83ebda224cSfranklin5   if (Na->left) {
849566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3));
859566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3,Na->left,v3));
869566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3,1.0,v2));
87ebda224cSfranklin5   } else {
889566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3));
89ebda224cSfranklin5   }
90ebda224cSfranklin5   PetscFunctionReturn(0);
91ebda224cSfranklin5 }
92ebda224cSfranklin5 
93ebda224cSfranklin5 PetscErrorCode MatMultHermitianTranspose_Normal(Mat N,Vec x,Vec y)
94ebda224cSfranklin5 {
95ebda224cSfranklin5   Mat_Normal     *Na = (Mat_Normal*)N->data;
96ebda224cSfranklin5   Vec            in;
97ebda224cSfranklin5 
98ebda224cSfranklin5   PetscFunctionBegin;
99ebda224cSfranklin5   in = x;
100ebda224cSfranklin5   if (Na->left) {
101ebda224cSfranklin5     if (!Na->leftwork) {
1029566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(Na->left,&Na->leftwork));
103ebda224cSfranklin5     }
1049566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in));
105ebda224cSfranklin5     in   = Na->leftwork;
106ebda224cSfranklin5   }
1079566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A,in,Na->w));
1089566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y));
109ebda224cSfranklin5   if (Na->right) {
1109566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(y,Na->right,y));
111ebda224cSfranklin5   }
1129566063dSJacob Faibussowitsch   PetscCall(VecScale(y,Na->scale));
113ebda224cSfranklin5   PetscFunctionReturn(0);
114ebda224cSfranklin5 }
115ebda224cSfranklin5 
116ebda224cSfranklin5 PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
117ebda224cSfranklin5 {
118ebda224cSfranklin5   Mat_Normal     *Na = (Mat_Normal*)N->data;
119ebda224cSfranklin5   Vec            in;
120ebda224cSfranklin5 
121ebda224cSfranklin5   PetscFunctionBegin;
122ebda224cSfranklin5   in = v1;
123ebda224cSfranklin5   if (Na->left) {
124ebda224cSfranklin5     if (!Na->leftwork) {
1259566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(Na->left,&Na->leftwork));
126ebda224cSfranklin5     }
1279566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in));
128ebda224cSfranklin5     in   = Na->leftwork;
129ebda224cSfranklin5   }
1309566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A,in,Na->w));
1319566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w,Na->scale));
132ebda224cSfranklin5   if (Na->right) {
1339566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3));
1349566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3,Na->right,v3));
1359566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3,1.0,v2));
136ebda224cSfranklin5   } else {
1379566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3));
138ebda224cSfranklin5   }
139ebda224cSfranklin5   PetscFunctionReturn(0);
140ebda224cSfranklin5 }
141ebda224cSfranklin5 
142ebda224cSfranklin5 PetscErrorCode MatDestroyHermitian_Normal(Mat N)
143ebda224cSfranklin5 {
144ebda224cSfranklin5   Mat_Normal     *Na = (Mat_Normal*)N->data;
145ebda224cSfranklin5 
146ebda224cSfranklin5   PetscFunctionBegin;
1479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
148*a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
1499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
1509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->left));
1519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->right));
1529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->leftwork));
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rightwork));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->data));
1559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMatHermitian_C",NULL));
156ebda224cSfranklin5   PetscFunctionReturn(0);
157ebda224cSfranklin5 }
158ebda224cSfranklin5 
159ebda224cSfranklin5 /*
160ebda224cSfranklin5       Slow, nonscalable version
161ebda224cSfranklin5 */
162ebda224cSfranklin5 PetscErrorCode MatGetDiagonalHermitian_Normal(Mat N,Vec v)
163ebda224cSfranklin5 {
164ebda224cSfranklin5   Mat_Normal        *Na = (Mat_Normal*)N->data;
165ebda224cSfranklin5   Mat               A   = Na->A;
166ebda224cSfranklin5   PetscInt          i,j,rstart,rend,nnz;
167ebda224cSfranklin5   const PetscInt    *cols;
168ebda224cSfranklin5   PetscScalar       *diag,*work,*values;
169ebda224cSfranklin5   const PetscScalar *mvalues;
170ebda224cSfranklin5 
171ebda224cSfranklin5   PetscFunctionBegin;
1729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work));
1739566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work,A->cmap->N));
1749566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
175ebda224cSfranklin5   for (i=rstart; i<rend; i++) {
1769566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,i,&nnz,&cols,&mvalues));
177ebda224cSfranklin5     for (j=0; j<nnz; j++) {
178a8f6171dSBarry Smith       work[cols[j]] += mvalues[j]*PetscConj(mvalues[j]);
179ebda224cSfranklin5     }
1809566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,i,&nnz,&cols,&mvalues));
181ebda224cSfranklin5   }
1821c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N)));
183ebda224cSfranklin5   rstart = N->cmap->rstart;
184ebda224cSfranklin5   rend   = N->cmap->rend;
1859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v,&values));
1869566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values,diag+rstart,rend-rstart));
1879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v,&values));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag,work));
1899566063dSJacob Faibussowitsch   PetscCall(VecScale(v,Na->scale));
190ebda224cSfranklin5   PetscFunctionReturn(0);
191ebda224cSfranklin5 }
192ebda224cSfranklin5 
193*a48507d3SJose E. Roman PetscErrorCode MatGetDiagonalBlockHermitian_Normal(Mat N,Mat *D)
194*a48507d3SJose E. Roman {
195*a48507d3SJose E. Roman   Mat_Normal *Na = (Mat_Normal*)N->data;
196*a48507d3SJose E. Roman   Mat        M,A = Na->A;
197*a48507d3SJose E. Roman 
198*a48507d3SJose E. Roman   PetscFunctionBegin;
199*a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A,&M));
200*a48507d3SJose E. Roman   PetscCall(MatCreateNormalHermitian(M,&Na->D));
201*a48507d3SJose E. Roman   *D = Na->D;
202*a48507d3SJose E. Roman   PetscFunctionReturn(0);
203*a48507d3SJose E. Roman }
204*a48507d3SJose E. Roman 
20565f45395SPierre Jolivet PetscErrorCode MatNormalGetMatHermitian_Normal(Mat A,Mat *M)
20665f45395SPierre Jolivet {
20765f45395SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal*)A->data;
20865f45395SPierre Jolivet 
20965f45395SPierre Jolivet   PetscFunctionBegin;
21065f45395SPierre Jolivet   *M = Aa->A;
21165f45395SPierre Jolivet   PetscFunctionReturn(0);
21265f45395SPierre Jolivet }
21365f45395SPierre Jolivet 
21465f45395SPierre Jolivet /*@
21565f45395SPierre Jolivet       MatNormalHermitianGetMat - Gets the Mat object stored inside a MATNORMALHERMITIAN
21665f45395SPierre Jolivet 
21765f45395SPierre Jolivet    Logically collective on Mat
21865f45395SPierre Jolivet 
21965f45395SPierre Jolivet    Input Parameter:
22065f45395SPierre Jolivet .   A  - the MATNORMALHERMITIAN matrix
22165f45395SPierre Jolivet 
22265f45395SPierre Jolivet    Output Parameter:
22365f45395SPierre Jolivet .   M - the matrix object stored inside A
22465f45395SPierre Jolivet 
22565f45395SPierre Jolivet    Level: intermediate
22665f45395SPierre Jolivet 
227db781477SPatrick Sanan .seealso: `MatCreateNormalHermitian()`
22865f45395SPierre Jolivet 
22965f45395SPierre Jolivet @*/
23065f45395SPierre Jolivet PetscErrorCode MatNormalHermitianGetMat(Mat A,Mat *M)
23165f45395SPierre Jolivet {
23265f45395SPierre Jolivet   PetscFunctionBegin;
23365f45395SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
23465f45395SPierre Jolivet   PetscValidType(A,1);
23565f45395SPierre Jolivet   PetscValidPointer(M,2);
236cac4c232SBarry Smith   PetscUseMethod(A,"MatNormalGetMatHermitian_C",(Mat,Mat*),(A,M));
23765f45395SPierre Jolivet   PetscFunctionReturn(0);
23865f45395SPierre Jolivet }
23965f45395SPierre Jolivet 
240ebda224cSfranklin5 /*@
241ebda224cSfranklin5       MatCreateNormalHermitian - Creates a new matrix object that behaves like (A*)'*A.
242ebda224cSfranklin5 
243ebda224cSfranklin5    Collective on Mat
244ebda224cSfranklin5 
245ebda224cSfranklin5    Input Parameter:
246ebda224cSfranklin5 .   A  - the (possibly rectangular complex) matrix
247ebda224cSfranklin5 
248ebda224cSfranklin5    Output Parameter:
249ebda224cSfranklin5 .   N - the matrix that represents (A*)'*A
250ebda224cSfranklin5 
251ebda224cSfranklin5    Level: intermediate
252ebda224cSfranklin5 
25395452b02SPatrick Sanan    Notes:
25495452b02SPatrick Sanan     The product (A*)'*A is NOT actually formed! Rather the new matrix
255ebda224cSfranklin5           object performs the matrix-vector product by first multiplying by
256ebda224cSfranklin5           A and then (A*)'
257ebda224cSfranklin5 @*/
258ebda224cSfranklin5 PetscErrorCode  MatCreateNormalHermitian(Mat A,Mat *N)
259ebda224cSfranklin5 {
260ebda224cSfranklin5   PetscInt       m,n;
261ebda224cSfranklin5   Mat_Normal     *Na;
2625cb049f5SJose E. Roman   VecType        vtype;
263ebda224cSfranklin5 
264ebda224cSfranklin5   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
2669566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),N));
2679566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE));
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*N,MATNORMALHERMITIAN));
2699566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap,&(*N)->rmap));
2709566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap,&(*N)->cmap));
271ebda224cSfranklin5 
2729566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(*N,&Na));
273ebda224cSfranklin5   (*N)->data = (void*) Na;
2749566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
275ebda224cSfranklin5   Na->A      = A;
276ebda224cSfranklin5   Na->scale  = 1.0;
277ebda224cSfranklin5 
2789566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,NULL,&Na->w));
279ebda224cSfranklin5 
280ebda224cSfranklin5   (*N)->ops->destroy          = MatDestroyHermitian_Normal;
281ebda224cSfranklin5   (*N)->ops->mult             = MatMultHermitian_Normal;
282ebda224cSfranklin5   (*N)->ops->multtranspose    = MatMultHermitianTranspose_Normal;
283ebda224cSfranklin5   (*N)->ops->multtransposeadd = MatMultHermitianTransposeAdd_Normal;
284ebda224cSfranklin5   (*N)->ops->multadd          = MatMultHermitianAdd_Normal;
285ebda224cSfranklin5   (*N)->ops->getdiagonal      = MatGetDiagonalHermitian_Normal;
286*a48507d3SJose E. Roman   (*N)->ops->getdiagonalblock = MatGetDiagonalBlockHermitian_Normal;
287ebda224cSfranklin5   (*N)->ops->scale            = MatScaleHermitian_Normal;
288ebda224cSfranklin5   (*N)->ops->diagonalscale    = MatDiagonalScaleHermitian_Normal;
289ebda224cSfranklin5   (*N)->assembled             = PETSC_TRUE;
2905cb049f5SJose E. Roman   (*N)->preallocated          = PETSC_TRUE;
29112ab1b64SPierre Jolivet 
2929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMatHermitian_C",MatNormalGetMatHermitian_Normal));
2939566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N,MAT_HERMITIAN,PETSC_TRUE));
2949566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A,&vtype));
2959566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N,vtype));
2965cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE)
2979566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N,A->boundtocpu));
2985cb049f5SJose E. Roman #endif
299ebda224cSfranklin5   PetscFunctionReturn(0);
300ebda224cSfranklin5 }
301