xref: /petsc/src/snes/tests/ex7.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Solves u`` + u^{2} = f with Newton-like methods. Using\n\
3*c4762a1bSJed Brown  matrix-free techniques with user-provided explicit preconditioner matrix.\n\n";
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #include <petscsnes.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
8*c4762a1bSJed Brown extern PetscErrorCode FormJacobianNoMatrix(SNES,Vec,Mat,Mat,void*);
9*c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
10*c4762a1bSJed Brown extern PetscErrorCode FormFunctioni(void *,PetscInt,Vec,PetscScalar *);
11*c4762a1bSJed Brown extern PetscErrorCode OtherFunctionForDifferencing(void*,Vec,Vec);
12*c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(SNES,Vec);
13*c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown typedef struct {
16*c4762a1bSJed Brown   PetscViewer viewer;
17*c4762a1bSJed Brown } MonitorCtx;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown typedef struct {
20*c4762a1bSJed Brown   PetscBool variant;
21*c4762a1bSJed Brown } AppCtx;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown int main(int argc,char **argv)
24*c4762a1bSJed Brown {
25*c4762a1bSJed Brown   SNES           snes;                 /* SNES context */
26*c4762a1bSJed Brown   SNESType       type = SNESNEWTONLS;        /* default nonlinear solution method */
27*c4762a1bSJed Brown   Vec            x,r,F,U;              /* vectors */
28*c4762a1bSJed Brown   Mat            J,B;                  /* Jacobian matrix-free, explicit preconditioner */
29*c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
30*c4762a1bSJed Brown   PetscScalar    h,xp = 0.0,v;
31*c4762a1bSJed Brown   PetscInt       its,n = 5,i;
32*c4762a1bSJed Brown   PetscErrorCode ierr;
33*c4762a1bSJed Brown   PetscBool      puremf = PETSC_FALSE;
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
36*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-variant",&user.variant);CHKERRQ(ierr);
38*c4762a1bSJed Brown   h    = 1.0/(n-1);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   /* Set up data structures */
41*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&x);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = VecDuplicate(x,&F);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = VecDuplicate(x,&U);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   /* create explict matrix preconditioner */
49*c4762a1bSJed Brown   ierr         = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&B);CHKERRQ(ierr);
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   /* Store right-hand-side of PDE and exact solution */
52*c4762a1bSJed Brown   for (i=0; i<n; i++) {
53*c4762a1bSJed Brown     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
54*c4762a1bSJed Brown     ierr = VecSetValues(F,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
55*c4762a1bSJed Brown     v    = xp*xp*xp;
56*c4762a1bSJed Brown     ierr = VecSetValues(U,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
57*c4762a1bSJed Brown     xp  += h;
58*c4762a1bSJed Brown   }
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /* Create nonlinear solver */
61*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = SNESSetType(snes,type);CHKERRQ(ierr);
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown   /* Set various routines and options */
65*c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,F);CHKERRQ(ierr);
66*c4762a1bSJed Brown   if (user.variant) {
67*c4762a1bSJed Brown     /* this approach is not normally needed, one should use the MatCreateSNESMF() below usually */
68*c4762a1bSJed Brown     ierr = MatCreateMFFD(PETSC_COMM_WORLD,n,n,n,n,&J);CHKERRQ(ierr);
69*c4762a1bSJed Brown     ierr = MatMFFDSetFunction(J,(PetscErrorCode (*)(void*, Vec, Vec))SNESComputeFunction,snes);CHKERRQ(ierr);
70*c4762a1bSJed Brown     ierr = MatMFFDSetFunctioni(J,FormFunctioni);CHKERRQ(ierr);
71*c4762a1bSJed Brown     /* Use the matrix free operator for both the Jacobian used to define the linear system and used to define the preconditioner */
72*c4762a1bSJed Brown     /* This tests MatGetDiagonal() for MATMFFD */
73*c4762a1bSJed Brown     ierr = PetscOptionsHasName(NULL,NULL,"-puremf",&puremf);CHKERRQ(ierr);
74*c4762a1bSJed Brown   } else {
75*c4762a1bSJed Brown     /* create matrix free matrix for Jacobian */
76*c4762a1bSJed Brown     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
77*c4762a1bSJed Brown     /* demonstrates differencing a different function than FormFunction() to apply a matrix operator */
78*c4762a1bSJed Brown     /* note we use the same context for this function as FormFunction, the F vector */
79*c4762a1bSJed Brown     ierr = MatMFFDSetFunction(J,OtherFunctionForDifferencing,F);CHKERRQ(ierr);
80*c4762a1bSJed Brown   }
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown   /* Set various routines and options */
83*c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,puremf ? J : B,puremf ? FormJacobianNoMatrix : FormJacobian,&user);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   /* Solve nonlinear system */
87*c4762a1bSJed Brown   ierr = FormInitialGuess(snes,x);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   /* Free data structures */
93*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);  ierr = VecDestroy(&r);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);  ierr = VecDestroy(&F);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);  ierr = MatDestroy(&B);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = PetscFinalize();
98*c4762a1bSJed Brown   return ierr;
99*c4762a1bSJed Brown }
100*c4762a1bSJed Brown /* --------------------  Evaluate Function F(x) --------------------- */
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown PetscErrorCode  FormFunction(SNES snes,Vec x,Vec f,void *dummy)
103*c4762a1bSJed Brown {
104*c4762a1bSJed Brown   const PetscScalar *xx,*FF;
105*c4762a1bSJed Brown   PetscScalar       *ff,d;
106*c4762a1bSJed Brown   PetscInt          i,n;
107*c4762a1bSJed Brown   PetscErrorCode    ierr;
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown   ierr  = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr  = VecGetArray(f,&ff);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr  = VecGetArrayRead((Vec) dummy,&FF);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr  = VecGetSize(x,&n);CHKERRQ(ierr);
113*c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
114*c4762a1bSJed Brown   ff[0] = xx[0];
115*c4762a1bSJed Brown   for (i=1; i<n-1; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
116*c4762a1bSJed Brown   ff[n-1] = xx[n-1] - 1.0;
117*c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr    = VecRestoreArray(f,&ff);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr    = VecRestoreArrayRead((Vec)dummy,&FF);CHKERRQ(ierr);
120*c4762a1bSJed Brown   return 0;
121*c4762a1bSJed Brown }
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown PetscErrorCode  FormFunctioni(void *dummy,PetscInt i,Vec x,PetscScalar *s)
124*c4762a1bSJed Brown {
125*c4762a1bSJed Brown   const PetscScalar *xx,*FF;
126*c4762a1bSJed Brown   PetscScalar       d;
127*c4762a1bSJed Brown   PetscInt          n;
128*c4762a1bSJed Brown   PetscErrorCode    ierr;
129*c4762a1bSJed Brown   SNES              snes = (SNES) dummy;
130*c4762a1bSJed Brown   Vec               F;
131*c4762a1bSJed Brown 
132*c4762a1bSJed Brown   ierr  = SNESGetFunction(snes,NULL,NULL,(void**)&F);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr  = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr  = VecGetArrayRead(F,&FF);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr  = VecGetSize(x,&n);CHKERRQ(ierr);
136*c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
137*c4762a1bSJed Brown   if (i == 0) {
138*c4762a1bSJed Brown     *s = xx[0];
139*c4762a1bSJed Brown   } else if (i == n-1) {
140*c4762a1bSJed Brown     *s = xx[n-1] - 1.0;
141*c4762a1bSJed Brown   } else {
142*c4762a1bSJed Brown     *s = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
143*c4762a1bSJed Brown   }
144*c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(F,&FF);CHKERRQ(ierr);
146*c4762a1bSJed Brown   return 0;
147*c4762a1bSJed Brown }
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown /*
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown    Example function that when differenced produces the same matrix free Jacobian as FormFunction()
152*c4762a1bSJed Brown    this is provided to show how a user can provide a different function
153*c4762a1bSJed Brown */
154*c4762a1bSJed Brown PetscErrorCode  OtherFunctionForDifferencing(void *dummy,Vec x,Vec f)
155*c4762a1bSJed Brown {
156*c4762a1bSJed Brown   PetscErrorCode ierr;
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown   ierr = FormFunction(NULL,x,f,dummy);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = VecShift(f,1.0);CHKERRQ(ierr);
160*c4762a1bSJed Brown   return 0;
161*c4762a1bSJed Brown }
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown /* --------------------  Form initial approximation ----------------- */
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown PetscErrorCode  FormInitialGuess(SNES snes,Vec x)
166*c4762a1bSJed Brown {
167*c4762a1bSJed Brown   PetscErrorCode ierr;
168*c4762a1bSJed Brown   PetscScalar    pfive = .50;
169*c4762a1bSJed Brown   ierr = VecSet(x,pfive);CHKERRQ(ierr);
170*c4762a1bSJed Brown   return 0;
171*c4762a1bSJed Brown }
172*c4762a1bSJed Brown /* --------------------  Evaluate Jacobian F'(x) -------------------- */
173*c4762a1bSJed Brown /*  Evaluates a matrix that is used to precondition the matrix-free
174*c4762a1bSJed Brown     jacobian. In this case, the explict preconditioner matrix is
175*c4762a1bSJed Brown     also EXACTLY the Jacobian. In general, it would be some lower
176*c4762a1bSJed Brown     order, simplified apprioximation */
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown PetscErrorCode  FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
179*c4762a1bSJed Brown {
180*c4762a1bSJed Brown   const PetscScalar *xx;
181*c4762a1bSJed Brown   PetscScalar       A[3],d;
182*c4762a1bSJed Brown   PetscInt          i,n,j[3];
183*c4762a1bSJed Brown   PetscErrorCode    ierr;
184*c4762a1bSJed Brown   AppCtx            *user = (AppCtx*) dummy;
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
188*c4762a1bSJed Brown   d    = (PetscReal)(n - 1); d = d*d;
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown   i    = 0; A[0] = 1.0;
191*c4762a1bSJed Brown   ierr = MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES);CHKERRQ(ierr);
192*c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
193*c4762a1bSJed Brown     j[0] = i - 1; j[1] = i;                   j[2] = i + 1;
194*c4762a1bSJed Brown     A[0] = d;     A[1] = -2.0*d + 2.0*xx[i];  A[2] = d;
195*c4762a1bSJed Brown     ierr = MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);
196*c4762a1bSJed Brown   }
197*c4762a1bSJed Brown   i     = n-1; A[0] = 1.0;
198*c4762a1bSJed Brown   ierr  = MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr  = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr  = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr  = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown   if (user->variant) {
204*c4762a1bSJed Brown     ierr = MatMFFDSetBase(jac,x,NULL);CHKERRQ(ierr);
205*c4762a1bSJed Brown   }
206*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
207*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
208*c4762a1bSJed Brown   return 0;
209*c4762a1bSJed Brown }
210*c4762a1bSJed Brown 
211*c4762a1bSJed Brown PetscErrorCode  FormJacobianNoMatrix(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
212*c4762a1bSJed Brown {
213*c4762a1bSJed Brown   PetscErrorCode    ierr;
214*c4762a1bSJed Brown   AppCtx            *user = (AppCtx*) dummy;
215*c4762a1bSJed Brown 
216*c4762a1bSJed Brown   if (user->variant) {
217*c4762a1bSJed Brown     ierr = MatMFFDSetBase(jac,x,NULL);CHKERRQ(ierr);
218*c4762a1bSJed Brown   }
219*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
221*c4762a1bSJed Brown   return 0;
222*c4762a1bSJed Brown }
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown /* --------------------  User-defined monitor ----------------------- */
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown PetscErrorCode  Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy)
227*c4762a1bSJed Brown {
228*c4762a1bSJed Brown   PetscErrorCode ierr;
229*c4762a1bSJed Brown   MonitorCtx     *monP = (MonitorCtx*) dummy;
230*c4762a1bSJed Brown   Vec            x;
231*c4762a1bSJed Brown   MPI_Comm       comm;
232*c4762a1bSJed Brown 
233*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = PetscFPrintf(comm,stdout,"iter = %D, SNES Function norm %g \n",its,(double)fnorm);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = VecView(x,monP->viewer);CHKERRQ(ierr);
237*c4762a1bSJed Brown   return 0;
238*c4762a1bSJed Brown }
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown 
241*c4762a1bSJed Brown /*TEST
242*c4762a1bSJed Brown 
243*c4762a1bSJed Brown    test:
244*c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short
245*c4762a1bSJed Brown 
246*c4762a1bSJed Brown    test:
247*c4762a1bSJed Brown       suffix: 2
248*c4762a1bSJed Brown       args: -variant -ksp_gmres_cgs_refinement_type refine_always  -snes_monitor_short
249*c4762a1bSJed Brown       output_file: output/ex7_1.out
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown    # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning
252*c4762a1bSJed Brown    test:
253*c4762a1bSJed Brown       suffix: 3
254*c4762a1bSJed Brown       args: -variant -pc_type jacobi -snes_view -ksp_monitor
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown    # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning
257*c4762a1bSJed Brown    test:
258*c4762a1bSJed Brown       suffix: 4
259*c4762a1bSJed Brown       args: -variant -pc_type jacobi -puremf  -snes_view -ksp_monitor
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown TEST*/
262