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