xref: /petsc/src/ts/tutorials/ex24.c (revision 3c633725528e547aaaa9b672a746f5d686a276e1)
1c4762a1bSJed Brown static char help[] = "Pseudotransient continuation to solve a many-variable system that comes from the 2 variable Rosenbrock function + trivial.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
6c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
7c4762a1bSJed Brown static PetscErrorCode MonitorObjective(TS,PetscInt,PetscReal,Vec,void*);
8c4762a1bSJed Brown 
9c4762a1bSJed Brown typedef struct {
10c4762a1bSJed Brown   PetscInt  n;
11c4762a1bSJed Brown   PetscBool monitor_short;
12c4762a1bSJed Brown } Ctx;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown int main(int argc,char **argv)
15c4762a1bSJed Brown {
16c4762a1bSJed Brown   TS             ts;            /* time integration context */
17c4762a1bSJed Brown   Vec            X;             /* solution, residual vectors */
18c4762a1bSJed Brown   Mat            J;             /* Jacobian matrix */
19c4762a1bSJed Brown   PetscErrorCode ierr;
20c4762a1bSJed Brown   PetscScalar    *x;
21c4762a1bSJed Brown   PetscReal      ftime;
22c4762a1bSJed Brown   PetscInt       i,steps,nits,lits;
23c4762a1bSJed Brown   PetscBool      view_final;
24c4762a1bSJed Brown   Ctx            ctx;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27c4762a1bSJed Brown   ctx.n = 3;
28c4762a1bSJed Brown   ierr  = PetscOptionsGetInt(NULL,NULL,"-n",&ctx.n,NULL);CHKERRQ(ierr);
29*3c633725SBarry Smith   PetscCheck(ctx.n >= 2,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The dimension specified with -n must be at least 2");
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   view_final = PETSC_FALSE;
32c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-view_final",&view_final,NULL);CHKERRQ(ierr);
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   ctx.monitor_short = PETSC_FALSE;
35c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor_short",&ctx.monitor_short,NULL);CHKERRQ(ierr);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown      Create Jacobian matrix data structure and state vector
39c4762a1bSJed Brown   */
40c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx.n,ctx.n);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = MatCreateVecs(J,&X,NULL);CHKERRQ(ierr);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* Create time integration context */
47c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
48c4762a1bSJed Brown   ierr = TSSetType(ts,TSPSEUDO);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,FormIFunction,&ctx);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&ctx);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,1000);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,1e-3);CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = TSMonitorSet(ts,MonitorObjective,&ctx,NULL);CHKERRQ(ierr);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Customize time integrator; set runtime options
58c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
63c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64c4762a1bSJed Brown   ierr = VecSet(X,0.0);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
66c4762a1bSJed Brown #if 1
67c4762a1bSJed Brown   x[0] = 5.;
68c4762a1bSJed Brown   x[1] = -5.;
69c4762a1bSJed Brown   for (i=2; i<ctx.n; i++) x[i] = 5.;
70c4762a1bSJed Brown #else
71c4762a1bSJed Brown   x[0] = 1.0;
72c4762a1bSJed Brown   x[1] = 15.0;
73c4762a1bSJed Brown   for (i=2; i<ctx.n; i++) x[i] = 10.0;
74c4762a1bSJed Brown #endif
75c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   ierr = TSSolve(ts,X);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = TSGetSNESIterations(ts,&nits);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = TSGetKSPIterations(ts,&lits);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Time integrator took (%D,%D,%D) iterations to reach final time %g\n",steps,nits,lits,(double)ftime);CHKERRQ(ierr);
83c4762a1bSJed Brown   if (view_final) {
84c4762a1bSJed Brown     ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
89c4762a1bSJed Brown      are no longer needed.
90c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = PetscFinalize();
96c4762a1bSJed Brown   return ierr;
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown static PetscErrorCode MonitorObjective(TS ts,PetscInt step,PetscReal t,Vec X,void *ictx)
100c4762a1bSJed Brown {
101c4762a1bSJed Brown   Ctx               *ctx = (Ctx*)ictx;
102c4762a1bSJed Brown   PetscErrorCode    ierr;
103c4762a1bSJed Brown   const PetscScalar *x;
104c4762a1bSJed Brown   PetscScalar       f;
105c4762a1bSJed Brown   PetscReal         dt,gnorm;
106c4762a1bSJed Brown   PetscInt          i,snesit,linit;
107c4762a1bSJed Brown   SNES              snes;
108c4762a1bSJed Brown   Vec               Xdot,F;
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   PetscFunctionBeginUser;
111c4762a1bSJed Brown   /* Compute objective functional */
112c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
113c4762a1bSJed Brown   f    = 0;
114c4762a1bSJed Brown   for (i=0; i<ctx->n-1; i++) f += PetscSqr(1. - x[i]) + 100. * PetscSqr(x[i+1] - PetscSqr(x[i]));
115c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Compute norm of gradient */
118c4762a1bSJed Brown   ierr = VecDuplicate(X,&Xdot);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = VecDuplicate(X,&F);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = FormIFunction(ts,t,X,Xdot,F,ictx);CHKERRQ(ierr);
122c4762a1bSJed Brown   ierr = VecNorm(F,NORM_2,&gnorm);CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = VecDestroy(&Xdot);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = VecDestroy(&F);CHKERRQ(ierr);
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&snesit);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = SNESGetLinearSolveIterations(snes,&linit);CHKERRQ(ierr);
130c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,
131c4762a1bSJed Brown                      (ctx->monitor_short
132c4762a1bSJed Brown                       ? "%3D t=%10.1e  dt=%10.1e  f=%10.1e  df=%10.1e  it=(%2D,%3D)\n"
133c4762a1bSJed Brown                       : "%3D t=%10.4e  dt=%10.4e  f=%10.4e  df=%10.4e  it=(%2D,%3D)\n"),
134c4762a1bSJed Brown                      step,(double)t,(double)dt,(double)PetscRealPart(f),(double)gnorm,snesit,linit);CHKERRQ(ierr);
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /* ------------------------------------------------------------------- */
139c4762a1bSJed Brown /*
140c4762a1bSJed Brown    FormIFunction - Evaluates nonlinear function, F(X,Xdot) = Xdot + grad(objective(X))
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    Input Parameters:
143c4762a1bSJed Brown +  ts   - the TS context
144c4762a1bSJed Brown .  t - time
145c4762a1bSJed Brown .  X    - input vector
146c4762a1bSJed Brown .  Xdot - time derivative
147c4762a1bSJed Brown -  ctx  - optional user-defined context
148c4762a1bSJed Brown 
149c4762a1bSJed Brown    Output Parameter:
150c4762a1bSJed Brown .  F - function vector
151c4762a1bSJed Brown  */
152c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ictx)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   PetscErrorCode    ierr;
155c4762a1bSJed Brown   const PetscScalar *x;
156c4762a1bSJed Brown   PetscScalar       *f;
157c4762a1bSJed Brown   PetscInt          i;
158c4762a1bSJed Brown   Ctx               *ctx = (Ctx*)ictx;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBeginUser;
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown     Get pointers to vector data.
163c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
164c4762a1bSJed Brown     the data array.  Otherwise, the routine is implementation dependent.
165c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
166c4762a1bSJed Brown     the array.
167c4762a1bSJed Brown   */
168c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = VecZeroEntries(F);CHKERRQ(ierr);
170c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* Compute gradient of objective */
173c4762a1bSJed Brown   for (i=0; i<ctx->n-1; i++) {
174c4762a1bSJed Brown     PetscScalar a,a0,a1;
175c4762a1bSJed Brown     a       = x[i+1] - PetscSqr(x[i]);
176c4762a1bSJed Brown     a0      = -2.*x[i];
177c4762a1bSJed Brown     a1      = 1.;
178c4762a1bSJed Brown     f[i]   += -2.*(1. - x[i]) + 200.*a*a0;
179c4762a1bSJed Brown     f[i+1] += 200.*a*a1;
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown   /* Restore vectors */
182c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = VecAXPY(F,1.0,Xdot);CHKERRQ(ierr);
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown /* ------------------------------------------------------------------- */
188c4762a1bSJed Brown /*
189c4762a1bSJed Brown    FormIJacobian - Evaluates Jacobian matrix.
190c4762a1bSJed Brown 
191c4762a1bSJed Brown    Input Parameters:
192c4762a1bSJed Brown +  ts - the TS context
193c4762a1bSJed Brown .  t - pseudo-time
194c4762a1bSJed Brown .  X - input vector
195c4762a1bSJed Brown .  Xdot - time derivative
196c4762a1bSJed Brown .  shift - multiplier for mass matrix
197c4762a1bSJed Brown .  dummy - user-defined context
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    Output Parameters:
200c4762a1bSJed Brown .  J - Jacobian matrix
201c4762a1bSJed Brown .  B - optionally different preconditioning matrix
202c4762a1bSJed Brown .  flag - flag indicating matrix structure
203c4762a1bSJed Brown */
204c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat J,Mat B,void *ictx)
205c4762a1bSJed Brown {
206c4762a1bSJed Brown   const PetscScalar *x;
207c4762a1bSJed Brown   PetscErrorCode    ierr;
208c4762a1bSJed Brown   PetscInt          i;
209c4762a1bSJed Brown   Ctx               *ctx = (Ctx*)ictx;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   PetscFunctionBeginUser;
212c4762a1bSJed Brown   ierr = MatZeroEntries(B);CHKERRQ(ierr);
213c4762a1bSJed Brown   /*
214c4762a1bSJed Brown      Get pointer to vector data
215c4762a1bSJed Brown   */
216c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /*
219c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
220c4762a1bSJed Brown   */
221c4762a1bSJed Brown   for (i=0; i<ctx->n-1; i++) {
222c4762a1bSJed Brown     PetscInt    rowcol[2];
223c4762a1bSJed Brown     PetscScalar v[2][2],a,a0,a1,a00,a01,a10,a11;
224c4762a1bSJed Brown     rowcol[0] = i;
225c4762a1bSJed Brown     rowcol[1] = i+1;
226c4762a1bSJed Brown     a         = x[i+1] - PetscSqr(x[i]);
227c4762a1bSJed Brown     a0        = -2.*x[i];
228c4762a1bSJed Brown     a00       = -2.;
229c4762a1bSJed Brown     a01       = 0.;
230c4762a1bSJed Brown     a1        = 1.;
231c4762a1bSJed Brown     a10       = 0.;
232c4762a1bSJed Brown     a11       = 0.;
233c4762a1bSJed Brown     v[0][0]   = 2. + 200.*(a*a00 + a0*a0);
234c4762a1bSJed Brown     v[0][1]   = 200.*(a*a01 + a1*a0);
235c4762a1bSJed Brown     v[1][0]   = 200.*(a*a10 + a0*a1);
236c4762a1bSJed Brown     v[1][1]   = 200.*(a*a11 + a1*a1);
237c4762a1bSJed Brown     ierr      = MatSetValues(B,2,rowcol,2,rowcol,&v[0][0],ADD_VALUES);CHKERRQ(ierr);
238c4762a1bSJed Brown   }
239c4762a1bSJed Brown   for (i=0; i<ctx->n; i++) {
240c4762a1bSJed Brown     ierr = MatSetValue(B,i,i,(PetscScalar)shift,ADD_VALUES);CHKERRQ(ierr);
241c4762a1bSJed Brown   }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /*
246c4762a1bSJed Brown      Assemble matrix
247c4762a1bSJed Brown   */
248c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
250c4762a1bSJed Brown   if (J != B) {
251c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown   PetscFunctionReturn(0);
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown /*TEST
258c4762a1bSJed Brown 
259c4762a1bSJed Brown     test:
260c4762a1bSJed Brown       requires: !single
261c4762a1bSJed Brown 
262c4762a1bSJed Brown     test:
263c4762a1bSJed Brown       args: -pc_type lu -ts_dt 1e-5 -ts_max_time 1e5 -n 50 -monitor_short -snes_max_it 5 -snes_type newtonls -ts_max_snes_failures -1
264c4762a1bSJed Brown       requires: !single
265c4762a1bSJed Brown       suffix: 2
266c4762a1bSJed Brown 
267c4762a1bSJed Brown TEST*/
268