Actual source code: ex36SE.c

  1: static char help[] = "Transistor amplifier (semi-explicit).\n";

  3: /*F
  4:   [I 0] [y'] = f(t,y,z)
  5:   [0 0] [z'] = g(t,y,z)
  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 IFunctionSemiExplicit(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] = -400 * PetscSinReal(200 * PETSC_PI * t) + 1000 * y[3] + ydot[0];
 46:   f[1] = 0.5 - 1 / (2. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.)) + (500 * y[1]) / 9. + ydot[1];
 47:   f[2] = -222.5522222222222 + 33 / (100. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.)) + (1000 * y[4]) / 27. + ydot[2];
 48:   f[3] = 0.0006666766666666667 - 1 / (1.e8 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.)) + PetscSinReal(200 * PETSC_PI * t) / 2500. + y[0] / 4500. - (11 * y[3]) / 9000.;
 49:   f[4] = 0.0006676566666666666 - 99 / (1.e8 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.)) + y[2] / 9000. - y[4] / 4500.;

 51:   PetscCall(VecRestoreArrayRead(Y, &y));
 52:   PetscCall(VecRestoreArrayRead(Ydot, &ydot));
 53:   PetscCall(VecRestoreArray(F, &f));
 54:   PetscFunctionReturn(PETSC_SUCCESS);
 55: }

 57: /*
 58:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 59: */
 60: static PetscErrorCode IJacobianSemiExplicit(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *ctx)
 61: {
 62:   PetscInt           rowcol[] = {0, 1, 2, 3, 4};
 63:   const PetscScalar *y, *ydot;
 64:   PetscScalar        J[5][5];

 66:   PetscFunctionBeginUser;
 67:   PetscCall(VecGetArrayRead(Y, &y));
 68:   PetscCall(VecGetArrayRead(Ydot, &ydot));

 70:   PetscCall(PetscMemzero(J, sizeof(J)));

 72:   J[0][0] = a;
 73:   J[0][3] = 1000;
 74:   J[1][0] = 250 / (13. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 75:   J[1][1] = 55.55555555555556 + a + 250 / (13. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 76:   J[1][3] = -250 / (13. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 77:   J[2][0] = -165 / (13. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 78:   J[2][1] = -165 / (13. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 79:   J[2][2] = a;
 80:   J[2][3] = 165 / (13. * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 81:   J[2][4] = 37.03703703703704;
 82:   J[3][0] = 0.00022222222222222223 + 1 / (2.6e6 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 83:   J[3][1] = 1 / (2.6e6 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 84:   J[3][3] = -0.0012222222222222222 - 1 / (2.6e6 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 85:   J[4][0] = 99 / (2.6e6 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 86:   J[4][1] = 99 / (2.6e6 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 87:   J[4][2] = 0.00011111111111111112;
 88:   J[4][3] = -99 / (2.6e6 * PetscExpReal((500 * (y[0] + y[1] - y[3])) / 13.));
 89:   J[4][4] = -0.00022222222222222223;

 91:   PetscCall(MatSetValues(B, 5, rowcol, 5, rowcol, &J[0][0], INSERT_VALUES));

 93:   PetscCall(VecRestoreArrayRead(Y, &y));
 94:   PetscCall(VecRestoreArrayRead(Ydot, &ydot));

 96:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 97:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
 98:   if (A != B) {
 99:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
100:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
101:   }
102:   PetscFunctionReturn(PETSC_SUCCESS);
103: }

105: int main(int argc, char **argv)
106: {
107:   TS           ts; /* ODE integrator */
108:   Vec          Y;  /* solution will be stored here */
109:   Mat          A;  /* Jacobian matrix */
110:   PetscMPIInt  size;
111:   PetscInt     n = 5;
112:   PetscScalar *y;

114:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115:      Initialize program
116:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117:   PetscFunctionBeginUser;
118:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
119:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
120:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:     Create necessary matrix and vectors
124:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
126:   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
127:   PetscCall(MatSetFromOptions(A));
128:   PetscCall(MatSetUp(A));

130:   PetscCall(MatCreateVecs(A, &Y, NULL));

132:   PetscCall(VecGetArray(Y, &y));
133:   y[0] = -3.0;
134:   y[1] = 3.0;
135:   y[2] = 6.0;
136:   y[3] = 0.0;
137:   y[4] = 6.0;
138:   PetscCall(VecRestoreArray(Y, &y));

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:      Create timestepping solver context
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
144:   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
145:   PetscCall(TSSetType(ts, TSARKIMEX));
146:   PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
147:   PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE));
148:   /*PetscCall(TSSetType(ts,TSROSW));*/
149:   PetscCall(TSSetIFunction(ts, NULL, IFunctionSemiExplicit, NULL));
150:   PetscCall(TSSetIJacobian(ts, A, A, IJacobianSemiExplicit, NULL));

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:      Set initial conditions
154:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155:   PetscCall(TSSetSolution(ts, Y));

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158:      Set solver options
159:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160:   PetscCall(TSSetMaxTime(ts, 0.15));
161:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
162:   PetscCall(TSSetTimeStep(ts, .001));
163:   PetscCall(TSSetFromOptions(ts));

165:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:      Do Time stepping
167:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168:   PetscCall(TSSolve(ts, Y));

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
172:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   PetscCall(MatDestroy(&A));
174:   PetscCall(VecDestroy(&Y));
175:   PetscCall(TSDestroy(&ts));
176:   PetscCall(PetscFinalize());
177:   return 0;
178: }