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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 26c4762a1bSJed Brown ctx.n = 3; 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&ctx.n,NULL)); 283c633725SBarry Smith PetscCheck(ctx.n >= 2,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The dimension specified with -n must be at least 2"); 29c4762a1bSJed Brown 30c4762a1bSJed Brown view_final = PETSC_FALSE; 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_final",&view_final,NULL)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown ctx.monitor_short = PETSC_FALSE; 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor_short",&ctx.monitor_short,NULL)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown Create Jacobian matrix data structure and state vector 38c4762a1bSJed Brown */ 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx.n,ctx.n)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(J)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(J)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(J,&X,NULL)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Create time integration context */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSPSEUDO)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,&ctx)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,&ctx)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,1000)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,1e-3)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,MonitorObjective,&ctx,NULL)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56c4762a1bSJed Brown Customize time integrator; set runtime options 57c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 585f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X,0.0)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(X,&x)); 65c4762a1bSJed Brown #if 1 66c4762a1bSJed Brown x[0] = 5.; 67c4762a1bSJed Brown x[1] = -5.; 68c4762a1bSJed Brown for (i=2; i<ctx.n; i++) x[i] = 5.; 69c4762a1bSJed Brown #else 70c4762a1bSJed Brown x[0] = 1.0; 71c4762a1bSJed Brown x[1] = 15.0; 72c4762a1bSJed Brown for (i=2; i<ctx.n; i++) x[i] = 10.0; 73c4762a1bSJed Brown #endif 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(X,&x)); 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,X)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSNESIterations(ts,&nits)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetKSPIterations(ts,&lits)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Time integrator took (%D,%D,%D) iterations to reach final time %g\n",steps,nits,lits,(double)ftime)); 82c4762a1bSJed Brown if (view_final) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 87c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 88c4762a1bSJed Brown are no longer needed. 89c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 90c4762a1bSJed Brown 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 94*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 95*b122ec5aSJacob Faibussowitsch return 0; 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown static PetscErrorCode MonitorObjective(TS ts,PetscInt step,PetscReal t,Vec X,void *ictx) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown Ctx *ctx = (Ctx*)ictx; 101c4762a1bSJed Brown PetscErrorCode ierr; 102c4762a1bSJed Brown const PetscScalar *x; 103c4762a1bSJed Brown PetscScalar f; 104c4762a1bSJed Brown PetscReal dt,gnorm; 105c4762a1bSJed Brown PetscInt i,snesit,linit; 106c4762a1bSJed Brown SNES snes; 107c4762a1bSJed Brown Vec Xdot,F; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBeginUser; 110c4762a1bSJed Brown /* Compute objective functional */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 112c4762a1bSJed Brown f = 0; 113c4762a1bSJed Brown for (i=0; i<ctx->n-1; i++) f += PetscSqr(1. - x[i]) + 100. * PetscSqr(x[i+1] - PetscSqr(x[i])); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Compute norm of gradient */ 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&Xdot)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&F)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Xdot)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(FormIFunction(ts,t,X,Xdot,F,ictx)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(F,NORM_2,&gnorm)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Xdot)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&F)); 124c4762a1bSJed Brown 1255f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSNES(ts,&snes)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&snesit)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetLinearSolveIterations(snes,&linit)); 129c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, 130c4762a1bSJed Brown (ctx->monitor_short 131c4762a1bSJed Brown ? "%3D t=%10.1e dt=%10.1e f=%10.1e df=%10.1e it=(%2D,%3D)\n" 132c4762a1bSJed Brown : "%3D t=%10.4e dt=%10.4e f=%10.4e df=%10.4e it=(%2D,%3D)\n"), 133c4762a1bSJed Brown step,(double)t,(double)dt,(double)PetscRealPart(f),(double)gnorm,snesit,linit);CHKERRQ(ierr); 134c4762a1bSJed Brown PetscFunctionReturn(0); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown FormIFunction - Evaluates nonlinear function, F(X,Xdot) = Xdot + grad(objective(X)) 140c4762a1bSJed Brown 141c4762a1bSJed Brown Input Parameters: 142c4762a1bSJed Brown + ts - the TS context 143c4762a1bSJed Brown . t - time 144c4762a1bSJed Brown . X - input vector 145c4762a1bSJed Brown . Xdot - time derivative 146c4762a1bSJed Brown - ctx - optional user-defined context 147c4762a1bSJed Brown 148c4762a1bSJed Brown Output Parameter: 149c4762a1bSJed Brown . F - function vector 150c4762a1bSJed Brown */ 151c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ictx) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown const PetscScalar *x; 154c4762a1bSJed Brown PetscScalar *f; 155c4762a1bSJed Brown PetscInt i; 156c4762a1bSJed Brown Ctx *ctx = (Ctx*)ictx; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBeginUser; 159c4762a1bSJed Brown /* 160c4762a1bSJed Brown Get pointers to vector data. 161c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 162c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 163c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 164c4762a1bSJed Brown the array. 165c4762a1bSJed Brown */ 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(F)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* Compute gradient of objective */ 171c4762a1bSJed Brown for (i=0; i<ctx->n-1; i++) { 172c4762a1bSJed Brown PetscScalar a,a0,a1; 173c4762a1bSJed Brown a = x[i+1] - PetscSqr(x[i]); 174c4762a1bSJed Brown a0 = -2.*x[i]; 175c4762a1bSJed Brown a1 = 1.; 176c4762a1bSJed Brown f[i] += -2.*(1. - x[i]) + 200.*a*a0; 177c4762a1bSJed Brown f[i+1] += 200.*a*a1; 178c4762a1bSJed Brown } 179c4762a1bSJed Brown /* Restore vectors */ 1805f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(F,1.0,Xdot)); 183c4762a1bSJed Brown PetscFunctionReturn(0); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown FormIJacobian - Evaluates Jacobian matrix. 188c4762a1bSJed Brown 189c4762a1bSJed Brown Input Parameters: 190c4762a1bSJed Brown + ts - the TS context 191c4762a1bSJed Brown . t - pseudo-time 192c4762a1bSJed Brown . X - input vector 193c4762a1bSJed Brown . Xdot - time derivative 194c4762a1bSJed Brown . shift - multiplier for mass matrix 195c4762a1bSJed Brown . dummy - user-defined context 196c4762a1bSJed Brown 197c4762a1bSJed Brown Output Parameters: 198c4762a1bSJed Brown . J - Jacobian matrix 199c4762a1bSJed Brown . B - optionally different preconditioning matrix 200c4762a1bSJed Brown . flag - flag indicating matrix structure 201c4762a1bSJed Brown */ 202c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat J,Mat B,void *ictx) 203c4762a1bSJed Brown { 204c4762a1bSJed Brown const PetscScalar *x; 205c4762a1bSJed Brown PetscInt i; 206c4762a1bSJed Brown Ctx *ctx = (Ctx*)ictx; 207c4762a1bSJed Brown 208c4762a1bSJed Brown PetscFunctionBeginUser; 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(B)); 210c4762a1bSJed Brown /* 211c4762a1bSJed Brown Get pointer to vector data 212c4762a1bSJed Brown */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* 216c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 217c4762a1bSJed Brown */ 218c4762a1bSJed Brown for (i=0; i<ctx->n-1; i++) { 219c4762a1bSJed Brown PetscInt rowcol[2]; 220c4762a1bSJed Brown PetscScalar v[2][2],a,a0,a1,a00,a01,a10,a11; 221c4762a1bSJed Brown rowcol[0] = i; 222c4762a1bSJed Brown rowcol[1] = i+1; 223c4762a1bSJed Brown a = x[i+1] - PetscSqr(x[i]); 224c4762a1bSJed Brown a0 = -2.*x[i]; 225c4762a1bSJed Brown a00 = -2.; 226c4762a1bSJed Brown a01 = 0.; 227c4762a1bSJed Brown a1 = 1.; 228c4762a1bSJed Brown a10 = 0.; 229c4762a1bSJed Brown a11 = 0.; 230c4762a1bSJed Brown v[0][0] = 2. + 200.*(a*a00 + a0*a0); 231c4762a1bSJed Brown v[0][1] = 200.*(a*a01 + a1*a0); 232c4762a1bSJed Brown v[1][0] = 200.*(a*a10 + a0*a1); 233c4762a1bSJed Brown v[1][1] = 200.*(a*a11 + a1*a1); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&v[0][0],ADD_VALUES)); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown for (i=0; i<ctx->n; i++) { 2375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(B,i,i,(PetscScalar)shift,ADD_VALUES)); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* 243c4762a1bSJed Brown Assemble matrix 244c4762a1bSJed Brown */ 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 247c4762a1bSJed Brown if (J != B) { 2485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown PetscFunctionReturn(0); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /*TEST 255c4762a1bSJed Brown 256c4762a1bSJed Brown test: 257c4762a1bSJed Brown requires: !single 258c4762a1bSJed Brown 259c4762a1bSJed Brown test: 260c4762a1bSJed 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 261c4762a1bSJed Brown requires: !single 262c4762a1bSJed Brown suffix: 2 263c4762a1bSJed Brown 264c4762a1bSJed Brown TEST*/ 265