Actual source code: ex23.c
1: /*$Id: ex23.c,v 1.6 2001/03/23 23:24:25 balay Exp $*/
3: static char help[] = "Solves PDE problem from ex22.cnn";
5: #include petscda.h
6: #include petscpf.h
7: #include petscsnes.h
9: /*
11: In this example the PDE is
12: Uxx = 2,
13: u(0) = .25
14: u(1) = 0
16: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
18: Use the usual centered finite differences.
20: Note we treat the problem as non-linear though it happens to be linear
22: See ex22.c for a design optimization problem
24: */
26: typedef struct {
27: PetscViewer u_viewer;
28: PetscViewer fu_viewer;
29: } UserCtx;
31: extern int FormFunction(SNES,Vec,Vec,void*);
33: int main(int argc,char **argv)
34: {
35: int ierr,N = 5;
36: UserCtx user;
37: DA da;
38: DMMG *dmmg;
40: PetscInitialize(&argc,&argv,PETSC_NULL,help);
41: PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);
43: /* Hardwire several options; can be changed at command line */
44: PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
45: PetscOptionsSetValue("-ksp_type","fgmres");
46: PetscOptionsSetValue("-ksp_max_it","5");
47: PetscOptionsSetValue("-pc_mg_type","full");
48: PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
49: PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
50: PetscOptionsSetValue("-mg_coarse_ksp_max_it","3");
51: PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
52: PetscOptionsSetValue("-snes_mf_type","wp");
53: PetscOptionsSetValue("-snes_mf_compute_norma","no");
54: PetscOptionsSetValue("-snes_mf_compute_normu","no");
55: PetscOptionsSetValue("-snes_eq_ls","basic");
56: /* PetscOptionsSetValue("-snes_eq_ls","basicnonorms"); */
57: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
58:
59: /* Create a global vector from a da arrays */
60: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,N,1,1,PETSC_NULL,&da);
62: /* create graphics windows */
63: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);
64: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - discretization of function",-1,-1,-1,-1,&user.fu_viewer);
66: /* create nonlinear multi-level solver */
67: DMMGCreate(PETSC_COMM_WORLD,2,&user,&dmmg);
68: DMMGSetUseMatrixFree(dmmg);
69: DMMGSetDM(dmmg,(DM)da);
70: DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
71: DMMGSolve(dmmg);
72: DMMGDestroy(dmmg);
74: DADestroy(da);
75: PetscViewerDestroy(user.u_viewer);
76: PetscViewerDestroy(user.fu_viewer);
78: PetscFinalize();
79: return 0;
80: }
81:
82: /*
83: This local function acts on the ghosted version of U (accessed via DAGetLocalVector())
84: BUT the global, nonghosted version of FU
86: */
87: int FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
88: {
89: DMMG dmmg = (DMMG)dummy;
90: int ierr,xs,xm,i,N;
91: Scalar *u,*fu,d,h;
92: Vec vu;
93: DA da = (DA) dmmg->dm;
96: DAGetLocalVector(da,&vu);
97: DAGlobalToLocalBegin(da,U,INSERT_VALUES,vu);
98: DAGlobalToLocalEnd(da,U,INSERT_VALUES,vu);
100: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
101: DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
102: DAVecGetArray(da,vu,(void**)&u);
103: DAVecGetArray(da,FU,(void**)&fu);
104: d = N-1.0;
105: h = 1.0/d;
107: for (i=xs; i<xs+xm; i++) {
108: if (i == 0) fu[0] = 2.0*d*(u[0] - .25) + h*u[0]*u[0];
109: else if (i == N-1) fu[N-1] = 2.0*d*u[N-1] + h*u[N-1]*u[N-1];
110: else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h) + h*u[i]*u[i];
111: }
113: DAVecRestoreArray(da,vu,(void**)&u);
114: DAVecRestoreArray(da,FU,(void**)&fu);
115: DARestoreLocalVector(da,&vu);
116: PetscLogFlops(9*N);
117: return(0);
118: }