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 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 PetscErrorCode ierr; 37 const PetscScalar *y,*ydot; 38 PetscScalar *f; 39 40 PetscFunctionBegin; 41 /* The next three lines allow us to access the entries of the vectors directly */ 42 ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 43 ierr = VecGetArrayRead(Ydot,&ydot);CHKERRQ(ierr); 44 ierr = VecGetArrayWrite(F,&f);CHKERRQ(ierr); 45 46 f[0]= PetscSinReal(200*PETSC_PI*t)/2500. - y[0]/1000. - ydot[0]/1.e6 + ydot[1]/1.e6; 47 f[1]=0.0006666766666666667 - PetscExpReal((500*(y[1] - y[2]))/13.)/1.e8 - y[1]/4500. + ydot[0]/1.e6 - ydot[1]/1.e6; 48 f[2]=-1.e-6 + PetscExpReal((500*(y[1] - y[2]))/13.)/1.e6 - y[2]/9000. - ydot[2]/500000.; 49 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; 50 f[4]=-y[4]/9000. + (3*ydot[3])/1.e6 - (3*ydot[4])/1.e6; 51 52 ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 53 ierr = VecRestoreArrayRead(Ydot,&ydot);CHKERRQ(ierr); 54 ierr = VecRestoreArrayWrite(F,&f);CHKERRQ(ierr); 55 PetscFunctionReturn(0); 56 } 57 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 PetscErrorCode ierr; 64 PetscInt rowcol[] = {0,1,2,3,4}; 65 const PetscScalar *y,*ydot; 66 PetscScalar J[5][5]; 67 68 PetscFunctionBegin; 69 ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 70 ierr = VecGetArrayRead(Ydot,&ydot);CHKERRQ(ierr); 71 72 ierr = PetscMemzero(J,sizeof(J));CHKERRQ(ierr); 73 74 J[0][0]=-0.001 - a/1.e6; 75 J[0][1]=a/1.e6; 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 88 89 ierr = MatSetValues(B,5,rowcol,5,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 90 91 ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 92 ierr = VecRestoreArrayRead(Ydot,&ydot);CHKERRQ(ierr); 93 94 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 95 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96 if (A != B) { 97 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99 } 100 PetscFunctionReturn(0); 101 } 102 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 PetscErrorCode ierr; 109 PetscMPIInt size; 110 PetscInt n = 5; 111 PetscScalar *y; 112 113 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114 Initialize program 115 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 116 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 117 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 118 if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 119 120 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121 Create necessary matrix and vectors 122 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 123 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 124 ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 125 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 126 ierr = MatSetUp(A);CHKERRQ(ierr); 127 128 ierr = MatCreateVecs(A,&Y,NULL);CHKERRQ(ierr); 129 130 ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 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 ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 137 138 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139 Create timestepping solver context 140 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 142 ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 143 ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 144 ierr = TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);CHKERRQ(ierr); 145 ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr); 146 /*ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr);*/ 147 ierr = TSSetIFunction(ts,NULL,IFunctionImplicit,NULL);CHKERRQ(ierr); 148 ierr = TSSetIJacobian(ts,A,A,IJacobianImplicit,NULL);CHKERRQ(ierr); 149 150 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151 Set initial conditions 152 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153 ierr = TSSetSolution(ts,Y);CHKERRQ(ierr); 154 155 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156 Set solver options 157 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158 ierr = TSSetMaxTime(ts,0.15);CHKERRQ(ierr); 159 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 160 ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr); 161 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 162 163 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 164 Do time stepping 165 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 166 ierr = TSSolve(ts,Y);CHKERRQ(ierr); 167 168 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169 Free work space. All PETSc objects should be destroyed when they are no longer needed. 170 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 171 ierr = MatDestroy(&A);CHKERRQ(ierr); 172 ierr = VecDestroy(&Y);CHKERRQ(ierr); 173 ierr = TSDestroy(&ts);CHKERRQ(ierr); 174 ierr = PetscFinalize(); 175 return ierr; 176 } 177 178 /*TEST 179 build: 180 requires: !single !complex 181 test: 182 183 TEST*/ 184