Actual source code: ex27.c
2: static char help[] = "Reads a PETSc matrix and vector from a file and solves the normal equations.\n\n";
3: /*T
4: Concepts: KSP^solving a linear system
5: Concepts: Normal equations
6: Processors: n
7: T*/
9: /*
10: Include "petscksp.h" so that we can use KSP solvers. Note that this file
11: automatically includes:
12: petsc.h - base PETSc routines petscvec.h - vectors
13: petscsys.h - system routines petscmat.h - matrices
14: petscis.h - index sets petscksp.h - Krylov subspace methods
15: petscviewer.h - viewers petscpc.h - preconditioners
16: */
17: #include petscksp.h
22: int main(int argc,char **args)
23: {
24: KSP ksp; /* linear solver context */
25: Mat A,N; /* matrix */
26: Vec x,b,u,Ab; /* approx solution, RHS, exact solution */
27: PetscViewer fd; /* viewer */
28: char file[PETSC_MAX_PATH_LEN]; /* input file name */
29: PetscErrorCode ierr,ierrp;
30: PetscInt its,n,m;
31: PetscReal norm;
32: PetscScalar zero = 0.0,none = -1.0;
34: PetscInitialize(&argc,&args,(char *)0,help);
37: /*
38: Determine files from which we read the linear system
39: (matrix and right-hand-side vector).
40: */
41: PetscOptionsGetString(PETSC_NULL,"-f",file,PETSC_MAX_PATH_LEN-1,PETSC_NULL);
43: /* -----------------------------------------------------------
44: Beginning of linear solver loop
45: ----------------------------------------------------------- */
46: /*
47: Loop through the linear solve 2 times.
48: - The intention here is to preload and solve a small system;
49: then load another (larger) system and solve it as well.
50: This process preloads the instructions with the smaller
51: system so that more accurate performance monitoring (via
52: -log_summary) can be done with the larger one (that actually
53: is the system of interest).
54: */
55: PreLoadBegin(PETSC_FALSE,"Load system");
57: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
58: Load system
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: /*
62: Open binary file. Note that we use PETSC_FILE_RDONLY to indicate
63: reading from this file.
64: */
65: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_FILE_RDONLY,&fd);
67: /*
68: Load the matrix and vector; then destroy the viewer.
69: */
70: MatLoad(fd,MATMPIAIJ,&A);
71: PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
72: ierrp = VecLoad(fd,PETSC_NULL,&b);
73: PetscPopErrorHandler();
74: MatGetLocalSize(A,&m,&n);
75: if (ierrp) { /* if file contains no RHS, then use a vector of all ones */
76: PetscScalar one = 1.0;
77: VecCreate(PETSC_COMM_WORLD,&b);
78: VecSetSizes(b,m,PETSC_DECIDE);
79: VecSetFromOptions(b);
80: VecSet(b,one);
81: }
82: PetscViewerDestroy(fd);
84: /*
85: If the loaded matrix is larger than the vector (due to being padded
86: to match the block size of the system), then create a new padded vector.
87: */
88: {
89: PetscInt j,mvec,start,end,idex;
90: Vec tmp;
91: PetscScalar *bold;
93: /* Create a new vector b by padding the old one */
94: VecCreate(PETSC_COMM_WORLD,&tmp);
95: VecSetSizes(tmp,m,PETSC_DECIDE);
96: VecSetFromOptions(tmp);
97: VecGetOwnershipRange(b,&start,&end);
98: VecGetLocalSize(b,&mvec);
99: VecGetArray(b,&bold);
100: for (j=0; j<mvec; j++) {
101: idex = start+j;
102: VecSetValues(tmp,1,&idex,bold+j,INSERT_VALUES);
103: }
104: VecRestoreArray(b,&bold);
105: VecDestroy(b);
106: VecAssemblyBegin(tmp);
107: VecAssemblyEnd(tmp);
108: b = tmp;
109: }
110: VecDuplicate(b,&u);
111: VecCreateMPI(PETSC_COMM_WORLD,n,PETSC_DECIDE,&x);
112: VecDuplicate(x,&Ab);
113: VecSet(x,zero);
115: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
116: Setup solve for system
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: /*
120: Conclude profiling last stage; begin profiling next stage.
121: */
122: PreLoadStage("KSPSetUp");
124: MatCreateNormal(A,&N);
125: MatMultTranspose(A,b,Ab);
127: /*
128: Create linear solver; set operators; set runtime options.
129: */
130: KSPCreate(PETSC_COMM_WORLD,&ksp);
132: KSPSetOperators(ksp,N,N,SAME_NONZERO_PATTERN);
133: KSPSetFromOptions(ksp);
135: /*
136: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
137: enable more precise profiling of setting up the preconditioner.
138: These calls are optional, since both will be called within
139: KSPSolve() if they haven't been called already.
140: */
141: KSPSetUp(ksp);
142: KSPSetUpOnBlocks(ksp);
144: /*
145: Solve system
146: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: /*
149: Begin profiling next stage
150: */
151: PreLoadStage("KSPSolve");
153: /*
154: Solve linear system
155: */
156: KSPSolve(ksp,Ab,x);
158: /*
159: Conclude profiling this stage
160: */
161: PreLoadStage("Cleanup");
163: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
164: Check error, print output, free data structures.
165: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167: /*
168: Check error
169: */
170: MatMult(A,x,u);
171: VecAXPY(u,none,b);
172: VecNorm(u,NORM_2,&norm);
173: KSPGetIterationNumber(ksp,&its);
174: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
175: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A\n",norm);
177: /*
178: Free work space. All PETSc objects should be destroyed when they
179: are no longer needed.
180: */
181: MatDestroy(A); VecDestroy(b);
182: MatDestroy(N); VecDestroy(Ab);
183: VecDestroy(u); VecDestroy(x);
184: KSPDestroy(ksp);
185: PreLoadEnd();
186: /* -----------------------------------------------------------
187: End of linear solver loop
188: ----------------------------------------------------------- */
190: PetscFinalize();
191: return 0;
192: }