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