xref: /petsc/src/mat/tests/ex90.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Tests MatPtAP() \n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*
7c4762a1bSJed Brown  * This is an extremely simple example to test MatPtAP. It is very useful when developing and debugging the code.
8c4762a1bSJed Brown  *
9c4762a1bSJed Brown  * A =
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    1   2   0   4
12c4762a1bSJed Brown    0   1   2   0
13c4762a1bSJed Brown    2   0   4   0
14c4762a1bSJed Brown    0   1   2   1
15c4762a1bSJed Brown  *
16c4762a1bSJed Brown  *
17c4762a1bSJed Brown  *
18c4762a1bSJed Brown  * P =
19c4762a1bSJed Brown 
20c4762a1bSJed Brown    1.00000   0.00000
21c4762a1bSJed Brown    0.30000   0.50000
22c4762a1bSJed Brown    0.00000   0.80000
23c4762a1bSJed Brown    0.90000   0.00000
24c4762a1bSJed Brown  *
25c4762a1bSJed Brown  *
26c4762a1bSJed Brown  *AP =
27c4762a1bSJed Brown  *
28c4762a1bSJed Brown  * 5.20000   1.00000
29c4762a1bSJed Brown    0.30000   2.10000
30c4762a1bSJed Brown    2.00000   3.20000
31c4762a1bSJed Brown    1.20000   2.10000
32c4762a1bSJed Brown  *
33c4762a1bSJed Brown  * PT =
34c4762a1bSJed Brown 
35c4762a1bSJed Brown    1.00000   0.30000   0.00000   0.90000
36c4762a1bSJed Brown    0.00000   0.50000   0.80000   0.00000
37c4762a1bSJed Brown 
38c4762a1bSJed Brown  *
39c4762a1bSJed Brown  * C =
40c4762a1bSJed Brown 
41c4762a1bSJed Brown    6.3700   3.5200
42c4762a1bSJed Brown    1.7500   3.6100
43c4762a1bSJed Brown  *
44c4762a1bSJed Brown  * */
45c4762a1bSJed Brown 
46c4762a1bSJed Brown int main(int argc,char **argv)
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   PetscErrorCode ierr;
49c4762a1bSJed Brown   Mat            A,P,PtAP;
50c4762a1bSJed Brown   PetscInt       i1[] = {0, 3, 5}, i2[] = {0,2,5};
51c4762a1bSJed Brown   PetscInt       j1[] = {0, 1, 3, 1, 2}, j2[] = {0, 2, 1, 2, 3};
52c4762a1bSJed Brown   PetscScalar    a1[] = {1, 2, 4, 1, 2}, a2[] = {2, 4, 1, 2, 1};
53c4762a1bSJed Brown   PetscInt       pi1[] = {0,1,3}, pi2[] = {0,1,2};
54c4762a1bSJed Brown   PetscInt       pj1[] = {0, 0, 1}, pj2[] = {1,0};
55c4762a1bSJed Brown   PetscScalar    pa1[] = {1, 0.3, 0.5}, pa2[] = {0.8, 0.9};
56c4762a1bSJed Brown   MPI_Comm       comm;
57c4762a1bSJed Brown   PetscMPIInt    rank,size;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
60c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
61ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
62ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
63*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 2,comm,PETSC_ERR_ARG_INCOMP,"You have to use two processor cores to run this example ");
64c4762a1bSJed Brown   ierr = MatCreateMPIAIJWithArrays(comm,2,2,PETSC_DETERMINE,PETSC_DETERMINE,rank? i2:i1,rank? j2:j1,rank? a2:a1,&A);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = MatCreateMPIAIJWithArrays(comm,2,1,PETSC_DETERMINE,PETSC_DETERMINE,rank? pi2:pi1,rank? pj2:pj1,rank? pa2:pa1,&P);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,1.1,&PtAP);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = MatView(P,NULL);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = MatView(PtAP,NULL);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,1.1,&PtAP);CHKERRQ(ierr);
71c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = MatView(P,NULL);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = MatView(PtAP,NULL);CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = MatDestroy(&P);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatDestroy(&PtAP);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = PetscFinalize();
78c4762a1bSJed Brown   return ierr;
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown /*TEST
82c4762a1bSJed Brown    test:
83c4762a1bSJed Brown      nsize: 2
84c4762a1bSJed Brown      args:   -matptap_via allatonce
85c4762a1bSJed Brown      output_file: output/ex90_1.out
86c4762a1bSJed Brown 
87c4762a1bSJed Brown    test:
88c4762a1bSJed Brown      nsize: 2
89c4762a1bSJed Brown      suffix: merged
90c4762a1bSJed Brown      args:   -matptap_via allatonce_merged
91c4762a1bSJed Brown      output_file: output/ex90_1.out
92c4762a1bSJed Brown 
93c4762a1bSJed Brown TEST*/
94