Actual source code: adapthist.c

  1: #include <petsc/private/tshistoryimpl.h>

  3: typedef struct {
  4:   TSHistory hist;
  5:   PetscBool bw;
  6: } TSAdapt_History;

  8: static PetscErrorCode TSAdaptChoose_History(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
  9: {
 10:   PetscInt        step;
 11:   TSAdapt_History *thadapt = (TSAdapt_History*)adapt->data;

 14:   TSGetStepNumber(ts,&step);
 15:   TSHistoryGetTimeStep(thadapt->hist,thadapt->bw,step+1,next_h);
 16:   *accept  = PETSC_TRUE;
 17:   *next_sc = 0;
 18:   *wlte    = -1;
 19:   *wltea   = -1;
 20:   *wlter   = -1;
 21:   return 0;
 22: }

 24: static PetscErrorCode TSAdaptReset_History(TSAdapt adapt)
 25: {
 26:   TSAdapt_History *thadapt = (TSAdapt_History*)adapt->data;

 28:   TSHistoryDestroy(&thadapt->hist);
 29:   return 0;
 30: }

 32: static PetscErrorCode TSAdaptDestroy_History(TSAdapt adapt)
 33: {
 34:   TSAdaptReset_History(adapt);
 35:   PetscFree(adapt->data);
 36:   return 0;
 37: }

 39: /* this is not public, as TSHistory is not a public object */
 40: PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt adapt, TSHistory hist, PetscBool backward)
 41: {
 42:   PetscReal      *hist_t;
 43:   PetscInt       n;
 44:   PetscBool      flg;

 48:   PetscObjectTypeCompare((PetscObject)adapt,TSADAPTHISTORY,&flg);
 49:   if (!flg) return 0;
 50:   TSHistoryGetHistory(hist,&n,(const PetscReal**)&hist_t,NULL,NULL);
 51:   TSAdaptHistorySetHistory(adapt,n,hist_t,backward);
 52:   return 0;
 53: }

 55: /*@
 56:    TSAdaptHistoryGetStep - Gets time and time step for a given step number in the history

 58:    Logically Collective on TSAdapt

 60:    Input Parameters:
 61: +  adapt    - the TSAdapt context
 62: -  step     - the step number

 64:    Output Parameters:
 65: +  t  - the time corresponding to the requested step (can be NULL)
 66: -  dt - the time step to be taken at the requested step (can be NULL)

 68:    Notes: The time history is internally copied, and the user can free the hist array. The user still needs to specify the termination of the solve via TSSetMaxSteps().

 70:    Level: advanced

 72: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptHistorySetTrajectory(), TSADAPTHISTORY
 73: @*/
 74: PetscErrorCode TSAdaptHistoryGetStep(TSAdapt adapt, PetscInt step, PetscReal *t, PetscReal *dt)
 75: {
 76:   TSAdapt_History *thadapt;
 77:   PetscBool       flg;

 81:   PetscObjectTypeCompare((PetscObject)adapt,TSADAPTHISTORY,&flg);
 83:   thadapt = (TSAdapt_History*)adapt->data;
 84:   TSHistoryGetTimeStep(thadapt->hist,thadapt->bw,step,dt);
 85:   TSHistoryGetTime(thadapt->hist,thadapt->bw,step,t);
 86:   return 0;
 87: }

 89: /*@
 90:    TSAdaptHistorySetHistory - Sets a time history in the adaptor

 92:    Logically Collective on TSAdapt

 94:    Input Parameters:
 95: +  adapt    - the TSAdapt context
 96: .  n        - size of the time history
 97: .  hist     - the time history
 98: -  backward - if the time history has to be followed backward

100:    Notes: The time history is internally copied, and the user can free the hist array. The user still needs to specify the termination of the solve via TSSetMaxSteps().

102:    Level: advanced

104: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptHistorySetTrajectory(), TSADAPTHISTORY
105: @*/
106: PetscErrorCode TSAdaptHistorySetHistory(TSAdapt adapt, PetscInt n, PetscReal hist[], PetscBool backward)
107: {
108:   TSAdapt_History *thadapt;
109:   PetscBool       flg;

115:   PetscObjectTypeCompare((PetscObject)adapt,TSADAPTHISTORY,&flg);
116:   if (!flg) return 0;
117:   thadapt = (TSAdapt_History*)adapt->data;
118:   TSHistoryDestroy(&thadapt->hist);
119:   TSHistoryCreate(PetscObjectComm((PetscObject)adapt),&thadapt->hist);
120:   TSHistorySetHistory(thadapt->hist,n,hist,NULL,PETSC_FALSE);
121:   thadapt->bw = backward;
122:   return 0;
123: }

125: /*@
126:    TSAdaptHistorySetTrajectory - Sets a time history in the adaptor from a given TSTrajectory

128:    Logically Collective on TSAdapt

130:    Input Parameters:
131: +  adapt    - the TSAdapt context
132: .  tj       - the TSTrajectory context
133: -  backward - if the time history has to be followed backward

135:    Notes: The time history is internally copied, and the user can destroy the TSTrajectory if not needed. The user still needs to specify the termination of the solve via TSSetMaxSteps().

137:    Level: advanced

139: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptHistorySetHistory(), TSADAPTHISTORY
140: @*/
141: PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt adapt, TSTrajectory tj, PetscBool backward)
142: {
143:   PetscBool       flg;

148:   PetscObjectTypeCompare((PetscObject)adapt,TSADAPTHISTORY,&flg);
149:   if (!flg) return 0;
150:   TSAdaptHistorySetTSHistory(adapt,tj->tsh,backward);
151:   return 0;
152: }

154: /*MC
155:    TSADAPTHISTORY - Time stepping controller that follows a given time history, used for Tangent Linear Model simulations

157:    Level: developer

159: .seealso: TS, TSAdapt, TSGetAdapt(), TSAdaptHistorySetHistory()
160: M*/
161: PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt adapt)
162: {
163:   TSAdapt_History *thadapt;

165:   PetscNew(&thadapt);
166:   adapt->matchstepfac[0] = PETSC_SMALL; /* prevent from accumulation errors */
167:   adapt->matchstepfac[1] = 0.0;         /* we will always match the final step, prevent TSAdaptChoose to mess with it */
168:   adapt->data            = thadapt;

170:   adapt->ops->choose  = TSAdaptChoose_History;
171:   adapt->ops->reset   = TSAdaptReset_History;
172:   adapt->ops->destroy = TSAdaptDestroy_History;
173:   return 0;
174: }