1c4762a1bSJed Brown static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, non-square\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
AssembleMatrix(MPI_Comm comm,Mat * A)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, 5, 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, 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
Compare2(Vec * X,const char * test)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
CheckMatrices(Mat A,Mat B,Vec left,Vec right,Vec X,Vec Y,Vec X1,Vec Y1)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));
62c4762a1bSJed Brown
639566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, ltmp[0]));
649566063dSJacob Faibussowitsch PetscCall(MatMult(B, X, ltmp[1]));
659566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp, "MatMult"));
66c4762a1bSJed Brown
679566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Y, rtmp[0]));
689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Y, rtmp[1]));
699566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp, "MatMultTranspose"));
70c4762a1bSJed Brown
719566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1, ltmp[0]));
729566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1, ltmp[1]));
739566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, X, ltmp[0], ltmp[0]));
749566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B, X, ltmp[1], ltmp[1]));
759566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp, "MatMultAdd v2==v3"));
76c4762a1bSJed Brown
779566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, X, Y1, ltmp[0]));
789566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B, X, Y1, ltmp[1]));
799566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp, "MatMultAdd v2!=v3"));
80c4762a1bSJed Brown
819566063dSJacob Faibussowitsch PetscCall(VecCopy(X1, rtmp[0]));
829566063dSJacob Faibussowitsch PetscCall(VecCopy(X1, rtmp[1]));
839566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, Y, rtmp[0], rtmp[0]));
849566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B, Y, rtmp[1], rtmp[1]));
859566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2==v3"));
86c4762a1bSJed Brown
879566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, Y, X1, rtmp[0]));
889566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B, Y, X1, rtmp[1]));
899566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2!=v3"));
90c4762a1bSJed Brown
919566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2, <mp));
929566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2, &rtmp));
93*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown
main(int argc,char * argv[])96d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
97d71ae5a4SJacob Faibussowitsch {
98c4762a1bSJed Brown Mat A, B, Asub, Bsub;
99c4762a1bSJed Brown PetscInt ms, idxrow[3], idxcol[4];
100c4762a1bSJed Brown Vec left, right, X, Y, X1, Y1;
101c4762a1bSJed Brown IS isrow, iscol;
102c4762a1bSJed Brown PetscBool random = PETSC_TRUE;
103c4762a1bSJed Brown
104327415f7SBarry Smith PetscFunctionBeginUser;
1059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1069566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &A));
1079566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &B));
1089566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRIX, NULL));
1099566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRICES, NULL));
1109566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &ms, NULL));
111c4762a1bSJed Brown
112c4762a1bSJed Brown idxrow[0] = ms + 1;
113c4762a1bSJed Brown idxrow[1] = ms + 2;
114c4762a1bSJed Brown idxrow[2] = ms + 4;
1159566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 3, idxrow, PETSC_USE_POINTER, &isrow));
116c4762a1bSJed Brown
117c4762a1bSJed Brown idxcol[0] = ms + 1;
118c4762a1bSJed Brown idxcol[1] = ms + 2;
119c4762a1bSJed Brown idxcol[2] = ms + 4;
120c4762a1bSJed Brown idxcol[3] = ms + 5;
1219566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 4, idxcol, PETSC_USE_POINTER, &iscol));
122c4762a1bSJed Brown
1239566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isrow, iscol, MAT_INITIAL_MATRIX, &Asub));
1249566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B, isrow, iscol, MAT_INITIAL_MATRIX, &Bsub));
125c4762a1bSJed Brown
1269566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Asub, &right, &left));
1279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &X));
1289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &X1));
1299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &Y));
1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &Y1));
131c4762a1bSJed Brown
1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-random", &random, NULL));
133c4762a1bSJed Brown if (random) {
1349566063dSJacob Faibussowitsch PetscCall(VecSetRandom(right, NULL));
1359566063dSJacob Faibussowitsch PetscCall(VecSetRandom(left, NULL));
1369566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X, NULL));
1379566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y, NULL));
1389566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X1, NULL));
1399566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y1, NULL));
140c4762a1bSJed Brown } else {
1419566063dSJacob Faibussowitsch PetscCall(VecSet(right, 1.0));
1429566063dSJacob Faibussowitsch PetscCall(VecSet(left, 2.0));
1439566063dSJacob Faibussowitsch PetscCall(VecSet(X, 3.0));
1449566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 4.0));
1459566063dSJacob Faibussowitsch PetscCall(VecSet(X1, 3.0));
1469566063dSJacob Faibussowitsch PetscCall(VecSet(Y1, 4.0));
147c4762a1bSJed Brown }
1489566063dSJacob Faibussowitsch PetscCall(CheckMatrices(Asub, Bsub, left, right, X, Y, X1, Y1));
1499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow));
1509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol));
1519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asub));
1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bsub));
1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left));
1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right));
1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1));
1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1));
1619566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
162b122ec5aSJacob Faibussowitsch return 0;
163c4762a1bSJed Brown }
164c4762a1bSJed Brown
165c4762a1bSJed Brown /*TEST
166c4762a1bSJed Brown
167c4762a1bSJed Brown test:
168c4762a1bSJed Brown nsize: 3
169c4762a1bSJed Brown
170c4762a1bSJed Brown TEST*/
171