xref: /petsc/src/snes/tests/ex5.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 /*T
6c4762a1bSJed Brown    Concepts: SNES^basic uniprocessor example
7c4762a1bSJed Brown    Processors: 1
8c4762a1bSJed Brown T*/
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
12c4762a1bSJed Brown    file automatically includes:
13c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
14c4762a1bSJed Brown      petscmat.h - matrices
15c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
16c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
17c4762a1bSJed Brown      petscksp.h   - linear solvers
18c4762a1bSJed Brown */
19c4762a1bSJed Brown 
20c4762a1bSJed Brown #include <petscsnes.h>
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown    User-defined routines
24c4762a1bSJed Brown */
25c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
26c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
27c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
28c4762a1bSJed Brown 
29c4762a1bSJed Brown int main(int argc,char **argv)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   SNES           snes;                   /* SNES context */
32c4762a1bSJed Brown   Vec            x,r,F,U;                /* vectors */
33c4762a1bSJed Brown   Mat            J;                      /* Jacobian matrix */
34c4762a1bSJed Brown   PetscInt       its,n = 5,i,maxit,maxf,lens[3] = {1,2,2};
35c4762a1bSJed Brown   PetscMPIInt    size;
36c4762a1bSJed Brown   PetscScalar    h,xp,v,none = -1.0;
37c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
38c4762a1bSJed Brown   KSP            ksp;
39c4762a1bSJed Brown   PC             pc;
40c4762a1bSJed Brown 
41*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
425f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
432c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!");
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
45c4762a1bSJed Brown   h    = 1.0/(n-1);
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48c4762a1bSJed Brown      Create nonlinear solver context
49c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCVPBJACOBI));
55c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
57c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58c4762a1bSJed Brown   /*
59c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
60c4762a1bSJed Brown   */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,n));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /*
69c4762a1bSJed Brown      Set function evaluation routine and vector
70c4762a1bSJed Brown   */
715f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormFunction,(void*)F));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
75c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76c4762a1bSJed Brown 
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetVariableBlockSizes(J,3,lens));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /*
84c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
85c4762a1bSJed Brown      routine. User can override with:
86c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
87c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
88c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
89c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
90c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
91c4762a1bSJed Brown                          products within Newton-Krylov method
92c4762a1bSJed Brown   */
93c4762a1bSJed Brown 
945f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
98c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
102c4762a1bSJed Brown   */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution"));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
108c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
109c4762a1bSJed Brown   */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*
113c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
114c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
115c4762a1bSJed Brown      the option -snes_view
116c4762a1bSJed Brown   */
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown      Initialize application:
122c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
123c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   xp = 0.0;
126c4762a1bSJed Brown   for (i=0; i<n; i++) {
127c4762a1bSJed Brown     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES));
129c4762a1bSJed Brown     v    = xp*xp*xp;
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(U,1,&i,&v,INSERT_VALUES));
131c4762a1bSJed Brown     xp  += h;
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
136c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
139c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
140c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
141c4762a1bSJed Brown      this vector to zero by calling VecSet().
142c4762a1bSJed Brown   */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %D\n\n",its));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown      Check solution and clean up
150c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /*
153c4762a1bSJed Brown      Check the error
154c4762a1bSJed Brown   */
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,none,U));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /*
160c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
161c4762a1bSJed Brown      are no longer needed.
162c4762a1bSJed Brown   */
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));  CHKERRQ(VecDestroy(&r));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));  CHKERRQ(VecDestroy(&F));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));  CHKERRQ(SNESDestroy(&snes));
166*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
167*b122ec5aSJacob Faibussowitsch   return 0;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown /* ------------------------------------------------------------------- */
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
172c4762a1bSJed Brown 
173c4762a1bSJed Brown    Input/Output Parameter:
174c4762a1bSJed Brown .  x - the solution vector
175c4762a1bSJed Brown */
176c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
177c4762a1bSJed Brown {
178c4762a1bSJed Brown   PetscScalar    pfive = .50;
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
180c4762a1bSJed Brown   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown /* ------------------------------------------------------------------- */
183c4762a1bSJed Brown /*
184c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Input Parameters:
187c4762a1bSJed Brown .  snes - the SNES context
188c4762a1bSJed Brown .  x - input vector
189c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
190c4762a1bSJed Brown 
191c4762a1bSJed Brown    Output Parameter:
192c4762a1bSJed Brown .  f - function vector
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    Note:
195c4762a1bSJed Brown    The user-defined context can contain any application-specific data
196c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
197c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
198c4762a1bSJed Brown    a vector containing the right-hand-side of the discretized PDE.
199c4762a1bSJed Brown  */
200c4762a1bSJed Brown 
201c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   Vec               g = (Vec)ctx;
204c4762a1bSJed Brown   const PetscScalar *xx,*gg;
205c4762a1bSJed Brown   PetscScalar       *ff,d;
206c4762a1bSJed Brown   PetscInt          i,n;
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /*
209c4762a1bSJed Brown      Get pointers to vector data.
210c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
211c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
212c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
213c4762a1bSJed Brown          the array.
214c4762a1bSJed Brown   */
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&ff));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(g,&gg));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /*
220c4762a1bSJed Brown      Compute function
221c4762a1bSJed Brown   */
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
223c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
224c4762a1bSJed Brown   ff[0] = xx[0];
225c4762a1bSJed 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];
226c4762a1bSJed Brown   ff[n-1] = xx[n-1] - 1.0;
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /*
229c4762a1bSJed Brown      Restore vectors
230c4762a1bSJed Brown   */
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&ff));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(g,&gg));
234c4762a1bSJed Brown   return 0;
235c4762a1bSJed Brown }
236c4762a1bSJed Brown /* ------------------------------------------------------------------- */
237c4762a1bSJed Brown /*
238c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    Input Parameters:
241c4762a1bSJed Brown .  snes - the SNES context
242c4762a1bSJed Brown .  x - input vector
243c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    Output Parameters:
246c4762a1bSJed Brown .  jac - Jacobian matrix
247c4762a1bSJed Brown .  B - optionally different preconditioning matrix
248c4762a1bSJed Brown 
249c4762a1bSJed Brown */
250c4762a1bSJed Brown 
251c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
252c4762a1bSJed Brown {
253c4762a1bSJed Brown   const PetscScalar *xx;
254c4762a1bSJed Brown   PetscScalar       A[3],d;
255c4762a1bSJed Brown   PetscInt          i,n,j[3];
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /*
258c4762a1bSJed Brown      Get pointer to vector data
259c4762a1bSJed Brown   */
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /*
263c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
264c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
265c4762a1bSJed Brown         row at once.
266c4762a1bSJed Brown   */
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
268c4762a1bSJed Brown   d    = (PetscReal)(n - 1); d = d*d;
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   /*
271c4762a1bSJed Brown      Interior grid points
272c4762a1bSJed Brown   */
273c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
274c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
275c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
2765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /*
280c4762a1bSJed Brown      Boundary points
281c4762a1bSJed Brown   */
282c4762a1bSJed Brown   i = 0;   A[0] = 1.0;
283c4762a1bSJed Brown 
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   i = n-1; A[0] = 1.0;
287c4762a1bSJed Brown 
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   /*
291c4762a1bSJed Brown      Restore vector
292c4762a1bSJed Brown   */
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /*
296c4762a1bSJed Brown      Assemble matrix
297c4762a1bSJed Brown   */
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
300c4762a1bSJed Brown   if (jac != B) {
3015f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
3025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
303c4762a1bSJed Brown   }
304c4762a1bSJed Brown   return 0;
305c4762a1bSJed Brown }
306c4762a1bSJed Brown 
307c4762a1bSJed Brown /*TEST
308c4762a1bSJed Brown 
309c4762a1bSJed Brown    test:
310c4762a1bSJed Brown       args: -snes_monitor_short -snes_view -ksp_monitor
311c4762a1bSJed Brown 
312c4762a1bSJed Brown    # this is just a test for SNESKSPTRASPOSEONLY and KSPSolveTranspose to behave properly
313c4762a1bSJed Brown    # the solution is wrong on purpose
314c4762a1bSJed Brown    test:
315c4762a1bSJed Brown       requires: !single !complex
316c4762a1bSJed Brown       suffix: transpose_only
317c4762a1bSJed 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
318c4762a1bSJed Brown 
319c4762a1bSJed Brown TEST*/
320