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