xref: /petsc/src/mat/tests/ex102.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateLRC()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Vec            x,b,c=NULL;
9d751fb92SStefano Zampini   Mat            A,U,V,LR,X,LRe;
10d751fb92SStefano Zampini   PetscInt       M = 5, N = 7;
11c4762a1bSJed Brown   PetscBool      flg;
12c4762a1bSJed Brown 
13*327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&M,NULL));
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19c4762a1bSJed Brown          Create the sparse matrix
20c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
229566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
239566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(A,"A_"));
249566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
259566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A,5,NULL));
269566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
289566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A,NULL));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31c4762a1bSJed Brown          Create the dense matrices
32c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
339566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&U));
349566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(U,PETSC_DECIDE,PETSC_DECIDE,M,3));
359566063dSJacob Faibussowitsch   PetscCall(MatSetType(U,MATDENSE));
369566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(U,"U_"));
379566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(U));
389566063dSJacob Faibussowitsch   PetscCall(MatSetUp(U));
399566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(U,NULL));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&V));
429566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(V,PETSC_DECIDE,PETSC_DECIDE,N,3));
439566063dSJacob Faibussowitsch   PetscCall(MatSetType(V,MATDENSE));
449566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(V,"V_"));
459566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(V));
469566063dSJacob Faibussowitsch   PetscCall(MatSetUp(V));
479566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(V,NULL));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50c4762a1bSJed Brown          Create a vector to hold the diagonal of C
51d751fb92SStefano Zampini          A sequential vector can be created as well on each process
52d751fb92SStefano Zampini          It is user responsibility to ensure the data in the vector
53d751fb92SStefano Zampini          is consistent across processors
54c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-use_c",&flg));
56c4762a1bSJed Brown   if (flg) {
579566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,3,&c));
589566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(c,NULL));
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62c4762a1bSJed Brown          Create low rank correction matrix
63c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-low_rank",&flg));
65c4762a1bSJed Brown   if (flg) {
66c4762a1bSJed Brown     /* create a low-rank matrix, with no A-matrix */
679566063dSJacob Faibussowitsch     PetscCall(MatCreateLRC(NULL,U,c,V,&LR));
689566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
69c4762a1bSJed Brown   } else {
709566063dSJacob Faibussowitsch     PetscCall(MatCreateLRC(A,U,c,V,&LR));
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74d751fb92SStefano Zampini          Create the low rank correction matrix explicitly to check for
75d751fb92SStefano Zampini          correctness
76d751fb92SStefano Zampini      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
779566063dSJacob Faibussowitsch   PetscCall(MatHermitianTranspose(V,MAT_INITIAL_MATRIX,&X));
789566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(X,c,NULL));
799566063dSJacob Faibussowitsch   PetscCall(MatMatMult(U,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&LRe));
809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
811baa6e33SBarry Smith   if (A) PetscCall(MatAYPX(LRe,1.0,A,DIFFERENT_NONZERO_PATTERN));
82d751fb92SStefano Zampini 
83d751fb92SStefano Zampini   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84c4762a1bSJed Brown          Create test vectors
85c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
869566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(LR,&x,&b));
879566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x,NULL));
889566063dSJacob Faibussowitsch   PetscCall(MatMult(LR,x,b));
899566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(LR,b,x));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
92d751fb92SStefano Zampini 
93d751fb92SStefano Zampini   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94d751fb92SStefano Zampini          Check correctness
95d751fb92SStefano Zampini      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
969566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(LR,LRe,10,&flg));
979566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error in MatMult\n"));
98d751fb92SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
999566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeEqual(LR,LRe,10,&flg));
1009566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error in MatMultTranspose\n"));
101d751fb92SStefano Zampini #endif
102d751fb92SStefano Zampini 
1039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&LRe));
1059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&U));
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&V));
1079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&c));
1089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&LR));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /*
111c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
112c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
113c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
114c4762a1bSJed Brown          options are chosen (e.g., -log_view).
115c4762a1bSJed Brown   */
1169566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
117b122ec5aSJacob Faibussowitsch   return 0;
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown /*TEST
121c4762a1bSJed Brown 
122d751fb92SStefano Zampini    testset:
123d751fb92SStefano Zampini       output_file: output/ex102_1.out
124d751fb92SStefano Zampini       nsize: {{1 2}}
125d751fb92SStefano Zampini       args: -low_rank {{0 1}} -use_c {{0 1}}
126c4762a1bSJed Brown       test:
127d751fb92SStefano Zampini         suffix: standard
128c4762a1bSJed Brown       test:
129d751fb92SStefano Zampini         suffix: cuda
130d751fb92SStefano Zampini         requires: cuda
131d751fb92SStefano Zampini         args: -A_mat_type aijcusparse -U_mat_type densecuda -V_mat_type densecuda
132c4762a1bSJed Brown 
133c4762a1bSJed Brown TEST*/
134