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