Actual source code: ex7.c
2: static char help[] = "Nonlinear, time-dependent PDE in 2d.\n";
4: /*
5: Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
6: Include "petscts.h" so that we can use SNES solvers. Note that this
7: file automatically includes:
8: petscsys.h - base PETSc routines petscvec.h - vectors
9: petscmat.h - matrices
10: petscis.h - index sets petscksp.h - Krylov subspace methods
11: petscviewer.h - viewers petscpc.h - preconditioners
12: petscksp.h - linear solvers
13: */
14: #include <petscdm.h>
15: #include <petscdmda.h>
16: #include <petscts.h>
18: /*
19: User-defined routines
20: */
21: extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec);
22: extern PetscErrorCode MyTSMonitor(TS,PetscInt,PetscReal,Vec,void*);
23: extern PetscErrorCode MySNESMonitor(SNES,PetscInt,PetscReal,PetscViewerAndFormat*);
25: int main(int argc,char **argv)
26: {
27: TS ts; /* time integrator */
28: SNES snes;
29: Vec x,r; /* solution, residual vectors */
30: DM da;
31: PetscViewerAndFormat *vf;
33: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: Initialize program
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36: PetscInitialize(&argc,&argv,(char*)0,help);
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Create distributed array (DMDA) to manage parallel grid and vectors
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
41: DMSetFromOptions(da);
42: DMSetUp(da);
44: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: Extract global vectors from DMDA; then duplicate for remaining
46: vectors that are the same types
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48: DMCreateGlobalVector(da,&x);
49: VecDuplicate(x,&r);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create timestepping solver context
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: TSCreate(PETSC_COMM_WORLD,&ts);
55: TSSetProblemType(ts,TS_NONLINEAR);
56: TSSetRHSFunction(ts,NULL,FormFunction,da);
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Create matrix data structure; set Jacobian evaluation routine
61: Set Jacobian matrix data structure and default Jacobian evaluation
62: routine. User can override with:
63: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
64: (unless user explicitly sets preconditioner)
65: -snes_mf_operator : form preconditioning matrix as set by the user,
66: but use matrix-free approx for Jacobian-vector
67: products within Newton-Krylov method
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: TSSetMaxTime(ts,1.0);
72: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
73: TSMonitorSet(ts,MyTSMonitor,PETSC_VIEWER_STDOUT_WORLD,NULL);
74: TSSetDM(ts,da);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Customize nonlinear solver
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: TSSetType(ts,TSBEULER);
79: TSGetSNES(ts,&snes);
80: PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
81: SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Set initial conditions
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: FormInitialSolution(da,x);
87: TSSetTimeStep(ts,.0001);
88: TSSetSolution(ts,x);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Set runtime options
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: TSSetFromOptions(ts);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Solve nonlinear system
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: TSSolve(ts,x);
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Free work space. All PETSc objects should be destroyed when they
102: are no longer needed.
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: VecDestroy(&x);
105: VecDestroy(&r);
106: TSDestroy(&ts);
107: DMDestroy(&da);
109: PetscFinalize();
110: return 0;
111: }
112: /* ------------------------------------------------------------------- */
113: /*
114: FormFunction - Evaluates nonlinear function, F(x).
116: Input Parameters:
117: . ts - the TS context
118: . X - input vector
119: . ptr - optional user-defined context, as set by SNESSetFunction()
121: Output Parameter:
122: . F - function vector
123: */
124: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
125: {
126: DM da;
127: PetscInt i,j,Mx,My,xs,ys,xm,ym;
128: PetscReal two = 2.0,hx,hy,sx,sy;
129: PetscScalar u,uxx,uyy,**x,**f;
130: Vec localX;
133: TSGetDM(ts,&da);
134: DMGetLocalVector(da,&localX);
135: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
137: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
138: hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
140: /*
141: Scatter ghost points to local vector,using the 2-step process
142: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
143: By placing code between these two statements, computations can be
144: done while messages are in transition.
145: */
146: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
147: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
149: /*
150: Get pointers to vector data
151: */
152: DMDAVecGetArrayRead(da,localX,&x);
153: DMDAVecGetArray(da,F,&f);
155: /*
156: Get local grid boundaries
157: */
158: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
160: /*
161: Compute function over the locally owned part of the grid
162: */
163: for (j=ys; j<ys+ym; j++) {
164: for (i=xs; i<xs+xm; i++) {
165: if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
166: f[j][i] = x[j][i];
167: continue;
168: }
169: u = x[j][i];
170: uxx = (two*u - x[j][i-1] - x[j][i+1])*sx;
171: uyy = (two*u - x[j-1][i] - x[j+1][i])*sy;
172: /* f[j][i] = -(uxx + uyy); */
173: f[j][i] = -u*(uxx + uyy) - (4.0 - 1.0)*((x[j][i+1] - x[j][i-1])*(x[j][i+1] - x[j][i-1])*.25*sx +
174: (x[j+1][i] - x[j-1][i])*(x[j+1][i] - x[j-1][i])*.25*sy);
175: }
176: }
178: /*
179: Restore vectors
180: */
181: DMDAVecRestoreArrayRead(da,localX,&x);
182: DMDAVecRestoreArray(da,F,&f);
183: DMRestoreLocalVector(da,&localX);
184: PetscLogFlops(11.0*ym*xm);
185: return 0;
186: }
188: /* ------------------------------------------------------------------- */
189: PetscErrorCode FormInitialSolution(DM da,Vec U)
190: {
191: PetscInt i,j,xs,ys,xm,ym,Mx,My;
192: PetscScalar **u;
193: PetscReal hx,hy,x,y,r;
196: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
198: hx = 1.0/(PetscReal)(Mx-1);
199: hy = 1.0/(PetscReal)(My-1);
201: /*
202: Get pointers to vector data
203: */
204: DMDAVecGetArray(da,U,&u);
206: /*
207: Get local grid boundaries
208: */
209: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
211: /*
212: Compute function over the locally owned part of the grid
213: */
214: for (j=ys; j<ys+ym; j++) {
215: y = j*hy;
216: for (i=xs; i<xs+xm; i++) {
217: x = i*hx;
218: r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
219: if (r < .125) u[j][i] = PetscExpReal(-30.0*r*r*r);
220: else u[j][i] = 0.0;
221: }
222: }
224: /*
225: Restore vectors
226: */
227: DMDAVecRestoreArray(da,U,&u);
228: return 0;
229: }
231: PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
232: {
233: PetscReal norm;
234: MPI_Comm comm;
237: if (step < 0) return 0; /* step of -1 indicates an interpolated solution */
238: VecNorm(v,NORM_2,&norm);
239: PetscObjectGetComm((PetscObject)ts,&comm);
240: PetscPrintf(comm,"timestep %D time %g norm %g\n",step,(double)ptime,(double)norm);
241: return 0;
242: }
244: /*
245: MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
246: Input Parameters:
247: snes - the SNES context
248: its - iteration number
249: fnorm - 2-norm function value (may be estimated)
250: ctx - optional user-defined context for private data for the
251: monitor routine, as set by SNESMonitorSet()
252: */
253: PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf)
254: {
256: SNESMonitorDefaultShort(snes,its,fnorm,vf);
257: return 0;
258: }
260: /*TEST
262: test:
263: args: -ts_max_steps 5
265: test:
266: suffix: 2
267: args: -ts_max_steps 5 -snes_mf_operator
269: test:
270: suffix: 3
271: args: -ts_max_steps 5 -snes_mf -pc_type none
273: TEST*/