1c4762a1bSJed Brown static char help[] = "Tests MatCreateLRC()\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Vec x, b, c = NULL;
8d751fb92SStefano Zampini Mat A, U, V, LR, X, LRe;
9d751fb92SStefano Zampini PetscInt M = 5, N = 7;
10c4762a1bSJed Brown PetscBool flg;
11c4762a1bSJed Brown
12327415f7SBarry Smith PetscFunctionBeginUser;
13c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &M, NULL));
159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL));
16c4762a1bSJed Brown
17c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18c4762a1bSJed Brown Create the sparse matrix
19c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
229566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "A_"));
239566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
249566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL));
259566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 5, NULL, 5, NULL));
269566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
279566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A, NULL));
28c4762a1bSJed Brown
29c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30c4762a1bSJed Brown Create the dense matrices
31c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
329566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &U));
339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(U, PETSC_DECIDE, PETSC_DECIDE, M, 3));
349566063dSJacob Faibussowitsch PetscCall(MatSetType(U, MATDENSE));
359566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(U, "U_"));
369566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(U));
379566063dSJacob Faibussowitsch PetscCall(MatSetUp(U));
389566063dSJacob Faibussowitsch PetscCall(MatSetRandom(U, NULL));
39c4762a1bSJed Brown
409566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &V));
419566063dSJacob Faibussowitsch PetscCall(MatSetSizes(V, PETSC_DECIDE, PETSC_DECIDE, N, 3));
429566063dSJacob Faibussowitsch PetscCall(MatSetType(V, MATDENSE));
439566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(V, "V_"));
449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(V));
459566063dSJacob Faibussowitsch PetscCall(MatSetUp(V));
469566063dSJacob Faibussowitsch PetscCall(MatSetRandom(V, NULL));
47c4762a1bSJed Brown
48c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49c4762a1bSJed Brown Create a vector to hold the diagonal of C
50d751fb92SStefano Zampini A sequential vector can be created as well on each process
51d751fb92SStefano Zampini It is user responsibility to ensure the data in the vector
52d751fb92SStefano Zampini is consistent across processors
53c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
549566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-use_c", &flg));
55c4762a1bSJed Brown if (flg) {
5677433607SBarry Smith PetscCall(VecCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, PETSC_DECIDE, 3, &c));
579566063dSJacob Faibussowitsch PetscCall(VecSetRandom(c, NULL));
58c4762a1bSJed Brown }
59c4762a1bSJed Brown
60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61c4762a1bSJed Brown Create low rank correction matrix
62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
639566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-low_rank", &flg));
64c4762a1bSJed Brown if (flg) {
65c4762a1bSJed Brown /* create a low-rank matrix, with no A-matrix */
669566063dSJacob Faibussowitsch PetscCall(MatCreateLRC(NULL, U, c, V, &LR));
679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
68c4762a1bSJed Brown } else {
699566063dSJacob Faibussowitsch PetscCall(MatCreateLRC(A, U, c, V, &LR));
70c4762a1bSJed Brown }
71c4762a1bSJed Brown
72c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73d751fb92SStefano Zampini Create the low rank correction matrix explicitly to check for
74d751fb92SStefano Zampini correctness
75d751fb92SStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
769566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(V, MAT_INITIAL_MATRIX, &X));
779566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(X, c, NULL));
78fb842aefSJose E. Roman PetscCall(MatMatMult(U, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &LRe));
799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X));
801baa6e33SBarry Smith if (A) PetscCall(MatAYPX(LRe, 1.0, A, DIFFERENT_NONZERO_PATTERN));
81d751fb92SStefano Zampini
82d751fb92SStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown Create test vectors
84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
859566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(LR, &x, &b));
869566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, NULL));
879566063dSJacob Faibussowitsch PetscCall(MatMult(LR, x, b));
889566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(LR, b, x));
899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
91d751fb92SStefano Zampini
92d751fb92SStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93d751fb92SStefano Zampini Check correctness
94d751fb92SStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
959566063dSJacob Faibussowitsch PetscCall(MatMultEqual(LR, LRe, 10, &flg));
969566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error in MatMult\n"));
97d751fb92SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
989566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeEqual(LR, LRe, 10, &flg));
999566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error in MatMultTranspose\n"));
100d751fb92SStefano Zampini #endif
101d751fb92SStefano Zampini
1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&LRe));
1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&U));
1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&V));
1069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c));
1079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&LR));
108c4762a1bSJed Brown
109c4762a1bSJed Brown /*
110c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine
111c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI
112c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime
113c4762a1bSJed Brown options are chosen (e.g., -log_view).
114c4762a1bSJed Brown */
1159566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
116b122ec5aSJacob Faibussowitsch return 0;
117c4762a1bSJed Brown }
118c4762a1bSJed Brown
119c4762a1bSJed Brown /*TEST
120c4762a1bSJed Brown
121d751fb92SStefano Zampini testset:
122*3886731fSPierre Jolivet output_file: output/empty.out
123d751fb92SStefano Zampini nsize: {{1 2}}
124d751fb92SStefano Zampini args: -low_rank {{0 1}} -use_c {{0 1}}
125c4762a1bSJed Brown test:
126d751fb92SStefano Zampini suffix: standard
127c4762a1bSJed Brown test:
128d751fb92SStefano Zampini suffix: cuda
129d751fb92SStefano Zampini requires: cuda
130d751fb92SStefano Zampini args: -A_mat_type aijcusparse -U_mat_type densecuda -V_mat_type densecuda
131c4762a1bSJed Brown
132c4762a1bSJed Brown TEST*/
133