xref: /petsc/src/mat/tests/ex177.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests various routines in MatKAIJ format.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #define IMAX 15
5c4762a1bSJed Brown 
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat          A, B, TA;
9c4762a1bSJed Brown   PetscScalar *S, *T;
10c4762a1bSJed Brown   PetscViewer  fd;
11c4762a1bSJed Brown   char         file[PETSC_MAX_PATH_LEN];
12c4762a1bSJed Brown   PetscInt     m, n, M, N, p = 1, q = 1, i, j;
13c4762a1bSJed Brown   PetscMPIInt  rank, size;
14c4762a1bSJed Brown   PetscBool    flg;
15c4762a1bSJed Brown 
16327415f7SBarry Smith   PetscFunctionBeginUser;
17c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
20c4762a1bSJed Brown 
21910cf402Sprj-   /* Load AIJ matrix A */
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
2328b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
269566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, fd));
279566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Get dof, then create S and T */
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-q", &q, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(p * q, &S, p * q, &T));
33c4762a1bSJed Brown   for (i = 0; i < p * q; i++) S[i] = 0;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   for (i = 0; i < p; i++) {
36c4762a1bSJed Brown     for (j = 0; j < q; j++) {
37910cf402Sprj-       /* Set some random non-zero values */
38c4762a1bSJed Brown       S[i + p * j] = ((PetscReal)((i + 1) * (j + 1))) / ((PetscReal)(p + q));
39c4762a1bSJed Brown       T[i + p * j] = ((PetscReal)((p - i) + j)) / ((PetscReal)(p * q));
40c4762a1bSJed Brown     }
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Test KAIJ when both S & T are not NULL */
44c4762a1bSJed Brown 
45910cf402Sprj-   /* Create KAIJ matrix TA */
469566063dSJacob Faibussowitsch   PetscCall(MatCreateKAIJ(A, p, q, S, T, &TA));
479566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
489566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatConvert(TA, MATAIJ, MAT_INITIAL_MATRIX, &B));
51c4762a1bSJed Brown 
52910cf402Sprj-   /* Test MatKAIJGetScaledIdentity() */
539566063dSJacob Faibussowitsch   PetscCall(MatKAIJGetScaledIdentity(TA, &flg));
5428b400f6SJacob Faibussowitsch   PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in Test 1: MatKAIJGetScaledIdentity()");
55c4762a1bSJed Brown   /* Test MatMult() */
569566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(TA, B, 10, &flg));
5728b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 1: MatMult() for KAIJ matrix");
58c4762a1bSJed Brown   /* Test MatMultAdd() */
599566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(TA, B, 10, &flg));
6028b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 1: MatMultAdd() for KAIJ matrix");
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&TA));
639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Test KAIJ when S is NULL */
66c4762a1bSJed Brown 
67910cf402Sprj-   /* Create KAIJ matrix TA */
689566063dSJacob Faibussowitsch   PetscCall(MatCreateKAIJ(A, p, q, NULL, T, &TA));
699566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(MatConvert(TA, MATAIJ, MAT_INITIAL_MATRIX, &B));
73c4762a1bSJed Brown 
74910cf402Sprj-   /* Test MatKAIJGetScaledIdentity() */
759566063dSJacob Faibussowitsch   PetscCall(MatKAIJGetScaledIdentity(TA, &flg));
7628b400f6SJacob Faibussowitsch   PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in Test 2: MatKAIJGetScaledIdentity()");
77c4762a1bSJed Brown   /* Test MatMult() */
789566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(TA, B, 10, &flg));
7928b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 2: MatMult() for KAIJ matrix");
80c4762a1bSJed Brown   /* Test MatMultAdd() */
819566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(TA, B, 10, &flg));
8228b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 2: MatMultAdd() for KAIJ matrix");
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&TA));
859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Test KAIJ when T is NULL */
88c4762a1bSJed Brown 
89910cf402Sprj-   /* Create KAIJ matrix TA */
909566063dSJacob Faibussowitsch   PetscCall(MatCreateKAIJ(A, p, q, S, NULL, &TA));
919566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
929566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(MatConvert(TA, MATAIJ, MAT_INITIAL_MATRIX, &B));
95c4762a1bSJed Brown 
96910cf402Sprj-   /* Test MatKAIJGetScaledIdentity() */
979566063dSJacob Faibussowitsch   PetscCall(MatKAIJGetScaledIdentity(TA, &flg));
9828b400f6SJacob Faibussowitsch   PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in Test 3: MatKAIJGetScaledIdentity()");
99c4762a1bSJed Brown   /* Test MatMult() */
1009566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(TA, B, 10, &flg));
10128b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 3: MatMult() for KAIJ matrix");
102c4762a1bSJed Brown   /* Test MatMultAdd() */
1039566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(TA, B, 10, &flg));
10428b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 3: MatMultAdd() for KAIJ matrix");
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&TA));
1079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
108c4762a1bSJed Brown 
10915229ffcSPierre Jolivet   /* Test KAIJ when T is an identity matrix */
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   if (p == q) {
112c4762a1bSJed Brown     for (i = 0; i < p; i++) {
113c4762a1bSJed Brown       for (j = 0; j < q; j++) {
114c4762a1bSJed Brown         if (i == j) T[i + j * p] = 1.0;
115c4762a1bSJed Brown         else T[i + j * p] = 0.0;
116c4762a1bSJed Brown       }
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown 
119910cf402Sprj-     /* Create KAIJ matrix TA */
1209566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(A, p, q, S, T, &TA));
1219566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, &n));
1229566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, &N));
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch     PetscCall(MatConvert(TA, MATAIJ, MAT_INITIAL_MATRIX, &B));
125c4762a1bSJed Brown 
126910cf402Sprj-     /* Test MatKAIJGetScaledIdentity() */
1279566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetScaledIdentity(TA, &flg));
12828b400f6SJacob Faibussowitsch     PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in Test 4: MatKAIJGetScaledIdentity()");
129c4762a1bSJed Brown     /* Test MatMult() */
1309566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(TA, B, 10, &flg));
13128b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 4: MatMult() for KAIJ matrix");
132c4762a1bSJed Brown     /* Test MatMultAdd() */
1339566063dSJacob Faibussowitsch     PetscCall(MatMultAddEqual(TA, B, 10, &flg));
13428b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error in Test 4: MatMultAdd() for KAIJ matrix");
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&TA));
1379566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
138910cf402Sprj- 
1399566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(A, p, q, NULL, T, &TA));
140910cf402Sprj-     /* Test MatKAIJGetScaledIdentity() */
1419566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetScaledIdentity(TA, &flg));
14228b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in Test 5: MatKAIJGetScaledIdentity()");
1439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&TA));
144910cf402Sprj- 
145910cf402Sprj-     for (i = 0; i < p; i++) {
146910cf402Sprj-       for (j = 0; j < q; j++) {
147910cf402Sprj-         if (i == j) S[i + j * p] = T[i + j * p] = 2.0;
148910cf402Sprj-         else S[i + j * p] = T[i + j * p] = 0.0;
149910cf402Sprj-       }
150910cf402Sprj-     }
1519566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(A, p, q, S, T, &TA));
152910cf402Sprj-     /* Test MatKAIJGetScaledIdentity() */
1539566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetScaledIdentity(TA, &flg));
15428b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in Test 6: MatKAIJGetScaledIdentity()");
1559566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&TA));
156c4762a1bSJed Brown   }
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* Done with all tests */
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch   PetscCall(PetscFree2(S, T));
1619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
163b122ec5aSJacob Faibussowitsch   return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown /*TEST
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   build:
169c4762a1bSJed Brown     requires: !complex
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   test:
172dfd57a17SPierre Jolivet     requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
173*3886731fSPierre Jolivet     output_file: output/empty.out
174c4762a1bSJed Brown     nsize: {{1 2 3 4}}
175c4762a1bSJed Brown     args: -f ${DATAFILESPATH}/matrices/small -p {{2 3 7}} -q {{3 7}} -viewer_binary_skip_info
176c4762a1bSJed Brown 
177c4762a1bSJed Brown TEST*/
178