Actual source code: ex36A.c
1: static char help[] = "Transistor amplifier (autonomous).\n";
3: /*F
4: M y'=f(y)
6: Useful options: -ts_monitor_lg_solution -ts_monitor_lg_timestep -lg_indicate_data_points 0
7: F*/
9: /*
10: Include "petscts.h" so that we can use TS solvers. Note that this
11: file automatically includes:
12: petscsys.h - base PETSc routines petscvec.h - vectors
13: petscmat.h - matrices
14: petscis.h - index sets petscksp.h - Krylov subspace methods
15: petscviewer.h - viewers petscpc.h - preconditioners
16: petscksp.h - linear solvers
17: */
18: #include <petscts.h>
20: FILE *gfilepointer_data, *gfilepointer_info;
22: /* Defines the source */
23: PetscErrorCode Ue(PetscScalar t, PetscScalar *U)
24: {
25: PetscFunctionBeginUser;
26: U = 0.4 * sin(200 * pi * t);
27: PetscFunctionReturn(PETSC_SUCCESS);
28: }
29: * /
31: /*
32: Defines the DAE passed to the time solver
33: */
34: static PetscErrorCode IFunctionImplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *ctx)
35: {
36: const PetscScalar *y, *ydot;
37: PetscScalar *f;
39: PetscFunctionBeginUser;
40: /* The next three lines allow us to access the entries of the vectors directly */
41: PetscCall(VecGetArrayRead(Y, &y));
42: PetscCall(VecGetArrayRead(Ydot, &ydot));
43: PetscCall(VecGetArray(F, &f));
45: f[0] = PetscSinReal(200 * PETSC_PI * y[5]) / 2500. - y[0] / 1000. - ydot[0] / 1.e6 + ydot[1] / 1.e6;
46: f[1] = 0.0006666766666666667 - PetscExpReal((500 * (y[1] - y[2])) / 13.) / 1.e8 - y[1] / 4500. + ydot[0] / 1.e6 - ydot[1] / 1.e6;
47: f[2] = -1.e-6 + PetscExpReal((500 * (y[1] - y[2])) / 13.) / 1.e6 - y[2] / 9000. - ydot[2] / 500000.;
48: f[3] = 0.0006676566666666666 - (99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 1.e8 - y[3] / 9000. - (3 * ydot[3]) / 1.e6 + (3 * ydot[4]) / 1.e6;
49: f[4] = -y[4] / 9000. + (3 * ydot[3]) / 1.e6 - (3 * ydot[4]) / 1.e6;
50: f[5] = -1 + ydot[5];
52: PetscCall(VecRestoreArrayRead(Y, &y));
53: PetscCall(VecRestoreArrayRead(Ydot, &ydot));
54: PetscCall(VecRestoreArray(F, &f));
55: PetscFunctionReturn(PETSC_SUCCESS);
56: }
58: /*
59: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
60: */
61: static PetscErrorCode IJacobianImplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *ctx)
62: {
63: PetscInt rowcol[] = {0, 1, 2, 3, 4, 5};
64: const PetscScalar *y, *ydot;
65: PetscScalar J[6][6];
67: PetscFunctionBeginUser;
68: PetscCall(VecGetArrayRead(Y, &y));
69: PetscCall(VecGetArrayRead(Ydot, &ydot));
71: PetscCall(PetscMemzero(J, sizeof(J)));
73: J[0][0] = -0.001 - a / 1.e6;
74: J[0][1] = a / 1.e6;
75: J[0][5] = (2 * PETSC_PI * PetscCosReal(200 * PETSC_PI * y[5])) / 25.;
76: J[1][0] = a / 1.e6;
77: J[1][1] = -0.00022222222222222223 - a / 1.e6 - PetscExpReal((500 * (y[1] - y[2])) / 13.) / 2.6e6;
78: J[1][2] = PetscExpReal((500 * (y[1] - y[2])) / 13.) / 2.6e6;
79: J[2][1] = PetscExpReal((500 * (y[1] - y[2])) / 13.) / 26000.;
80: J[2][2] = -0.00011111111111111112 - a / 500000. - PetscExpReal((500 * (y[1] - y[2])) / 13.) / 26000.;
81: J[3][1] = (-99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 2.6e6;
82: J[3][2] = (99 * PetscExpReal((500 * (y[1] - y[2])) / 13.)) / 2.6e6;
83: J[3][3] = -0.00011111111111111112 - (3 * a) / 1.e6;
84: J[3][4] = (3 * a) / 1.e6;
85: J[4][3] = (3 * a) / 1.e6;
86: J[4][4] = -0.00011111111111111112 - (3 * a) / 1.e6;
87: J[5][5] = a;
89: PetscCall(MatSetValues(B, 6, rowcol, 6, rowcol, &J[0][0], INSERT_VALUES));
91: PetscCall(VecRestoreArrayRead(Y, &y));
92: PetscCall(VecRestoreArrayRead(Ydot, &ydot));
94: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
95: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
96: if (A != B) {
97: PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
98: PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
99: }
100: PetscFunctionReturn(PETSC_SUCCESS);
101: }
103: int main(int argc, char **argv)
104: {
105: TS ts; /* ODE integrator */
106: Vec Y; /* solution will be stored here */
107: Mat A; /* Jacobian matrix */
108: PetscMPIInt size;
109: PetscInt n = 6;
110: PetscScalar *y;
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Initialize program
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: PetscFunctionBeginUser;
116: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
117: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
118: PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Create necessary matrix and vectors
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
124: PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
125: PetscCall(MatSetFromOptions(A));
126: PetscCall(MatSetUp(A));
128: PetscCall(MatCreateVecs(A, &Y, NULL));
130: PetscCall(VecGetArray(Y, &y));
131: y[0] = 0.0;
132: y[1] = 3.0;
133: y[2] = y[1];
134: y[3] = 6.0;
135: y[4] = 0.0;
136: y[5] = 0.0;
137: PetscCall(VecRestoreArray(Y, &y));
139: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: Create timestepping solver context
141: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
143: PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
144: PetscCall(TSSetType(ts, TSARKIMEX));
145: PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
146: PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE));
147: /*PetscCall(TSSetType(ts,TSROSW));*/
148: PetscCall(TSSetIFunction(ts, NULL, IFunctionImplicit, NULL));
149: PetscCall(TSSetIJacobian(ts, A, A, IJacobianImplicit, NULL));
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Set initial conditions
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: PetscCall(TSSetSolution(ts, Y));
156: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: Set solver options
158: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159: PetscCall(TSSetMaxTime(ts, 0.15));
160: PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
161: PetscCall(TSSetTimeStep(ts, .001));
162: PetscCall(TSSetFromOptions(ts));
164: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165: Do Time stepping
166: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167: PetscCall(TSSolve(ts, Y));
169: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170: Free work space. All PETSc objects should be destroyed when they are no longer needed.
171: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172: PetscCall(MatDestroy(&A));
173: PetscCall(VecDestroy(&Y));
174: PetscCall(TSDestroy(&ts));
175: PetscCall(PetscFinalize());
176: return 0;
177: }