[petsc-users] Implementing TS routine
Zhang, Hong
hongzhang at anl.gov
Thu Jul 4 23:18:20 CDT 2019
On Jul 3, 2019, at 3:10 PM, Manuel Valera <mvalera-w at sdsu.edu<mailto:mvalera-w at sdsu.edu>> wrote:
Thanks Zhang for your answer,
I ended up getting a compiling and running TS routine... that does not give me the answers, and i am not sure what i am doing wrong,
My TS code so far looks like this:
(...initialization...)
call TSCreate(PETSC_COMM_WORLD,ts,ierr)
call TSSetProblemType(ts,TS_NONLINEAR,ierr)
call TSSetSolution(ts,gTemperature,ierr)
call TSSetRHSFunction(ts,PETSC_NULL_VEC,FormRHSFunction,PETSC_NULL_VEC,ierr)
The second last argument should be the user context. If you pass NULL, the ctx variable in your FormRHSFunction will be NULL.
call TSSetType(ts,TSSSP,ierr)
call TSSetTimeStep(ts,dt,ierr)
call TSSetDM(ts,daScalars,ierr)
call TSSetMaxTime(ts,TotalTime,ierr)
call TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr)
call TSSetFromOptions(ts,ierr)
call TSSetUp(ts,ierr)
(... then inside main calculations loop i have...)
call TSSolve(ts,gTemperature,ierr)
(...and my RHSFunction looks like this...)
subroutine FormRHSFunction(ts,t,t_loc,rhs_loc,ctx,ierr)
use petscvec
use petscts
use petscdmda
use petscdm
USE utils
USE dmdaobjs
USE probsize
USE modelparams
implicit none
TS :: ts
PetscReal :: t
Vec :: ctx,t_loc,rhs_loc
PetscErrorCode :: ierr
call TSGetDM(ts,daScalars,ierr)
call DMGetLocalVector(daScalars,t_loc,ierr);
t_loc is the input, it should be not modified.
call DMGlobalToLocalBegin(daScalars,gTemperature,INSERT_VALUES,t_loc,ierr);CHKERRQ(ierr)
call DMGlobalToLocalEnd(daScalars,gTemperature,INSERT_VALUES,t_loc,ierr);CHKERRQ(ierr)
FormRHSFunction is supposed to implement rhs_loc = f(t,t_loc). So you want to scatter ghost points for t_loc, not gTemperature.
call DMDASolveTExplicit(3)
call DMGetLocalVector(daScalars,rhs_loc,ierr);
call DMGlobalToLocalBegin(daScalars,TRHSS_ts,INSERT_VALUES,rhs_loc,ierr);CHKERRQ(ierr)
call DMGlobalToLocalEnd(daScalars,TRHSS_ts,INSERT_VALUES,rhs_loc,ierr);CHKERRQ(ierr)
There is no need to scatter ghost points for rhs_loc.
call DMRestoreLocalVector(daScalars,t_loc,ierr);CHKERRQ(ierr)
call DMRestoreLocalVector(daScalars,rhs_loc,ierr);CHKERRQ(ierr)
end subroutine FormRHSFunction
Where DMDASolveTExplicit(3) is the old function to calculate time integration with runge kutta, modified to only generate the f(t,u) which in this case is rhs_loc,
I still have several doubts:
* Will this explicit RHS calculation work with TS routines? my time integration is explicit so far and it would involve a fair deal of extra work to make it implicit.
For explicit time integration, one needs to provide only RHSFunction.
* This 't' parameter in the RHSFunction is controlled by PETSC? i am not providing anything for it directly, where is it coming from?
It is controlled by PETSc. If your problem is autonomous (the RHS does not depend on t), it can be simply ignored.
* Do i need to provide a Jacobian or the TS routine will try to guess one? This is related to the first point where my time scheme being explicit does not use a jacobian.
For explicit time integration, Jacobian is not needed.
Hong
Thanks, any help is appreciated, if you see any errors or need more information please let me know,
Regards,
Manuel
On Wed, Jun 26, 2019 at 9:54 PM Zhang, Hong <hongzhang at anl.gov<mailto:hongzhang at anl.gov>> wrote:
On Jun 26, 2019, at 4:17 PM, Manuel Valera via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:
Hi PETSc,
I am trying to implement the Time stepping routines in my model, i have a working runge-kutta time scheme that goes to the following steps:
* Interpolate u,v,w to the time advancing variable position.
* Calculate nonlinear coefficients and advect velocities with a forward-backward shock capturing scheme.
* Calculate the variable laplacian
* Sum terms to create RHS (nonlinear advected velocities + laplacian)
* Finally, the runge kutta integration is done in a typical way that looks like:
newtemp(t+1) = prevtemp(t) + dt*RHS
So my questions are:
* I think my problem is nonlinear, but is being made linearized by the advecting scheme, is this correct? this is to know if i should use the linear or nonlinear routines for TS.
TSComputeRHSFunctionLinear is just a convenience function for linear ODEs in the form udot=Au. Using it won’t buy you much. So for TS starters, it is fine to assume your problem is nonlinear and think of the form udot=f(t,u) where f is the RHS function.
* How do i know what are the appropriate routines i should be using here? so far i think i should use the following:
call TSCreate(PETSC_COMM_WORLD,ts,ierr)
call TSSetProblemType(ts,TS_LINEAR,ierr)
call TSSetTimeStep(ts,dt,ierr)
call TSSetFromOptions(ts,ierr)
call TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,NULL,ierr)
call TSSolve(ts,loctemperature,ierr)
Should i use call TSSetRHSJacobian for the temperature jacobian in this case?
I would suggest to write your own RHSFunction and set it to TS with TSSetRHSFunction().
I am using https://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex4.c.html as a general guide, is there a more appropriate example?
ex4 is a good example. In addition, ex9 uses finite volume method with slope limiters to solve a variety of problems such as advection equation, burgers equation and shallow water equation. It might be an overkill, but it seems to be close to the problem you are trying to solve. Note that you might want to use the SSP methods (-ts_type ssp) instead of the classic Runge-Kutta methods for problems with shocks.
Hong (Mr.)
The dt value and total timesteps are controlled by the model,
Thanks for your help,
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190705/f89dd97e/attachment-0001.html>
More information about the petsc-users
mailing list