Actual source code: ex41.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Parallel bouncing ball example to test TS event feature.\n";

  4: /*
  5:   The dynamics of the bouncing ball is described by the ODE
  6:                   u1_t = u2
  7:                   u2_t = -9.8
  8:   
  9:   Each processor is assigned one ball.
 10:   
 11:   The event function routine checks for the ball hitting the
 12:   ground (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by
 13:   a factor of 0.9 and its height set to 1.0*rank.
 14: */

 16: #include <petscts.h>

 20: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx)
 21: {
 23:   PetscScalar    *u;

 26:   /* Event for ball height */
 27:   VecGetArray(U,&u);
 28:   fvalue[0] = u[0];
 29:   VecRestoreArray(U,&u);
 30:   return(0);
 31: }

 35: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,void* ctx)
 36: {
 38:   PetscScalar    *u;
 39:   PetscMPIInt    rank;

 42:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 43:   if (nevents) {
 44:     VecGetArray(U,&u);
 45:     PetscPrintf(PETSC_COMM_SELF,"Processor [%d]: Ball hit the ground at t = %f seconds\n",rank,t);
 46:     /* Set new initial conditions with .9 attenuation */
 47:     u[0] = 1.0*rank;
 48:     u[1] = -0.9*u[1];
 49:     VecRestoreArray(U,&u);
 50:   }
 51:   TSSetSolution(ts,U);
 52:   return(0);
 53: }

 57: /*
 58:      Defines the ODE passed to the ODE solver
 59: */
 60: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
 61: {
 63:   PetscScalar    *u,*udot,*f;

 66:   /*  The next three lines allow us to access the entries of the vectors directly */
 67:   VecGetArray(U,&u);
 68:   VecGetArray(Udot,&udot);
 69:   VecGetArray(F,&f);

 71:   f[0] = udot[0] - u[1];
 72:   f[1] = udot[1] + 9.8;

 74:   VecRestoreArray(U,&u);
 75:   VecRestoreArray(Udot,&udot);
 76:   VecRestoreArray(F,&f);
 77:   return(0);
 78: }

 82: /*
 83:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 84: */
 85: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
 86: {
 88:   PetscInt       rowcol[2];
 89:   PetscScalar    *u,*udot,J[2][2];
 90:   PetscInt        rstart;

 93:   VecGetArray(U,&u);
 94:   VecGetArray(Udot,&udot);

 96:   MatGetOwnershipRange(A,&rstart,NULL);
 97:   rowcol[0] = rstart; rowcol[1] = rstart+1;

 99:   J[0][0] = a;                       J[0][1] = -1;
100:   J[1][0] = 0.0;                     J[1][1] = a;
101:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);

103:   VecRestoreArray(U,&u);
104:   VecRestoreArray(Udot,&udot);

106:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
107:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
108:   if (A != B) {
109:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
110:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
111:   }
112:   return(0);
113: }
114: 
117: int main(int argc,char **argv)
118: {
119:   TS             ts;            /* ODE integrator */
120:   Vec            U;             /* solution will be stored here */
121:   Mat            A;             /* Jacobian matrix */
123:   PetscMPIInt    rank;
124:   PetscInt       n = 2,direction=-1;
125:   PetscScalar    *u;
126:   PetscBool      terminate=PETSC_FALSE;

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:      Initialize program
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:   PetscInitialize(&argc,&argv,(char*)0,help);
132:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

134:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135:     Create necessary matrix and vectors
136:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137:   MatCreate(PETSC_COMM_WORLD,&A);
138:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
139:   MatSetType(A,MATDENSE);
140:   MatSetFromOptions(A);
141:   MatSetUp(A);

143:   MatGetVecs(A,&U,NULL);

145:   VecGetArray(U,&u);
146:   u[0] = 1.0*rank;
147:   u[1] = 20.0;
148:   VecRestoreArray(U,&u);

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:      Create timestepping solver context
152:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153:   TSCreate(PETSC_COMM_WORLD,&ts);
154:   TSSetProblemType(ts,TS_NONLINEAR);
155:   TSSetType(ts,TSROSW);
156:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,NULL);
157:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,NULL);

159:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160:      Set initial conditions
161:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162:   TSSetSolution(ts,U);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Set solver options
166:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167:   TSSetDuration(ts,1000,30.0);
168:   TSSetInitialTimeStep(ts,0.0,0.1);
169:   TSSetFromOptions(ts);
170: 
171:   TSSetEventMonitor(ts,1,&direction,&terminate,EventFunction,PostEventFunction,NULL);

173:   TSAdapt adapt;
174:   TSGetAdapt(ts,&adapt);
175:   /* The adapative time step controller could take very large timesteps resulting in 
176:      the same event occuring multiple times in the same interval. A max. step 
177:      limit is enforced here to avoid this issue. 
178:   */
179:   TSAdaptSetStepLimits(adapt,0.0,0.5);
180:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181:      Run timestepping solver
182:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183:     TSSolve(ts,U);

185:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
187:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188: 
189:   MatDestroy(&A);
190:   VecDestroy(&U);
191:   TSDestroy(&ts);
192: 
193:   PetscFinalize();
194:   return(0);
195: }