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