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 (%D,%D,%D) 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 PetscErrorCode ierr; 102 const PetscScalar *x; 103 PetscScalar f; 104 PetscReal dt,gnorm; 105 PetscInt i,snesit,linit; 106 SNES snes; 107 Vec Xdot,F; 108 109 PetscFunctionBeginUser; 110 /* Compute objective functional */ 111 PetscCall(VecGetArrayRead(X,&x)); 112 f = 0; 113 for (i=0; i<ctx->n-1; i++) f += PetscSqr(1. - x[i]) + 100. * PetscSqr(x[i+1] - PetscSqr(x[i])); 114 PetscCall(VecRestoreArrayRead(X,&x)); 115 116 /* Compute norm of gradient */ 117 PetscCall(VecDuplicate(X,&Xdot)); 118 PetscCall(VecDuplicate(X,&F)); 119 PetscCall(VecZeroEntries(Xdot)); 120 PetscCall(FormIFunction(ts,t,X,Xdot,F,ictx)); 121 PetscCall(VecNorm(F,NORM_2,&gnorm)); 122 PetscCall(VecDestroy(&Xdot)); 123 PetscCall(VecDestroy(&F)); 124 125 PetscCall(TSGetTimeStep(ts,&dt)); 126 PetscCall(TSGetSNES(ts,&snes)); 127 PetscCall(SNESGetIterationNumber(snes,&snesit)); 128 PetscCall(SNESGetLinearSolveIterations(snes,&linit)); 129 ierr = PetscPrintf(PETSC_COMM_WORLD, 130 (ctx->monitor_short 131 ? "%3D t=%10.1e dt=%10.1e f=%10.1e df=%10.1e it=(%2D,%3D)\n" 132 : "%3D t=%10.4e dt=%10.4e f=%10.4e df=%10.4e it=(%2D,%3D)\n"), 133 step,(double)t,(double)dt,(double)PetscRealPart(f),(double)gnorm,snesit,linit);PetscCall(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 /* ------------------------------------------------------------------- */ 138 /* 139 FormIFunction - Evaluates nonlinear function, F(X,Xdot) = Xdot + grad(objective(X)) 140 141 Input Parameters: 142 + ts - the TS context 143 . t - time 144 . X - input vector 145 . Xdot - time derivative 146 - ctx - optional user-defined context 147 148 Output Parameter: 149 . F - function vector 150 */ 151 static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ictx) 152 { 153 const PetscScalar *x; 154 PetscScalar *f; 155 PetscInt i; 156 Ctx *ctx = (Ctx*)ictx; 157 158 PetscFunctionBeginUser; 159 /* 160 Get pointers to vector data. 161 - For default PETSc vectors, VecGetArray() returns a pointer to 162 the data array. Otherwise, the routine is implementation dependent. 163 - You MUST call VecRestoreArray() when you no longer need access to 164 the array. 165 */ 166 PetscCall(VecGetArrayRead(X,&x)); 167 PetscCall(VecZeroEntries(F)); 168 PetscCall(VecGetArray(F,&f)); 169 170 /* Compute gradient of objective */ 171 for (i=0; i<ctx->n-1; i++) { 172 PetscScalar a,a0,a1; 173 a = x[i+1] - PetscSqr(x[i]); 174 a0 = -2.*x[i]; 175 a1 = 1.; 176 f[i] += -2.*(1. - x[i]) + 200.*a*a0; 177 f[i+1] += 200.*a*a1; 178 } 179 /* Restore vectors */ 180 PetscCall(VecRestoreArrayRead(X,&x)); 181 PetscCall(VecRestoreArray(F,&f)); 182 PetscCall(VecAXPY(F,1.0,Xdot)); 183 PetscFunctionReturn(0); 184 } 185 /* ------------------------------------------------------------------- */ 186 /* 187 FormIJacobian - Evaluates Jacobian matrix. 188 189 Input Parameters: 190 + ts - the TS context 191 . t - pseudo-time 192 . X - input vector 193 . Xdot - time derivative 194 . shift - multiplier for mass matrix 195 . dummy - user-defined context 196 197 Output Parameters: 198 . J - Jacobian matrix 199 . B - optionally different preconditioning matrix 200 . flag - flag indicating matrix structure 201 */ 202 static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat J,Mat B,void *ictx) 203 { 204 const PetscScalar *x; 205 PetscInt i; 206 Ctx *ctx = (Ctx*)ictx; 207 208 PetscFunctionBeginUser; 209 PetscCall(MatZeroEntries(B)); 210 /* 211 Get pointer to vector data 212 */ 213 PetscCall(VecGetArrayRead(X,&x)); 214 215 /* 216 Compute Jacobian entries and insert into matrix. 217 */ 218 for (i=0; i<ctx->n-1; i++) { 219 PetscInt rowcol[2]; 220 PetscScalar v[2][2],a,a0,a1,a00,a01,a10,a11; 221 rowcol[0] = i; 222 rowcol[1] = i+1; 223 a = x[i+1] - PetscSqr(x[i]); 224 a0 = -2.*x[i]; 225 a00 = -2.; 226 a01 = 0.; 227 a1 = 1.; 228 a10 = 0.; 229 a11 = 0.; 230 v[0][0] = 2. + 200.*(a*a00 + a0*a0); 231 v[0][1] = 200.*(a*a01 + a1*a0); 232 v[1][0] = 200.*(a*a10 + a0*a1); 233 v[1][1] = 200.*(a*a11 + a1*a1); 234 PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&v[0][0],ADD_VALUES)); 235 } 236 for (i=0; i<ctx->n; i++) { 237 PetscCall(MatSetValue(B,i,i,(PetscScalar)shift,ADD_VALUES)); 238 } 239 240 PetscCall(VecRestoreArrayRead(X,&x)); 241 242 /* 243 Assemble matrix 244 */ 245 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 246 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 247 if (J != B) { 248 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 249 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 250 } 251 PetscFunctionReturn(0); 252 } 253 254 /*TEST 255 256 test: 257 requires: !single 258 259 test: 260 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 261 requires: !single 262 suffix: 2 263 264 TEST*/ 265