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