xref: /petsc/src/tao/matrix/adamat.c (revision 2286efddd54511ab18e8e2adb1e023c4bf8f0b92)
17fb60732SBarry Smith #include <petsc/private/matimpl.h> /*I  "petscmat.h"  I*/
2e0877f53SBarry Smith 
33b242c63SJacob Faibussowitsch static PetscErrorCode MatCreateADA(Mat, Vec, Vec, Mat *);
40ebee16dSLisandro Dalcin 
5e0877f53SBarry Smith typedef struct {
6e0877f53SBarry Smith   Mat      A;
7e0877f53SBarry Smith   Vec      D1;
8e0877f53SBarry Smith   Vec      D2;
9e0877f53SBarry Smith   Vec      W;
10e0877f53SBarry Smith   Vec      W2;
11e0877f53SBarry Smith   Vec      ADADiag;
12e0877f53SBarry Smith   PetscInt GotDiag;
13e0877f53SBarry Smith } _p_TaoMatADACtx;
14e0877f53SBarry Smith typedef _p_TaoMatADACtx *TaoMatADACtx;
15e0877f53SBarry Smith 
MatMult_ADA(Mat mat,Vec a,Vec y)16d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ADA(Mat mat, Vec a, Vec y)
17d71ae5a4SJacob Faibussowitsch {
18e0877f53SBarry Smith   TaoMatADACtx ctx;
19e0877f53SBarry Smith   PetscReal    one = 1.0;
20e0877f53SBarry Smith 
21e0877f53SBarry Smith   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
239566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->A, a, ctx->W));
241baa6e33SBarry Smith   if (ctx->D1) PetscCall(VecPointwiseMult(ctx->W, ctx->D1, ctx->W));
259566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(ctx->A, ctx->W, y));
26e0877f53SBarry Smith   if (ctx->D2) {
279566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(ctx->W2, ctx->D2, a));
289566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, one, ctx->W2));
29e0877f53SBarry Smith   }
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31e0877f53SBarry Smith }
32e0877f53SBarry Smith 
MatMultTranspose_ADA(Mat mat,Vec a,Vec y)33d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_ADA(Mat mat, Vec a, Vec y)
34d71ae5a4SJacob Faibussowitsch {
35e0877f53SBarry Smith   PetscFunctionBegin;
369566063dSJacob Faibussowitsch   PetscCall(MatMult_ADA(mat, a, y));
373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38e0877f53SBarry Smith }
39e0877f53SBarry Smith 
MatDiagonalSet_ADA(Mat M,Vec D,InsertMode mode)40d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDiagonalSet_ADA(Mat M, Vec D, InsertMode mode)
41d71ae5a4SJacob Faibussowitsch {
42e0877f53SBarry Smith   TaoMatADACtx ctx;
43e0877f53SBarry Smith   PetscReal    zero = 0.0, one = 1.0;
44e0877f53SBarry Smith 
45e0877f53SBarry Smith   PetscFunctionBegin;
463c859ba3SBarry Smith   PetscCheck(mode != INSERT_VALUES, PetscObjectComm((PetscObject)M), PETSC_ERR_SUP, "Cannot insert diagonal entries of this matrix type, can only add");
479566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M, &ctx));
48e0877f53SBarry Smith   if (!ctx->D2) {
499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(D, &ctx->D2));
509566063dSJacob Faibussowitsch     PetscCall(VecSet(ctx->D2, zero));
51e0877f53SBarry Smith   }
529566063dSJacob Faibussowitsch   PetscCall(VecAXPY(ctx->D2, one, D));
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54e0877f53SBarry Smith }
55e0877f53SBarry Smith 
MatDestroy_ADA(Mat mat)56d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ADA(Mat mat)
57d71ae5a4SJacob Faibussowitsch {
58e0877f53SBarry Smith   TaoMatADACtx ctx;
59e0877f53SBarry Smith 
60e0877f53SBarry Smith   PetscFunctionBegin;
619566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->W));
639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->W2));
649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->ADADiag));
659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->A));
669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->D1));
679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->D2));
689566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70e0877f53SBarry Smith }
71e0877f53SBarry Smith 
MatView_ADA(Mat mat,PetscViewer viewer)72d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_ADA(Mat mat, PetscViewer viewer)
73d71ae5a4SJacob Faibussowitsch {
74e0877f53SBarry Smith   PetscFunctionBegin;
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76e0877f53SBarry Smith }
77e0877f53SBarry Smith 
MatShift_ADA(Mat Y,PetscReal a)78d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_ADA(Mat Y, PetscReal a)
79d71ae5a4SJacob Faibussowitsch {
80e0877f53SBarry Smith   TaoMatADACtx ctx;
81e0877f53SBarry Smith 
82e0877f53SBarry Smith   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(Y, &ctx));
849566063dSJacob Faibussowitsch   PetscCall(VecShift(ctx->D2, a));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86e0877f53SBarry Smith }
87e0877f53SBarry Smith 
MatDuplicate_ADA(Mat mat,MatDuplicateOption op,Mat * M)88d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_ADA(Mat mat, MatDuplicateOption op, Mat *M)
89d71ae5a4SJacob Faibussowitsch {
90e0877f53SBarry Smith   TaoMatADACtx ctx;
91e0877f53SBarry Smith   Mat          A2;
92e0877f53SBarry Smith   Vec          D1b = NULL, D2b;
93e0877f53SBarry Smith 
94e0877f53SBarry Smith   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
969566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ctx->A, op, &A2));
97e0877f53SBarry Smith   if (ctx->D1) {
989566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ctx->D1, &D1b));
999566063dSJacob Faibussowitsch     PetscCall(VecCopy(ctx->D1, D1b));
100e0877f53SBarry Smith   }
1019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ctx->D2, &D2b));
1029566063dSJacob Faibussowitsch   PetscCall(VecCopy(ctx->D2, D2b));
1039566063dSJacob Faibussowitsch   PetscCall(MatCreateADA(A2, D1b, D2b, M));
1041baa6e33SBarry Smith   if (ctx->D1) PetscCall(PetscObjectDereference((PetscObject)D1b));
1059566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)D2b));
1069566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)A2));
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108e0877f53SBarry Smith }
109e0877f53SBarry Smith 
MatEqual_ADA(Mat A,Mat B,PetscBool * flg)110d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatEqual_ADA(Mat A, Mat B, PetscBool *flg)
111d71ae5a4SJacob Faibussowitsch {
112e0877f53SBarry Smith   TaoMatADACtx ctx1, ctx2;
113e0877f53SBarry Smith 
114e0877f53SBarry Smith   PetscFunctionBegin;
1159566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx1));
1169566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(B, &ctx2));
1179566063dSJacob Faibussowitsch   PetscCall(VecEqual(ctx1->D2, ctx2->D2, flg));
1181baa6e33SBarry Smith   if (*flg == PETSC_TRUE) PetscCall(VecEqual(ctx1->D1, ctx2->D1, flg));
1191baa6e33SBarry Smith   if (*flg == PETSC_TRUE) PetscCall(MatEqual(ctx1->A, ctx2->A, flg));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121e0877f53SBarry Smith }
122e0877f53SBarry Smith 
MatScale_ADA(Mat mat,PetscReal a)123d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_ADA(Mat mat, PetscReal a)
124d71ae5a4SJacob Faibussowitsch {
125e0877f53SBarry Smith   TaoMatADACtx ctx;
126e0877f53SBarry Smith 
127e0877f53SBarry Smith   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1299566063dSJacob Faibussowitsch   PetscCall(VecScale(ctx->D1, a));
1301baa6e33SBarry Smith   if (ctx->D2) PetscCall(VecScale(ctx->D2, a));
1313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
132e0877f53SBarry Smith }
133e0877f53SBarry Smith 
MatTranspose_ADA(Mat mat,MatReuse reuse,Mat * B)134d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_ADA(Mat mat, MatReuse reuse, Mat *B)
135d71ae5a4SJacob Faibussowitsch {
136e0877f53SBarry Smith   TaoMatADACtx ctx;
137e0877f53SBarry Smith 
138e0877f53SBarry Smith   PetscFunctionBegin;
1397fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(mat, *B));
1409566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
141cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
1429566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, B));
143cf37664fSBarry Smith   } else if (reuse == MAT_REUSE_MATRIX) {
1449566063dSJacob Faibussowitsch     PetscCall(MatCopy(mat, *B, SAME_NONZERO_PATTERN));
145cf37664fSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Does not support inplace transpose");
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147e0877f53SBarry Smith }
148e0877f53SBarry Smith 
MatADAComputeDiagonal(Mat mat)149d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatADAComputeDiagonal(Mat mat)
150d71ae5a4SJacob Faibussowitsch {
151e0877f53SBarry Smith   PetscInt     i, m, n, low, high;
152e0877f53SBarry Smith   PetscScalar *dtemp, *dptr;
153e0877f53SBarry Smith   TaoMatADACtx ctx;
154e0877f53SBarry Smith 
155e0877f53SBarry Smith   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1579566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(mat, &low, &high));
1589566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &m, &n));
159e0877f53SBarry Smith 
1609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &dtemp));
161e0877f53SBarry Smith   for (i = 0; i < n; i++) {
1629566063dSJacob Faibussowitsch     PetscCall(MatGetColumnVector(ctx->A, ctx->W, i));
1639566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(ctx->W, ctx->W, ctx->W));
1649566063dSJacob Faibussowitsch     PetscCall(VecDotBegin(ctx->D1, ctx->W, dtemp + i));
165e0877f53SBarry Smith   }
16648a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(VecDotEnd(ctx->D1, ctx->W, dtemp + i));
167e0877f53SBarry Smith 
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->ADADiag, &dptr));
169ad540459SPierre Jolivet   for (i = low; i < high; i++) dptr[i - low] = dtemp[i];
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->ADADiag, &dptr));
1719566063dSJacob Faibussowitsch   PetscCall(PetscFree(dtemp));
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173e0877f53SBarry Smith }
174e0877f53SBarry Smith 
MatGetDiagonal_ADA(Mat mat,Vec v)175d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_ADA(Mat mat, Vec v)
176d71ae5a4SJacob Faibussowitsch {
177e0877f53SBarry Smith   PetscReal    one = 1.0;
178e0877f53SBarry Smith   TaoMatADACtx ctx;
179e0877f53SBarry Smith 
180e0877f53SBarry Smith   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1829566063dSJacob Faibussowitsch   PetscCall(MatADAComputeDiagonal(mat));
1839566063dSJacob Faibussowitsch   PetscCall(VecCopy(ctx->ADADiag, v));
1841baa6e33SBarry Smith   if (ctx->D2) PetscCall(VecAXPY(v, one, ctx->D2));
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
186e0877f53SBarry Smith }
187e0877f53SBarry Smith 
MatCreateSubMatrix_ADA(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat * newmat)188d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_ADA(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat)
189d71ae5a4SJacob Faibussowitsch {
190e0877f53SBarry Smith   PetscInt     low, high;
191e0877f53SBarry Smith   IS           ISrow;
192e0877f53SBarry Smith   Vec          D1, D2;
193e0877f53SBarry Smith   Mat          Atemp;
194e0877f53SBarry Smith   TaoMatADACtx ctx;
195e0877f53SBarry Smith   PetscBool    isequal;
196e0877f53SBarry Smith 
197e0877f53SBarry Smith   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(ISEqual(isrow, iscol, &isequal));
1993c859ba3SBarry Smith   PetscCheck(isequal, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for identical column and row indices");
2009566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
201e0877f53SBarry Smith 
2029566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(ctx->A, &low, &high));
2039566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)mat), high - low, low, 1, &ISrow));
2049566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(ctx->A, ISrow, iscol, cll, &Atemp));
2059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ISrow));
206e0877f53SBarry Smith 
207e0877f53SBarry Smith   if (ctx->D1) {
2089566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ctx->D1, &D1));
2099566063dSJacob Faibussowitsch     PetscCall(VecCopy(ctx->D1, D1));
210e0877f53SBarry Smith   } else {
211e0877f53SBarry Smith     D1 = NULL;
212e0877f53SBarry Smith   }
213e0877f53SBarry Smith 
214e0877f53SBarry Smith   if (ctx->D2) {
215e0877f53SBarry Smith     Vec D2sub;
216e0877f53SBarry Smith 
2179566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(ctx->D2, isrow, &D2sub));
2189566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(D2sub, &D2));
2199566063dSJacob Faibussowitsch     PetscCall(VecCopy(D2sub, D2));
2209566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(ctx->D2, isrow, &D2sub));
221e0877f53SBarry Smith   } else {
222e0877f53SBarry Smith     D2 = NULL;
223e0877f53SBarry Smith   }
224e0877f53SBarry Smith 
2259566063dSJacob Faibussowitsch   PetscCall(MatCreateADA(Atemp, D1, D2, newmat));
2269566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(*newmat, &ctx));
2279566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)Atemp));
2281baa6e33SBarry Smith   if (ctx->D1) PetscCall(PetscObjectDereference((PetscObject)D1));
2291baa6e33SBarry Smith   if (ctx->D2) PetscCall(PetscObjectDereference((PetscObject)D2));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231e0877f53SBarry Smith }
232e0877f53SBarry Smith 
MatCreateSubMatrices_ADA(Mat A,PetscInt n,IS * irow,IS * icol,MatReuse scall,Mat ** B)233d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_ADA(Mat A, PetscInt n, IS *irow, IS *icol, MatReuse scall, Mat **B)
234d71ae5a4SJacob Faibussowitsch {
235e0877f53SBarry Smith   PetscInt i;
236e0877f53SBarry Smith 
237e0877f53SBarry Smith   PetscFunctionBegin;
23848a46eb9SPierre Jolivet   if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B));
23948a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_ADA(A, irow[i], icol[i], scall, &(*B)[i]));
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241e0877f53SBarry Smith }
242e0877f53SBarry Smith 
MatGetColumnVector_ADA(Mat mat,Vec Y,PetscInt col)243d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetColumnVector_ADA(Mat mat, Vec Y, PetscInt col)
244d71ae5a4SJacob Faibussowitsch {
245e0877f53SBarry Smith   PetscInt    low, high;
246e0877f53SBarry Smith   PetscScalar zero = 0.0, one = 1.0;
247e0877f53SBarry Smith 
248e0877f53SBarry Smith   PetscFunctionBegin;
2499566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, zero));
2509566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(Y, &low, &high));
25148a46eb9SPierre Jolivet   if (col >= low && col < high) PetscCall(VecSetValue(Y, col, one, INSERT_VALUES));
2529566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y));
2539566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));
2549566063dSJacob Faibussowitsch   PetscCall(MatMult_ADA(mat, Y, Y));
2553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
256e0877f53SBarry Smith }
257e0877f53SBarry Smith 
MatConvert_ADA(Mat mat,MatType newtype,Mat * NewMat)258d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_ADA(Mat mat, MatType newtype, Mat *NewMat)
259d71ae5a4SJacob Faibussowitsch {
260e0877f53SBarry Smith   PetscMPIInt  size;
261e0877f53SBarry Smith   PetscBool    sametype, issame, isdense, isseqdense;
262e0877f53SBarry Smith   TaoMatADACtx ctx;
263e0877f53SBarry Smith 
264e0877f53SBarry Smith   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
267e0877f53SBarry Smith 
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, newtype, &sametype));
2699566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSAME, &issame));
2709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIDENSE, &isdense));
2719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQDENSE, &isseqdense));
272e0877f53SBarry Smith 
273e0877f53SBarry Smith   if (sametype || issame) {
2749566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, NewMat));
275e0877f53SBarry Smith   } else if (isdense) {
276e0877f53SBarry Smith     PetscInt           i, j, low, high, m, n, M, N;
277e0877f53SBarry Smith     const PetscScalar *dptr;
278e0877f53SBarry Smith     Vec                X;
279e0877f53SBarry Smith 
2809566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ctx->D2, &X));
2819566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mat, &M, &N));
2829566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mat, &m, &n));
2839566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)mat), m, m, N, N, NULL, NewMat));
2849566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(*NewMat, &low, &high));
285e0877f53SBarry Smith     for (i = 0; i < M; i++) {
2869566063dSJacob Faibussowitsch       PetscCall(MatGetColumnVector_ADA(mat, X, i));
2879566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &dptr));
28848a46eb9SPierre Jolivet       for (j = 0; j < high - low; j++) PetscCall(MatSetValue(*NewMat, low + j, i, dptr[j], INSERT_VALUES));
2899566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &dptr));
290e0877f53SBarry Smith     }
2919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(*NewMat, MAT_FINAL_ASSEMBLY));
2929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(*NewMat, MAT_FINAL_ASSEMBLY));
2939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&X));
294e0877f53SBarry Smith   } else if (isseqdense && size == 1) {
295e0877f53SBarry Smith     PetscInt           i, j, low, high, m, n, M, N;
296e0877f53SBarry Smith     const PetscScalar *dptr;
297e0877f53SBarry Smith     Vec                X;
298e0877f53SBarry Smith 
2999566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ctx->D2, &X));
3009566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mat, &M, &N));
3019566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mat, &m, &n));
3029566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PetscObjectComm((PetscObject)mat), N, N, NULL, NewMat));
3039566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(*NewMat, &low, &high));
304e0877f53SBarry Smith     for (i = 0; i < M; i++) {
3059566063dSJacob Faibussowitsch       PetscCall(MatGetColumnVector_ADA(mat, X, i));
3069566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &dptr));
30748a46eb9SPierre Jolivet       for (j = 0; j < high - low; j++) PetscCall(MatSetValue(*NewMat, low + j, i, dptr[j], INSERT_VALUES));
3089566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &dptr));
309e0877f53SBarry Smith     }
3109566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(*NewMat, MAT_FINAL_ASSEMBLY));
3119566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(*NewMat, MAT_FINAL_ASSEMBLY));
3129566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&X));
313691b26d3SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No support to convert objects to that type");
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
315e0877f53SBarry Smith }
316e0877f53SBarry Smith 
MatNorm_ADA(Mat mat,NormType type,PetscReal * norm)317d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_ADA(Mat mat, NormType type, PetscReal *norm)
318d71ae5a4SJacob Faibussowitsch {
319e0877f53SBarry Smith   TaoMatADACtx ctx;
320e0877f53SBarry Smith 
321e0877f53SBarry Smith   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
323e0877f53SBarry Smith   if (type == NORM_FROBENIUS) {
324e0877f53SBarry Smith     *norm = 1.0;
325e0877f53SBarry Smith   } else if (type == NORM_1 || type == NORM_INFINITY) {
326e0877f53SBarry Smith     *norm = 1.0;
327e0877f53SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No two norm");
3283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
329e0877f53SBarry Smith }
330a7e14dcfSSatish Balay 
3313b242c63SJacob Faibussowitsch /*
332a7e14dcfSSatish Balay    MatCreateADA - Creates a matrix M=A^T D1 A + D2 where D1, D2 are diagonal
333a7e14dcfSSatish Balay 
334c3339decSBarry Smith    Collective
335a7e14dcfSSatish Balay 
336a7e14dcfSSatish Balay    Input Parameters:
337a7e14dcfSSatish Balay +  mat - matrix of arbitrary type
3380c0fd78eSBarry Smith .  d1 - A vector defining a diagonal matrix
3390c0fd78eSBarry Smith -  d2 - A vector defining a diagonal matrix
340a7e14dcfSSatish Balay 
3412fe279fdSBarry Smith    Output Parameter:
342a7e14dcfSSatish Balay .  J - New matrix whose operations are defined in terms of mat, D1, and D2.
343a7e14dcfSSatish Balay 
344c3339decSBarry Smith    Level: developer
345c3339decSBarry Smith 
346c3339decSBarry Smith    Note:
347a7e14dcfSSatish Balay    The user provides the input data and is responsible for destroying
34827430b45SBarry Smith    this data after matrix `J` has been destroyed.
349a7e14dcfSSatish Balay 
350c3339decSBarry Smith .seealso: `Mat`, `MatCreate()`
3513b242c63SJacob Faibussowitsch */
MatCreateADA(Mat mat,Vec d1,Vec d2,Mat * J)352ba38deedSJacob Faibussowitsch static PetscErrorCode MatCreateADA(Mat mat, Vec d1, Vec d2, Mat *J)
353d71ae5a4SJacob Faibussowitsch {
354367daffbSBarry Smith   MPI_Comm     comm = PetscObjectComm((PetscObject)mat);
355a7e14dcfSSatish Balay   TaoMatADACtx ctx;
356a7e14dcfSSatish Balay   PetscInt     nloc, n;
357a7e14dcfSSatish Balay 
358a7e14dcfSSatish Balay   PetscFunctionBegin;
3599566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
360a7e14dcfSSatish Balay   ctx->A  = mat;
361a7e14dcfSSatish Balay   ctx->D1 = d1;
362a7e14dcfSSatish Balay   ctx->D2 = d2;
363a7e14dcfSSatish Balay   if (d1) {
3649566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(d1, &ctx->W));
3659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)d1));
366a7e14dcfSSatish Balay   } else {
3670e660641SBarry Smith     ctx->W = NULL;
368a7e14dcfSSatish Balay   }
369a7e14dcfSSatish Balay   if (d2) {
3709566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(d2, &ctx->W2));
3719566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(d2, &ctx->ADADiag));
3729566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)d2));
373a7e14dcfSSatish Balay   } else {
3740e660641SBarry Smith     ctx->W2      = NULL;
3750e660641SBarry Smith     ctx->ADADiag = NULL;
376a7e14dcfSSatish Balay   }
377a7e14dcfSSatish Balay 
378a7e14dcfSSatish Balay   ctx->GotDiag = 0;
3799566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
380a7e14dcfSSatish Balay 
3819566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(d2, &nloc));
3829566063dSJacob Faibussowitsch   PetscCall(VecGetSize(d2, &n));
383a7e14dcfSSatish Balay 
3849566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(comm, nloc, nloc, n, n, ctx, J));
3859566063dSJacob Faibussowitsch   PetscCall(MatShellSetManageScalingShifts(*J));
386*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (PetscErrorCodeFn *)MatMult_ADA));
387*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_ADA));
388*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (PetscErrorCodeFn *)MatView_ADA));
389*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_ADA));
390*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_DIAGONAL_SET, (PetscErrorCodeFn *)MatDiagonalSet_ADA));
391*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_SHIFT, (PetscErrorCodeFn *)MatShift_ADA));
392*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_EQUAL, (PetscErrorCodeFn *)MatEqual_ADA));
393*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_SCALE, (PetscErrorCodeFn *)MatScale_ADA));
394*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_TRANSPOSE, (PetscErrorCodeFn *)MatTranspose_ADA));
395*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_ADA));
396*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRICES, (PetscErrorCodeFn *)MatCreateSubMatrices_ADA));
397*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_NORM, (PetscErrorCodeFn *)MatNorm_ADA));
398*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_DUPLICATE, (PetscErrorCodeFn *)MatDuplicate_ADA));
399*57d50842SBarry Smith   PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRIX, (PetscErrorCodeFn *)MatCreateSubMatrix_ADA));
400a7e14dcfSSatish Balay 
4019566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*J, MAT_SYMMETRIC, PETSC_TRUE));
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
403a7e14dcfSSatish Balay }
404