Actual source code: ex7.c
1: /*$Id: ex7.c,v 1.3 2001/04/14 03:58:23 bsmith Exp $*/
3: /* Program usage: mpirun -np <procs> ex5 [-help] [all PETSc options] */
5: static char help[] = "Nonlinear, time-dependent PDE in 2d.n";
8: /*
9: Include "petscda.h" so that we can use distributed arrays (DAs).
10: Include "petscts.h" so that we can use SNES solvers. Note that this
11: file 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: petscsles.h - linear solvers
17: */
18: #include petscda.h
19: #include petscts.h
22: /*
23: User-defined routines
24: */
25: extern int FormFunction(TS,double,Vec,Vec,void*),FormInitialSolution(DA,Vec);
26: extern int Monitor(TS,int,double,Vec,void*);
28: int main(int argc,char **argv)
29: {
30: TS ts; /* nonlinear solver */
31: Vec x,r; /* solution, residual vectors */
32: Mat J; /* Jacobian matrix */
33: int steps; /* iterations for convergence */
34: int ierr;
35: DA da;
36: MatFDColoring matfdcoloring;
37: ISColoring iscoloring;
38: double ftime;
40: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: Initialize program
42: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: PetscInitialize(&argc,&argv,(char *)0,help);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Create distributed array (DA) to manage parallel grid and vectors
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,
51: 1,1,PETSC_NULL,PETSC_NULL,&da);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Extract global vectors from DA; then duplicate for remaining
55: vectors that are the same types
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: DACreateGlobalVector(da,&x);
58: VecDuplicate(x,&r);
60: TSCreate(PETSC_COMM_WORLD,TS_NONLINEAR,&ts);
61: TSSetRHSFunction(ts,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: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: DAGetColoring(da,IS_COLORING_GLOBAL,MATMPIAIJ,&iscoloring,&J);
77: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
78: ISColoringDestroy(iscoloring);
79: MatFDColoringSetFunction(matfdcoloring,(int (*)(void))FormFunction,da);
80: MatFDColoringSetFromOptions(matfdcoloring);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Customize nonlinear solver; set runtime options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: TSSetInitialTimeStep(ts,0.0,.0001);
86: TSSetType(ts,TS_BEULER);
87: TSSetDuration(ts,100,1.0);
88: TSSetFromOptions(ts);
89: TSSetMonitor(ts,Monitor,0,0);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Set initial conditions
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: FormInitialSolution(da,x);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Solve nonlinear system
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: TSSetSolution(ts,x);
100: TSStep(ts,&steps,&ftime);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Free work space. All PETSc objects should be destroyed when they
105: are no longer needed.
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: MatDestroy(J);
109: MatFDColoringDestroy(matfdcoloring);
110: VecDestroy(x);
111: VecDestroy(r);
112: TSDestroy(ts);
113: DADestroy(da);
114: PetscFinalize();
116: return(0);
117: }
118: /* ------------------------------------------------------------------- */
119: /*
120: FormFunction - Evaluates nonlinear function, F(x).
122: Input Parameters:
123: . ts - the TS context
124: . X - input vector
125: . ptr - optional user-defined context, as set by SNESSetFunction()
127: Output Parameter:
128: . F - function vector
129: */
130: int FormFunction(TS ts,double time,Vec X,Vec F,void *ptr)
131: {
132: DA da = (DA)ptr;
133: int ierr,i,j,Mx,My,xs,ys,xm,ym;
134: double two = 2.0,hx,hy,hxdhy,hydhx,sx,sy;
135: Scalar u,uxx,uyy,**x,**f;
136: Vec localX;
139: DAGetLocalVector(da,&localX);
140: DAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
141: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
143: hx = 1.0/(double)(Mx-1); sx = 1.0/(hx*hx);
144: hy = 1.0/(double)(My-1); sy = 1.0/(hy*hy);
145: hxdhy = hx/hy;
146: hydhx = hy/hx;
148: /*
149: Scatter ghost points to local vector,using the 2-step process
150: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
151: By placing code between these two statements, computations can be
152: done while messages are in transition.
153: */
154: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
155: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
157: /*
158: Get pointers to vector data
159: */
160: DAVecGetArray(da,localX,(void**)&x);
161: DAVecGetArray(da,F,(void**)&f);
163: /*
164: Get local grid boundaries
165: */
166: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
168: /*
169: Compute function over the locally owned part of the grid
170: */
171: for (j=ys; j<ys+ym; j++) {
172: for (i=xs; i<xs+xm; i++) {
173: if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
174: f[j][i] = x[j][i];
175: continue;
176: }
177: u = x[j][i];
178: uxx = (two*u - x[j][i-1] - x[j][i+1])*sx;
179: uyy = (two*u - x[j-1][i] - x[j+1][i])*sy;
180: /* f[j][i] = -(uxx + uyy); */
181: 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 +
182: (x[j+1][i] - x[j-1][i])*(x[j+1][i] - x[j-1][i])*.25*sy);
183: }
184: }
186: /*
187: Restore vectors
188: */
189: DAVecRestoreArray(da,localX,(void**)&x);
190: DAVecRestoreArray(da,F,(void**)&f);
191: DARestoreLocalVector(da,&localX);
192: PetscLogFlops(11*ym*xm);
193: return(0);
194: }
196: /* ------------------------------------------------------------------- */
197: int FormInitialSolution(DA da,Vec U)
198: {
199: int ierr,i,j,xs,ys,xm,ym,Mx,My;
200: Scalar **u;
201: double hx,hy,x,y,r;
204: DAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
205: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
207: hx = 1.0/(double)(Mx-1);
208: hy = 1.0/(double)(My-1);
210: /*
211: Get pointers to vector data
212: */
213: DAVecGetArray(da,U,(void**)&u);
215: /*
216: Get local grid boundaries
217: */
218: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
220: /*
221: Compute function over the locally owned part of the grid
222: */
223: for (j=ys; j<ys+ym; j++) {
224: y = j*hy;
225: for (i=xs; i<xs+xm; i++) {
226: x = i*hx;
227: r = PetscSqrtScalar((x-.5)*(x-.5) + (y-.5)*(y-.5));
228: if (r < .125) {
229: u[j][i] = PetscExpScalar(-30.0*r*r*r);
230: } else {
231: u[j][i] = 0.0;
232: }
233: }
234: }
236: /*
237: Restore vectors
238: */
239: DAVecRestoreArray(da,U,(void**)&u);
240: return(0);
241: }
243: int Monitor(TS ts,int step,double ptime,Vec v,void *ctx)
244: {
245: int ierr;
246: double norm;
247: MPI_Comm comm;
250: VecNorm(v,NORM_2,&norm);
251: PetscObjectGetComm((PetscObject)ts,&comm);
252: PetscPrintf(comm,"timestep %d time %g norm %gn",step,ptime,norm);
253: return(0);
254: }