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 14c4762a1bSJed Brown int main(int argc,char **argv) 15c4762a1bSJed Brown { 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 25*327415f7SBarry Smith PetscFunctionBeginUser; 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,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 97c4762a1bSJed Brown static PetscErrorCode MonitorObjective(TS ts,PetscInt step,PetscReal t,Vec X,void *ictx) 98c4762a1bSJed Brown { 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)); 12763a3b9bcSJacob Faibussowitsch 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" 12863a3b9bcSJacob Faibussowitsch : "%3" PetscInt_FMT " t=%10.4e dt=%10.4e f=%10.4e df=%10.4e it=(%2" PetscInt_FMT ",%3" PetscInt_FMT ")\n"), 129d0609cedSBarry Smith step,(double)t,(double)dt,(double)PetscRealPart(f),(double)gnorm,snesit,linit)); 130c4762a1bSJed Brown PetscFunctionReturn(0); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 134c4762a1bSJed Brown /* 135c4762a1bSJed Brown FormIFunction - Evaluates nonlinear function, F(X,Xdot) = Xdot + grad(objective(X)) 136c4762a1bSJed Brown 137c4762a1bSJed Brown Input Parameters: 138c4762a1bSJed Brown + ts - the TS context 139c4762a1bSJed Brown . t - time 140c4762a1bSJed Brown . X - input vector 141c4762a1bSJed Brown . Xdot - time derivative 142c4762a1bSJed Brown - ctx - optional user-defined context 143c4762a1bSJed Brown 144c4762a1bSJed Brown Output Parameter: 145c4762a1bSJed Brown . F - function vector 146c4762a1bSJed Brown */ 147c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ictx) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown const PetscScalar *x; 150c4762a1bSJed Brown PetscScalar *f; 151c4762a1bSJed Brown PetscInt i; 152c4762a1bSJed Brown Ctx *ctx = (Ctx*)ictx; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBeginUser; 155c4762a1bSJed Brown /* 156c4762a1bSJed Brown Get pointers to vector data. 157c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 158c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 159c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 160c4762a1bSJed Brown the array. 161c4762a1bSJed Brown */ 1629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 1639566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1649566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Compute gradient of objective */ 167c4762a1bSJed Brown for (i=0; i<ctx->n-1; i++) { 168c4762a1bSJed Brown PetscScalar a,a0,a1; 169c4762a1bSJed Brown a = x[i+1] - PetscSqr(x[i]); 170c4762a1bSJed Brown a0 = -2.*x[i]; 171c4762a1bSJed Brown a1 = 1.; 172c4762a1bSJed Brown f[i] += -2.*(1. - x[i]) + 200.*a*a0; 173c4762a1bSJed Brown f[i+1] += 200.*a*a1; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown /* Restore vectors */ 1769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 1779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 1789566063dSJacob Faibussowitsch PetscCall(VecAXPY(F,1.0,Xdot)); 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 182c4762a1bSJed Brown /* 183c4762a1bSJed Brown FormIJacobian - Evaluates Jacobian matrix. 184c4762a1bSJed Brown 185c4762a1bSJed Brown Input Parameters: 186c4762a1bSJed Brown + ts - the TS context 187c4762a1bSJed Brown . t - pseudo-time 188c4762a1bSJed Brown . X - input vector 189c4762a1bSJed Brown . Xdot - time derivative 190c4762a1bSJed Brown . shift - multiplier for mass matrix 191c4762a1bSJed Brown . dummy - user-defined context 192c4762a1bSJed Brown 193c4762a1bSJed Brown Output Parameters: 194c4762a1bSJed Brown . J - Jacobian matrix 195c4762a1bSJed Brown . B - optionally different preconditioning matrix 196c4762a1bSJed Brown . flag - flag indicating matrix structure 197c4762a1bSJed Brown */ 198c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat J,Mat B,void *ictx) 199c4762a1bSJed Brown { 200c4762a1bSJed Brown const PetscScalar *x; 201c4762a1bSJed Brown PetscInt i; 202c4762a1bSJed Brown Ctx *ctx = (Ctx*)ictx; 203c4762a1bSJed Brown 204c4762a1bSJed Brown PetscFunctionBeginUser; 2059566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B)); 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Get pointer to vector data 208c4762a1bSJed Brown */ 2099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown for (i=0; i<ctx->n-1; i++) { 215c4762a1bSJed Brown PetscInt rowcol[2]; 216c4762a1bSJed Brown PetscScalar v[2][2],a,a0,a1,a00,a01,a10,a11; 217c4762a1bSJed Brown rowcol[0] = i; 218c4762a1bSJed Brown rowcol[1] = i+1; 219c4762a1bSJed Brown a = x[i+1] - PetscSqr(x[i]); 220c4762a1bSJed Brown a0 = -2.*x[i]; 221c4762a1bSJed Brown a00 = -2.; 222c4762a1bSJed Brown a01 = 0.; 223c4762a1bSJed Brown a1 = 1.; 224c4762a1bSJed Brown a10 = 0.; 225c4762a1bSJed Brown a11 = 0.; 226c4762a1bSJed Brown v[0][0] = 2. + 200.*(a*a00 + a0*a0); 227c4762a1bSJed Brown v[0][1] = 200.*(a*a01 + a1*a0); 228c4762a1bSJed Brown v[1][0] = 200.*(a*a10 + a0*a1); 229c4762a1bSJed Brown v[1][1] = 200.*(a*a11 + a1*a1); 2309566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&v[0][0],ADD_VALUES)); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown for (i=0; i<ctx->n; i++) { 2339566063dSJacob Faibussowitsch PetscCall(MatSetValue(B,i,i,(PetscScalar)shift,ADD_VALUES)); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 2369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Assemble matrix 240c4762a1bSJed Brown */ 2419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 243c4762a1bSJed Brown if (J != B) { 2449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown PetscFunctionReturn(0); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown /*TEST 251c4762a1bSJed Brown 252c4762a1bSJed Brown test: 253c4762a1bSJed Brown requires: !single 254c4762a1bSJed Brown 255c4762a1bSJed Brown test: 256c4762a1bSJed 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 257c4762a1bSJed Brown requires: !single 258c4762a1bSJed Brown suffix: 2 259c4762a1bSJed Brown 260c4762a1bSJed Brown TEST*/ 261