xref: /petsc/src/mat/tests/ex53.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ format.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #define IMAX 15
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat                A, B, C, At, Bt;
8c4762a1bSJed Brown   PetscViewer        fd;
9c4762a1bSJed Brown   char               file[PETSC_MAX_PATH_LEN];
10c4762a1bSJed Brown   PetscRandom        rand;
11c4762a1bSJed Brown   Vec                xx, yy, s1, s2;
12c4762a1bSJed Brown   PetscReal          s1norm, s2norm, rnorm, tol = 1.e-10;
13c4762a1bSJed Brown   PetscInt           rstart, rend, rows[2], cols[2], m, n, i, j, M, N, ct, row, ncols1, ncols2, bs;
14c4762a1bSJed Brown   PetscMPIInt        rank, size;
15c4762a1bSJed Brown   const PetscInt    *cols1, *cols2;
16c4762a1bSJed Brown   PetscScalar        vals1[4], vals2[4], v;
17c4762a1bSJed Brown   const PetscScalar *v1, *v2;
18c4762a1bSJed Brown   PetscBool          flg;
19c4762a1bSJed Brown 
20327415f7SBarry Smith   PetscFunctionBeginUser;
21c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* Check out if MatLoad() works */
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
2728b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Input file not specified");
289566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
309566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATBAIJ));
319566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, fd));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
379566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
389566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
399566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &xx));
409566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(xx, m, PETSC_DECIDE));
419566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(xx));
429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(xx, &s1));
439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(xx, &s2));
449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(xx, &yy));
459566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Test MatNorm() */
489566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_FROBENIUS, &s1norm));
499566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_FROBENIUS, &s2norm));
50c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm - s1norm) / s2norm;
5148a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
529566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_INFINITY, &s1norm));
539566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_INFINITY, &s2norm));
54c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm - s1norm) / s2norm;
5548a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
569566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_1, &s1norm));
579566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_1, &s2norm));
58c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm - s1norm) / s2norm;
5948a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Test MatMult() */
62c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
639566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rand));
649566063dSJacob Faibussowitsch     PetscCall(MatMult(A, xx, s1));
659566063dSJacob Faibussowitsch     PetscCall(MatMult(B, xx, s2));
669566063dSJacob Faibussowitsch     PetscCall(VecAXPY(s2, -1.0, s1));
679566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &rnorm));
6848a46eb9SPierre Jolivet     if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatMult - Norm2=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)rnorm, bs));
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* test MatMultAdd() */
72c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
739566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rand));
749566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(yy, rand));
759566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A, xx, yy, s1));
769566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(B, xx, yy, s2));
779566063dSJacob Faibussowitsch     PetscCall(VecAXPY(s2, -1.0, s1));
789566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &rnorm));
7948a46eb9SPierre Jolivet     if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatMultAdd - Norm2=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)rnorm, bs));
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Test MatMultTranspose() */
83c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
849566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rand));
859566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, xx, s1));
869566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(B, xx, s2));
879566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
889566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
89c4762a1bSJed Brown     rnorm = s2norm - s1norm;
9048a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatMultTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
93c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
949566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rand));
959566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(yy, rand));
969566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(A, xx, yy, s1));
979566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(B, xx, yy, s2));
989566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
999566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
100c4762a1bSJed Brown     rnorm = s2norm - s1norm;
10148a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error: MatMultTransposeAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Check MatGetValues() */
1059566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
1069566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
109c4762a1bSJed Brown     /* Create random row numbers ad col numbers */
1109566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &v));
111c4762a1bSJed Brown     cols[0] = (int)(PetscRealPart(v) * N);
1129566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &v));
113c4762a1bSJed Brown     cols[1] = (int)(PetscRealPart(v) * N);
1149566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &v));
115c4762a1bSJed Brown     rows[0] = rstart + (int)(PetscRealPart(v) * m);
1169566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &v));
117c4762a1bSJed Brown     rows[1] = rstart + (int)(PetscRealPart(v) * m);
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch     PetscCall(MatGetValues(A, 2, rows, 2, cols, vals1));
1209566063dSJacob Faibussowitsch     PetscCall(MatGetValues(B, 2, rows, 2, cols, vals2));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown     for (j = 0; j < 4; j++) {
123c4762a1bSJed Brown       if (vals1[j] != vals2[j]) {
1249566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]: Error: MatGetValues rstart = %2" PetscInt_FMT "  row = %2" PetscInt_FMT " col = %2" PetscInt_FMT " val1 = %e val2 = %e bs = %" PetscInt_FMT "\n", rank, rstart, rows[j / 2], cols[j % 2], (double)PetscRealPart(vals1[j]), (double)PetscRealPart(vals2[j]), bs));
125c4762a1bSJed Brown       }
126c4762a1bSJed Brown     }
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Test MatGetRow()/ MatRestoreRow() */
130c4762a1bSJed Brown   for (ct = 0; ct < 100; ct++) {
1319566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &v));
1328fffc762SJacob Faibussowitsch     row = rstart + (PetscInt)(PetscRealPart(v) * m);
1339566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, row, &ncols1, &cols1, &v1));
1349566063dSJacob Faibussowitsch     PetscCall(MatGetRow(B, row, &ncols2, &cols2, &v2));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown     for (i = 0, j = 0; i < ncols1 && j < ncols2; j++) {
137c4762a1bSJed Brown       while (cols2[j] != cols1[i]) i++;
13808401ef6SPierre Jolivet       PetscCheck(v1[i] == v2[j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetRow() failed - vals incorrect.");
139c4762a1bSJed Brown     }
14008401ef6SPierre Jolivet     PetscCheck(j >= ncols2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetRow() failed - cols incorrect");
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, row, &ncols1, &cols1, &v1));
1439566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(B, row, &ncols2, &cols2, &v2));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Test MatConvert() */
1479566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSAME, MAT_INITIAL_MATRIX, &C));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* See if MatMult Says both are same */
150c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
1519566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rand));
1529566063dSJacob Faibussowitsch     PetscCall(MatMult(A, xx, s1));
1539566063dSJacob Faibussowitsch     PetscCall(MatMult(C, xx, s2));
1549566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
1559566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
156c4762a1bSJed Brown     rnorm = s2norm - s1norm;
15748a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error in MatConvert: MatMult - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
158c4762a1bSJed Brown   }
1599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Test MatTranspose() */
1629566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &At));
1639566063dSJacob Faibussowitsch   PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &Bt));
164c4762a1bSJed Brown   for (i = 0; i < IMAX; i++) {
1659566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rand));
1669566063dSJacob Faibussowitsch     PetscCall(MatMult(At, xx, s1));
1679566063dSJacob Faibussowitsch     PetscCall(MatMult(Bt, xx, s2));
1689566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
1699566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
170c4762a1bSJed Brown     rnorm = s2norm - s1norm;
17148a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Error in MatConvert:MatMult - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", rank, (double)s1norm, (double)s2norm, bs));
172c4762a1bSJed Brown   }
1739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&At));
1749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bt));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xx));
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yy));
1809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
1819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
1829566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
184b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*TEST
188c4762a1bSJed Brown 
189c4762a1bSJed Brown    build:
190c4762a1bSJed Brown       requires: !complex
191c4762a1bSJed Brown 
192c4762a1bSJed Brown    test:
193dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
194c4762a1bSJed Brown       nsize: 3
195c4762a1bSJed Brown       args: -matload_block_size 1 -f ${DATAFILESPATH}/matrices/small
196*3886731fSPierre Jolivet       output_file: output/empty.out
197c4762a1bSJed Brown 
198c4762a1bSJed Brown    test:
199c4762a1bSJed Brown       suffix: 2
200dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
201c4762a1bSJed Brown       nsize: 3
202c4762a1bSJed Brown       args: -matload_block_size 2 -f ${DATAFILESPATH}/matrices/small
203*3886731fSPierre Jolivet       output_file: output/empty.out
204c4762a1bSJed Brown 
205c4762a1bSJed Brown    test:
206c4762a1bSJed Brown       suffix: 3
207dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
208c4762a1bSJed Brown       nsize: 3
209c4762a1bSJed Brown       args: -matload_block_size 4 -f ${DATAFILESPATH}/matrices/small
210*3886731fSPierre Jolivet       output_file: output/empty.out
211c4762a1bSJed Brown 
212c4762a1bSJed Brown    test:
213c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
214c4762a1bSJed Brown       suffix: 4
215dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
216c4762a1bSJed Brown       nsize: 3
217c4762a1bSJed Brown       args: -matload_block_size 5 -f ${DATAFILESPATH}/matrices/small
218*3886731fSPierre Jolivet       output_file: output/empty.out
219c4762a1bSJed Brown 
220c4762a1bSJed Brown    test:
221c4762a1bSJed Brown       suffix: 5
222dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
223c4762a1bSJed Brown       nsize: 3
224c4762a1bSJed Brown       args: -matload_block_size 6 -f ${DATAFILESPATH}/matrices/small
225*3886731fSPierre Jolivet       output_file: output/empty.out
226c4762a1bSJed Brown 
227c4762a1bSJed Brown    test:
228c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
229c4762a1bSJed Brown       suffix: 6
230dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
231c4762a1bSJed Brown       nsize: 3
232c4762a1bSJed Brown       args: -matload_block_size 7 -f ${DATAFILESPATH}/matrices/small
233*3886731fSPierre Jolivet       output_file: output/empty.out
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    test:
236c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
237c4762a1bSJed Brown       suffix: 7
238dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
239c4762a1bSJed Brown       nsize: 3
240c4762a1bSJed Brown       args: -matload_block_size 8 -f ${DATAFILESPATH}/matrices/small
241*3886731fSPierre Jolivet       output_file: output/empty.out
242c4762a1bSJed Brown 
243c4762a1bSJed Brown    test:
244c4762a1bSJed Brown       suffix: 8
245dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
246c4762a1bSJed Brown       nsize: 4
247c4762a1bSJed Brown       args: -matload_block_size 3 -f ${DATAFILESPATH}/matrices/small
248*3886731fSPierre Jolivet       output_file: output/empty.out
249c4762a1bSJed Brown 
250c4762a1bSJed Brown TEST*/
251