Actual source code: ex102.c
petsc-3.6.1 2015-08-06
2: static char help[] = "Tests MatCreateLRC()\n\n";
4: /*T
5: Concepts: Low rank correction
7: Processors: n
8: T*/
10: #include <petscmat.h>
14: int main(int argc,char **args)
15: {
16: Vec x,b;
17: Mat A,U,V,LR;
18: PetscInt i,j,Ii,J,Istart,Iend,m = 8,n = 7,rstart,rend;
20: PetscBool flg;
21: PetscScalar *u,a;
23: PetscInitialize(&argc,&args,(char*)0,help);
24: PetscOptionsGetInt(NULL,"-m",&m,NULL);
25: PetscOptionsGetInt(NULL,"-n",&n,NULL);
27: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28: Create the sparse matrix
29: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
30: MatCreate(PETSC_COMM_WORLD,&A);
31: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
32: MatSetFromOptions(A);
33: MatSetUp(A);
34: MatGetOwnershipRange(A,&Istart,&Iend);
35: for (Ii=Istart; Ii<Iend; Ii++) {
36: a = -1.0; i = Ii/n; j = Ii - i*n;
37: if (i>0) {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
38: if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
39: if (j>0) {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
40: if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
41: a = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&a,INSERT_VALUES);
42: }
43: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
44: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Create the dense matrices
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: MatCreate(PETSC_COMM_WORLD,&U);
50: MatSetSizes(U,PETSC_DECIDE,PETSC_DECIDE,m*n,3);
51: MatSetType(U,MATDENSE);
52: MatSetUp(U);
53: MatGetOwnershipRange(U,&rstart,&rend);
54: MatDenseGetArray(U,&u);
55: for (i=rstart; i<rend; i++) {
56: u[i-rstart] = (PetscReal)i;
57: u[i+rend-2*rstart] = (PetscReal)1000*i;
58: u[i+2*rend-3*rstart] = (PetscReal)100000*i;
59: }
60: MatDenseRestoreArray(U,&u);
61: MatAssemblyBegin(U,MAT_FINAL_ASSEMBLY);
62: MatAssemblyEnd(U,MAT_FINAL_ASSEMBLY);
65: MatCreate(PETSC_COMM_WORLD,&V);
66: MatSetSizes(V,PETSC_DECIDE,PETSC_DECIDE,m*n,3);
67: MatSetType(V,MATDENSE);
68: MatSetUp(V);
69: MatGetOwnershipRange(U,&rstart,&rend);
70: MatDenseGetArray(V,&u);
71: for (i=rstart; i<rend; i++) {
72: u[i-rstart] = (PetscReal)i;
73: u[i+rend-2*rstart] = (PetscReal)1.2*i;
74: u[i+2*rend-3*rstart] = (PetscReal)1.67*i+2;
75: }
76: MatDenseRestoreArray(V,&u);
77: MatAssemblyBegin(V,MAT_FINAL_ASSEMBLY);
78: MatAssemblyEnd(V,MAT_FINAL_ASSEMBLY);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Create low rank created matrix
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: MatCreateLRC(A,U,V,&LR);
83: MatSetUp(LR);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Create test vectors
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: VecCreate(PETSC_COMM_WORLD,&x);
89: VecSetSizes(x,PETSC_DECIDE,m*n);
90: VecSetFromOptions(x);
91: VecDuplicate(x,&b);
92: VecGetOwnershipRange(x,&rstart,&rend);
93: VecGetArray(x,&u);
94: for (i=rstart; i<rend; i++) u[i-rstart] = (PetscScalar)i;
95: VecRestoreArray(x,&u);
97: MatMult(LR,x,b);
98: /*
99: View the product if desired
100: */
101: PetscOptionsHasName(NULL,"-view_product",&flg);
102: if (flg) {VecView(b,PETSC_VIEWER_STDOUT_WORLD);}
104: VecDestroy(&x);
105: VecDestroy(&b);
106: /* you can destroy the matrices in any order you like */
107: MatDestroy(&A);
108: MatDestroy(&U);
109: MatDestroy(&V);
110: MatDestroy(&LR);
112: /*
113: Always call PetscFinalize() before exiting a program. This routine
114: - finalizes the PETSc libraries as well as MPI
115: - provides summary and diagnostic information if certain runtime
116: options are chosen (e.g., -log_summary).
117: */
118: PetscFinalize();
119: return 0;
120: }