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