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> 1721052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h> 18474dd773SHong Zhang 1902555c94SHong Zhang static PetscErrorCode TSReset_RK_MultirateNonsplit(TS ts) 20e5bffa22SHong Zhang { 21e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 22e5bffa22SHong Zhang RKTableau tab = rk->tableau; 23e5bffa22SHong Zhang PetscErrorCode ierr; 24e5bffa22SHong Zhang 25e5bffa22SHong Zhang PetscFunctionBegin; 26e5bffa22SHong Zhang ierr = VecDestroy(&rk->X0);CHKERRQ(ierr); 27e5bffa22SHong Zhang ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 28e5bffa22SHong Zhang PetscFunctionReturn(0); 29e5bffa22SHong Zhang } 30e5bffa22SHong Zhang 3102555c94SHong Zhang static PetscErrorCode TSInterpolate_RK_MultirateNonsplit(TS ts,PetscReal itime,Vec X) 32e5bffa22SHong Zhang { 33e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 34e5bffa22SHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 3563a6f1b4SHong Zhang PetscReal h = ts->time_step; 36e5bffa22SHong Zhang PetscReal tt,t; 37e5bffa22SHong Zhang PetscScalar *b; 38e5bffa22SHong Zhang const PetscReal *B = rk->tableau->binterp; 39e5bffa22SHong Zhang PetscErrorCode ierr; 40e5bffa22SHong Zhang 41e5bffa22SHong Zhang PetscFunctionBegin; 42e5bffa22SHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 4363a6f1b4SHong Zhang t = (itime - rk->ptime)/h; 44e5bffa22SHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 45e5bffa22SHong Zhang for (i=0; i<s; i++) b[i] = 0; 46e5bffa22SHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 47e5bffa22SHong Zhang for (i=0; i<s; i++) { 48e5bffa22SHong Zhang b[i] += h * B[i*p+j] * tt; 49e5bffa22SHong Zhang } 50e5bffa22SHong Zhang } 51e5bffa22SHong Zhang ierr = VecCopy(rk->X0,X);CHKERRQ(ierr); 52e5bffa22SHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 53e5bffa22SHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 54e5bffa22SHong Zhang PetscFunctionReturn(0); 55e5bffa22SHong Zhang } 56e5bffa22SHong Zhang 5702555c94SHong Zhang static PetscErrorCode TSStepRefine_RK_MultirateNonsplit(TS ts) 5863a6f1b4SHong Zhang { 5963a6f1b4SHong Zhang TS previousts,subts; 6063a6f1b4SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 6163a6f1b4SHong Zhang RKTableau tab = rk->tableau; 6263a6f1b4SHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 6363a6f1b4SHong Zhang Vec vec_fast,subvec_fast; 6463a6f1b4SHong Zhang const PetscInt s = tab->s; 6563a6f1b4SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 6663a6f1b4SHong Zhang PetscScalar *w = rk->work; 6763a6f1b4SHong Zhang PetscInt i,j,k; 6863a6f1b4SHong Zhang PetscReal t = ts->ptime,h = ts->time_step; 6963a6f1b4SHong Zhang PetscErrorCode ierr; 7063a6f1b4SHong Zhang 71*362febeeSStefano Zampini PetscFunctionBegin; 7263a6f1b4SHong Zhang ierr = VecDuplicate(ts->vec_sol,&vec_fast);CHKERRQ(ierr); 7363a6f1b4SHong Zhang previousts = rk->subts_current; 7463a6f1b4SHong Zhang ierr = TSRHSSplitGetSubTS(rk->subts_current,"fast",&subts);CHKERRQ(ierr); 7563a6f1b4SHong Zhang ierr = TSRHSSplitGetSubTS(subts,"fast",&subts);CHKERRQ(ierr); 7663a6f1b4SHong Zhang for (k=0; k<rk->dtratio; k++) { 7763a6f1b4SHong Zhang for (i=0; i<s; i++) { 7802555c94SHong Zhang ierr = TSInterpolate_RK_MultirateNonsplit(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr); 7963a6f1b4SHong Zhang for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 8063a6f1b4SHong 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 */ 8163a6f1b4SHong Zhang ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr); 8263a6f1b4SHong Zhang ierr = VecMAXPY(vec_fast,i,w,YdotRHS);CHKERRQ(ierr); 8363a6f1b4SHong Zhang /* update the fast components in the stage value */ 8463a6f1b4SHong Zhang ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 8563a6f1b4SHong Zhang ierr = VecISCopy(Y[i],rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr); 8663a6f1b4SHong Zhang ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 8763a6f1b4SHong Zhang /* compute the stage RHS */ 8863a6f1b4SHong Zhang ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 8963a6f1b4SHong Zhang } 9063a6f1b4SHong Zhang ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr); 9163a6f1b4SHong Zhang ierr = TSEvaluateStep(ts,tab->order,vec_fast,NULL);CHKERRQ(ierr); 9263a6f1b4SHong Zhang /* update the fast components in the solution */ 9363a6f1b4SHong Zhang ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 9463a6f1b4SHong Zhang ierr = VecISCopy(ts->vec_sol,rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr); 9563a6f1b4SHong Zhang ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr); 9663a6f1b4SHong Zhang 9763a6f1b4SHong Zhang if (subts) { 9863a6f1b4SHong Zhang Vec *YdotRHS_copy; 9963a6f1b4SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,s,&YdotRHS_copy);CHKERRQ(ierr); 10063a6f1b4SHong Zhang rk->subts_current = rk->subts_fast; 10163a6f1b4SHong Zhang ts->ptime = t+k*h/rk->dtratio; 10263a6f1b4SHong Zhang ts->time_step = h/rk->dtratio; 10363a6f1b4SHong Zhang ierr = TSRHSSplitGetIS(rk->subts_current,"fast",&rk->is_fast);CHKERRQ(ierr); 10463a6f1b4SHong Zhang for (i=0; i<s; i++) { 10563a6f1b4SHong Zhang ierr = VecCopy(rk->YdotRHS_slow[i],YdotRHS_copy[i]);CHKERRQ(ierr); 10663a6f1b4SHong Zhang ierr = VecCopy(YdotRHS[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr); 10763a6f1b4SHong Zhang } 10863a6f1b4SHong Zhang 10902555c94SHong Zhang ierr = TSStepRefine_RK_MultirateNonsplit(ts);CHKERRQ(ierr); 11063a6f1b4SHong Zhang 11163a6f1b4SHong Zhang rk->subts_current = previousts; 11263a6f1b4SHong Zhang ts->ptime = t; 11363a6f1b4SHong Zhang ts->time_step = h; 11463a6f1b4SHong Zhang ierr = TSRHSSplitGetIS(previousts,"fast",&rk->is_fast);CHKERRQ(ierr); 11563a6f1b4SHong Zhang for (i=0; i<s; i++) { 11663a6f1b4SHong Zhang ierr = VecCopy(YdotRHS_copy[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr); 11763a6f1b4SHong Zhang } 11863a6f1b4SHong Zhang ierr = VecDestroyVecs(s,&YdotRHS_copy);CHKERRQ(ierr); 11963a6f1b4SHong Zhang } 12063a6f1b4SHong Zhang } 12163a6f1b4SHong Zhang ierr = VecDestroy(&vec_fast);CHKERRQ(ierr); 12263a6f1b4SHong Zhang PetscFunctionReturn(0); 12363a6f1b4SHong Zhang } 12463a6f1b4SHong Zhang 12502555c94SHong Zhang static PetscErrorCode TSStep_RK_MultirateNonsplit(TS ts) 126e5bffa22SHong Zhang { 127e5bffa22SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 128e5bffa22SHong Zhang RKTableau tab = rk->tableau; 129e5bffa22SHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow; 130e5bffa22SHong Zhang Vec stage_slow,sol_slow; /* vectors store the slow components */ 131e5bffa22SHong Zhang Vec subvec_slow; /* sub vector to store the slow components */ 132e5bffa22SHong Zhang IS is_slow = rk->is_slow; 133e5bffa22SHong Zhang const PetscInt s = tab->s; 134e5bffa22SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 135e5bffa22SHong Zhang PetscScalar *w = rk->work; 13663a6f1b4SHong Zhang PetscInt i,j,dtratio = rk->dtratio; 137e5bffa22SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 138e5bffa22SHong Zhang PetscErrorCode ierr; 139e5bffa22SHong Zhang 140e5bffa22SHong Zhang PetscFunctionBegin; 141e5bffa22SHong Zhang rk->status = TS_STEP_INCOMPLETE; 142e5bffa22SHong Zhang ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr); 143e5bffa22SHong Zhang ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr); 144e5bffa22SHong Zhang ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr); 145e5bffa22SHong Zhang for (i=0; i<s; i++) { 146e5bffa22SHong Zhang rk->stage_time = t + h*c[i]; 147e5bffa22SHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 148e5bffa22SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 149e5bffa22SHong Zhang for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 150e5bffa22SHong Zhang ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr); 151e5bffa22SHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr); 152e5bffa22SHong Zhang /* compute the stage RHS */ 153e5bffa22SHong Zhang ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 154e5bffa22SHong Zhang } 155e5bffa22SHong Zhang /* update the slow components in the solution */ 156e5bffa22SHong Zhang rk->YdotRHS = YdotRHS_slow; 157e5bffa22SHong Zhang rk->dtratio = 1; 158e5bffa22SHong Zhang ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr); 159e5bffa22SHong Zhang rk->dtratio = dtratio; 160e5bffa22SHong Zhang rk->YdotRHS = YdotRHS; 161e5bffa22SHong Zhang /* update the slow components in the solution */ 162e5bffa22SHong Zhang ierr = VecGetSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr); 163e5bffa22SHong Zhang ierr = VecISCopy(ts->vec_sol,is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr); 164e5bffa22SHong Zhang ierr = VecRestoreSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr); 165e5bffa22SHong Zhang 16663a6f1b4SHong Zhang rk->subts_current = ts; 16763a6f1b4SHong Zhang rk->ptime = t; 16863a6f1b4SHong Zhang rk->time_step = h; 16902555c94SHong Zhang ierr = TSStepRefine_RK_MultirateNonsplit(ts);CHKERRQ(ierr); 17063a6f1b4SHong Zhang 171e5bffa22SHong Zhang ts->ptime = t + ts->time_step; 172e5bffa22SHong Zhang ts->time_step = next_time_step; 173e5bffa22SHong Zhang rk->status = TS_STEP_COMPLETE; 17463a6f1b4SHong Zhang 175e5bffa22SHong Zhang /* free memory of work vectors */ 176e5bffa22SHong Zhang ierr = VecDestroy(&stage_slow);CHKERRQ(ierr); 177e5bffa22SHong Zhang ierr = VecDestroy(&sol_slow);CHKERRQ(ierr); 178e5bffa22SHong Zhang PetscFunctionReturn(0); 179e5bffa22SHong Zhang } 180e5bffa22SHong Zhang 18102555c94SHong Zhang static PetscErrorCode TSSetUp_RK_MultirateNonsplit(TS ts) 1820fe4d17eSHong Zhang { 1830fe4d17eSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1840fe4d17eSHong Zhang RKTableau tab = rk->tableau; 1850fe4d17eSHong Zhang PetscErrorCode ierr; 1860fe4d17eSHong Zhang 1870fe4d17eSHong Zhang PetscFunctionBegin; 1880fe4d17eSHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr); 1890fe4d17eSHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr); 1900fe4d17eSHong 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"); 1910fe4d17eSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr); 1920fe4d17eSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr); 1930fe4d17eSHong 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"); 1940fe4d17eSHong Zhang ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr); 1950fe4d17eSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 1960fe4d17eSHong Zhang rk->subts_current = rk->subts_fast; 1970fe4d17eSHong Zhang 19802555c94SHong Zhang ts->ops->step = TSStep_RK_MultirateNonsplit; 19902555c94SHong Zhang ts->ops->interpolate = TSInterpolate_RK_MultirateNonsplit; 2000fe4d17eSHong Zhang PetscFunctionReturn(0); 2010fe4d17eSHong Zhang } 2020fe4d17eSHong Zhang 20363a6f1b4SHong Zhang /* 20463a6f1b4SHong Zhang Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest. 20563a6f1b4SHong Zhang */ 20663a6f1b4SHong Zhang static PetscErrorCode TSCopyDM(TS tssrc,TS tsdest) 20763a6f1b4SHong Zhang { 20863a6f1b4SHong Zhang DM newdm,dmsrc,dmdest; 20963a6f1b4SHong Zhang PetscErrorCode ierr; 21063a6f1b4SHong Zhang 21163a6f1b4SHong Zhang PetscFunctionBegin; 21263a6f1b4SHong Zhang ierr = TSGetDM(tssrc,&dmsrc);CHKERRQ(ierr); 21363a6f1b4SHong Zhang ierr = DMClone(dmsrc,&newdm);CHKERRQ(ierr); 21463a6f1b4SHong Zhang ierr = TSGetDM(tsdest,&dmdest);CHKERRQ(ierr); 21563a6f1b4SHong Zhang ierr = DMCopyDMTS(dmdest,newdm);CHKERRQ(ierr); 21663a6f1b4SHong Zhang ierr = DMCopyDMSNES(dmdest,newdm);CHKERRQ(ierr); 21763a6f1b4SHong Zhang ierr = TSSetDM(tsdest,newdm);CHKERRQ(ierr); 21863a6f1b4SHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 21963a6f1b4SHong Zhang PetscFunctionReturn(0); 22063a6f1b4SHong Zhang } 22163a6f1b4SHong Zhang 22202555c94SHong Zhang static PetscErrorCode TSReset_RK_MultirateSplit(TS ts) 223474dd773SHong Zhang { 224474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 225474dd773SHong Zhang PetscErrorCode ierr; 226474dd773SHong Zhang 227474dd773SHong Zhang PetscFunctionBegin; 228bb42530cSHong Zhang if (rk->subts_slow) { 229bb42530cSHong Zhang ierr = PetscFree(rk->subts_slow->data);CHKERRQ(ierr); 230bb42530cSHong Zhang rk->subts_slow = NULL; 231bb42530cSHong Zhang } 232bb42530cSHong Zhang if (rk->subts_fast) { 23363a6f1b4SHong Zhang ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr); 23463a6f1b4SHong Zhang ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr); 235bb42530cSHong Zhang ierr = VecDestroy(&rk->X0);CHKERRQ(ierr); 236bf0cca7dSHong Zhang ierr = TSReset_RK_MultirateSplit(rk->subts_fast);CHKERRQ(ierr); 237bb42530cSHong Zhang ierr = PetscFree(rk->subts_fast->data);CHKERRQ(ierr); 238bb42530cSHong Zhang rk->subts_fast = NULL; 239bb42530cSHong Zhang } 240474dd773SHong Zhang PetscFunctionReturn(0); 241474dd773SHong Zhang } 242474dd773SHong Zhang 24302555c94SHong Zhang static PetscErrorCode TSInterpolate_RK_MultirateSplit(TS ts,PetscReal itime,Vec X) 244474dd773SHong Zhang { 245474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 24663a6f1b4SHong Zhang Vec Xslow; 247474dd773SHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 248474dd773SHong Zhang PetscReal h; 249474dd773SHong Zhang PetscReal tt,t; 250474dd773SHong Zhang PetscScalar *b; 251474dd773SHong Zhang const PetscReal *B = rk->tableau->binterp; 252474dd773SHong Zhang PetscErrorCode ierr; 253474dd773SHong Zhang 254474dd773SHong Zhang PetscFunctionBegin; 255474dd773SHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 256474dd773SHong Zhang 257474dd773SHong Zhang switch (rk->status) { 258474dd773SHong Zhang case TS_STEP_INCOMPLETE: 259474dd773SHong Zhang case TS_STEP_PENDING: 260474dd773SHong Zhang h = ts->time_step; 261474dd773SHong Zhang t = (itime - ts->ptime)/h; 262474dd773SHong Zhang break; 263474dd773SHong Zhang case TS_STEP_COMPLETE: 264474dd773SHong Zhang h = ts->ptime - ts->ptime_prev; 265474dd773SHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 266474dd773SHong Zhang break; 267474dd773SHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 268474dd773SHong Zhang } 269474dd773SHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 270474dd773SHong Zhang for (i=0; i<s; i++) b[i] = 0; 271474dd773SHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 272474dd773SHong Zhang for (i=0; i<s; i++) { 273474dd773SHong Zhang b[i] += h * B[i*p+j] * tt; 274474dd773SHong Zhang } 275474dd773SHong Zhang } 27663a6f1b4SHong Zhang for (i=0; i<s; i++) { 27763a6f1b4SHong Zhang ierr = VecGetSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr); 27863a6f1b4SHong Zhang } 27963a6f1b4SHong Zhang ierr = VecGetSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr); 28063a6f1b4SHong Zhang ierr = VecISCopy(rk->X0,rk->is_slow,SCATTER_REVERSE,Xslow);CHKERRQ(ierr); 28163a6f1b4SHong Zhang ierr = VecMAXPY(Xslow,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 28263a6f1b4SHong Zhang ierr = VecRestoreSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr); 28363a6f1b4SHong Zhang for (i=0; i<s; i++) { 28463a6f1b4SHong Zhang ierr = VecRestoreSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr); 28563a6f1b4SHong Zhang } 286474dd773SHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 287474dd773SHong Zhang PetscFunctionReturn(0); 288474dd773SHong Zhang } 289474dd773SHong Zhang 290474dd773SHong Zhang /* 291474dd773SHong Zhang This is for partitioned RHS multirate RK method 292474dd773SHong Zhang The step completion formula is 293474dd773SHong Zhang 294474dd773SHong Zhang x1 = x0 + h b^T YdotRHS 295474dd773SHong Zhang 296474dd773SHong Zhang */ 29702555c94SHong Zhang static PetscErrorCode TSEvaluateStep_RK_MultirateSplit(TS ts,PetscInt order,Vec X,PetscBool *done) 298474dd773SHong Zhang { 299474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 300474dd773SHong Zhang RKTableau tab = rk->tableau; 301474dd773SHong Zhang Vec Xslow,Xfast; /* subvectors of X which store slow components and fast components respectively */ 302474dd773SHong Zhang PetscScalar *w = rk->work; 303474dd773SHong Zhang PetscReal h = ts->time_step; 304474dd773SHong Zhang PetscInt s = tab->s,j; 305474dd773SHong Zhang PetscErrorCode ierr; 306474dd773SHong Zhang 307474dd773SHong Zhang PetscFunctionBegin; 308474dd773SHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 309474dd773SHong Zhang if (rk->slow) { 310474dd773SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]; 311474dd773SHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr); 31263a6f1b4SHong Zhang ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr); 313feb237baSPierre Jolivet ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr); 314474dd773SHong Zhang } else { 315474dd773SHong Zhang for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j]; 316474dd773SHong Zhang ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr); 31763a6f1b4SHong Zhang ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr); 318474dd773SHong Zhang ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr); 319474dd773SHong Zhang } 320474dd773SHong Zhang PetscFunctionReturn(0); 321474dd773SHong Zhang } 322474dd773SHong Zhang 32302555c94SHong Zhang static PetscErrorCode TSStepRefine_RK_MultirateSplit(TS ts) 32463a6f1b4SHong Zhang { 32563a6f1b4SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 32663a6f1b4SHong Zhang TS subts_fast = rk->subts_fast,currentlevelts; 32763a6f1b4SHong Zhang TS_RK *subrk_fast = (TS_RK*)subts_fast->data; 32863a6f1b4SHong Zhang RKTableau tab = rk->tableau; 32963a6f1b4SHong Zhang Vec *Y = rk->Y; 33063a6f1b4SHong Zhang Vec *YdotRHS = rk->YdotRHS,*YdotRHS_fast = rk->YdotRHS_fast; 33163a6f1b4SHong Zhang Vec Yfast,Xfast; 33263a6f1b4SHong Zhang const PetscInt s = tab->s; 33363a6f1b4SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 33463a6f1b4SHong Zhang PetscScalar *w = rk->work; 33563a6f1b4SHong Zhang PetscInt i,j,k; 33663a6f1b4SHong Zhang PetscReal t = ts->ptime,h = ts->time_step; 33763a6f1b4SHong Zhang PetscErrorCode ierr; 33863a6f1b4SHong Zhang 339*362febeeSStefano Zampini PetscFunctionBegin; 34063a6f1b4SHong Zhang for (k=0; k<rk->dtratio; k++) { 34163a6f1b4SHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 34263a6f1b4SHong Zhang for (i=0; i<s; i++) { 34363a6f1b4SHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 34463a6f1b4SHong Zhang } 34563a6f1b4SHong Zhang /* propogate fast component using small time steps */ 34663a6f1b4SHong Zhang for (i=0; i<s; i++) { 34763a6f1b4SHong Zhang /* stage value for slow components */ 34802555c94SHong Zhang ierr = TSInterpolate_RK_MultirateSplit(rk->ts_root,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr); 34963a6f1b4SHong Zhang currentlevelts = rk->ts_root; 35063a6f1b4SHong Zhang while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */ 35163a6f1b4SHong Zhang currentlevelts = ((TS_RK*)currentlevelts->data)->subts_fast; 35202555c94SHong Zhang ierr = TSInterpolate_RK_MultirateSplit(currentlevelts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr); 35363a6f1b4SHong Zhang } 35463a6f1b4SHong Zhang for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 35563a6f1b4SHong Zhang subrk_fast->stage_time = t + h/rk->dtratio*c[i]; 35663a6f1b4SHong Zhang ierr = TSPreStage(subts_fast,subrk_fast->stage_time);CHKERRQ(ierr); 35763a6f1b4SHong Zhang /* stage value for fast components */ 35863a6f1b4SHong Zhang ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 35963a6f1b4SHong Zhang ierr = VecCopy(Xfast,Yfast);CHKERRQ(ierr); 36063a6f1b4SHong Zhang ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 36163a6f1b4SHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 36263a6f1b4SHong Zhang ierr = TSPostStage(subts_fast,subrk_fast->stage_time,i,Y);CHKERRQ(ierr); 36363a6f1b4SHong Zhang /* compute the stage RHS for fast components */ 36463a6f1b4SHong Zhang ierr = TSComputeRHSFunction(subts_fast,t+k*h*rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 36563a6f1b4SHong Zhang } 36663a6f1b4SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 36763a6f1b4SHong Zhang /* update the value of fast components using fast time step */ 36863a6f1b4SHong Zhang rk->slow = PETSC_FALSE; 36902555c94SHong Zhang ierr = TSEvaluateStep_RK_MultirateSplit(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 37063a6f1b4SHong Zhang for (i=0; i<s; i++) { 37163a6f1b4SHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 37263a6f1b4SHong Zhang } 37363a6f1b4SHong Zhang 37463a6f1b4SHong Zhang if (subrk_fast->subts_fast) { 37563a6f1b4SHong Zhang subts_fast->ptime = t+k*h/rk->dtratio; 37663a6f1b4SHong Zhang subts_fast->time_step = h/rk->dtratio; 37702555c94SHong Zhang ierr = TSStepRefine_RK_MultirateSplit(subts_fast);CHKERRQ(ierr); 37863a6f1b4SHong Zhang } 37963a6f1b4SHong Zhang /* update the fast components of the solution */ 38063a6f1b4SHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 38163a6f1b4SHong Zhang ierr = VecISCopy(rk->X0,rk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr); 38263a6f1b4SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr); 38363a6f1b4SHong Zhang } 38463a6f1b4SHong Zhang PetscFunctionReturn(0); 38563a6f1b4SHong Zhang } 38663a6f1b4SHong Zhang 38702555c94SHong Zhang static PetscErrorCode TSStep_RK_MultirateSplit(TS ts) 388474dd773SHong Zhang { 389474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 390474dd773SHong Zhang RKTableau tab = rk->tableau; 391474dd773SHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 39263a6f1b4SHong Zhang Vec *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow; 393474dd773SHong Zhang Vec Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively */ 394474dd773SHong Zhang const PetscInt s = tab->s; 395474dd773SHong Zhang const PetscReal *A = tab->A,*c = tab->c; 396474dd773SHong Zhang PetscScalar *w = rk->work; 39763a6f1b4SHong Zhang PetscInt i,j; 398474dd773SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 399474dd773SHong Zhang PetscErrorCode ierr; 400474dd773SHong Zhang 401474dd773SHong Zhang PetscFunctionBegin; 402474dd773SHong Zhang rk->status = TS_STEP_INCOMPLETE; 403474dd773SHong Zhang for (i=0; i<s; i++) { 404474dd773SHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr); 405474dd773SHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 406474dd773SHong Zhang } 407474dd773SHong Zhang ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr); 408474dd773SHong Zhang /* propogate both slow and fast components using large time steps */ 409474dd773SHong Zhang for (i=0; i<s; i++) { 410474dd773SHong Zhang rk->stage_time = t + h*c[i]; 411474dd773SHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 412474dd773SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 413474dd773SHong Zhang ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 414bb42530cSHong Zhang ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 415474dd773SHong Zhang for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 416474dd773SHong Zhang ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 417bb42530cSHong Zhang ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr); 418474dd773SHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 419bb42530cSHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 420474dd773SHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr); 421474dd773SHong Zhang ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 422474dd773SHong Zhang ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 423474dd773SHong Zhang } 424474dd773SHong Zhang rk->slow = PETSC_TRUE; 42563a6f1b4SHong Zhang /* update the slow components of the solution using slow time step */ 42602555c94SHong Zhang ierr = TSEvaluateStep_RK_MultirateSplit(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 42763a6f1b4SHong Zhang /* YdotRHS will be used for interpolation during refinement */ 428474dd773SHong Zhang for (i=0; i<s; i++) { 429474dd773SHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr); 430474dd773SHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 431474dd773SHong Zhang } 43263a6f1b4SHong Zhang 43302555c94SHong Zhang ierr = TSStepRefine_RK_MultirateSplit(ts);CHKERRQ(ierr); 43463a6f1b4SHong Zhang 435bb42530cSHong Zhang ts->ptime = t + ts->time_step; 436474dd773SHong Zhang ts->time_step = next_time_step; 437474dd773SHong Zhang rk->status = TS_STEP_COMPLETE; 438474dd773SHong Zhang PetscFunctionReturn(0); 439474dd773SHong Zhang } 440474dd773SHong Zhang 44102555c94SHong Zhang static PetscErrorCode TSSetUp_RK_MultirateSplit(TS ts) 4420fe4d17eSHong Zhang { 4430fe4d17eSHong Zhang TS_RK *rk = (TS_RK*)ts->data,*nextlevelrk,*currentlevelrk; 4440fe4d17eSHong Zhang TS nextlevelts; 4450fe4d17eSHong Zhang Vec X0; 4460fe4d17eSHong Zhang PetscErrorCode ierr; 4470fe4d17eSHong Zhang 4480fe4d17eSHong Zhang PetscFunctionBegin; 4490fe4d17eSHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr); 4500fe4d17eSHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr); 4510fe4d17eSHong 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"); 4520fe4d17eSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr); 4530fe4d17eSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr); 4540fe4d17eSHong 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"); 4550fe4d17eSHong Zhang 4560fe4d17eSHong Zhang ierr = VecDuplicate(ts->vec_sol,&X0);CHKERRQ(ierr); 4570fe4d17eSHong Zhang /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */ 4580fe4d17eSHong Zhang currentlevelrk = rk; 4590fe4d17eSHong Zhang while (currentlevelrk->subts_fast) { 4600fe4d17eSHong Zhang ierr = PetscMalloc1(rk->tableau->s,¤tlevelrk->YdotRHS_fast);CHKERRQ(ierr); 4610fe4d17eSHong Zhang ierr = PetscMalloc1(rk->tableau->s,¤tlevelrk->YdotRHS_slow);CHKERRQ(ierr); 4620fe4d17eSHong Zhang ierr = PetscObjectReference((PetscObject)X0);CHKERRQ(ierr); 4630fe4d17eSHong Zhang currentlevelrk->X0 = X0; 4640fe4d17eSHong Zhang currentlevelrk->ts_root = ts; 4650fe4d17eSHong Zhang 4660fe4d17eSHong Zhang /* set up the ts for the slow part */ 4670fe4d17eSHong Zhang nextlevelts = currentlevelrk->subts_slow; 4680fe4d17eSHong Zhang ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr); 4690fe4d17eSHong Zhang nextlevelrk->tableau = rk->tableau; 4700fe4d17eSHong Zhang nextlevelrk->work = rk->work; 4710fe4d17eSHong Zhang nextlevelrk->Y = rk->Y; 4720fe4d17eSHong Zhang nextlevelrk->YdotRHS = rk->YdotRHS; 4730fe4d17eSHong Zhang nextlevelts->data = (void*)nextlevelrk; 4740fe4d17eSHong Zhang ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr); 4750fe4d17eSHong Zhang ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr); 4760fe4d17eSHong Zhang 4770fe4d17eSHong Zhang /* set up the ts for the fast part */ 4780fe4d17eSHong Zhang nextlevelts = currentlevelrk->subts_fast; 4790fe4d17eSHong Zhang ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr); 4800fe4d17eSHong Zhang nextlevelrk->tableau = rk->tableau; 4810fe4d17eSHong Zhang nextlevelrk->work = rk->work; 4820fe4d17eSHong Zhang nextlevelrk->Y = rk->Y; 4830fe4d17eSHong Zhang nextlevelrk->YdotRHS = rk->YdotRHS; 4840fe4d17eSHong Zhang nextlevelrk->dtratio = rk->dtratio; 4850fe4d17eSHong Zhang ierr = TSRHSSplitGetIS(nextlevelts,"slow",&nextlevelrk->is_slow);CHKERRQ(ierr); 4860fe4d17eSHong Zhang ierr = TSRHSSplitGetSubTS(nextlevelts,"slow",&nextlevelrk->subts_slow);CHKERRQ(ierr); 4870fe4d17eSHong Zhang ierr = TSRHSSplitGetIS(nextlevelts,"fast",&nextlevelrk->is_fast);CHKERRQ(ierr); 4880fe4d17eSHong Zhang ierr = TSRHSSplitGetSubTS(nextlevelts,"fast",&nextlevelrk->subts_fast);CHKERRQ(ierr); 4890fe4d17eSHong Zhang nextlevelts->data = (void*)nextlevelrk; 4900fe4d17eSHong Zhang ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr); 4910fe4d17eSHong Zhang ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr); 4920fe4d17eSHong Zhang 4930fe4d17eSHong Zhang currentlevelrk = nextlevelrk; 4940fe4d17eSHong Zhang } 4950fe4d17eSHong Zhang ierr = VecDestroy(&X0);CHKERRQ(ierr); 4960fe4d17eSHong Zhang 49702555c94SHong Zhang ts->ops->step = TSStep_RK_MultirateSplit; 49802555c94SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_RK_MultirateSplit; 49902555c94SHong Zhang ts->ops->interpolate = TSInterpolate_RK_MultirateSplit; 5000fe4d17eSHong Zhang PetscFunctionReturn(0); 5010fe4d17eSHong Zhang } 5020fe4d17eSHong Zhang 50321052c0cSHong Zhang PetscErrorCode TSRKSetMultirate_RK(TS ts,PetscBool use_multirate) 504474dd773SHong Zhang { 505474dd773SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 506474dd773SHong Zhang PetscErrorCode ierr; 507474dd773SHong Zhang 508474dd773SHong Zhang PetscFunctionBegin; 509474dd773SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5100fe4d17eSHong Zhang rk->use_multirate = use_multirate; 5110fe4d17eSHong Zhang if (use_multirate) { 512474dd773SHong Zhang rk->dtratio = 2; 51302555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateSplit_C",TSSetUp_RK_MultirateSplit);CHKERRQ(ierr); 51402555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateSplit_C",TSReset_RK_MultirateSplit);CHKERRQ(ierr); 51502555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateNonsplit_C",TSSetUp_RK_MultirateNonsplit);CHKERRQ(ierr); 51602555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateNonsplit_C",TSReset_RK_MultirateNonsplit);CHKERRQ(ierr); 5170fe4d17eSHong Zhang } else { 5180fe4d17eSHong Zhang rk->dtratio = 0; 51902555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateSplit_C",NULL);CHKERRQ(ierr); 52002555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateSplit_C",NULL);CHKERRQ(ierr); 52102555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateNonsplit_C",NULL);CHKERRQ(ierr); 52202555c94SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateNonsplit_C",NULL);CHKERRQ(ierr); 523474dd773SHong Zhang } 524474dd773SHong Zhang PetscFunctionReturn(0); 525474dd773SHong Zhang } 526474dd773SHong Zhang 52721052c0cSHong Zhang PetscErrorCode TSRKGetMultirate_RK(TS ts,PetscBool *use_multirate) 5280fe4d17eSHong Zhang { 5290fe4d17eSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5300fe4d17eSHong Zhang 5310fe4d17eSHong Zhang PetscFunctionBegin; 5320fe4d17eSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5330fe4d17eSHong Zhang *use_multirate = rk->use_multirate; 5340fe4d17eSHong Zhang PetscFunctionReturn(0); 5350fe4d17eSHong Zhang } 536