17fb60732SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 2e0877f53SBarry Smith 30ebee16dSLisandro Dalcin PETSC_INTERN 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 169371c9d4SSatish Balay static PetscErrorCode MatMult_ADA(Mat mat, Vec a, Vec y) { 17e0877f53SBarry Smith TaoMatADACtx ctx; 18e0877f53SBarry Smith PetscReal one = 1.0; 19e0877f53SBarry Smith 20e0877f53SBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 229566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->A, a, ctx->W)); 231baa6e33SBarry Smith if (ctx->D1) PetscCall(VecPointwiseMult(ctx->W, ctx->D1, ctx->W)); 249566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->A, ctx->W, y)); 25e0877f53SBarry Smith if (ctx->D2) { 269566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ctx->W2, ctx->D2, a)); 279566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, one, ctx->W2)); 28e0877f53SBarry Smith } 29e0877f53SBarry Smith PetscFunctionReturn(0); 30e0877f53SBarry Smith } 31e0877f53SBarry Smith 329371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_ADA(Mat mat, Vec a, Vec y) { 33e0877f53SBarry Smith PetscFunctionBegin; 349566063dSJacob Faibussowitsch PetscCall(MatMult_ADA(mat, a, y)); 35e0877f53SBarry Smith PetscFunctionReturn(0); 36e0877f53SBarry Smith } 37e0877f53SBarry Smith 389371c9d4SSatish Balay static PetscErrorCode MatDiagonalSet_ADA(Mat M, Vec D, InsertMode mode) { 39e0877f53SBarry Smith TaoMatADACtx ctx; 40e0877f53SBarry Smith PetscReal zero = 0.0, one = 1.0; 41e0877f53SBarry Smith 42e0877f53SBarry Smith PetscFunctionBegin; 433c859ba3SBarry Smith PetscCheck(mode != INSERT_VALUES, PetscObjectComm((PetscObject)M), PETSC_ERR_SUP, "Cannot insert diagonal entries of this matrix type, can only add"); 449566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, &ctx)); 45e0877f53SBarry Smith if (!ctx->D2) { 469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(D, &ctx->D2)); 479566063dSJacob Faibussowitsch PetscCall(VecSet(ctx->D2, zero)); 48e0877f53SBarry Smith } 499566063dSJacob Faibussowitsch PetscCall(VecAXPY(ctx->D2, one, D)); 50e0877f53SBarry Smith PetscFunctionReturn(0); 51e0877f53SBarry Smith } 52e0877f53SBarry Smith 539371c9d4SSatish Balay static PetscErrorCode MatDestroy_ADA(Mat mat) { 54e0877f53SBarry Smith TaoMatADACtx ctx; 55e0877f53SBarry Smith 56e0877f53SBarry Smith PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->W)); 599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->W2)); 609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->ADADiag)); 619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->A)); 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->D1)); 639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->D2)); 649566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 65e0877f53SBarry Smith PetscFunctionReturn(0); 66e0877f53SBarry Smith } 67e0877f53SBarry Smith 689371c9d4SSatish Balay static PetscErrorCode MatView_ADA(Mat mat, PetscViewer viewer) { 69e0877f53SBarry Smith PetscFunctionBegin; 70e0877f53SBarry Smith PetscFunctionReturn(0); 71e0877f53SBarry Smith } 72e0877f53SBarry Smith 739371c9d4SSatish Balay static PetscErrorCode MatShift_ADA(Mat Y, PetscReal a) { 74e0877f53SBarry Smith TaoMatADACtx ctx; 75e0877f53SBarry Smith 76e0877f53SBarry Smith PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(Y, &ctx)); 789566063dSJacob Faibussowitsch PetscCall(VecShift(ctx->D2, a)); 79e0877f53SBarry Smith PetscFunctionReturn(0); 80e0877f53SBarry Smith } 81e0877f53SBarry Smith 829371c9d4SSatish Balay static PetscErrorCode MatDuplicate_ADA(Mat mat, MatDuplicateOption op, Mat *M) { 83e0877f53SBarry Smith TaoMatADACtx ctx; 84e0877f53SBarry Smith Mat A2; 85e0877f53SBarry Smith Vec D1b = NULL, D2b; 86e0877f53SBarry Smith 87e0877f53SBarry Smith PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 899566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ctx->A, op, &A2)); 90e0877f53SBarry Smith if (ctx->D1) { 919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->D1, &D1b)); 929566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->D1, D1b)); 93e0877f53SBarry Smith } 949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->D2, &D2b)); 959566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->D2, D2b)); 969566063dSJacob Faibussowitsch PetscCall(MatCreateADA(A2, D1b, D2b, M)); 971baa6e33SBarry Smith if (ctx->D1) PetscCall(PetscObjectDereference((PetscObject)D1b)); 989566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)D2b)); 999566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)A2)); 100e0877f53SBarry Smith PetscFunctionReturn(0); 101e0877f53SBarry Smith } 102e0877f53SBarry Smith 1039371c9d4SSatish Balay static PetscErrorCode MatEqual_ADA(Mat A, Mat B, PetscBool *flg) { 104e0877f53SBarry Smith TaoMatADACtx ctx1, ctx2; 105e0877f53SBarry Smith 106e0877f53SBarry Smith PetscFunctionBegin; 1079566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx1)); 1089566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(B, &ctx2)); 1099566063dSJacob Faibussowitsch PetscCall(VecEqual(ctx1->D2, ctx2->D2, flg)); 1101baa6e33SBarry Smith if (*flg == PETSC_TRUE) PetscCall(VecEqual(ctx1->D1, ctx2->D1, flg)); 1111baa6e33SBarry Smith if (*flg == PETSC_TRUE) PetscCall(MatEqual(ctx1->A, ctx2->A, flg)); 112e0877f53SBarry Smith PetscFunctionReturn(0); 113e0877f53SBarry Smith } 114e0877f53SBarry Smith 1159371c9d4SSatish Balay static PetscErrorCode MatScale_ADA(Mat mat, PetscReal a) { 116e0877f53SBarry Smith TaoMatADACtx ctx; 117e0877f53SBarry Smith 118e0877f53SBarry Smith PetscFunctionBegin; 1199566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1209566063dSJacob Faibussowitsch PetscCall(VecScale(ctx->D1, a)); 1211baa6e33SBarry Smith if (ctx->D2) PetscCall(VecScale(ctx->D2, a)); 122e0877f53SBarry Smith PetscFunctionReturn(0); 123e0877f53SBarry Smith } 124e0877f53SBarry Smith 1259371c9d4SSatish Balay static PetscErrorCode MatTranspose_ADA(Mat mat, MatReuse reuse, Mat *B) { 126e0877f53SBarry Smith TaoMatADACtx ctx; 127e0877f53SBarry Smith 128e0877f53SBarry Smith PetscFunctionBegin; 1297fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(mat, *B)); 1309566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 131cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 1329566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, B)); 133cf37664fSBarry Smith } else if (reuse == MAT_REUSE_MATRIX) { 1349566063dSJacob Faibussowitsch PetscCall(MatCopy(mat, *B, SAME_NONZERO_PATTERN)); 135cf37664fSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Does not support inplace transpose"); 136e0877f53SBarry Smith PetscFunctionReturn(0); 137e0877f53SBarry Smith } 138e0877f53SBarry Smith 1399371c9d4SSatish Balay static PetscErrorCode MatADAComputeDiagonal(Mat mat) { 140e0877f53SBarry Smith PetscInt i, m, n, low, high; 141e0877f53SBarry Smith PetscScalar *dtemp, *dptr; 142e0877f53SBarry Smith TaoMatADACtx ctx; 143e0877f53SBarry Smith 144e0877f53SBarry Smith PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &low, &high)); 1479566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &m, &n)); 148e0877f53SBarry Smith 1499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dtemp)); 150e0877f53SBarry Smith for (i = 0; i < n; i++) { 1519566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(ctx->A, ctx->W, i)); 1529566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ctx->W, ctx->W, ctx->W)); 1539566063dSJacob Faibussowitsch PetscCall(VecDotBegin(ctx->D1, ctx->W, dtemp + i)); 154e0877f53SBarry Smith } 155*48a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(VecDotEnd(ctx->D1, ctx->W, dtemp + i)); 156e0877f53SBarry Smith 1579566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->ADADiag, &dptr)); 1589371c9d4SSatish Balay for (i = low; i < high; i++) { dptr[i - low] = dtemp[i]; } 1599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->ADADiag, &dptr)); 1609566063dSJacob Faibussowitsch PetscCall(PetscFree(dtemp)); 161e0877f53SBarry Smith PetscFunctionReturn(0); 162e0877f53SBarry Smith } 163e0877f53SBarry Smith 1649371c9d4SSatish Balay static PetscErrorCode MatGetDiagonal_ADA(Mat mat, Vec v) { 165e0877f53SBarry Smith PetscReal one = 1.0; 166e0877f53SBarry Smith TaoMatADACtx ctx; 167e0877f53SBarry Smith 168e0877f53SBarry Smith PetscFunctionBegin; 1699566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1709566063dSJacob Faibussowitsch PetscCall(MatADAComputeDiagonal(mat)); 1719566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->ADADiag, v)); 1721baa6e33SBarry Smith if (ctx->D2) PetscCall(VecAXPY(v, one, ctx->D2)); 173e0877f53SBarry Smith PetscFunctionReturn(0); 174e0877f53SBarry Smith } 175e0877f53SBarry Smith 1769371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrix_ADA(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) { 177e0877f53SBarry Smith PetscInt low, high; 178e0877f53SBarry Smith IS ISrow; 179e0877f53SBarry Smith Vec D1, D2; 180e0877f53SBarry Smith Mat Atemp; 181e0877f53SBarry Smith TaoMatADACtx ctx; 182e0877f53SBarry Smith PetscBool isequal; 183e0877f53SBarry Smith 184e0877f53SBarry Smith PetscFunctionBegin; 1859566063dSJacob Faibussowitsch PetscCall(ISEqual(isrow, iscol, &isequal)); 1863c859ba3SBarry Smith PetscCheck(isequal, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for identical column and row indices"); 1879566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 188e0877f53SBarry Smith 1899566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ctx->A, &low, &high)); 1909566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)mat), high - low, low, 1, &ISrow)); 1919566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(ctx->A, ISrow, iscol, cll, &Atemp)); 1929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ISrow)); 193e0877f53SBarry Smith 194e0877f53SBarry Smith if (ctx->D1) { 1959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->D1, &D1)); 1969566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->D1, D1)); 197e0877f53SBarry Smith } else { 198e0877f53SBarry Smith D1 = NULL; 199e0877f53SBarry Smith } 200e0877f53SBarry Smith 201e0877f53SBarry Smith if (ctx->D2) { 202e0877f53SBarry Smith Vec D2sub; 203e0877f53SBarry Smith 2049566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->D2, isrow, &D2sub)); 2059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(D2sub, &D2)); 2069566063dSJacob Faibussowitsch PetscCall(VecCopy(D2sub, D2)); 2079566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->D2, isrow, &D2sub)); 208e0877f53SBarry Smith } else { 209e0877f53SBarry Smith D2 = NULL; 210e0877f53SBarry Smith } 211e0877f53SBarry Smith 2129566063dSJacob Faibussowitsch PetscCall(MatCreateADA(Atemp, D1, D2, newmat)); 2139566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(*newmat, &ctx)); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)Atemp)); 2151baa6e33SBarry Smith if (ctx->D1) PetscCall(PetscObjectDereference((PetscObject)D1)); 2161baa6e33SBarry Smith if (ctx->D2) PetscCall(PetscObjectDereference((PetscObject)D2)); 217e0877f53SBarry Smith PetscFunctionReturn(0); 218e0877f53SBarry Smith } 219e0877f53SBarry Smith 2209371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrices_ADA(Mat A, PetscInt n, IS *irow, IS *icol, MatReuse scall, Mat **B) { 221e0877f53SBarry Smith PetscInt i; 222e0877f53SBarry Smith 223e0877f53SBarry Smith PetscFunctionBegin; 224*48a46eb9SPierre Jolivet if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B)); 225*48a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_ADA(A, irow[i], icol[i], scall, &(*B)[i])); 226e0877f53SBarry Smith PetscFunctionReturn(0); 227e0877f53SBarry Smith } 228e0877f53SBarry Smith 2299371c9d4SSatish Balay static PetscErrorCode MatGetColumnVector_ADA(Mat mat, Vec Y, PetscInt col) { 230e0877f53SBarry Smith PetscInt low, high; 231e0877f53SBarry Smith PetscScalar zero = 0.0, one = 1.0; 232e0877f53SBarry Smith 233e0877f53SBarry Smith PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(VecSet(Y, zero)); 2359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(Y, &low, &high)); 236*48a46eb9SPierre Jolivet if (col >= low && col < high) PetscCall(VecSetValue(Y, col, one, INSERT_VALUES)); 2379566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y)); 2389566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y)); 2399566063dSJacob Faibussowitsch PetscCall(MatMult_ADA(mat, Y, Y)); 240e0877f53SBarry Smith PetscFunctionReturn(0); 241e0877f53SBarry Smith } 242e0877f53SBarry Smith 2439371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_ADA(Mat mat, MatType newtype, Mat *NewMat) { 244e0877f53SBarry Smith PetscMPIInt size; 245e0877f53SBarry Smith PetscBool sametype, issame, isdense, isseqdense; 246e0877f53SBarry Smith TaoMatADACtx ctx; 247e0877f53SBarry Smith 248e0877f53SBarry Smith PetscFunctionBegin; 2499566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 251e0877f53SBarry Smith 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, newtype, &sametype)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSAME, &issame)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIDENSE, &isdense)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQDENSE, &isseqdense)); 256e0877f53SBarry Smith 257e0877f53SBarry Smith if (sametype || issame) { 2589566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, NewMat)); 259e0877f53SBarry Smith } else if (isdense) { 260e0877f53SBarry Smith PetscInt i, j, low, high, m, n, M, N; 261e0877f53SBarry Smith const PetscScalar *dptr; 262e0877f53SBarry Smith Vec X; 263e0877f53SBarry Smith 2649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->D2, &X)); 2659566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &M, &N)); 2669566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &m, &n)); 2679566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)mat), m, m, N, N, NULL, NewMat)); 2689566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*NewMat, &low, &high)); 269e0877f53SBarry Smith for (i = 0; i < M; i++) { 2709566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector_ADA(mat, X, i)); 2719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &dptr)); 272*48a46eb9SPierre Jolivet for (j = 0; j < high - low; j++) PetscCall(MatSetValue(*NewMat, low + j, i, dptr[j], INSERT_VALUES)); 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &dptr)); 274e0877f53SBarry Smith } 2759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*NewMat, MAT_FINAL_ASSEMBLY)); 2769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*NewMat, MAT_FINAL_ASSEMBLY)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 278e0877f53SBarry Smith } else if (isseqdense && size == 1) { 279e0877f53SBarry Smith PetscInt i, j, low, high, m, n, M, N; 280e0877f53SBarry Smith const PetscScalar *dptr; 281e0877f53SBarry Smith Vec X; 282e0877f53SBarry Smith 2839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->D2, &X)); 2849566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &M, &N)); 2859566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &m, &n)); 2869566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PetscObjectComm((PetscObject)mat), N, N, NULL, NewMat)); 2879566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*NewMat, &low, &high)); 288e0877f53SBarry Smith for (i = 0; i < M; i++) { 2899566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector_ADA(mat, X, i)); 2909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &dptr)); 291*48a46eb9SPierre Jolivet for (j = 0; j < high - low; j++) PetscCall(MatSetValue(*NewMat, low + j, i, dptr[j], INSERT_VALUES)); 2929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &dptr)); 293e0877f53SBarry Smith } 2949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*NewMat, MAT_FINAL_ASSEMBLY)); 2959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*NewMat, MAT_FINAL_ASSEMBLY)); 2969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 297691b26d3SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No support to convert objects to that type"); 298e0877f53SBarry Smith PetscFunctionReturn(0); 299e0877f53SBarry Smith } 300e0877f53SBarry Smith 3019371c9d4SSatish Balay static PetscErrorCode MatNorm_ADA(Mat mat, NormType type, PetscReal *norm) { 302e0877f53SBarry Smith TaoMatADACtx ctx; 303e0877f53SBarry Smith 304e0877f53SBarry Smith PetscFunctionBegin; 3059566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 306e0877f53SBarry Smith if (type == NORM_FROBENIUS) { 307e0877f53SBarry Smith *norm = 1.0; 308e0877f53SBarry Smith } else if (type == NORM_1 || type == NORM_INFINITY) { 309e0877f53SBarry Smith *norm = 1.0; 310e0877f53SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No two norm"); 311e0877f53SBarry Smith PetscFunctionReturn(0); 312e0877f53SBarry Smith } 313a7e14dcfSSatish Balay 314a7e14dcfSSatish Balay /*@C 315a7e14dcfSSatish Balay MatCreateADA - Creates a matrix M=A^T D1 A + D2 where D1, D2 are diagonal 316a7e14dcfSSatish Balay 317a7e14dcfSSatish Balay Collective on matrix 318a7e14dcfSSatish Balay 319a7e14dcfSSatish Balay Input Parameters: 320a7e14dcfSSatish Balay + mat - matrix of arbitrary type 3210c0fd78eSBarry Smith . d1 - A vector defining a diagonal matrix 3220c0fd78eSBarry Smith - d2 - A vector defining a diagonal matrix 323a7e14dcfSSatish Balay 324a7e14dcfSSatish Balay Output Parameters: 325a7e14dcfSSatish Balay . J - New matrix whose operations are defined in terms of mat, D1, and D2. 326a7e14dcfSSatish Balay 327a7e14dcfSSatish Balay Notes: 328a7e14dcfSSatish Balay The user provides the input data and is responsible for destroying 329a7e14dcfSSatish Balay this data after matrix J has been destroyed. 330a7e14dcfSSatish Balay 331a7e14dcfSSatish Balay Level: developer 332a7e14dcfSSatish Balay 333db781477SPatrick Sanan .seealso: `MatCreate()` 334a7e14dcfSSatish Balay @*/ 3359371c9d4SSatish Balay PetscErrorCode MatCreateADA(Mat mat, Vec d1, Vec d2, Mat *J) { 336367daffbSBarry Smith MPI_Comm comm = PetscObjectComm((PetscObject)mat); 337a7e14dcfSSatish Balay TaoMatADACtx ctx; 338a7e14dcfSSatish Balay PetscInt nloc, n; 339a7e14dcfSSatish Balay 340a7e14dcfSSatish Balay PetscFunctionBegin; 3419566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 342a7e14dcfSSatish Balay ctx->A = mat; 343a7e14dcfSSatish Balay ctx->D1 = d1; 344a7e14dcfSSatish Balay ctx->D2 = d2; 345a7e14dcfSSatish Balay if (d1) { 3469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(d1, &ctx->W)); 3479566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)d1)); 348a7e14dcfSSatish Balay } else { 3490e660641SBarry Smith ctx->W = NULL; 350a7e14dcfSSatish Balay } 351a7e14dcfSSatish Balay if (d2) { 3529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(d2, &ctx->W2)); 3539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(d2, &ctx->ADADiag)); 3549566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)d2)); 355a7e14dcfSSatish Balay } else { 3560e660641SBarry Smith ctx->W2 = NULL; 3570e660641SBarry Smith ctx->ADADiag = NULL; 358a7e14dcfSSatish Balay } 359a7e14dcfSSatish Balay 360a7e14dcfSSatish Balay ctx->GotDiag = 0; 3619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 362a7e14dcfSSatish Balay 3639566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(d2, &nloc)); 3649566063dSJacob Faibussowitsch PetscCall(VecGetSize(d2, &n)); 365a7e14dcfSSatish Balay 3669566063dSJacob Faibussowitsch PetscCall(MatCreateShell(comm, nloc, nloc, n, n, ctx, J)); 3679566063dSJacob Faibussowitsch PetscCall(MatShellSetManageScalingShifts(*J)); 3689566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))MatMult_ADA)); 3699566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))MatDestroy_ADA)); 3709566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (void (*)(void))MatView_ADA)); 3719566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_ADA)); 3729566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_DIAGONAL_SET, (void (*)(void))MatDiagonalSet_ADA)); 3739566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_SHIFT, (void (*)(void))MatShift_ADA)); 3749566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_EQUAL, (void (*)(void))MatEqual_ADA)); 3759566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_SCALE, (void (*)(void))MatScale_ADA)); 3769566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_ADA)); 3779566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_ADA)); 3789566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRICES, (void (*)(void))MatCreateSubMatrices_ADA)); 3799566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_NORM, (void (*)(void))MatNorm_ADA)); 3809566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_ADA)); 3819566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRIX, (void (*)(void))MatCreateSubMatrix_ADA)); 382a7e14dcfSSatish Balay 3839566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)(*J), (PetscObject)ctx->W)); 3849566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)(*J))); 385a7e14dcfSSatish Balay 3869566063dSJacob Faibussowitsch PetscCall(MatSetOption(*J, MAT_SYMMETRIC, PETSC_TRUE)); 387a7e14dcfSSatish Balay PetscFunctionReturn(0); 388a7e14dcfSSatish Balay } 389