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