Actual source code: ex1.c
1: /*$Id: ex1.c,v 1.26 2001/08/07 03:04:16 balay Exp $*/
3: static char help[] = "Newton's method to solve a two-variable system, sequentially.nn";
5: /*T
6: Concepts: SNES^basic uniprocessor example
7: Processors: 1
8: T*/
10: /*
11: Include "petscsnes.h" so that we can use SNES solvers. Note that this
12: file automatically includes:
13: petsc.h - base PETSc routines petscvec.h - vectors
14: petscsys.h - system routines petscmat.h - matrices
15: petscis.h - index sets petscksp.h - Krylov subspace methods
16: petscviewer.h - viewers petscpc.h - preconditioners
17: petscsles.h - linear solvers
18: */
19: #include petscsnes.h
21: /*
22: User-defined routines
23: */
24: extern int FormJacobian1(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
25: extern int FormFunction1(SNES,Vec,Vec,void*);
26: extern int FormJacobian2(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
27: extern int FormFunction2(SNES,Vec,Vec,void*);
29: int main(int argc,char **argv)
30: {
31: SNES snes; /* nonlinear solver context */
32: SLES sles; /* linear solver context */
33: PC pc; /* preconditioner context */
34: KSP ksp; /* Krylov subspace method context */
35: Vec x,r; /* solution, residual vectors */
36: Mat J; /* Jacobian matrix */
37: int ierr,its,size;
38: PetscScalar pfive = .5,*xx;
39: PetscTruth flg;
41: PetscInitialize(&argc,&argv,(char *)0,help);
42: MPI_Comm_size(PETSC_COMM_WORLD,&size);
43: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Create nonlinear solver context
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: SNESCreate(PETSC_COMM_WORLD,&snes);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create matrix and vector data structures; set corresponding routines
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /*
56: Create vectors for solution and nonlinear function
57: */
58: VecCreateSeq(PETSC_COMM_SELF,2,&x);
59: VecDuplicate(x,&r);
61: /*
62: Create Jacobian matrix data structure
63: */
64: MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,2,2,&J);
65: MatSetFromOptions(J);
67: PetscOptionsHasName(PETSC_NULL,"-hard",&flg);
68: if (!flg) {
69: /*
70: Set function evaluation routine and vector.
71: */
72: SNESSetFunction(snes,r,FormFunction1,PETSC_NULL);
74: /*
75: Set Jacobian matrix data structure and Jacobian evaluation routine
76: */
77: SNESSetJacobian(snes,J,J,FormJacobian1,PETSC_NULL);
78: } else {
79: SNESSetFunction(snes,r,FormFunction2,PETSC_NULL);
80: SNESSetJacobian(snes,J,J,FormJacobian2,PETSC_NULL);
81: }
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Customize nonlinear solver; set runtime options
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: /*
88: Set linear solver defaults for this problem. By extracting the
89: SLES, KSP, and PC contexts from the SNES context, we can then
90: directly call any SLES, KSP, and PC routines to set various options.
91: */
92: SNESGetSLES(snes,&sles);
93: SLESGetKSP(sles,&ksp);
94: SLESGetPC(sles,&pc);
95: PCSetType(pc,PCNONE);
96: KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);
98: /*
99: Set SNES/SLES/KSP/PC runtime options, e.g.,
100: -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
101: These options will override those specified above as long as
102: SNESSetFromOptions() is called _after_ any other customization
103: routines.
104: */
105: SNESSetFromOptions(snes);
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Evaluate initial guess; then solve nonlinear system
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: if (!flg) {
111: VecSet(&pfive,x);
112: } else {
113: VecGetArray(x,&xx);
114: xx[0] = 2.0; xx[1] = 3.0;
115: VecRestoreArray(x,&xx);
116: }
117: /*
118: Note: The user should initialize the vector, x, with the initial guess
119: for the nonlinear solver prior to calling SNESSolve(). In particular,
120: to employ an initial guess of zero, the user should explicitly set
121: this vector to zero by calling VecSet().
122: */
124: SNESSolve(snes,x,&its);
125: if (flg) {
126: Vec f;
127: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
128: SNESGetFunction(snes,&f,0,0);
129: VecView(r,PETSC_VIEWER_STDOUT_WORLD);
130: }
132: PetscPrintf(PETSC_COMM_SELF,"number of Newton iterations = %dnn",its);
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Free work space. All PETSc objects should be destroyed when they
136: are no longer needed.
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: VecDestroy(x); VecDestroy(r);
140: MatDestroy(J); SNESDestroy(snes);
142: PetscFinalize();
143: return 0;
144: }
145: /* ------------------------------------------------------------------- */
146: /*
147: FormFunction1 - Evaluates nonlinear function, F(x).
149: Input Parameters:
150: . snes - the SNES context
151: . x - input vector
152: . dummy - optional user-defined context (not used here)
154: Output Parameter:
155: . f - function vector
156: */
157: int FormFunction1(SNES snes,Vec x,Vec f,void *dummy)
158: {
159: int ierr;
160: PetscScalar *xx,*ff;
162: /*
163: Get pointers to vector data.
164: - For default PETSc vectors, VecGetArray() returns a pointer to
165: the data array. Otherwise, the routine is implementation dependent.
166: - You MUST call VecRestoreArray() when you no longer need access to
167: the array.
168: */
169: VecGetArray(x,&xx);
170: VecGetArray(f,&ff);
172: /*
173: Compute function
174: */
175: ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0;
176: ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0;
178: /*
179: Restore vectors
180: */
181: VecRestoreArray(x,&xx);
182: VecRestoreArray(f,&ff);
184: return 0;
185: }
186: /* ------------------------------------------------------------------- */
187: /*
188: FormJacobian1 - Evaluates Jacobian matrix.
190: Input Parameters:
191: . snes - the SNES context
192: . x - input vector
193: . dummy - optional user-defined context (not used here)
195: Output Parameters:
196: . jac - Jacobian matrix
197: . B - optionally different preconditioning matrix
198: . flag - flag indicating matrix structure
199: */
200: int FormJacobian1(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
201: {
202: PetscScalar *xx,A[4];
203: int ierr,idx[2] = {0,1};
205: /*
206: Get pointer to vector data
207: */
208: VecGetArray(x,&xx);
210: /*
211: Compute Jacobian entries and insert into matrix.
212: - Since this is such a small problem, we set all entries for
213: the matrix at once.
214: */
215: A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0];
216: A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1];
217: MatSetValues(*jac,2,idx,2,idx,A,INSERT_VALUES);
218: *flag = SAME_NONZERO_PATTERN;
220: /*
221: Restore vector
222: */
223: VecRestoreArray(x,&xx);
225: /*
226: Assemble matrix
227: */
228: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
229: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
231: return 0;
232: }
235: /* ------------------------------------------------------------------- */
236: int FormFunction2(SNES snes,Vec x,Vec f,void *dummy)
237: {
238: int ierr;
239: PetscScalar *xx,*ff;
241: /*
242: Get pointers to vector data.
243: - For default PETSc vectors, VecGetArray() returns a pointer to
244: the data array. Otherwise, the routine is implementation dependent.
245: - You MUST call VecRestoreArray() when you no longer need access to
246: the array.
247: */
248: VecGetArray(x,&xx);
249: VecGetArray(f,&ff);
251: /*
252: Compute function
253: */
254: ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0];
255: ff[1] = xx[1];
257: /*
258: Restore vectors
259: */
260: VecRestoreArray(x,&xx);
261: VecRestoreArray(f,&ff);
263: return 0;
264: }
265: /* ------------------------------------------------------------------- */
266: int FormJacobian2(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
267: {
268: PetscScalar *xx,A[4];
269: int ierr,idx[2] = {0,1};
271: /*
272: Get pointer to vector data
273: */
274: VecGetArray(x,&xx);
276: /*
277: Compute Jacobian entries and insert into matrix.
278: - Since this is such a small problem, we set all entries for
279: the matrix at once.
280: */
281: A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0;
282: A[2] = 0.0; A[3] = 1.0;
283: MatSetValues(*jac,2,idx,2,idx,A,INSERT_VALUES);
284: *flag = SAME_NONZERO_PATTERN;
286: /*
287: Restore vector
288: */
289: VecRestoreArray(x,&xx);
291: /*
292: Assemble matrix
293: */
294: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
295: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
297: return 0;
298: }