xref: /petsc/src/mat/impls/normal/normm.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
127d8b95cSPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/
2c8a8475eSBarry Smith 
3c8a8475eSBarry Smith typedef struct {
47cfd0b8cSBarry Smith   Mat A;
5a48507d3SJose E. Roman   Mat D; /* local submatrix for diagonal part */
627d8b95cSPierre Jolivet   Vec w;
7c8a8475eSBarry Smith } Mat_Normal;
8c8a8475eSBarry Smith 
MatIncreaseOverlap_Normal(Mat A,PetscInt is_max,IS is[],PetscInt ov)966976f2fSJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Normal(Mat A, PetscInt is_max, IS is[], PetscInt ov)
10d71ae5a4SJacob Faibussowitsch {
1127d8b95cSPierre Jolivet   Mat_Normal *a;
12fa80d070SPierre Jolivet   Mat         pattern;
13fa80d070SPierre Jolivet 
14fa80d070SPierre Jolivet   PetscFunctionBegin;
1508401ef6SPierre Jolivet   PetscCheck(ov >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified");
1627d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
179566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, a->A, NULL, &pattern));
189566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(pattern, MATPRODUCT_AtB));
199566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(pattern));
209566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(pattern));
219566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(pattern, is_max, is, ov));
229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pattern));
233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24fa80d070SPierre Jolivet }
25fa80d070SPierre Jolivet 
MatCreateSubMatrices_Normal(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat * submat[])2666976f2fSJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Normal(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
27d71ae5a4SJacob Faibussowitsch {
2827d8b95cSPierre Jolivet   Mat_Normal *a;
2927d8b95cSPierre Jolivet   Mat         B, *suba;
30fa80d070SPierre Jolivet   IS         *row;
31b9c875b8SPierre Jolivet   PetscScalar shift, scale;
32fa80d070SPierre Jolivet   PetscInt    M;
33fa80d070SPierre Jolivet 
34fa80d070SPierre Jolivet   PetscFunctionBegin;
3527d8b95cSPierre Jolivet   PetscCheck(irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
36b9c875b8SPierre Jolivet   PetscCall(MatShellGetScalingShifts(mat, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
3727d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(mat, &a));
3827d8b95cSPierre Jolivet   B = a->A;
3948a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
409566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &row));
429566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
439566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row[0]));
44fa80d070SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
459566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
46fa80d070SPierre Jolivet   for (M = 0; M < n; ++M) {
479566063dSJacob Faibussowitsch     PetscCall(MatCreateNormal(suba[M], *submat + M));
48b9c875b8SPierre Jolivet     PetscCall(MatShift((*submat)[M], shift));
49b9c875b8SPierre Jolivet     PetscCall(MatScale((*submat)[M], scale));
50fa80d070SPierre Jolivet   }
519566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row[0]));
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(row));
539566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(n, &suba));
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55fa80d070SPierre Jolivet }
56fa80d070SPierre Jolivet 
MatPermute_Normal(Mat A,IS rowp,IS colp,Mat * B)5766976f2fSJacob Faibussowitsch static PetscErrorCode MatPermute_Normal(Mat A, IS rowp, IS colp, Mat *B)
58d71ae5a4SJacob Faibussowitsch {
5927d8b95cSPierre Jolivet   Mat_Normal *a;
6027d8b95cSPierre Jolivet   Mat         C, Aa;
61fa80d070SPierre Jolivet   IS          row;
62b9c875b8SPierre Jolivet   PetscScalar shift, scale;
63fa80d070SPierre Jolivet 
64fa80d070SPierre Jolivet   PetscFunctionBegin;
6508401ef6SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
66b9c875b8SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
6727d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
6827d8b95cSPierre Jolivet   Aa = a->A;
699566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
709566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row));
719566063dSJacob Faibussowitsch   PetscCall(MatPermute(Aa, row, colp, &C));
729566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row));
739566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
75b9c875b8SPierre Jolivet   PetscCall(MatShift(*B, shift));
76b9c875b8SPierre Jolivet   PetscCall(MatScale(*B, scale));
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78fa80d070SPierre Jolivet }
79fa80d070SPierre Jolivet 
MatDuplicate_Normal(Mat A,MatDuplicateOption op,Mat * B)8066976f2fSJacob Faibussowitsch static PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B)
81d71ae5a4SJacob Faibussowitsch {
8227d8b95cSPierre Jolivet   Mat_Normal *a;
83fa80d070SPierre Jolivet   Mat         C;
84fa80d070SPierre Jolivet 
85fa80d070SPierre Jolivet   PetscFunctionBegin;
8627d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
879566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(a->A, op, &C));
889566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
9027d8b95cSPierre Jolivet   if (op == MAT_COPY_VALUES) PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92fa80d070SPierre Jolivet }
93fa80d070SPierre Jolivet 
MatCopy_Normal(Mat A,Mat B,MatStructure str)9466976f2fSJacob Faibussowitsch static PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str)
95d71ae5a4SJacob Faibussowitsch {
9627d8b95cSPierre Jolivet   Mat_Normal *a, *b;
97fa80d070SPierre Jolivet 
98fa80d070SPierre Jolivet   PetscFunctionBegin;
9927d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
10027d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(B, &b));
1019566063dSJacob Faibussowitsch   PetscCall(MatCopy(a->A, b->A, str));
1023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103fa80d070SPierre Jolivet }
104fa80d070SPierre Jolivet 
MatMult_Normal(Mat N,Vec x,Vec y)10566976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y)
106d71ae5a4SJacob Faibussowitsch {
10727d8b95cSPierre Jolivet   Mat_Normal *Na;
108c8a8475eSBarry Smith 
109c8a8475eSBarry Smith   PetscFunctionBegin;
11027d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
11127d8b95cSPierre Jolivet   PetscCall(MatMult(Na->A, x, Na->w));
1129566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A, Na->w, y));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114db090513SMatthew Knepley }
115db090513SMatthew Knepley 
MatDestroy_Normal(Mat N)11666976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_Normal(Mat N)
117d71ae5a4SJacob Faibussowitsch {
11827d8b95cSPierre Jolivet   Mat_Normal *Na;
119c8a8475eSBarry Smith 
120c8a8475eSBarry Smith   PetscFunctionBegin;
12127d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
123a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
1249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
12527d8b95cSPierre Jolivet   PetscCall(PetscFree(Na));
1269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL));
1279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_seqaij_C", NULL));
1289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_mpiaij_C", NULL));
1292f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
1302f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_hypre_C", NULL));
1312f0a5c5aSJose E. Roman #endif
1329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_seqdense_C", NULL));
1339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_mpidense_C", NULL));
13427d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatShellSetContext_C", NULL));
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136c8a8475eSBarry Smith }
137c8a8475eSBarry Smith 
13817a6fd94SBarry Smith /*
13917a6fd94SBarry Smith       Slow, nonscalable version
14017a6fd94SBarry Smith */
MatGetDiagonal_Normal(Mat N,Vec v)14166976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Normal(Mat N, Vec v)
142d71ae5a4SJacob Faibussowitsch {
14327d8b95cSPierre Jolivet   Mat_Normal        *Na;
14427d8b95cSPierre Jolivet   Mat                A;
145b24ad042SBarry Smith   PetscInt           i, j, rstart, rend, nnz;
146b24ad042SBarry Smith   const PetscInt    *cols;
147e91c04dfSPierre Jolivet   PetscScalar       *work, *values;
148b3cc6726SBarry Smith   const PetscScalar *mvalues;
14917a6fd94SBarry Smith 
15017a6fd94SBarry Smith   PetscFunctionBegin;
15127d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
15227d8b95cSPierre Jolivet   A = Na->A;
153e91c04dfSPierre Jolivet   PetscCall(PetscMalloc1(A->cmap->N, &work));
1549566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
1559566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
15617a6fd94SBarry Smith   for (i = rstart; i < rend; i++) {
1579566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
158ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * mvalues[j];
1599566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
16017a6fd94SBarry Smith   }
161e91c04dfSPierre Jolivet   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, work, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
162d0f46423SBarry Smith   rstart = N->cmap->rstart;
163d0f46423SBarry Smith   rend   = N->cmap->rend;
1649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
165e91c04dfSPierre Jolivet   PetscCall(PetscArraycpy(values, work + rstart, rend - rstart));
1669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
167e91c04dfSPierre Jolivet   PetscCall(PetscFree(work));
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16917a6fd94SBarry Smith }
170c8a8475eSBarry Smith 
MatGetDiagonalBlock_Normal(Mat N,Mat * D)17166976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D)
172d71ae5a4SJacob Faibussowitsch {
17327d8b95cSPierre Jolivet   Mat_Normal *Na;
17427d8b95cSPierre Jolivet   Mat         M, A;
175a48507d3SJose E. Roman 
176a48507d3SJose E. Roman   PetscFunctionBegin;
17727d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
17827d8b95cSPierre Jolivet   A = Na->A;
179a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
180a48507d3SJose E. Roman   PetscCall(MatCreateNormal(M, &Na->D));
181a48507d3SJose E. Roman   *D = Na->D;
1823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
183a48507d3SJose E. Roman }
184a48507d3SJose E. Roman 
MatNormalGetMat_Normal(Mat A,Mat * M)18566976f2fSJacob Faibussowitsch static PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M)
186d71ae5a4SJacob Faibussowitsch {
18727d8b95cSPierre Jolivet   Mat_Normal *Aa;
188fa80d070SPierre Jolivet 
189fa80d070SPierre Jolivet   PetscFunctionBegin;
19027d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &Aa));
191fa80d070SPierre Jolivet   *M = Aa->A;
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
193fa80d070SPierre Jolivet }
194fa80d070SPierre Jolivet 
195fa80d070SPierre Jolivet /*@
19611a5261eSBarry Smith   MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL`
197fa80d070SPierre Jolivet 
1982ef1f0ffSBarry Smith   Logically Collective
199fa80d070SPierre Jolivet 
200fa80d070SPierre Jolivet   Input Parameter:
20111a5261eSBarry Smith . A - the `MATNORMAL` matrix
202fa80d070SPierre Jolivet 
203fa80d070SPierre Jolivet   Output Parameter:
2042ef1f0ffSBarry Smith . M - the matrix object stored inside `A`
205fa80d070SPierre Jolivet 
206fa80d070SPierre Jolivet   Level: intermediate
207fa80d070SPierre Jolivet 
2081cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()`
209fa80d070SPierre Jolivet @*/
MatNormalGetMat(Mat A,Mat * M)210d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat(Mat A, Mat *M)
211d71ae5a4SJacob Faibussowitsch {
212fa80d070SPierre Jolivet   PetscFunctionBegin;
213fa80d070SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
214fa80d070SPierre Jolivet   PetscValidType(A, 1);
2154f572ea9SToby Isaac   PetscAssertPointer(M, 2);
216cac4c232SBarry Smith   PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M));
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218fa80d070SPierre Jolivet }
219fa80d070SPierre Jolivet 
MatConvert_Normal_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)22066976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
221d71ae5a4SJacob Faibussowitsch {
22227d8b95cSPierre Jolivet   Mat_Normal *Aa;
223fa80d070SPierre Jolivet   Mat         B;
224353cda34SPierre Jolivet   Vec         left, right, dshift;
225353cda34SPierre Jolivet   PetscScalar scale, shift;
226fa80d070SPierre Jolivet   PetscInt    m, n, M, N;
227fa80d070SPierre Jolivet 
228fa80d070SPierre Jolivet   PetscFunctionBegin;
22927d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &Aa));
230353cda34SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, &dshift, &left, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
2319566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
2329566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
233fa80d070SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
234fa80d070SPierre Jolivet     B = *newmat;
2359566063dSJacob Faibussowitsch     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
236fa80d070SPierre Jolivet   } else {
2379566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
2389566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
2399566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
2409566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
2419566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
242fa80d070SPierre Jolivet   }
2439566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(B));
244ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
245ac530a7eSPierre Jolivet   else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
2469566063dSJacob Faibussowitsch   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
247353cda34SPierre Jolivet   PetscCall(MatDiagonalScale(*newmat, left, right));
248353cda34SPierre Jolivet   PetscCall(MatScale(*newmat, scale));
249353cda34SPierre Jolivet   PetscCall(MatShift(*newmat, shift));
250353cda34SPierre Jolivet   if (dshift) PetscCall(MatDiagonalSet(*newmat, dshift, ADD_VALUES));
2513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
252fa80d070SPierre Jolivet }
253fa80d070SPierre Jolivet 
2542f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
MatConvert_Normal_HYPRE(Mat A,MatType type,MatReuse reuse,Mat * B)25566976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
2562f0a5c5aSJose E. Roman {
2572f0a5c5aSJose E. Roman   PetscFunctionBegin;
2582f0a5c5aSJose E. Roman   if (reuse == MAT_INITIAL_MATRIX) {
2592f0a5c5aSJose E. Roman     PetscCall(MatConvert(A, MATAIJ, reuse, B));
2602f0a5c5aSJose E. Roman     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
2612f0a5c5aSJose E. Roman   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2632f0a5c5aSJose E. Roman }
2642f0a5c5aSJose E. Roman #endif
2652f0a5c5aSJose E. Roman 
266fa80d070SPierre Jolivet typedef struct {
267fa80d070SPierre Jolivet   Mat work[2];
268fa80d070SPierre Jolivet } Normal_Dense;
269fa80d070SPierre Jolivet 
MatProductNumeric_Normal_Dense(Mat C)27066976f2fSJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Normal_Dense(Mat C)
271d71ae5a4SJacob Faibussowitsch {
272fa80d070SPierre Jolivet   Mat           A, B;
273fa80d070SPierre Jolivet   Normal_Dense *contents;
274fa80d070SPierre Jolivet   Mat_Normal   *a;
275b9c875b8SPierre Jolivet   Vec           right;
276b9c875b8SPierre Jolivet   PetscScalar  *array, scale;
277fa80d070SPierre Jolivet 
278fa80d070SPierre Jolivet   PetscFunctionBegin;
2798f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
280fa80d070SPierre Jolivet   A = C->product->A;
281fa80d070SPierre Jolivet   B = C->product->B;
28227d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
283fa80d070SPierre Jolivet   contents = (Normal_Dense *)C->product->data;
28428b400f6SJacob Faibussowitsch   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
28505694d27SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
286b9c875b8SPierre Jolivet   if (right) {
2879566063dSJacob Faibussowitsch     PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN));
288b9c875b8SPierre Jolivet     PetscCall(MatDiagonalScale(C, right, NULL));
289fa80d070SPierre Jolivet   }
2909566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[0]));
2919566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
2929566063dSJacob Faibussowitsch   PetscCall(MatDensePlaceArray(contents->work[1], array));
2939566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[1]));
2949566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
2959566063dSJacob Faibussowitsch   PetscCall(MatDenseResetArray(contents->work[1]));
2969566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
2979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
299b9c875b8SPierre Jolivet   PetscCall(MatScale(C, scale));
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301fa80d070SPierre Jolivet }
302fa80d070SPierre Jolivet 
MatNormal_DenseDestroy(PetscCtxRt ctx)303*2a8381b2SBarry Smith static PetscErrorCode MatNormal_DenseDestroy(PetscCtxRt ctx)
304d71ae5a4SJacob Faibussowitsch {
305cc1eb50dSBarry Smith   Normal_Dense *contents = *(Normal_Dense **)ctx;
306fa80d070SPierre Jolivet 
307fa80d070SPierre Jolivet   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work));
3099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work + 1));
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(contents));
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
312fa80d070SPierre Jolivet }
313fa80d070SPierre Jolivet 
MatProductSymbolic_Normal_Dense(Mat C)31466976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C)
315d71ae5a4SJacob Faibussowitsch {
316fa80d070SPierre Jolivet   Mat           A, B;
317fa80d070SPierre Jolivet   Normal_Dense *contents = NULL;
318fa80d070SPierre Jolivet   Mat_Normal   *a;
319b9c875b8SPierre Jolivet   Vec           right;
320b9c875b8SPierre Jolivet   PetscScalar  *array, scale;
321fa80d070SPierre Jolivet   PetscInt      n, N, m, M;
322fa80d070SPierre Jolivet 
323fa80d070SPierre Jolivet   PetscFunctionBegin;
3248f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
32528b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
326fa80d070SPierre Jolivet   A = C->product->A;
327fa80d070SPierre Jolivet   B = C->product->B;
32805694d27SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
32927d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
3309566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &m, &n));
3319566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C, &M, &N));
332fa80d070SPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
3339566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(B, NULL, &n));
3349566063dSJacob Faibussowitsch     PetscCall(MatGetSize(B, NULL, &N));
3359566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
3369566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, NULL));
3379566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, m, n, M, N));
338fa80d070SPierre Jolivet   }
3399566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
3409566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
3419566063dSJacob Faibussowitsch   PetscCall(PetscNew(&contents));
342fa80d070SPierre Jolivet   C->product->data    = contents;
343fa80d070SPierre Jolivet   C->product->destroy = MatNormal_DenseDestroy;
344b9c875b8SPierre Jolivet   if (right) PetscCall(MatProductCreate(a->A, C, NULL, contents->work));
345b9c875b8SPierre Jolivet   else PetscCall(MatProductCreate(a->A, B, NULL, contents->work));
3469566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB));
3479566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[0]));
3489566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[0]));
3499566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1));
3509566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB));
3519566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[1]));
3529566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[1]));
3539566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
3549566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array));
3559566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array));
3569566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
357fa80d070SPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
359fa80d070SPierre Jolivet }
360fa80d070SPierre Jolivet 
MatProductSetFromOptions_Normal_Dense(Mat C)36166976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C)
362d71ae5a4SJacob Faibussowitsch {
363fa80d070SPierre Jolivet   Mat_Product *product = C->product;
364fa80d070SPierre Jolivet 
365fa80d070SPierre Jolivet   PetscFunctionBegin;
366794904c8SPierre Jolivet   if (product->type == MATPRODUCT_AB) C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
3673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
368fa80d070SPierre Jolivet }
369fa80d070SPierre Jolivet 
3702ef1f0ffSBarry Smith /*MC
3712ef1f0ffSBarry Smith   MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A
3722ef1f0ffSBarry Smith 
3732ef1f0ffSBarry Smith   Level: intermediate
3742ef1f0ffSBarry Smith 
37527d8b95cSPierre Jolivet   Developer Notes:
37627d8b95cSPierre Jolivet   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
37727d8b95cSPierre Jolivet 
37827d8b95cSPierre Jolivet   Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage
37927d8b95cSPierre Jolivet 
3801cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
3812ef1f0ffSBarry Smith M*/
3822ef1f0ffSBarry Smith 
383c8a8475eSBarry Smith /*@
384789736e1SBarry Smith   MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like $A^T A$.
385c8a8475eSBarry Smith 
386c3339decSBarry Smith   Collective
387c8a8475eSBarry Smith 
388c8a8475eSBarry Smith   Input Parameter:
389c8a8475eSBarry Smith . A - the (possibly rectangular) matrix
390c8a8475eSBarry Smith 
391c8a8475eSBarry Smith   Output Parameter:
392789736e1SBarry Smith . N - the matrix that represents $A^T A $
393c8a8475eSBarry Smith 
394c30e7cdbSSatish Balay   Level: intermediate
395c30e7cdbSSatish Balay 
39695452b02SPatrick Sanan   Notes:
397789736e1SBarry Smith   The product $A^T A$ is NOT actually formed! Rather the new matrix
39811a5261eSBarry Smith   object performs the matrix-vector product, `MatMult()`, by first multiplying by
399789736e1SBarry Smith   $A$ and then $A^T$
400789736e1SBarry Smith 
401789736e1SBarry Smith   If `MatGetFactor()` is called on this matrix with `MAT_FACTOR_QR` then the inner matrix `A` is used for the factorization
40211a5261eSBarry Smith 
4031cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
404c8a8475eSBarry Smith @*/
MatCreateNormal(Mat A,Mat * N)405d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormal(Mat A, Mat *N)
406d71ae5a4SJacob Faibussowitsch {
407c8a8475eSBarry Smith   Mat_Normal *Na;
4089ca0e1b6SStefano Zampini   VecType     vtype;
409c8a8475eSBarry Smith 
410c8a8475eSBarry Smith   PetscFunctionBegin;
4119566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
4129566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
4139566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
41427d8b95cSPierre Jolivet   PetscCall(MatSetType(*N, MATSHELL));
4154dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
41627d8b95cSPierre Jolivet   PetscCall(MatShellSetContext(*N, Na));
4179566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
418c3122656SLisandro Dalcin   Na->A = A;
4199566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
4202205254eSKarl Rupp 
421b4903688SPierre Jolivet   PetscCall(MatSetBlockSize(*N, A->cmap->bs));
42257d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_Normal));
42357d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Normal));
42457d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMult_Normal));
42557d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (PetscErrorCodeFn *)MatDuplicate_Normal));
42657d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_Normal));
42757d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL_BLOCK, (PetscErrorCodeFn *)MatGetDiagonalBlock_Normal));
42857d50842SBarry Smith   PetscCall(MatShellSetOperation(*N, MATOP_COPY, (PetscErrorCodeFn *)MatCopy_Normal));
429fa80d070SPierre Jolivet   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
430fa80d070SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
431fa80d070SPierre Jolivet   (*N)->ops->permute           = MatPermute_Normal;
4329ca0e1b6SStefano Zampini 
43327d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalGetMat_Normal));
43427d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ));
43527d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ));
4362f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
43727d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE));
4382f0a5c5aSJose E. Roman #endif
43927d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense));
44027d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense));
44127d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable));
44227d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
44327d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
4449566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE));
4459566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
4469566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
4479ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
4489566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
4499ca0e1b6SStefano Zampini #endif
45027d8b95cSPierre Jolivet   PetscCall(MatSetUp(*N));
45127d8b95cSPierre Jolivet   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL));
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
453c8a8475eSBarry Smith }
454