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