xref: /petsc/src/mat/tests/ex97.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, non-square\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
59371c9d4SSatish Balay static PetscErrorCode AssembleMatrix(MPI_Comm comm, Mat *A) {
6c4762a1bSJed Brown   Mat      B;
7c4762a1bSJed Brown   PetscInt i, ms, me;
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &B));
119566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, 5, 6, PETSC_DETERMINE, PETSC_DETERMINE));
129566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
139566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
149566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &ms, &me));
15*48a46eb9SPierre Jolivet   for (i = ms; i < me; i++) PetscCall(MatSetValue(B, i, i, 1.0 * i, INSERT_VALUES));
169566063dSJacob Faibussowitsch   PetscCall(MatSetValue(B, me - 1, me, me * me, INSERT_VALUES));
179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
19c4762a1bSJed Brown   *A = B;
20c4762a1bSJed Brown   PetscFunctionReturn(0);
21c4762a1bSJed Brown }
22c4762a1bSJed Brown 
239371c9d4SSatish Balay static PetscErrorCode Compare2(Vec *X, const char *test) {
24c4762a1bSJed Brown   PetscReal norm;
25c4762a1bSJed Brown   Vec       Y;
26c4762a1bSJed Brown   PetscInt  verbose = 0;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   PetscFunctionBegin;
299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X[0], &Y));
309566063dSJacob Faibussowitsch   PetscCall(VecCopy(X[0], Y));
319566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y, -1.0, X[1]));
329566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_INFINITY, &norm));
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-verbose", &verbose, NULL));
35c4762a1bSJed Brown   if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) {
369566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%30s: norm difference < sqrt(eps_machine)\n", test));
37c4762a1bSJed Brown   } else {
389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%30s: norm difference %g\n", test, (double)norm));
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown   if (verbose > 1) {
419566063dSJacob Faibussowitsch     PetscCall(VecView(X[0], PETSC_VIEWER_STDOUT_WORLD));
429566063dSJacob Faibussowitsch     PetscCall(VecView(X[1], PETSC_VIEWER_STDOUT_WORLD));
439566063dSJacob Faibussowitsch     PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD));
44c4762a1bSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
46c4762a1bSJed Brown   PetscFunctionReturn(0);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
499371c9d4SSatish Balay static PetscErrorCode CheckMatrices(Mat A, Mat B, Vec left, Vec right, Vec X, Vec Y, Vec X1, Vec Y1) {
50c4762a1bSJed Brown   Vec *ltmp, *rtmp;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(right, 2, &rtmp));
549566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(left, 2, &ltmp));
559566063dSJacob Faibussowitsch   PetscCall(MatScale(A, PETSC_PI));
569566063dSJacob Faibussowitsch   PetscCall(MatScale(B, PETSC_PI));
579566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A, left, right));
589566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(B, left, right));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(MatMult(A, X, ltmp[0]));
619566063dSJacob Faibussowitsch   PetscCall(MatMult(B, X, ltmp[1]));
629566063dSJacob Faibussowitsch   PetscCall(Compare2(ltmp, "MatMult"));
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, Y, rtmp[0]));
659566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(B, Y, rtmp[1]));
669566063dSJacob Faibussowitsch   PetscCall(Compare2(rtmp, "MatMultTranspose"));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch   PetscCall(VecCopy(Y1, ltmp[0]));
699566063dSJacob Faibussowitsch   PetscCall(VecCopy(Y1, ltmp[1]));
709566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A, X, ltmp[0], ltmp[0]));
719566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(B, X, ltmp[1], ltmp[1]));
729566063dSJacob Faibussowitsch   PetscCall(Compare2(ltmp, "MatMultAdd v2==v3"));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A, X, Y1, ltmp[0]));
759566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(B, X, Y1, ltmp[1]));
769566063dSJacob Faibussowitsch   PetscCall(Compare2(ltmp, "MatMultAdd v2!=v3"));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, rtmp[0]));
799566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, rtmp[1]));
809566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A, Y, rtmp[0], rtmp[0]));
819566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(B, Y, rtmp[1], rtmp[1]));
829566063dSJacob Faibussowitsch   PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2==v3"));
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A, Y, X1, rtmp[0]));
859566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(B, Y, X1, rtmp[1]));
869566063dSJacob Faibussowitsch   PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2!=v3"));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(2, &ltmp));
899566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(2, &rtmp));
90c4762a1bSJed Brown   PetscFunctionReturn(0);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
939371c9d4SSatish Balay int main(int argc, char *argv[]) {
94c4762a1bSJed Brown   Mat       A, B, Asub, Bsub;
95c4762a1bSJed Brown   PetscInt  ms, idxrow[3], idxcol[4];
96c4762a1bSJed Brown   Vec       left, right, X, Y, X1, Y1;
97c4762a1bSJed Brown   IS        isrow, iscol;
98c4762a1bSJed Brown   PetscBool random = PETSC_TRUE;
99c4762a1bSJed Brown 
100327415f7SBarry Smith   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1029566063dSJacob Faibussowitsch   PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &A));
1039566063dSJacob Faibussowitsch   PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &B));
1049566063dSJacob Faibussowitsch   PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRIX, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRICES, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &ms, NULL));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   idxrow[0] = ms + 1;
109c4762a1bSJed Brown   idxrow[1] = ms + 2;
110c4762a1bSJed Brown   idxrow[2] = ms + 4;
1119566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 3, idxrow, PETSC_USE_POINTER, &isrow));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   idxcol[0] = ms + 1;
114c4762a1bSJed Brown   idxcol[1] = ms + 2;
115c4762a1bSJed Brown   idxcol[2] = ms + 4;
116c4762a1bSJed Brown   idxcol[3] = ms + 5;
1179566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 4, idxcol, PETSC_USE_POINTER, &iscol));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isrow, iscol, MAT_INITIAL_MATRIX, &Asub));
1209566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(B, isrow, iscol, MAT_INITIAL_MATRIX, &Bsub));
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Asub, &right, &left));
1239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(right, &X));
1249566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(right, &X1));
1259566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(left, &Y));
1269566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(left, &Y1));
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-random", &random, NULL));
129c4762a1bSJed Brown   if (random) {
1309566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(right, NULL));
1319566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(left, NULL));
1329566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X, NULL));
1339566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Y, NULL));
1349566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X1, NULL));
1359566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Y1, NULL));
136c4762a1bSJed Brown   } else {
1379566063dSJacob Faibussowitsch     PetscCall(VecSet(right, 1.0));
1389566063dSJacob Faibussowitsch     PetscCall(VecSet(left, 2.0));
1399566063dSJacob Faibussowitsch     PetscCall(VecSet(X, 3.0));
1409566063dSJacob Faibussowitsch     PetscCall(VecSet(Y, 4.0));
1419566063dSJacob Faibussowitsch     PetscCall(VecSet(X1, 3.0));
1429566063dSJacob Faibussowitsch     PetscCall(VecSet(Y1, 4.0));
143c4762a1bSJed Brown   }
1449566063dSJacob Faibussowitsch   PetscCall(CheckMatrices(Asub, Bsub, left, right, X, Y, X1, Y1));
1459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
1469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
1479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asub));
1509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bsub));
1519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
1529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
1559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y1));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
158b122ec5aSJacob Faibussowitsch   return 0;
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /*TEST
162c4762a1bSJed Brown 
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown       nsize: 3
165c4762a1bSJed Brown 
166c4762a1bSJed Brown TEST*/
167