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