Actual source code: ex102.c
2: /* Program usage: mpirun -np <procs> ex2 [-help] [all PETSc options] */
4: static char help[] = "Tests MatCreateLRC()\n\n";
6: /*T
7: Concepts: Low rank correction
9: Processors: n
10: T*/
12: #include petscmat.h
16: int main(int argc,char **args)
17: {
18: Vec x,b;
19: Mat A,U,V,LR;
20: PetscInt i,j,Ii,J,Istart,Iend,m = 8,n = 7,rstart,rend;
22: PetscTruth flg;
23: PetscScalar *u,a;
25: PetscInitialize(&argc,&args,(char *)0,help);
26: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
27: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Create the sparse matrix
31: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32: MatCreate(PETSC_COMM_WORLD,&A);
33: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
34: MatSetFromOptions(A);
35: MatGetOwnershipRange(A,&Istart,&Iend);
36: for (Ii=Istart; Ii<Iend; Ii++) {
37: a = -1.0; i = Ii/n; j = Ii - i*n;
38: if (i>0) {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
39: if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
40: if (j>0) {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
41: if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&a,INSERT_VALUES);}
42: a = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&a,INSERT_VALUES);
43: }
44: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
45: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Create the dense matrices
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: MatCreate(PETSC_COMM_WORLD,&U);
51: MatSetSizes(U,PETSC_DECIDE,PETSC_DECIDE,m*n,3);
52: MatSetType(U,MATDENSE);
53: MatGetOwnershipRange(U,&rstart,&rend);
54: MatGetArray(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: MatRestoreArray(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: MatGetOwnershipRange(U,&rstart,&rend);
69: MatGetArray(V,&u);
70: for (i=rstart; i<rend; i++) {
71: u[i-rstart] = (PetscReal)i;
72: u[i+rend-2*rstart] = (PetscReal)1.2*i;
73: u[i+2*rend-3*rstart] = (PetscReal)1.67*i+2;
74: }
75: MatRestoreArray(V,&u);
76: MatAssemblyBegin(V,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(V,MAT_FINAL_ASSEMBLY);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Create low rank created matrix
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: MatCreateLRC(A,U,V,&LR);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Create test vectors
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: VecCreate(PETSC_COMM_WORLD,&x);
87: VecSetSizes(x,PETSC_DECIDE,m*n);
88: VecSetFromOptions(x);
89: VecDuplicate(x,&b);
90: VecGetOwnershipRange(x,&rstart,&rend);
91: VecGetArray(x,&u);
92: for (i=rstart; i<rend; i++) u[i-rstart] = (PetscScalar)i;
93: VecRestoreArray(x,&u);
95: MatMult(LR,x,b);
96: /*
97: View the product if desired
98: */
99: PetscOptionsHasName(PETSC_NULL,"-view_product",&flg);
100: if (flg) {VecView(b,PETSC_VIEWER_STDOUT_WORLD);}
102: VecDestroy(x);
103: VecDestroy(b);
104: /* you can destroy the matrices in any order you like */
105: MatDestroy(A);
106: MatDestroy(U);
107: MatDestroy(V);
108: MatDestroy(LR);
110: /*
111: Always call PetscFinalize() before exiting a program. This routine
112: - finalizes the PETSc libraries as well as MPI
113: - provides summary and diagnostic information if certain runtime
114: options are chosen (e.g., -log_summary).
115: */
116: PetscFinalize();
117: return 0;
118: }