1474dd773SHong Zhang /* 2474dd773SHong Zhang Code for time stepping with the multi-rate Runge-Kutta method 3474dd773SHong Zhang 4474dd773SHong Zhang Notes: 5474dd773SHong Zhang 1) The general system is written as 6474dd773SHong Zhang Udot = F(t,U) for the nonsplit version of multi-rate RK, 7474dd773SHong Zhang user should give the indexes for both slow and fast components; 8474dd773SHong Zhang 2) The general system is written as 9474dd773SHong Zhang Usdot = Fs(t,Us,Uf) 10474dd773SHong Zhang Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits, 11474dd773SHong Zhang user should partioned RHS by themselves and also provide the indexes for both slow and fast components. 12474dd773SHong Zhang */ 13474dd773SHong Zhang 14474dd773SHong Zhang #include <petsc/private/tsimpl.h> 15474dd773SHong Zhang #include <petscdm.h> 16474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h> 17474dd773SHong Zhang 18474dd773SHong Zhang static TSRKType TSMRKDefault = TSRK2A; 19474dd773SHong Zhang 20e5bffa22SHong Zhang static PetscErrorCode TSSetUp_MRKNONSPLIT(TS ts) 21e5bffa22SHong Zhang { 22e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 23e5bffa22SHong Zhang RKTableau tab = rk->tableau; 24e5bffa22SHong Zhang PetscErrorCode ierr; 25e5bffa22SHong Zhang 26e5bffa22SHong Zhang PetscFunctionBegin; 27e5bffa22SHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr); 28e5bffa22SHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr); 29e5bffa22SHong Zhang if (!rk->is_slow || !rk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use multirate RK"); 30e5bffa22SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr); 31e5bffa22SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr); 32e5bffa22SHong Zhang if (!rk->subts_slow || !rk->subts_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS"); 33e5bffa22SHong Zhang ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr); 34e5bffa22SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 35*63a6f1b4SHong Zhang rk->subts_current = rk->subts_fast; 36e5bffa22SHong Zhang PetscFunctionReturn(0); 37e5bffa22SHong Zhang } 38e5bffa22SHong Zhang 39e5bffa22SHong Zhang static PetscErrorCode TSReset_MRKNONSPLIT(TS ts) 40e5bffa22SHong Zhang { 41e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 42e5bffa22SHong Zhang RKTableau tab = rk->tableau; 43e5bffa22SHong Zhang PetscErrorCode ierr; 44e5bffa22SHong Zhang 45e5bffa22SHong Zhang PetscFunctionBegin; 46e5bffa22SHong Zhang ierr = VecDestroy(&rk->X0);CHKERRQ(ierr); 47e5bffa22SHong Zhang ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 48e5bffa22SHong Zhang PetscFunctionReturn(0); 49e5bffa22SHong Zhang } 50e5bffa22SHong Zhang 51e5bffa22SHong Zhang static PetscErrorCode TSInterpolate_MRKNONSPLIT(TS ts,PetscReal itime,Vec X) 52e5bffa22SHong Zhang { 53e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 54e5bffa22SHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 55*63a6f1b4SHong Zhang PetscReal h = ts->time_step; 56e5bffa22SHong Zhang PetscReal tt,t; 57e5bffa22SHong Zhang PetscScalar *b; 58e5bffa22SHong Zhang const PetscReal *B = rk->tableau->binterp; 59e5bffa22SHong Zhang PetscErrorCode ierr; 60e5bffa22SHong Zhang 61e5bffa22SHong Zhang PetscFunctionBegin; 62e5bffa22SHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 63*63a6f1b4SHong Zhang t = (itime - rk->ptime)/h; 64e5bffa22SHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 65e5bffa22SHong Zhang for (i=0; i<s; i++) b[i] = 0; 66e5bffa22SHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 67e5bffa22SHong Zhang for (i=0; i<s; i++) { 68e5bffa22SHong Zhang b[i] += h * B[i*p+j] * tt; 69e5bffa22SHong Zhang } 70e5bffa22SHong Zhang } 71e5bffa22SHong Zhang ierr = VecCopy(rk->X0,X);CHKERRQ(ierr); 72e5bffa22SHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 73e5bffa22SHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 74e5bffa22SHong Zhang PetscFunctionReturn(0); 75e5bffa22SHong Zhang } 76e5bffa22SHong Zhang 77*63a6f1b4SHong Zhang static PetscErrorCode TSStepRefine_MRKNONSPLIT(TS ts) 78*63a6f1b4SHong Zhang { 79*63a6f1b4SHong Zhang TS previousts,subts; 80*63a6f1b4SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 81*63a6f1b4SHong Zhang RKTableau tab = rk->tableau; 82*63a6f1b4SHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 83*63a6f1b4SHong Zhang Vec vec_fast,subvec_fast; 84*63a6f1b4SHong Zhang const PetscInt s = tab->s; 85*63a6f1b4SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 86*63a6f1b4SHong Zhang PetscScalar *w = rk->work; 87*63a6f1b4SHong Zhang PetscInt i,j,k; 88*63a6f1b4SHong Zhang PetscReal t = ts->ptime,h = ts->time_step; 89*63a6f1b4SHong Zhang PetscErrorCode ierr; 90*63a6f1b4SHong Zhang 91*63a6f1b4SHong Zhang ierr = VecDuplicate(ts->vec_sol,&vec_fast);CHKERRQ(ierr); 92*63a6f1b4SHong Zhang previousts = rk->subts_current; 93*63a6f1b4SHong Zhang ierr = TSRHSSplitGetSubTS(rk->subts_current,"fast",&subts);CHKERRQ(ierr); 94*63a6f1b4SHong Zhang ierr = TSRHSSplitGetSubTS(subts,"fast",&subts);CHKERRQ(ierr); 95*63a6f1b4SHong Zhang for (k=0; k<rk->dtratio; k++) { 96*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 97*63a6f1b4SHong Zhang ierr = TSInterpolate_MRKNONSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr); 98*63a6f1b4SHong Zhang for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 99*63a6f1b4SHong Zhang /* update the fast components in the stage value, the slow components will be ignored, so we do not care the slow part in vec_fast */ 100*63a6f1b4SHong Zhang ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr); 101*63a6f1b4SHong Zhang ierr = VecMAXPY(vec_fast,i,w,YdotRHS);CHKERRQ(ierr); 102*63a6f1b4SHong Zhang /* update the fast components in the stage value */ 103*63a6f1b4SHong Zhang ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 104*63a6f1b4SHong Zhang ierr = VecISCopy(Y[i],rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr); 105*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 106*63a6f1b4SHong Zhang /* compute the stage RHS */ 107*63a6f1b4SHong Zhang ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 108*63a6f1b4SHong Zhang } 109*63a6f1b4SHong Zhang ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr); 110*63a6f1b4SHong Zhang ierr = TSEvaluateStep(ts,tab->order,vec_fast,NULL);CHKERRQ(ierr); 111*63a6f1b4SHong Zhang /* update the fast components in the solution */ 112*63a6f1b4SHong Zhang ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 113*63a6f1b4SHong Zhang ierr = VecISCopy(ts->vec_sol,rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr); 114*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 115*63a6f1b4SHong Zhang 116*63a6f1b4SHong Zhang if (subts) { 117*63a6f1b4SHong Zhang Vec *YdotRHS_copy; 118*63a6f1b4SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,s,&YdotRHS_copy);CHKERRQ(ierr); 119*63a6f1b4SHong Zhang rk->subts_current = rk->subts_fast; 120*63a6f1b4SHong Zhang ts->ptime = t+k*h/rk->dtratio; 121*63a6f1b4SHong Zhang ts->time_step = h/rk->dtratio; 122*63a6f1b4SHong Zhang ierr = TSRHSSplitGetIS(rk->subts_current,"fast",&rk->is_fast);CHKERRQ(ierr); 123*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 124*63a6f1b4SHong Zhang ierr = VecCopy(rk->YdotRHS_slow[i],YdotRHS_copy[i]);CHKERRQ(ierr); 125*63a6f1b4SHong Zhang ierr = VecCopy(YdotRHS[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr); 126*63a6f1b4SHong Zhang } 127*63a6f1b4SHong Zhang 128*63a6f1b4SHong Zhang ierr = TSStepRefine_MRKNONSPLIT(ts);CHKERRQ(ierr); 129*63a6f1b4SHong Zhang 130*63a6f1b4SHong Zhang rk->subts_current = previousts; 131*63a6f1b4SHong Zhang ts->ptime = t; 132*63a6f1b4SHong Zhang ts->time_step = h; 133*63a6f1b4SHong Zhang ierr = TSRHSSplitGetIS(previousts,"fast",&rk->is_fast);CHKERRQ(ierr); 134*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 135*63a6f1b4SHong Zhang ierr = VecCopy(YdotRHS_copy[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr); 136*63a6f1b4SHong Zhang } 137*63a6f1b4SHong Zhang ierr = VecDestroyVecs(s,&YdotRHS_copy);CHKERRQ(ierr); 138*63a6f1b4SHong Zhang } 139*63a6f1b4SHong Zhang } 140*63a6f1b4SHong Zhang ierr = VecDestroy(&vec_fast);CHKERRQ(ierr); 141*63a6f1b4SHong Zhang PetscFunctionReturn(0); 142*63a6f1b4SHong Zhang } 143*63a6f1b4SHong Zhang 144e5bffa22SHong Zhang static PetscErrorCode TSStep_MRKNONSPLIT(TS ts) 145e5bffa22SHong Zhang { 146e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 147e5bffa22SHong Zhang RKTableau tab = rk->tableau; 148e5bffa22SHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow; 149e5bffa22SHong Zhang Vec stage_slow,sol_slow; /* vectors store the slow components */ 150e5bffa22SHong Zhang Vec subvec_slow; /* sub vector to store the slow components */ 151e5bffa22SHong Zhang IS is_slow = rk->is_slow; 152e5bffa22SHong Zhang const PetscInt s = tab->s; 153e5bffa22SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 154e5bffa22SHong Zhang PetscScalar *w = rk->work; 155*63a6f1b4SHong Zhang PetscInt i,j,dtratio = rk->dtratio; 156e5bffa22SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 157e5bffa22SHong Zhang PetscErrorCode ierr; 158e5bffa22SHong Zhang 159e5bffa22SHong Zhang PetscFunctionBegin; 160e5bffa22SHong Zhang rk->status = TS_STEP_INCOMPLETE; 161e5bffa22SHong Zhang ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr); 162e5bffa22SHong Zhang ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr); 163e5bffa22SHong Zhang ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr); 164e5bffa22SHong Zhang for (i=0; i<s; i++) { 165e5bffa22SHong Zhang rk->stage_time = t + h*c[i]; 166e5bffa22SHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 167e5bffa22SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 168e5bffa22SHong Zhang for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 169e5bffa22SHong Zhang ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr); 170e5bffa22SHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 171e5bffa22SHong Zhang /* compute the stage RHS */ 172e5bffa22SHong Zhang ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 173e5bffa22SHong Zhang } 174e5bffa22SHong Zhang /* update the slow components in the solution */ 175e5bffa22SHong Zhang rk->YdotRHS = YdotRHS_slow; 176e5bffa22SHong Zhang rk->dtratio = 1; 177e5bffa22SHong Zhang ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr); 178e5bffa22SHong Zhang rk->dtratio = dtratio; 179e5bffa22SHong Zhang rk->YdotRHS = YdotRHS; 180e5bffa22SHong Zhang /* update the slow components in the solution */ 181e5bffa22SHong Zhang ierr = VecGetSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr); 182e5bffa22SHong Zhang ierr = VecISCopy(ts->vec_sol,is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr); 183e5bffa22SHong Zhang ierr = VecRestoreSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr); 184e5bffa22SHong Zhang 185*63a6f1b4SHong Zhang rk->subts_current = ts; 186*63a6f1b4SHong Zhang rk->ptime = t; 187*63a6f1b4SHong Zhang rk->time_step = h; 188*63a6f1b4SHong Zhang ierr = TSStepRefine_MRKNONSPLIT(ts);CHKERRQ(ierr); 189*63a6f1b4SHong Zhang 190e5bffa22SHong Zhang ts->ptime = t + ts->time_step; 191e5bffa22SHong Zhang ts->time_step = next_time_step; 192e5bffa22SHong Zhang rk->status = TS_STEP_COMPLETE; 193*63a6f1b4SHong Zhang 194e5bffa22SHong Zhang /* free memory of work vectors */ 195e5bffa22SHong Zhang ierr = VecDestroy(&stage_slow);CHKERRQ(ierr); 196e5bffa22SHong Zhang ierr = VecDestroy(&sol_slow);CHKERRQ(ierr); 197e5bffa22SHong Zhang PetscFunctionReturn(0); 198e5bffa22SHong Zhang } 199e5bffa22SHong Zhang 200*63a6f1b4SHong Zhang /* 201*63a6f1b4SHong Zhang Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest. 202*63a6f1b4SHong Zhang */ 203*63a6f1b4SHong Zhang static PetscErrorCode TSCopyDM(TS tssrc,TS tsdest) 204*63a6f1b4SHong Zhang { 205*63a6f1b4SHong Zhang DM newdm,dmsrc,dmdest; 206*63a6f1b4SHong Zhang PetscErrorCode ierr; 207*63a6f1b4SHong Zhang 208*63a6f1b4SHong Zhang PetscFunctionBegin; 209*63a6f1b4SHong Zhang ierr = TSGetDM(tssrc,&dmsrc);CHKERRQ(ierr); 210*63a6f1b4SHong Zhang ierr = DMClone(dmsrc,&newdm);CHKERRQ(ierr); 211*63a6f1b4SHong Zhang ierr = TSGetDM(tsdest,&dmdest);CHKERRQ(ierr); 212*63a6f1b4SHong Zhang ierr = DMCopyDMTS(dmdest,newdm);CHKERRQ(ierr); 213*63a6f1b4SHong Zhang ierr = DMCopyDMSNES(dmdest,newdm);CHKERRQ(ierr); 214*63a6f1b4SHong Zhang ierr = TSSetDM(tsdest,newdm);CHKERRQ(ierr); 215*63a6f1b4SHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 216*63a6f1b4SHong Zhang PetscFunctionReturn(0); 217*63a6f1b4SHong Zhang } 218*63a6f1b4SHong Zhang 219474dd773SHong Zhang static PetscErrorCode TSSetUp_MRKSPLIT(TS ts) 220474dd773SHong Zhang { 221bb42530cSHong Zhang TS_RK *rk = (TS_RK*)ts->data,*nextlevelrk,*currentlevelrk; 222bb42530cSHong Zhang TS nextlevelts; 223*63a6f1b4SHong Zhang Vec X0; 224474dd773SHong Zhang PetscErrorCode ierr; 225474dd773SHong Zhang 226474dd773SHong Zhang PetscFunctionBegin; 227474dd773SHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr); 228474dd773SHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr); 229474dd773SHong Zhang if (!rk->is_slow || !rk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use -ts_type bsi"); 230474dd773SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr); 231474dd773SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr); 232474dd773SHong Zhang if (!rk->subts_slow || !rk->subts_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS"); 233474dd773SHong Zhang 234*63a6f1b4SHong Zhang ierr = VecDuplicate(ts->vec_sol,&X0);CHKERRQ(ierr); 235bb42530cSHong Zhang /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */ 236bb42530cSHong Zhang currentlevelrk = rk; 237bb42530cSHong Zhang while (currentlevelrk->subts_fast) { 238*63a6f1b4SHong Zhang ierr = PetscMalloc1(rk->tableau->s,¤tlevelrk->YdotRHS_fast);CHKERRQ(ierr); 239*63a6f1b4SHong Zhang ierr = PetscMalloc1(rk->tableau->s,¤tlevelrk->YdotRHS_slow);CHKERRQ(ierr); 240*63a6f1b4SHong Zhang ierr = PetscObjectReference((PetscObject)X0);CHKERRQ(ierr); 241*63a6f1b4SHong Zhang currentlevelrk->X0 = X0; 242*63a6f1b4SHong Zhang currentlevelrk->ts_root = ts; 243*63a6f1b4SHong Zhang 244bb42530cSHong Zhang /* set up the ts for the slow part */ 245bb42530cSHong Zhang nextlevelts = currentlevelrk->subts_slow; 246bb42530cSHong Zhang ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr); 247bb42530cSHong Zhang nextlevelrk->tableau = rk->tableau; 248bb42530cSHong Zhang nextlevelrk->work = rk->work; 249bb42530cSHong Zhang nextlevelrk->Y = rk->Y; 250*63a6f1b4SHong Zhang nextlevelrk->YdotRHS = rk->YdotRHS; 251bb42530cSHong Zhang nextlevelts->data = (void*)nextlevelrk; 252*63a6f1b4SHong Zhang ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr); 253bb42530cSHong Zhang ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr); 254bb42530cSHong Zhang 255bb42530cSHong Zhang /* set up the ts for the fast part */ 256bb42530cSHong Zhang nextlevelts = currentlevelrk->subts_fast; 257bb42530cSHong Zhang ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr); 258bb42530cSHong Zhang nextlevelrk->tableau = rk->tableau; 259bb42530cSHong Zhang nextlevelrk->work = rk->work; 260bb42530cSHong Zhang nextlevelrk->Y = rk->Y; 261*63a6f1b4SHong Zhang nextlevelrk->YdotRHS = rk->YdotRHS; 262bb42530cSHong Zhang nextlevelrk->dtratio = rk->dtratio; 263bb42530cSHong Zhang ierr = TSRHSSplitGetIS(nextlevelts,"slow",&nextlevelrk->is_slow);CHKERRQ(ierr); 264bb42530cSHong Zhang ierr = TSRHSSplitGetSubTS(nextlevelts,"slow",&nextlevelrk->subts_slow);CHKERRQ(ierr); 265bb42530cSHong Zhang ierr = TSRHSSplitGetIS(nextlevelts,"fast",&nextlevelrk->is_fast);CHKERRQ(ierr); 266bb42530cSHong Zhang ierr = TSRHSSplitGetSubTS(nextlevelts,"fast",&nextlevelrk->subts_fast);CHKERRQ(ierr); 267bb42530cSHong Zhang nextlevelts->data = (void*)nextlevelrk; 268*63a6f1b4SHong Zhang ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr); 269bb42530cSHong Zhang ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr); 270bb42530cSHong Zhang 271bb42530cSHong Zhang currentlevelrk = nextlevelrk; 272bb42530cSHong Zhang } 273*63a6f1b4SHong Zhang ierr = VecDestroy(&X0);CHKERRQ(ierr); 274474dd773SHong Zhang PetscFunctionReturn(0); 275474dd773SHong Zhang } 276474dd773SHong Zhang 277474dd773SHong Zhang static PetscErrorCode TSReset_MRKSPLIT(TS ts) 278474dd773SHong Zhang { 279474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 280474dd773SHong Zhang PetscErrorCode ierr; 281474dd773SHong Zhang 282474dd773SHong Zhang PetscFunctionBegin; 283bb42530cSHong Zhang if (rk->subts_slow) { 284bb42530cSHong Zhang ierr = PetscFree(rk->subts_slow->data);CHKERRQ(ierr); 285bb42530cSHong Zhang rk->subts_slow = NULL; 286bb42530cSHong Zhang } 287bb42530cSHong Zhang if (rk->subts_fast) { 288*63a6f1b4SHong Zhang ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr); 289*63a6f1b4SHong Zhang ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr); 290bb42530cSHong Zhang ierr = VecDestroy(&rk->X0);CHKERRQ(ierr); 291bb42530cSHong Zhang ierr = TSReset_MRKSPLIT(rk->subts_fast); 292bb42530cSHong Zhang ierr = PetscFree(rk->subts_fast->data);CHKERRQ(ierr); 293bb42530cSHong Zhang rk->subts_fast = NULL; 294bb42530cSHong Zhang } 295474dd773SHong Zhang PetscFunctionReturn(0); 296474dd773SHong Zhang } 297474dd773SHong Zhang 298474dd773SHong Zhang static PetscErrorCode TSInterpolate_MRKSPLIT(TS ts,PetscReal itime,Vec X) 299474dd773SHong Zhang { 300474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 301*63a6f1b4SHong Zhang Vec Xslow; 302474dd773SHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 303474dd773SHong Zhang PetscReal h; 304474dd773SHong Zhang PetscReal tt,t; 305474dd773SHong Zhang PetscScalar *b; 306474dd773SHong Zhang const PetscReal *B = rk->tableau->binterp; 307474dd773SHong Zhang PetscErrorCode ierr; 308474dd773SHong Zhang 309474dd773SHong Zhang PetscFunctionBegin; 310474dd773SHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 311474dd773SHong Zhang 312474dd773SHong Zhang switch (rk->status) { 313474dd773SHong Zhang case TS_STEP_INCOMPLETE: 314474dd773SHong Zhang case TS_STEP_PENDING: 315474dd773SHong Zhang h = ts->time_step; 316474dd773SHong Zhang t = (itime - ts->ptime)/h; 317474dd773SHong Zhang break; 318474dd773SHong Zhang case TS_STEP_COMPLETE: 319474dd773SHong Zhang h = ts->ptime - ts->ptime_prev; 320474dd773SHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 321474dd773SHong Zhang break; 322474dd773SHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 323474dd773SHong Zhang } 324474dd773SHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 325474dd773SHong Zhang for (i=0; i<s; i++) b[i] = 0; 326474dd773SHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 327474dd773SHong Zhang for (i=0; i<s; i++) { 328474dd773SHong Zhang b[i] += h * B[i*p+j] * tt; 329474dd773SHong Zhang } 330474dd773SHong Zhang } 331*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 332*63a6f1b4SHong Zhang ierr = VecGetSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr); 333*63a6f1b4SHong Zhang } 334*63a6f1b4SHong Zhang ierr = VecGetSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr); 335*63a6f1b4SHong Zhang ierr = VecISCopy(rk->X0,rk->is_slow,SCATTER_REVERSE,Xslow);CHKERRQ(ierr); 336*63a6f1b4SHong Zhang ierr = VecMAXPY(Xslow,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 337*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr); 338*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 339*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr); 340*63a6f1b4SHong Zhang } 341474dd773SHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 342474dd773SHong Zhang PetscFunctionReturn(0); 343474dd773SHong Zhang } 344474dd773SHong Zhang 345474dd773SHong Zhang /* 346474dd773SHong Zhang This is for partitioned RHS multirate RK method 347474dd773SHong Zhang The step completion formula is 348474dd773SHong Zhang 349474dd773SHong Zhang x1 = x0 + h b^T YdotRHS 350474dd773SHong Zhang 351474dd773SHong Zhang */ 352474dd773SHong Zhang static PetscErrorCode TSEvaluateStep_MRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done) 353474dd773SHong Zhang { 354474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 355474dd773SHong Zhang RKTableau tab = rk->tableau; 356474dd773SHong Zhang Vec Xslow,Xfast; /* subvectors of X which store slow components and fast components respectively */ 357474dd773SHong Zhang PetscScalar *w = rk->work; 358474dd773SHong Zhang PetscReal h = ts->time_step; 359474dd773SHong Zhang PetscInt s = tab->s,j; 360474dd773SHong Zhang PetscErrorCode ierr; 361474dd773SHong Zhang 362474dd773SHong Zhang PetscFunctionBegin; 363474dd773SHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 364474dd773SHong Zhang if (rk->slow) { 365474dd773SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]; 366474dd773SHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr); 367*63a6f1b4SHong Zhang ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr); 368474dd773SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);; 369474dd773SHong Zhang } else { 370474dd773SHong Zhang for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j]; 371474dd773SHong Zhang ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr); 372*63a6f1b4SHong Zhang ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr); 373474dd773SHong Zhang ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr); 374474dd773SHong Zhang } 375474dd773SHong Zhang PetscFunctionReturn(0); 376474dd773SHong Zhang } 377474dd773SHong Zhang 378*63a6f1b4SHong Zhang static PetscErrorCode TSStepRefine_MRKSPLIT(TS ts) 379*63a6f1b4SHong Zhang { 380*63a6f1b4SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 381*63a6f1b4SHong Zhang TS subts_fast = rk->subts_fast,currentlevelts; 382*63a6f1b4SHong Zhang TS_RK *subrk_fast = (TS_RK*)subts_fast->data; 383*63a6f1b4SHong Zhang RKTableau tab = rk->tableau; 384*63a6f1b4SHong Zhang Vec *Y = rk->Y; 385*63a6f1b4SHong Zhang Vec *YdotRHS = rk->YdotRHS,*YdotRHS_fast = rk->YdotRHS_fast; 386*63a6f1b4SHong Zhang Vec Yfast,Xfast; 387*63a6f1b4SHong Zhang const PetscInt s = tab->s; 388*63a6f1b4SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 389*63a6f1b4SHong Zhang PetscScalar *w = rk->work; 390*63a6f1b4SHong Zhang PetscInt i,j,k; 391*63a6f1b4SHong Zhang PetscReal t = ts->ptime,h = ts->time_step; 392*63a6f1b4SHong Zhang PetscErrorCode ierr; 393*63a6f1b4SHong Zhang 394*63a6f1b4SHong Zhang for (k=0; k<rk->dtratio; k++) { 395*63a6f1b4SHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 396*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 397*63a6f1b4SHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 398*63a6f1b4SHong Zhang } 399*63a6f1b4SHong Zhang /* propogate fast component using small time steps */ 400*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 401*63a6f1b4SHong Zhang /* stage value for slow components */ 402*63a6f1b4SHong Zhang ierr = TSInterpolate_MRKSPLIT(rk->ts_root,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr); 403*63a6f1b4SHong Zhang currentlevelts = rk->ts_root; 404*63a6f1b4SHong Zhang while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */ 405*63a6f1b4SHong Zhang currentlevelts = ((TS_RK*)currentlevelts->data)->subts_fast; 406*63a6f1b4SHong Zhang ierr = TSInterpolate_MRKSPLIT(currentlevelts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr); 407*63a6f1b4SHong Zhang } 408*63a6f1b4SHong Zhang for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 409*63a6f1b4SHong Zhang subrk_fast->stage_time = t + h/rk->dtratio*c[i]; 410*63a6f1b4SHong Zhang ierr = TSPreStage(subts_fast,subrk_fast->stage_time);CHKERRQ(ierr); 411*63a6f1b4SHong Zhang /* stage value for fast components */ 412*63a6f1b4SHong Zhang ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 413*63a6f1b4SHong Zhang ierr = VecCopy(Xfast,Yfast);CHKERRQ(ierr); 414*63a6f1b4SHong Zhang ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 415*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 416*63a6f1b4SHong Zhang ierr = TSPostStage(subts_fast,subrk_fast->stage_time,i,Y);CHKERRQ(ierr); 417*63a6f1b4SHong Zhang /* compute the stage RHS for fast components */ 418*63a6f1b4SHong Zhang ierr = TSComputeRHSFunction(subts_fast,t+k*h*rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 419*63a6f1b4SHong Zhang } 420*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 421*63a6f1b4SHong Zhang /* update the value of fast components using fast time step */ 422*63a6f1b4SHong Zhang rk->slow = PETSC_FALSE; 423*63a6f1b4SHong Zhang ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 424*63a6f1b4SHong Zhang for (i=0; i<s; i++) { 425*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 426*63a6f1b4SHong Zhang } 427*63a6f1b4SHong Zhang 428*63a6f1b4SHong Zhang if (subrk_fast->subts_fast) { 429*63a6f1b4SHong Zhang subts_fast->ptime = t+k*h/rk->dtratio; 430*63a6f1b4SHong Zhang subts_fast->time_step = h/rk->dtratio; 431*63a6f1b4SHong Zhang ierr = TSStepRefine_MRKSPLIT(subts_fast);CHKERRQ(ierr); 432*63a6f1b4SHong Zhang } 433*63a6f1b4SHong Zhang /* update the fast components of the solution */ 434*63a6f1b4SHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 435*63a6f1b4SHong Zhang ierr = VecISCopy(rk->X0,rk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr); 436*63a6f1b4SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 437*63a6f1b4SHong Zhang } 438*63a6f1b4SHong Zhang PetscFunctionReturn(0); 439*63a6f1b4SHong Zhang } 440*63a6f1b4SHong Zhang 441474dd773SHong Zhang static PetscErrorCode TSStep_MRKSPLIT(TS ts) 442474dd773SHong Zhang { 443474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 444474dd773SHong Zhang RKTableau tab = rk->tableau; 445474dd773SHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 446*63a6f1b4SHong Zhang Vec *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow; 447474dd773SHong Zhang Vec Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively */ 448474dd773SHong Zhang const PetscInt s = tab->s; 449474dd773SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 450474dd773SHong Zhang PetscScalar *w = rk->work; 451*63a6f1b4SHong Zhang PetscInt i,j; 452474dd773SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 453474dd773SHong Zhang PetscErrorCode ierr; 454474dd773SHong Zhang 455474dd773SHong Zhang PetscFunctionBegin; 456474dd773SHong Zhang rk->status = TS_STEP_INCOMPLETE; 457474dd773SHong Zhang for (i=0; i<s; i++) { 458474dd773SHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr); 459474dd773SHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 460474dd773SHong Zhang } 461474dd773SHong Zhang ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr); 462474dd773SHong Zhang /* propogate both slow and fast components using large time steps */ 463474dd773SHong Zhang for (i=0; i<s; i++) { 464474dd773SHong Zhang rk->stage_time = t + h*c[i]; 465474dd773SHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 466474dd773SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 467474dd773SHong Zhang ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 468bb42530cSHong Zhang ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 469474dd773SHong Zhang for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 470474dd773SHong Zhang ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 471bb42530cSHong Zhang ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr); 472474dd773SHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 473bb42530cSHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 474474dd773SHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 475474dd773SHong Zhang ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 476474dd773SHong Zhang ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 477474dd773SHong Zhang } 478474dd773SHong Zhang rk->slow = PETSC_TRUE; 479*63a6f1b4SHong Zhang /* update the slow components of the solution using slow time step */ 480474dd773SHong Zhang ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 481*63a6f1b4SHong Zhang /* YdotRHS will be used for interpolation during refinement */ 482474dd773SHong Zhang for (i=0; i<s; i++) { 483474dd773SHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr); 484474dd773SHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 485474dd773SHong Zhang } 486*63a6f1b4SHong Zhang 487*63a6f1b4SHong Zhang ierr = TSStepRefine_MRKSPLIT(ts);CHKERRQ(ierr); 488*63a6f1b4SHong Zhang 489bb42530cSHong Zhang ts->ptime = t + ts->time_step; 490474dd773SHong Zhang ts->time_step = next_time_step; 491474dd773SHong Zhang rk->status = TS_STEP_COMPLETE; 492474dd773SHong Zhang PetscFunctionReturn(0); 493474dd773SHong Zhang } 494474dd773SHong Zhang 495474dd773SHong Zhang /*@C 496474dd773SHong Zhang TSRKSetMultirateType - Set the type of RK Multirate scheme 497474dd773SHong Zhang 498474dd773SHong Zhang Logically collective 499474dd773SHong Zhang 500474dd773SHong Zhang Input Parameter: 501474dd773SHong Zhang + ts - timestepping context 502474dd773SHong Zhang - mrktype - type of MRK-scheme 503474dd773SHong Zhang 504474dd773SHong Zhang Options Database: 505474dd773SHong Zhang . -ts_rk_multiarte_type - <none,nonsplit,split> 506474dd773SHong Zhang 507474dd773SHong Zhang Level: intermediate 508474dd773SHong Zhang @*/ 509474dd773SHong Zhang PetscErrorCode TSRKSetMultirateType(TS ts, TSMRKType mrktype) 510474dd773SHong Zhang { 511474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 512474dd773SHong Zhang PetscErrorCode ierr; 513474dd773SHong Zhang 514474dd773SHong Zhang PetscFunctionBegin; 515474dd773SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 516474dd773SHong Zhang switch(mrktype){ 517474dd773SHong Zhang case TSMRKNONE: 518474dd773SHong Zhang break; 519474dd773SHong Zhang case TSMRKNONSPLIT: 520474dd773SHong Zhang ts->ops->step = TSStep_MRKNONSPLIT; 521474dd773SHong Zhang ts->ops->interpolate = TSInterpolate_MRKNONSPLIT; 522474dd773SHong Zhang rk->dtratio = 2; 523474dd773SHong Zhang ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr); 524474dd773SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",TSSetUp_MRKNONSPLIT);CHKERRQ(ierr); 525474dd773SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",TSReset_MRKNONSPLIT);CHKERRQ(ierr); 526474dd773SHong Zhang break; 527474dd773SHong Zhang case TSMRKSPLIT: 528474dd773SHong Zhang ts->ops->step = TSStep_MRKSPLIT; 529474dd773SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MRKSPLIT; 530474dd773SHong Zhang ts->ops->interpolate = TSInterpolate_MRKSPLIT; 531474dd773SHong Zhang rk->dtratio = 2; 532474dd773SHong Zhang ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr); 533474dd773SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",TSSetUp_MRKSPLIT);CHKERRQ(ierr); 534474dd773SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",TSReset_MRKSPLIT);CHKERRQ(ierr); 535474dd773SHong Zhang break; 536474dd773SHong Zhang default : 537474dd773SHong Zhang SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mrktype); 538474dd773SHong Zhang } 539474dd773SHong Zhang PetscFunctionReturn(0); 540474dd773SHong Zhang } 541474dd773SHong Zhang 542