xref: /petsc/src/ts/tutorials/ex24.c (revision 2a1887a77e7b2c6e00dd0ba96d1387c839460237)
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 
main(int argc,char ** argv)14d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown   TS           ts; /* time integration context */
17c4762a1bSJed Brown   Vec          X;  /* solution, residual vectors */
18c4762a1bSJed Brown   Mat          J;  /* Jacobian matrix */
19c4762a1bSJed Brown   PetscScalar *x;
20c4762a1bSJed Brown   PetscReal    ftime;
21c4762a1bSJed Brown   PetscInt     i, steps, nits, lits;
22c4762a1bSJed Brown   PetscBool    view_final;
23c4762a1bSJed Brown   Ctx          ctx;
24c4762a1bSJed Brown 
25327415f7SBarry Smith   PetscFunctionBeginUser;
26c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
27c4762a1bSJed Brown   ctx.n = 3;
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &ctx.n, NULL));
293c633725SBarry 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;
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_final", &view_final, NULL));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   ctx.monitor_short = PETSC_FALSE;
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor_short", &ctx.monitor_short, NULL));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown      Create Jacobian matrix data structure and state vector
39c4762a1bSJed Brown   */
409566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, ctx.n, ctx.n));
429566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
439566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
449566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(J, &X, NULL));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* Create time integration context */
479566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
489566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSPSEUDO));
499566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &ctx));
509566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &ctx));
519566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1000));
529566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
539566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1e-3));
549566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MonitorObjective, &ctx, NULL));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Customize time integrator; set runtime options
58c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
599566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
63c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
649566063dSJacob Faibussowitsch   PetscCall(VecSet(X, 0.0));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
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
759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
789566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
799566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
809566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &nits));
819566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &lits));
8263a3b9bcSJacob Faibussowitsch   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));
831baa6e33SBarry Smith   if (view_final) PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
87c4762a1bSJed Brown      are no longer needed.
88c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
929566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
939566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
94b122ec5aSJacob Faibussowitsch   return 0;
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
MonitorObjective(TS ts,PetscInt step,PetscReal t,Vec X,void * ictx)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode MonitorObjective(TS ts, PetscInt step, PetscReal t, Vec X, void *ictx)
98d71ae5a4SJacob Faibussowitsch {
99c4762a1bSJed Brown   Ctx               *ctx = (Ctx *)ictx;
100c4762a1bSJed Brown   const PetscScalar *x;
101c4762a1bSJed Brown   PetscScalar        f;
102c4762a1bSJed Brown   PetscReal          dt, gnorm;
103c4762a1bSJed Brown   PetscInt           i, snesit, linit;
104c4762a1bSJed Brown   SNES               snes;
105c4762a1bSJed Brown   Vec                Xdot, F;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBeginUser;
108c4762a1bSJed Brown   /* Compute objective functional */
1099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
110c4762a1bSJed Brown   f = 0;
111c4762a1bSJed Brown   for (i = 0; i < ctx->n - 1; i++) f += PetscSqr(1. - x[i]) + 100. * PetscSqr(x[i + 1] - PetscSqr(x[i]));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Compute norm of gradient */
1159566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Xdot));
1169566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &F));
1179566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Xdot));
1189566063dSJacob Faibussowitsch   PetscCall(FormIFunction(ts, t, X, Xdot, F, ictx));
1199566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &gnorm));
1209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xdot));
1219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F));
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1249566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1259566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &snesit));
1269566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes, &linit));
12757508eceSPierre Jolivet   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" : "%3" PetscInt_FMT " t=%10.4e  dt=%10.4e  f=%10.4e  df=%10.4e  it=(%2" PetscInt_FMT ",%3" PetscInt_FMT ")\n", step, (double)t, (double)dt, (double)PetscRealPart(f), (double)gnorm, snesit, linit));
1283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /* ------------------------------------------------------------------- */
132c4762a1bSJed Brown /*
133c4762a1bSJed Brown    FormIFunction - Evaluates nonlinear function, F(X,Xdot) = Xdot + grad(objective(X))
134c4762a1bSJed Brown 
135c4762a1bSJed Brown    Input Parameters:
136c4762a1bSJed Brown +  ts   - the TS context
137c4762a1bSJed Brown .  t - time
138c4762a1bSJed Brown .  X    - input vector
139c4762a1bSJed Brown .  Xdot - time derivative
140c4762a1bSJed Brown -  ctx  - optional user-defined context
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    Output Parameter:
143c4762a1bSJed Brown .  F - function vector
144c4762a1bSJed Brown  */
FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void * ictx)145d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ictx)
146d71ae5a4SJacob Faibussowitsch {
147c4762a1bSJed Brown   const PetscScalar *x;
148c4762a1bSJed Brown   PetscScalar       *f;
149c4762a1bSJed Brown   PetscInt           i;
150c4762a1bSJed Brown   Ctx               *ctx = (Ctx *)ictx;
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   PetscFunctionBeginUser;
153c4762a1bSJed Brown   /*
154c4762a1bSJed Brown     Get pointers to vector data.
155c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
156c4762a1bSJed Brown     the data array.  Otherwise, the routine is implementation dependent.
157c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
158c4762a1bSJed Brown     the array.
159c4762a1bSJed Brown   */
1609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1619566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
1629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* Compute gradient of objective */
165c4762a1bSJed Brown   for (i = 0; i < ctx->n - 1; i++) {
166c4762a1bSJed Brown     PetscScalar a, a0, a1;
167c4762a1bSJed Brown     a  = x[i + 1] - PetscSqr(x[i]);
168c4762a1bSJed Brown     a0 = -2. * x[i];
169c4762a1bSJed Brown     a1 = 1.;
170c4762a1bSJed Brown     f[i] += -2. * (1. - x[i]) + 200. * a * a0;
171c4762a1bSJed Brown     f[i + 1] += 200. * a * a1;
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown   /* Restore vectors */
1749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1769566063dSJacob Faibussowitsch   PetscCall(VecAXPY(F, 1.0, Xdot));
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown /* ------------------------------------------------------------------- */
180c4762a1bSJed Brown /*
181c4762a1bSJed Brown    FormIJacobian - Evaluates Jacobian matrix.
182c4762a1bSJed Brown 
183c4762a1bSJed Brown    Input Parameters:
184c4762a1bSJed Brown +  ts - the TS context
185c4762a1bSJed Brown .  t - pseudo-time
186c4762a1bSJed Brown .  X - input vector
187c4762a1bSJed Brown .  Xdot - time derivative
188c4762a1bSJed Brown .  shift - multiplier for mass matrix
189c4762a1bSJed Brown .  dummy - user-defined context
190c4762a1bSJed Brown 
191c4762a1bSJed Brown    Output Parameters:
192c4762a1bSJed Brown .  J - Jacobian matrix
1937addb90fSBarry Smith .  B - optionally different matrix used to construct the preconditioner
194c4762a1bSJed Brown */
FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat J,Mat B,void * ictx)195d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal shift, Mat J, Mat B, void *ictx)
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown   const PetscScalar *x;
198c4762a1bSJed Brown   PetscInt           i;
199c4762a1bSJed Brown   Ctx               *ctx = (Ctx *)ictx;
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   PetscFunctionBeginUser;
2029566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
203c4762a1bSJed Brown   /*
204c4762a1bSJed Brown      Get pointer to vector data
205c4762a1bSJed Brown   */
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /*
209c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
210c4762a1bSJed Brown   */
211c4762a1bSJed Brown   for (i = 0; i < ctx->n - 1; i++) {
212c4762a1bSJed Brown     PetscInt    rowcol[2];
213c4762a1bSJed Brown     PetscScalar v[2][2], a, a0, a1, a00, a01, a10, a11;
214c4762a1bSJed Brown     rowcol[0] = i;
215c4762a1bSJed Brown     rowcol[1] = i + 1;
216c4762a1bSJed Brown     a         = x[i + 1] - PetscSqr(x[i]);
217c4762a1bSJed Brown     a0        = -2. * x[i];
218c4762a1bSJed Brown     a00       = -2.;
219c4762a1bSJed Brown     a01       = 0.;
220c4762a1bSJed Brown     a1        = 1.;
221c4762a1bSJed Brown     a10       = 0.;
222c4762a1bSJed Brown     a11       = 0.;
223c4762a1bSJed Brown     v[0][0]   = 2. + 200. * (a * a00 + a0 * a0);
224c4762a1bSJed Brown     v[0][1]   = 200. * (a * a01 + a1 * a0);
225c4762a1bSJed Brown     v[1][0]   = 200. * (a * a10 + a0 * a1);
226c4762a1bSJed Brown     v[1][1]   = 200. * (a * a11 + a1 * a1);
2279566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &v[0][0], ADD_VALUES));
228c4762a1bSJed Brown   }
22948a46eb9SPierre Jolivet   for (i = 0; i < ctx->n; i++) PetscCall(MatSetValue(B, i, i, (PetscScalar)shift, ADD_VALUES));
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /*
234c4762a1bSJed Brown      Assemble matrix
235c4762a1bSJed Brown   */
2369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
238c4762a1bSJed Brown   if (J != B) {
2399566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2409566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
241c4762a1bSJed Brown   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243c4762a1bSJed Brown }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown /*TEST
246c4762a1bSJed Brown 
247c4762a1bSJed Brown     test:
248c4762a1bSJed Brown       requires: !single
249c4762a1bSJed Brown 
250c4762a1bSJed Brown     test:
251*188af4bfSBarry Smith       args: -pc_type lu -ts_time_step 1e-5 -ts_max_time 1e5 -n 50 -monitor_short -snes_max_it 5 -snes_type newtonls -ts_max_snes_failures unlimited
252c4762a1bSJed Brown       requires: !single
253c4762a1bSJed Brown       suffix: 2
254c4762a1bSJed Brown 
255c4762a1bSJed Brown TEST*/
256