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