PetscErrorCode TSEventHandler(TS ts) { PetscErrorCode ierr; TSEvent event; PetscReal t; Vec U; PetscInt i; PetscReal dt,dt_min; PetscInt rollback=0,in[2],out[2]; PetscInt fvalue_sign,fvalueprev_sign; PetscFunctionBegin; PetscValidHeaderSpecific(ts,TS_CLASSID,1); if (!ts->event) PetscFunctionReturn(0); event = ts->event; ierr = TSGetTime(ts,&t);CHKERRQ(ierr); ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); if (event->status == TSEVENT_NONE) { if (ts->steps == 1) /* After first successful step */ event->timestep_orig = ts->ptime - ts->ptime_prev; event->timestep_prev = dt; } if (event->status == TSEVENT_RESET_NEXTSTEP) { /* Restore time step */ dt = PetscMin(event->timestep_orig,event->timestep_prev); ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); event->status = TSEVENT_NONE; } if (event->status == TSEVENT_NONE) { event->ptime_end = t; } ierr = VecLockPush(U);CHKERRQ(ierr); ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); ierr = VecLockPop(U);CHKERRQ(ierr); for (i=0; i < event->nevents; i++) { if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { event->status = TSEVENT_ZERO; event->fvalue_right[i] = event->fvalue[i]; continue; } fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { switch (event->direction[i]) { case -1: if (fvalue_sign < 0) { rollback = 1; /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->fvalue_right[i] = event->fvalue[i]; event->side[i] = 1; if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; event->status = TSEVENT_LOCATED_INTERVAL; } break; case 1: if (fvalue_sign > 0) { rollback = 1; /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->fvalue_right[i] = event->fvalue[i]; event->side[i] = 1; if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; event->status = TSEVENT_LOCATED_INTERVAL; } break; case 0: rollback = 1; /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->fvalue_right[i] = event->fvalue[i]; event->side[i] = 1; if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; event->status = TSEVENT_LOCATED_INTERVAL; break; } } } in[0] = event->status; in[1] = rollback; ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); event->status = (TSEventStatus)out[0]; rollback = out[1]; if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; event->nevents_zero = 0; if (event->status == TSEVENT_ZERO) { for (i=0; i < event->nevents; i++) { if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { event->events_zero[event->nevents_zero++] = i; if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);CHKERRQ(ierr); } event->zerocrossing[i] = PETSC_FALSE; } event->side[i] = 0; } ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr); dt = event->ptime_end - t; if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */ dt = event->timestep_prev; event->status = TSEVENT_NONE; } ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); event->iterctr = 0; PetscFunctionReturn(0); } if (event->status == TSEVENT_LOCATED_INTERVAL) { ierr = TSRollBack(ts);CHKERRQ(ierr); ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); event->status = TSEVENT_PROCESSING; event->ptime_right = t; } else { for (i=0; i < event->nevents; i++) { if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) { /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); event->side[i] = -1; } event->fvalue_prev[i] = event->fvalue[i]; } if (event->monitor && event->status == TSEVENT_PROCESSING) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Stepping forward as no event detected in interval [%g - %g]\n",event->iterctr,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->ptime_prev = t; } if (event->status == TSEVENT_PROCESSING) event->iterctr++; ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr); PetscFunctionReturn(0); }
static PetscErrorCode TSStep_Theta(TS ts) { TS_Theta *th = (TS_Theta*)ts->data; PetscInt its,lits,reject,next_scheme; PetscReal next_time_step; TSAdapt adapt; PetscBool stageok,accept = PETSC_TRUE; PetscErrorCode ierr; PetscFunctionBegin; th->status = TS_STEP_INCOMPLETE; ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) { PetscReal shift = 1./(th->Theta*ts->time_step); th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; ierr = TSPreStep(ts);CHKERRQ(ierr); ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); } if (th->extrapolate) { ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); } else { ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); } ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); ts->snes_its += its; ts->ksp_its += lits; ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr); ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr); if (!stageok) {accept = PETSC_FALSE; goto reject_step;} ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr); th->status = TS_STEP_PENDING; /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr); ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); if (!accept) { /* Roll back the current step */ ts->ptime += next_time_step; /* This will be undone in rollback */ th->status = TS_STEP_INCOMPLETE; ierr = TSRollBack(ts);CHKERRQ(ierr); goto reject_step; } /* ignore next_scheme for now */ ts->ptime += ts->time_step; ts->time_step = next_time_step; ts->steps++; th->status = TS_STEP_COMPLETE; break; reject_step: if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) { ts->reason = TS_DIVERGED_STEP_REJECTED; ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); } continue; } PetscFunctionReturn(0); }