xref: /petsc/src/mat/tests/ex102.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MatCreateLRC()\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*T
5*c4762a1bSJed Brown    Concepts: Low rank correction
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown    Processors: n
8*c4762a1bSJed Brown T*/
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown #include <petscmat.h>
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown int main(int argc,char **args)
13*c4762a1bSJed Brown {
14*c4762a1bSJed Brown   Vec            x,b,c=NULL;
15*c4762a1bSJed Brown   Mat            A,U,V,LR;
16*c4762a1bSJed Brown   PetscInt       i,j,Ii,J,Istart,Iend,m = 8,n = 7,rstart,rend;
17*c4762a1bSJed Brown   PetscErrorCode ierr;
18*c4762a1bSJed Brown   PetscBool      flg;
19*c4762a1bSJed Brown   PetscScalar    *u,a;
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
22*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26*c4762a1bSJed Brown          Create the sparse matrix
27*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
28*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
33*c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
34*c4762a1bSJed Brown     a = -1.0; i = Ii/n; j = Ii - i*n;
35*c4762a1bSJed Brown     if (i>0)   {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);CHKERRQ(ierr);}
36*c4762a1bSJed Brown     if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);CHKERRQ(ierr);}
37*c4762a1bSJed Brown     if (j>0)   {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);CHKERRQ(ierr);}
38*c4762a1bSJed Brown     if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);CHKERRQ(ierr);}
39*c4762a1bSJed Brown     a = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&a,INSERT_VALUES);CHKERRQ(ierr);
40*c4762a1bSJed Brown   }
41*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45*c4762a1bSJed Brown          Create the dense matrices
46*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&U);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = MatSetSizes(U,PETSC_DECIDE,PETSC_DECIDE,m*n,3);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = MatSetType(U,MATDENSE);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = MatSetUp(U);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(U,&rstart,&rend);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = MatDenseGetArray(U,&u);CHKERRQ(ierr);
53*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
54*c4762a1bSJed Brown     u[i-rstart]          = (PetscReal)i;
55*c4762a1bSJed Brown     u[i+rend-2*rstart]   = (PetscReal)1000*i;
56*c4762a1bSJed Brown     u[i+2*rend-3*rstart] = (PetscReal)100000*i;
57*c4762a1bSJed Brown   }
58*c4762a1bSJed Brown   ierr = MatDenseRestoreArray(U,&u);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatAssemblyBegin(U,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = MatAssemblyEnd(U,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&V);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = MatSetSizes(V,PETSC_DECIDE,PETSC_DECIDE,m*n,3);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = MatSetType(V,MATDENSE);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = MatSetUp(V);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(U,&rstart,&rend);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = MatDenseGetArray(V,&u);CHKERRQ(ierr);
68*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
69*c4762a1bSJed Brown     u[i-rstart]          = (PetscReal)i;
70*c4762a1bSJed Brown     u[i+rend-2*rstart]   = (PetscReal)1.2*i;
71*c4762a1bSJed Brown     u[i+2*rend-3*rstart] = (PetscReal)1.67*i+2;
72*c4762a1bSJed Brown   }
73*c4762a1bSJed Brown   ierr = MatDenseRestoreArray(V,&u);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = MatAssemblyBegin(V,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = MatAssemblyEnd(V,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78*c4762a1bSJed Brown          Create a vector to hold the diagonal of C
79*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-use_c",&flg);CHKERRQ(ierr);
81*c4762a1bSJed Brown   if (flg) {
82*c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,3,&c);CHKERRQ(ierr);
83*c4762a1bSJed Brown     ierr = VecGetArray(c,&u);CHKERRQ(ierr);
84*c4762a1bSJed Brown     u[0] = 2.0;
85*c4762a1bSJed Brown     u[1] = -1.0;
86*c4762a1bSJed Brown     u[2] = 1.0;
87*c4762a1bSJed Brown     ierr = VecRestoreArray(c,&u);CHKERRQ(ierr);
88*c4762a1bSJed Brown   }
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91*c4762a1bSJed Brown          Create low rank correction matrix
92*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-low_rank",&flg);CHKERRQ(ierr);
94*c4762a1bSJed Brown   if (flg) {
95*c4762a1bSJed Brown     /* create a low-rank matrix, with no A-matrix */
96*c4762a1bSJed Brown     ierr = MatCreateLRC(NULL,U,c,V,&LR);CHKERRQ(ierr);
97*c4762a1bSJed Brown   } else {
98*c4762a1bSJed Brown     ierr = MatCreateLRC(A,U,c,V,&LR);CHKERRQ(ierr);
99*c4762a1bSJed Brown   }
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102*c4762a1bSJed Brown          Create test vectors
103*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = VecSetSizes(x,PETSC_DECIDE,m*n);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = VecDuplicate(x,&b);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(x,&rstart,&rend);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = VecGetArray(x,&u);CHKERRQ(ierr);
110*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) u[i-rstart] = (PetscScalar)i;
111*c4762a1bSJed Brown   ierr = VecRestoreArray(x,&u);CHKERRQ(ierr);
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown   ierr = MatMult(LR,x,b);CHKERRQ(ierr);
114*c4762a1bSJed Brown   /*
115*c4762a1bSJed Brown      View the product if desired
116*c4762a1bSJed Brown   */
117*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-view_product",&flg);CHKERRQ(ierr);
118*c4762a1bSJed Brown   if (flg) {ierr = VecView(b,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
122*c4762a1bSJed Brown   /* you can destroy the matrices in any order you like */
123*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = MatDestroy(&U);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = MatDestroy(&V);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = VecDestroy(&c);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = MatDestroy(&LR);CHKERRQ(ierr);
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown   /*
130*c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
131*c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
132*c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
133*c4762a1bSJed Brown          options are chosen (e.g., -log_view).
134*c4762a1bSJed Brown   */
135*c4762a1bSJed Brown   ierr = PetscFinalize();
136*c4762a1bSJed Brown   return ierr;
137*c4762a1bSJed Brown }
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown /*TEST
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown    test:
143*c4762a1bSJed Brown       suffix: 1
144*c4762a1bSJed Brown       nsize: 2
145*c4762a1bSJed Brown       args: -view_product
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown    test:
148*c4762a1bSJed Brown       suffix: 2
149*c4762a1bSJed Brown       nsize: 2
150*c4762a1bSJed Brown       args: -low_rank -view_product
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown    test:
153*c4762a1bSJed Brown       suffix: 3
154*c4762a1bSJed Brown       nsize: 2
155*c4762a1bSJed Brown       args: -use_c -view_product
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown TEST*/
158