1 2 static char help[] = "Transistor amplifier.\n"; 3 4 /*F 5 ` This example illustrates the implementation of an implicit DAE index-1 of form M y'=f(t,y) with singular mass matrix, where 6 7 [ -C1 C1 ] 8 [ C1 -C1 ] 9 M =[ -C2 ]; Ck = k * 1e-06 10 [ -C3 C3] 11 [ C3 -C3] 12 13 [ -(U(t) - y[0])/1000 ] 14 [ -6/R + y[1]/4500 + 0.01 * h(y[1]-y[2]) ] 15 f(t,y)= [ y[2]/R - h(y[1]-y[2]) ] 16 [ (y[3]-6)/9000 + 0.99 * h([y1]-y[2]) ] 17 [ y[4]/9000 ] 18 19 U(t) = 0.4 * Sin(200 Pi t); h[V] = 1e-06 * Exp(V/0.026 - 1) ` 20 21 Useful options: -ts_monitor_lg_solution -ts_monitor_lg_timestep -lg_indicate_data_points 0 22 F*/ 23 24 /* 25 Include "petscts.h" so that we can use TS solvers. Note that this 26 file automatically includes: 27 petscsys.h - base PETSc routines petscvec.h - vectors 28 petscmat.h - matrices 29 petscis.h - index sets petscksp.h - Krylov subspace methods 30 petscviewer.h - viewers petscpc.h - preconditioners 31 petscksp.h - linear solvers 32 */ 33 #include <petscts.h> 34 35 FILE *gfilepointer_data, *gfilepointer_info; 36 37 /* Defines the source */ 38 PetscErrorCode Ue(PetscScalar t, PetscScalar *U) { 39 PetscFunctionBeginUser; 40 *U = 0.4 * PetscSinReal(200 * PETSC_PI * t); 41 PetscFunctionReturn(0); 42 } 43 44 /* 45 Defines the DAE passed to the time solver 46 */ 47 static PetscErrorCode IFunctionImplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *ctx) { 48 const PetscScalar *y, *ydot; 49 PetscScalar *f; 50 51 PetscFunctionBeginUser; 52 /* The next three lines allow us to access the entries of the vectors directly */ 53 PetscCall(VecGetArrayRead(Y, &y)); 54 PetscCall(VecGetArrayRead(Ydot, &ydot)); 55 PetscCall(VecGetArrayWrite(F, &f)); 56 57 f[0] = ydot[0] / 1.e6 - ydot[1] / 1.e6 - PetscSinReal(200 * PETSC_PI * t) / 2500. + y[0] / 1000.; 58 f[1] = -ydot[0] / 1.e6 + ydot[1] / 1.e6 - 0.0006666766666666667 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 1.e8 + y[1] / 4500.; 59 f[2] = ydot[2] / 500000. + 1.e-6 - PetscExpReal((500 * (y[1] - y[2])) / 13.) / 1.e6 + y[2] / 9000.; 60 f[3] = (3 * ydot[3]) / 1.e6 - (3 * ydot[4]) / 1.e6 - 0.0006676566666666666 + (99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 1.e8 + y[3] / 9000.; 61 f[4] = (3 * ydot[4]) / 1.e6 - (3 * ydot[3]) / 1.e6 + y[4] / 9000.; 62 63 PetscCall(VecRestoreArrayRead(Y, &y)); 64 PetscCall(VecRestoreArrayRead(Ydot, &ydot)); 65 PetscCall(VecRestoreArrayWrite(F, &f)); 66 PetscFunctionReturn(0); 67 } 68 69 /* 70 Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 71 */ 72 static PetscErrorCode IJacobianImplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *ctx) { 73 PetscInt rowcol[] = {0, 1, 2, 3, 4}; 74 const PetscScalar *y, *ydot; 75 PetscScalar J[5][5]; 76 77 PetscFunctionBeginUser; 78 PetscCall(VecGetArrayRead(Y, &y)); 79 PetscCall(VecGetArrayRead(Ydot, &ydot)); 80 81 PetscCall(PetscMemzero(J, sizeof(J))); 82 83 J[0][0] = a / 1.e6 + 0.001; 84 J[0][1] = -a / 1.e6; 85 J[1][0] = -a / 1.e6; 86 J[1][1] = a / 1.e6 + 0.00022222222222222223 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 2.6e6; 87 J[1][2] = -PetscExpReal((500 * (y[1] - y[2])) / 13.) / 2.6e6; 88 J[2][1] = -PetscExpReal((500 * (y[1] - y[2])) / 13.) / 26000.; 89 J[2][2] = a / 500000 + 0.00011111111111111112 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 26000.; 90 J[3][1] = (99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 2.6e6; 91 J[3][2] = (-99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 2.6e6; 92 J[3][3] = (3 * a) / 1.e6 + 0.00011111111111111112; 93 J[3][4] = -(3 * a) / 1.e6; 94 J[4][3] = -(3 * a) / 1.e6; 95 J[4][4] = (3 * a) / 1.e6 + 0.00011111111111111112; 96 97 PetscCall(MatSetValues(B, 5, rowcol, 5, rowcol, &J[0][0], INSERT_VALUES)); 98 99 PetscCall(VecRestoreArrayRead(Y, &y)); 100 PetscCall(VecRestoreArrayRead(Ydot, &ydot)); 101 102 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 103 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 104 if (A != B) { 105 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 106 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 107 } 108 PetscFunctionReturn(0); 109 } 110 111 int main(int argc, char **argv) { 112 TS ts; /* ODE integrator */ 113 Vec Y; /* solution will be stored here */ 114 Mat A; /* Jacobian matrix */ 115 PetscMPIInt size; 116 PetscInt n = 5; 117 PetscScalar *y; 118 119 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 120 Initialize program 121 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 122 PetscFunctionBeginUser; 123 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 124 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 125 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 126 127 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128 Create necessary matrix and vectors 129 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 131 PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 132 PetscCall(MatSetFromOptions(A)); 133 PetscCall(MatSetUp(A)); 134 135 PetscCall(MatCreateVecs(A, &Y, NULL)); 136 137 PetscCall(VecGetArray(Y, &y)); 138 y[0] = 0.0; 139 y[1] = 3.0; 140 y[2] = y[1]; 141 y[3] = 6.0; 142 y[4] = 0.0; 143 PetscCall(VecRestoreArray(Y, &y)); 144 145 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146 Create timestepping solver context 147 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 149 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 150 PetscCall(TSSetType(ts, TSARKIMEX)); 151 /* Must use ARKIMEX with fully implicit stages since mass matrix is not the indentity */ 152 PetscCall(TSARKIMEXSetType(ts, TSARKIMEXPRSSP2)); 153 PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1)); 154 /*PetscCall(TSSetType(ts,TSROSW));*/ 155 PetscCall(TSSetIFunction(ts, NULL, IFunctionImplicit, NULL)); 156 PetscCall(TSSetIJacobian(ts, A, A, IJacobianImplicit, NULL)); 157 158 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159 Set initial conditions 160 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 161 PetscCall(TSSetSolution(ts, Y)); 162 163 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 164 Set solver options 165 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 166 PetscCall(TSSetMaxTime(ts, 0.15)); 167 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 168 PetscCall(TSSetTimeStep(ts, .001)); 169 PetscCall(TSSetFromOptions(ts)); 170 171 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172 Do time stepping 173 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 174 PetscCall(TSSolve(ts, Y)); 175 176 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 177 Free work space. All PETSc objects should be destroyed when they are no longer needed. 178 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 179 PetscCall(MatDestroy(&A)); 180 PetscCall(VecDestroy(&Y)); 181 PetscCall(TSDestroy(&ts)); 182 PetscCall(PetscFinalize()); 183 return 0; 184 } 185 186 /*TEST 187 build: 188 requires: !single !complex 189 test: 190 args: -ts_monitor 191 192 TEST*/ 193