1c4762a1bSJed Brown static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, square, shifted (copied from ex97)\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode AssembleMatrix(MPI_Comm comm, Mat *A) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown Mat B; 8c4762a1bSJed Brown PetscInt i, ms, me; 9c4762a1bSJed Brown 10c4762a1bSJed Brown PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &B)); 129566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, 6, 6, PETSC_DETERMINE, PETSC_DETERMINE)); 139566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 149566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B, &ms, &me)); 1648a46eb9SPierre Jolivet for (i = ms; i < me; i++) PetscCall(MatSetValue(B, i, i, 1.0 * i, INSERT_VALUES)); 179566063dSJacob Faibussowitsch PetscCall(MatSetValue(B, me - 1, me - 1, me * me, INSERT_VALUES)); 189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 20c4762a1bSJed Brown *A = B; 21*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24d71ae5a4SJacob Faibussowitsch static PetscErrorCode Compare2(Vec *X, const char *test) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown PetscReal norm; 27c4762a1bSJed Brown Vec Y; 28c4762a1bSJed Brown PetscInt verbose = 0; 29c4762a1bSJed Brown 30c4762a1bSJed Brown PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X[0], &Y)); 329566063dSJacob Faibussowitsch PetscCall(VecCopy(X[0], Y)); 339566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1.0, X[1])); 349566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_INFINITY, &norm)); 35c4762a1bSJed Brown 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-verbose", &verbose, NULL)); 37c4762a1bSJed Brown if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) { 389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%30s: norm difference < sqrt(eps_machine)\n", test)); 39c4762a1bSJed Brown } else { 409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%30s: norm difference %g\n", test, (double)norm)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown if (verbose > 1) { 439566063dSJacob Faibussowitsch PetscCall(VecView(X[0], PETSC_VIEWER_STDOUT_WORLD)); 449566063dSJacob Faibussowitsch PetscCall(VecView(X[1], PETSC_VIEWER_STDOUT_WORLD)); 459566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD)); 46c4762a1bSJed Brown } 479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 48*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckMatrices(Mat A, Mat B, Vec left, Vec right, Vec X, Vec Y, Vec X1, Vec Y1) 52d71ae5a4SJacob Faibussowitsch { 53c4762a1bSJed Brown Vec *ltmp, *rtmp; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(right, 2, &rtmp)); 579566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(left, 2, <mp)); 589566063dSJacob Faibussowitsch PetscCall(MatScale(A, PETSC_PI)); 599566063dSJacob Faibussowitsch PetscCall(MatScale(B, PETSC_PI)); 609566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, left, right)); 619566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B, left, right)); 629566063dSJacob Faibussowitsch PetscCall(MatShift(A, PETSC_PI)); 639566063dSJacob Faibussowitsch PetscCall(MatShift(B, PETSC_PI)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, ltmp[0])); 669566063dSJacob Faibussowitsch PetscCall(MatMult(B, X, ltmp[1])); 679566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp, "MatMult")); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Y, rtmp[0])); 709566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Y, rtmp[1])); 719566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp, "MatMultTranspose")); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1, ltmp[0])); 749566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1, ltmp[1])); 759566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, X, ltmp[0], ltmp[0])); 769566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B, X, ltmp[1], ltmp[1])); 779566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp, "MatMultAdd v2==v3")); 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, X, Y1, ltmp[0])); 809566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B, X, Y1, ltmp[1])); 819566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp, "MatMultAdd v2!=v3")); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(VecCopy(X1, rtmp[0])); 849566063dSJacob Faibussowitsch PetscCall(VecCopy(X1, rtmp[1])); 859566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, Y, rtmp[0], rtmp[0])); 869566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B, Y, rtmp[1], rtmp[1])); 879566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2==v3")); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, Y, X1, rtmp[0])); 909566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B, Y, X1, rtmp[1])); 919566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2!=v3")); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2, <mp)); 949566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2, &rtmp)); 95*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 99d71ae5a4SJacob Faibussowitsch { 100c4762a1bSJed Brown Mat A, B, Asub, Bsub; 101c4762a1bSJed Brown PetscInt ms, idxrow[3], idxcol[3]; 102c4762a1bSJed Brown Vec left, right, X, Y, X1, Y1; 103c4762a1bSJed Brown IS isrow, iscol; 104c4762a1bSJed Brown PetscBool random = PETSC_TRUE; 105c4762a1bSJed Brown 106327415f7SBarry Smith PetscFunctionBeginUser; 1079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1089566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &A)); 1099566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &B)); 1109566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRIX, NULL)); 1119566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRICES, NULL)); 1129566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &ms, NULL)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown idxrow[0] = ms + 1; 115c4762a1bSJed Brown idxrow[1] = ms + 2; 116c4762a1bSJed Brown idxrow[2] = ms + 4; 1179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 3, idxrow, PETSC_USE_POINTER, &isrow)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown idxcol[0] = ms + 1; 120c4762a1bSJed Brown idxcol[1] = ms + 2; 121c4762a1bSJed Brown idxcol[2] = ms + 4; 1229566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 3, idxcol, PETSC_USE_POINTER, &iscol)); 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isrow, iscol, MAT_INITIAL_MATRIX, &Asub)); 1259566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B, isrow, iscol, MAT_INITIAL_MATRIX, &Bsub)); 126c4762a1bSJed Brown 1279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Asub, &right, &left)); 1289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &X)); 1299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &X1)); 1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &Y)); 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &Y1)); 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-random", &random, NULL)); 134c4762a1bSJed Brown if (random) { 1359566063dSJacob Faibussowitsch PetscCall(VecSetRandom(right, NULL)); 1369566063dSJacob Faibussowitsch PetscCall(VecSetRandom(left, NULL)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X, NULL)); 1389566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y, NULL)); 1399566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X1, NULL)); 1409566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y1, NULL)); 141c4762a1bSJed Brown } else { 1429566063dSJacob Faibussowitsch PetscCall(VecSet(right, 1.0)); 1439566063dSJacob Faibussowitsch PetscCall(VecSet(left, 2.0)); 1449566063dSJacob Faibussowitsch PetscCall(VecSet(X, 3.0)); 1459566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 4.0)); 1469566063dSJacob Faibussowitsch PetscCall(VecSet(X1, 3.0)); 1479566063dSJacob Faibussowitsch PetscCall(VecSet(Y1, 4.0)); 148c4762a1bSJed Brown } 1499566063dSJacob Faibussowitsch PetscCall(CheckMatrices(Asub, Bsub, left, right, X, Y, X1, Y1)); 1509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 1519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asub)); 1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bsub)); 1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 163b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /*TEST 167c4762a1bSJed Brown 168c4762a1bSJed Brown test: 169c4762a1bSJed Brown nsize: 3 170c4762a1bSJed Brown 171c4762a1bSJed Brown TEST*/ 172