1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown /*F 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown \begin{eqnarray} 7*c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) \\ 8*c4762a1bSJed Brown \frac{d \theta}{dt} = \omega - \omega_s 9*c4762a1bSJed Brown \end{eqnarray} 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown F*/ 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown /* 14*c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 15*c4762a1bSJed Brown file automatically includes: 16*c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 17*c4762a1bSJed Brown petscmat.h - matrices 18*c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 19*c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 20*c4762a1bSJed Brown petscksp.h - linear solvers 21*c4762a1bSJed Brown */ 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown #include <petscts.h> 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown typedef struct { 26*c4762a1bSJed Brown PetscScalar H,omega_s,E,V,X; 27*c4762a1bSJed Brown PetscRandom rand; 28*c4762a1bSJed Brown } AppCtx; 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown /* 31*c4762a1bSJed Brown Defines the ODE passed to the ODE solver 32*c4762a1bSJed Brown */ 33*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx) 34*c4762a1bSJed Brown { 35*c4762a1bSJed Brown PetscErrorCode ierr; 36*c4762a1bSJed Brown PetscScalar *f,r; 37*c4762a1bSJed Brown const PetscScalar *u,*udot; 38*c4762a1bSJed Brown static PetscScalar R = .4; 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown PetscFunctionBegin; 41*c4762a1bSJed Brown ierr = PetscRandomGetValue(ctx->rand,&r);CHKERRQ(ierr); 42*c4762a1bSJed Brown if (r > .9) R = .5; 43*c4762a1bSJed Brown if (r < .1) R = .4; 44*c4762a1bSJed Brown R = .4; 45*c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 46*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 49*c4762a1bSJed Brown f[0] = 2.0*ctx->H*udot[0]/ctx->omega_s + ctx->E*ctx->V*PetscSinScalar(u[1])/ctx->X - R; 50*c4762a1bSJed Brown f[1] = udot[1] - u[0] + ctx->omega_s; 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 55*c4762a1bSJed Brown PetscFunctionReturn(0); 56*c4762a1bSJed Brown } 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* 59*c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 60*c4762a1bSJed Brown */ 61*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx) 62*c4762a1bSJed Brown { 63*c4762a1bSJed Brown PetscErrorCode ierr; 64*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 65*c4762a1bSJed Brown PetscScalar J[2][2]; 66*c4762a1bSJed Brown const PetscScalar *u,*udot; 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown PetscFunctionBegin; 69*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 71*c4762a1bSJed Brown J[0][0] = 2.0*ctx->H*a/ctx->omega_s; J[0][1] = -ctx->E*ctx->V*PetscCosScalar(u[1])/ctx->X; 72*c4762a1bSJed Brown J[1][0] = -1.0; J[1][1] = a; 73*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 79*c4762a1bSJed Brown if (A != B) { 80*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82*c4762a1bSJed Brown } 83*c4762a1bSJed Brown PetscFunctionReturn(0); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown int main(int argc,char **argv) 87*c4762a1bSJed Brown { 88*c4762a1bSJed Brown TS ts; /* ODE integrator */ 89*c4762a1bSJed Brown Vec U; /* solution will be stored here */ 90*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 91*c4762a1bSJed Brown PetscErrorCode ierr; 92*c4762a1bSJed Brown PetscMPIInt size; 93*c4762a1bSJed Brown PetscInt n = 2; 94*c4762a1bSJed Brown AppCtx ctx; 95*c4762a1bSJed Brown PetscScalar *u; 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 98*c4762a1bSJed Brown Initialize program 99*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 100*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 101*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 102*c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105*c4762a1bSJed Brown Create necessary matrix and vectors 106*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 107*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 115*c4762a1bSJed Brown Set runtime options 116*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 117*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Reaction options","");CHKERRQ(ierr); 118*c4762a1bSJed Brown { 119*c4762a1bSJed Brown ctx.omega_s = 1.0; 120*c4762a1bSJed Brown ierr = PetscOptionsScalar("-omega_s","","",ctx.omega_s,&ctx.omega_s,NULL);CHKERRQ(ierr); 121*c4762a1bSJed Brown ctx.H = 1.0; 122*c4762a1bSJed Brown ierr = PetscOptionsScalar("-H","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 123*c4762a1bSJed Brown ctx.E = 1.0; 124*c4762a1bSJed Brown ierr = PetscOptionsScalar("-E","","",ctx.E,&ctx.E,NULL);CHKERRQ(ierr); 125*c4762a1bSJed Brown ctx.V = 1.0; 126*c4762a1bSJed Brown ierr = PetscOptionsScalar("-V","","",ctx.V,&ctx.V,NULL);CHKERRQ(ierr); 127*c4762a1bSJed Brown ctx.X = 1.0; 128*c4762a1bSJed Brown ierr = PetscOptionsScalar("-X","","",ctx.X,&ctx.X,NULL);CHKERRQ(ierr); 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 131*c4762a1bSJed Brown u[0] = 1; 132*c4762a1bSJed Brown u[1] = .7; 133*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = PetscOptionsGetVec(NULL,NULL,"-initial",U,NULL);CHKERRQ(ierr); 135*c4762a1bSJed Brown } 136*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&ctx.rand);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(ctx.rand);CHKERRQ(ierr); 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142*c4762a1bSJed Brown Create timestepping solver context 143*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151*c4762a1bSJed Brown Set initial conditions 152*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153*c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 154*c4762a1bSJed Brown 155*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156*c4762a1bSJed Brown Set solver options 157*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,2000.0);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 164*c4762a1bSJed Brown Solve nonlinear system 165*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 166*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 167*c4762a1bSJed Brown 168*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 170*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 171*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = PetscRandomDestroy(&ctx.rand);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = PetscFinalize(); 176*c4762a1bSJed Brown return ierr; 177*c4762a1bSJed Brown } 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown /*TEST 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown build: 183*c4762a1bSJed Brown requires: !complex !single 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown test: 186*c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_max_steps 10 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown test: 189*c4762a1bSJed Brown suffix: 2 190*c4762a1bSJed Brown args: -ts_max_steps 10 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown TEST*/ 193