xref: /petsc/src/ts/tutorials/ex36.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 {
40   PetscFunctionBeginUser;
41   *U = 0.4 * PetscSinReal(200 * PETSC_PI * t);
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46      Defines the DAE passed to the time solver
47 */
48 static PetscErrorCode IFunctionImplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *ctx)
49 {
50   const PetscScalar *y, *ydot;
51   PetscScalar       *f;
52 
53   PetscFunctionBeginUser;
54   /*  The next three lines allow us to access the entries of the vectors directly */
55   PetscCall(VecGetArrayRead(Y, &y));
56   PetscCall(VecGetArrayRead(Ydot, &ydot));
57   PetscCall(VecGetArrayWrite(F, &f));
58 
59   f[0] = ydot[0] / 1.e6 - ydot[1] / 1.e6 - PetscSinReal(200 * PETSC_PI * t) / 2500. + y[0] / 1000.;
60   f[1] = -ydot[0] / 1.e6 + ydot[1] / 1.e6 - 0.0006666766666666667 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 1.e8 + y[1] / 4500.;
61   f[2] = ydot[2] / 500000. + 1.e-6 - PetscExpReal((500 * (y[1] - y[2])) / 13.) / 1.e6 + y[2] / 9000.;
62   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.;
63   f[4] = (3 * ydot[4]) / 1.e6 - (3 * ydot[3]) / 1.e6 + y[4] / 9000.;
64 
65   PetscCall(VecRestoreArrayRead(Y, &y));
66   PetscCall(VecRestoreArrayRead(Ydot, &ydot));
67   PetscCall(VecRestoreArrayWrite(F, &f));
68   PetscFunctionReturn(0);
69 }
70 
71 /*
72      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
73 */
74 static PetscErrorCode IJacobianImplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *ctx)
75 {
76   PetscInt           rowcol[] = {0, 1, 2, 3, 4};
77   const PetscScalar *y, *ydot;
78   PetscScalar        J[5][5];
79 
80   PetscFunctionBeginUser;
81   PetscCall(VecGetArrayRead(Y, &y));
82   PetscCall(VecGetArrayRead(Ydot, &ydot));
83 
84   PetscCall(PetscMemzero(J, sizeof(J)));
85 
86   J[0][0] = a / 1.e6 + 0.001;
87   J[0][1] = -a / 1.e6;
88   J[1][0] = -a / 1.e6;
89   J[1][1] = a / 1.e6 + 0.00022222222222222223 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 2.6e6;
90   J[1][2] = -PetscExpReal((500 * (y[1] - y[2])) / 13.) / 2.6e6;
91   J[2][1] = -PetscExpReal((500 * (y[1] - y[2])) / 13.) / 26000.;
92   J[2][2] = a / 500000 + 0.00011111111111111112 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 26000.;
93   J[3][1] = (99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 2.6e6;
94   J[3][2] = (-99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 2.6e6;
95   J[3][3] = (3 * a) / 1.e6 + 0.00011111111111111112;
96   J[3][4] = -(3 * a) / 1.e6;
97   J[4][3] = -(3 * a) / 1.e6;
98   J[4][4] = (3 * a) / 1.e6 + 0.00011111111111111112;
99 
100   PetscCall(MatSetValues(B, 5, rowcol, 5, rowcol, &J[0][0], INSERT_VALUES));
101 
102   PetscCall(VecRestoreArrayRead(Y, &y));
103   PetscCall(VecRestoreArrayRead(Ydot, &ydot));
104 
105   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
106   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
107   if (A != B) {
108     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
109     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
110   }
111   PetscFunctionReturn(0);
112 }
113 
114 int main(int argc, char **argv)
115 {
116   TS           ts; /* ODE integrator */
117   Vec          Y;  /* solution will be stored here */
118   Mat          A;  /* Jacobian matrix */
119   PetscMPIInt  size;
120   PetscInt     n = 5;
121   PetscScalar *y;
122 
123   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124      Initialize program
125      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126   PetscFunctionBeginUser;
127   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
128   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
129   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
130 
131   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132     Create necessary matrix and vectors
133     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
135   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
136   PetscCall(MatSetFromOptions(A));
137   PetscCall(MatSetUp(A));
138 
139   PetscCall(MatCreateVecs(A, &Y, NULL));
140 
141   PetscCall(VecGetArray(Y, &y));
142   y[0] = 0.0;
143   y[1] = 3.0;
144   y[2] = y[1];
145   y[3] = 6.0;
146   y[4] = 0.0;
147   PetscCall(VecRestoreArray(Y, &y));
148 
149   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150      Create timestepping solver context
151      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
153   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
154   PetscCall(TSSetType(ts, TSARKIMEX));
155   /* Must use ARKIMEX with fully implicit stages since mass matrix is not the indentity */
156   PetscCall(TSARKIMEXSetType(ts, TSARKIMEXPRSSP2));
157   PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
158   /*PetscCall(TSSetType(ts,TSROSW));*/
159   PetscCall(TSSetIFunction(ts, NULL, IFunctionImplicit, NULL));
160   PetscCall(TSSetIJacobian(ts, A, A, IJacobianImplicit, NULL));
161 
162   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163      Set initial conditions
164    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165   PetscCall(TSSetSolution(ts, Y));
166 
167   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168      Set solver options
169    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170   PetscCall(TSSetMaxTime(ts, 0.15));
171   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
172   PetscCall(TSSetTimeStep(ts, .001));
173   PetscCall(TSSetFromOptions(ts));
174 
175   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176      Do time stepping
177      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178   PetscCall(TSSolve(ts, Y));
179 
180   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
182    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183   PetscCall(MatDestroy(&A));
184   PetscCall(VecDestroy(&Y));
185   PetscCall(TSDestroy(&ts));
186   PetscCall(PetscFinalize());
187   return 0;
188 }
189 
190 /*TEST
191     build:
192       requires: !single !complex
193     test:
194       args: -ts_monitor
195 
196 TEST*/
197