xref: /petsc/src/snes/tests/ex5.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3c4762a1bSJed Brown This example tests PCVPBJacobiSetBlocks().\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
7c4762a1bSJed Brown    file automatically includes:
8c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
9c4762a1bSJed Brown      petscmat.h - matrices
10c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
11c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
12c4762a1bSJed Brown      petscksp.h   - linear solvers
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petscsnes.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /*
18c4762a1bSJed Brown    User-defined routines
19c4762a1bSJed Brown */
20c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
21c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
22c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
23c4762a1bSJed Brown 
24c4762a1bSJed Brown int main(int argc,char **argv)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   SNES           snes;                   /* SNES context */
27c4762a1bSJed Brown   Vec            x,r,F,U;                /* vectors */
28c4762a1bSJed Brown   Mat            J;                      /* Jacobian matrix */
29c4762a1bSJed Brown   PetscInt       its,n = 5,i,maxit,maxf,lens[3] = {1,2,2};
30c4762a1bSJed Brown   PetscMPIInt    size;
31c4762a1bSJed Brown   PetscScalar    h,xp,v,none = -1.0;
32c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
33c4762a1bSJed Brown   KSP            ksp;
34c4762a1bSJed Brown   PC             pc;
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
38be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
40c4762a1bSJed Brown   h    = 1.0/(n-1);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43c4762a1bSJed Brown      Create nonlinear solver context
44c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
479566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes,&ksp));
489566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp,&pc));
499566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCVPBJACOBI));
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
52c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53c4762a1bSJed Brown   /*
54c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
55c4762a1bSJed Brown   */
569566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
579566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,PETSC_DECIDE,n));
589566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&F));
619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&U));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /*
64c4762a1bSJed Brown      Set function evaluation routine and vector
65c4762a1bSJed Brown   */
669566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)F));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
70c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
739566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
749566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
759566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
769566063dSJacob Faibussowitsch   PetscCall(MatSetVariableBlockSizes(J,3,lens));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /*
79c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
80c4762a1bSJed Brown      routine. User can override with:
81c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
82c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
83c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
84c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
85c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
86c4762a1bSJed Brown                          products within Newton-Krylov method
87c4762a1bSJed Brown   */
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
93c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /*
96c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
97c4762a1bSJed Brown   */
989566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution"));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
103c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
104c4762a1bSJed Brown   */
1059566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /*
108c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
109c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
110c4762a1bSJed Brown      the option -snes_view
111c4762a1bSJed Brown   */
1129566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
113*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Initialize application:
117c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
118c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   xp = 0.0;
121c4762a1bSJed Brown   for (i=0; i<n; i++) {
122c4762a1bSJed Brown     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
1239566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F,1,&i,&v,INSERT_VALUES));
124c4762a1bSJed Brown     v    = xp*xp*xp;
1259566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U,1,&i,&v,INSERT_VALUES));
126c4762a1bSJed Brown     xp  += h;
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
131c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4762a1bSJed Brown   /*
133c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
134c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
135c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
136c4762a1bSJed Brown      this vector to zero by calling VecSet().
137c4762a1bSJed Brown   */
1389566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1399566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,x));
1409566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
141*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %" PetscInt_FMT "\n\n",its));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown      Check solution and clean up
145c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown      Check the error
149c4762a1bSJed Brown   */
1509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x,none,U));
1519566063dSJacob Faibussowitsch   PetscCall(VecNorm(x,NORM_2,&norm));
152*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %" PetscInt_FMT "\n",(double)norm,its));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /*
155c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
156c4762a1bSJed Brown      are no longer needed.
157c4762a1bSJed Brown   */
1589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));  PetscCall(VecDestroy(&r));
1599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));  PetscCall(VecDestroy(&F));
1609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));  PetscCall(SNESDestroy(&snes));
1619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
162b122ec5aSJacob Faibussowitsch   return 0;
163c4762a1bSJed Brown }
164c4762a1bSJed Brown /* ------------------------------------------------------------------- */
165c4762a1bSJed Brown /*
166c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
167c4762a1bSJed Brown 
168c4762a1bSJed Brown    Input/Output Parameter:
169c4762a1bSJed Brown .  x - the solution vector
170c4762a1bSJed Brown */
171c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
172c4762a1bSJed Brown {
173c4762a1bSJed Brown   PetscScalar    pfive = .50;
1749566063dSJacob Faibussowitsch   PetscCall(VecSet(x,pfive));
175c4762a1bSJed Brown   return 0;
176c4762a1bSJed Brown }
177c4762a1bSJed Brown /* ------------------------------------------------------------------- */
178c4762a1bSJed Brown /*
179c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
180c4762a1bSJed Brown 
181c4762a1bSJed Brown    Input Parameters:
182c4762a1bSJed Brown .  snes - the SNES context
183c4762a1bSJed Brown .  x - input vector
184c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Output Parameter:
187c4762a1bSJed Brown .  f - function vector
188c4762a1bSJed Brown 
189c4762a1bSJed Brown    Note:
190c4762a1bSJed Brown    The user-defined context can contain any application-specific data
191c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
192c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
193c4762a1bSJed Brown    a vector containing the right-hand-side of the discretized PDE.
194c4762a1bSJed Brown  */
195c4762a1bSJed Brown 
196c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
197c4762a1bSJed Brown {
198c4762a1bSJed Brown   Vec               g = (Vec)ctx;
199c4762a1bSJed Brown   const PetscScalar *xx,*gg;
200c4762a1bSJed Brown   PetscScalar       *ff,d;
201c4762a1bSJed Brown   PetscInt          i,n;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /*
204c4762a1bSJed Brown      Get pointers to vector data.
205c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
206c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
207c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
208c4762a1bSJed Brown          the array.
209c4762a1bSJed Brown   */
2109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
2119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f,&ff));
2129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g,&gg));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /*
215c4762a1bSJed Brown      Compute function
216c4762a1bSJed Brown   */
2179566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&n));
218c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
219c4762a1bSJed Brown   ff[0] = xx[0];
220c4762a1bSJed 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] - gg[i];
221c4762a1bSJed Brown   ff[n-1] = xx[n-1] - 1.0;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /*
224c4762a1bSJed Brown      Restore vectors
225c4762a1bSJed Brown   */
2269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
2279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f,&ff));
2289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g,&gg));
229c4762a1bSJed Brown   return 0;
230c4762a1bSJed Brown }
231c4762a1bSJed Brown /* ------------------------------------------------------------------- */
232c4762a1bSJed Brown /*
233c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    Input Parameters:
236c4762a1bSJed Brown .  snes - the SNES context
237c4762a1bSJed Brown .  x - input vector
238c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    Output Parameters:
241c4762a1bSJed Brown .  jac - Jacobian matrix
242c4762a1bSJed Brown .  B - optionally different preconditioning matrix
243c4762a1bSJed Brown 
244c4762a1bSJed Brown */
245c4762a1bSJed Brown 
246c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
247c4762a1bSJed Brown {
248c4762a1bSJed Brown   const PetscScalar *xx;
249c4762a1bSJed Brown   PetscScalar       A[3],d;
250c4762a1bSJed Brown   PetscInt          i,n,j[3];
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /*
253c4762a1bSJed Brown      Get pointer to vector data
254c4762a1bSJed Brown   */
2559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /*
258c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
259c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
260c4762a1bSJed Brown         row at once.
261c4762a1bSJed Brown   */
2629566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&n));
263c4762a1bSJed Brown   d    = (PetscReal)(n - 1); d = d*d;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /*
266c4762a1bSJed Brown      Interior grid points
267c4762a1bSJed Brown   */
268c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
269c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
270c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
2719566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /*
275c4762a1bSJed Brown      Boundary points
276c4762a1bSJed Brown   */
277c4762a1bSJed Brown   i = 0;   A[0] = 1.0;
278c4762a1bSJed Brown 
2799566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   i = n-1; A[0] = 1.0;
282c4762a1bSJed Brown 
2839566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /*
286c4762a1bSJed Brown      Restore vector
287c4762a1bSJed Brown   */
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   /*
291c4762a1bSJed Brown      Assemble matrix
292c4762a1bSJed Brown   */
2939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
295c4762a1bSJed Brown   if (jac != B) {
2969566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
2979566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
298c4762a1bSJed Brown   }
299c4762a1bSJed Brown   return 0;
300c4762a1bSJed Brown }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown /*TEST
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown       args: -snes_monitor_short -snes_view -ksp_monitor
306c4762a1bSJed Brown 
307c4762a1bSJed Brown    # this is just a test for SNESKSPTRASPOSEONLY and KSPSolveTranspose to behave properly
308c4762a1bSJed Brown    # the solution is wrong on purpose
309c4762a1bSJed Brown    test:
310c4762a1bSJed Brown       requires: !single !complex
311c4762a1bSJed Brown       suffix: transpose_only
312c4762a1bSJed Brown       args: -snes_monitor_short -snes_view -ksp_monitor -snes_type ksptransposeonly -pc_type ilu -snes_test_jacobian -snes_test_jacobian_view -ksp_view_rhs -ksp_view_solution -ksp_view_mat_explicit -ksp_view_preconditioned_operator_explicit
313c4762a1bSJed Brown 
314c4762a1bSJed Brown TEST*/
315