<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class="">
<br class="">
<div><br class="">
<blockquote type="cite" class="">
<div class="">On Jul 3, 2019, at 3:10 PM, Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu" class="">mvalera-w@sdsu.edu</a>> wrote:</div>
<br class="Apple-interchange-newline">
<div class="">
<div dir="ltr" class="">Thanks Zhang for your answer,
<div class=""><br class="">
</div>
<div class="">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,</div>
<div class=""><br class="">
</div>
<div class="">My TS code so far looks like this:</div>
<div class=""><br class="">
</div>
<div class="">(...initialization...)</div>
<div class=""><br class="">
</div>
<div class=""><font face="courier new, monospace" class="">call TSCreate(PETSC_COMM_WORLD,ts,ierr)<br class="">
call TSSetProblemType(ts,TS_NONLINEAR,ierr)<br class="">
call TSSetSolution(ts,gTemperature,ierr)<br class="">
call TSSetRHSFunction(ts,PETSC_NULL_VEC,FormRHSFunction,PETSC_NULL_VEC,ierr)<br class="">
</font></div>
</div>
</div>
</blockquote>
<div><br class="">
</div>
<div>The second last argument should be the user context. If you pass NULL, the ctx variable in your FormRHSFunction will be NULL.  </div>
<br class="">
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""><font face="courier new, monospace" class="">call TSSetType(ts,TSSSP,ierr)<br class="">
call TSSetTimeStep(ts,dt,ierr)<br class="">
call TSSetDM(ts,daScalars,ierr)<br class="">
call TSSetMaxTime(ts,TotalTime,ierr)<br class="">
call TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr)<br class="">
call TSSetFromOptions(ts,ierr)<br class="">
call TSSetUp(ts,ierr)</font><br class="">
</div>
<div class=""><br class="">
</div>
<div class="">(... then inside main calculations loop i have...)</div>
<div class=""><br class="">
</div>
<div class=""><font face="courier new, monospace" class="">call TSSolve(ts,gTemperature,ierr)<br class="">
</font></div>
<div class=""><br class="">
</div>
<div class="">(...and my RHSFunction looks like this...)</div>
<div class=""><br class="">
</div>
<div class=""><font face="courier new, monospace" class="">subroutine FormRHSFunction(ts,t,t_loc,rhs_loc,ctx,ierr)<br class="">
<br class="">
  use petscvec<br class="">
  use petscts<br class="">
  use petscdmda<br class="">
  use petscdm<br class="">
<br class="">
  USE utils<br class="">
  USE dmdaobjs<br class="">
  USE probsize<br class="">
  USE modelparams<br class="">
<br class="">
  implicit none<br class="">
<br class="">
  TS            :: ts<br class="">
  PetscReal     :: t<br class="">
  Vec           :: ctx,t_loc,rhs_loc<br class="">
  PetscErrorCode    :: ierr<br class="">
<br class="">
<br class="">
<br class="">
        call TSGetDM(ts,daScalars,ierr)<br class="">
<br class="">
        call DMGetLocalVector(daScalars,t_loc,ierr);<br class="">
</font></div>
</div>
</div>
</blockquote>
<div><br class="">
</div>
t_loc is the input, it should be not modified.</div>
<div><br class="">
</div>
<div>
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""><font face="courier new, monospace" class="">        call DMGlobalToLocalBegin(daScalars,gTemperature,INSERT_VALUES,t_loc,ierr);CHKERRQ(ierr)<br class="">
        call DMGlobalToLocalEnd(daScalars,gTemperature,INSERT_VALUES,t_loc,ierr);CHKERRQ(ierr)<br class="">
</font></div>
</div>
</div>
</blockquote>
<div><br class="">
</div>
<div>FormRHSFunction is supposed to implement rhs_loc = f(t,t_loc). So you want to scatter ghost points for t_loc, not gTemperature.</div>
<br class="">
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""><font face="courier new, monospace" class="">        call DMDASolveTExplicit(3)<br class="">
<br class="">
        call DMGetLocalVector(daScalars,rhs_loc,ierr);<br class="">
        call DMGlobalToLocalBegin(daScalars,TRHSS_ts,INSERT_VALUES,rhs_loc,ierr);CHKERRQ(ierr)<br class="">
        call DMGlobalToLocalEnd(daScalars,TRHSS_ts,INSERT_VALUES,rhs_loc,ierr);CHKERRQ(ierr)<br class="">
</font></div>
</div>
</div>
</blockquote>
<div><br class="">
</div>
<div>There is no need to scatter ghost points for rhs_loc.</div>
<div><br class="">
</div>
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""><font face="courier new, monospace" class=""><br class="">
<br class="">
        call DMRestoreLocalVector(daScalars,t_loc,ierr);CHKERRQ(ierr)<br class="">
        call DMRestoreLocalVector(daScalars,rhs_loc,ierr);CHKERRQ(ierr)<br class="">
<br class="">
end subroutine FormRHSFunction</font><br class="">
</div>
<div class=""><br class="">
</div>
<div class="">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,</div>
<div class=""><br class="">
</div>
<div class="">I still have several doubts:</div>
<div class=""><br class="">
</div>
<div class="">
<ul class="">
<li class="">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. </li></ul>
</div>
</div>
</div>
</blockquote>
<div>For explicit time integration, one needs to provide only RHSFunction.</div>
<br class="">
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class="">
<ul class="">
<li class="">This 't' parameter in the RHSFunction is controlled by PETSC? i am not providing anything for it directly, where is it coming from?</li></ul>
</div>
</div>
</div>
</blockquote>
<div><br class="">
</div>
<div>It is controlled by PETSc. If your problem is autonomous (the RHS does not depend on t), it can be simply ignored.</div>
<div> </div>
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class="">
<ul class="">
<li class="">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.</li></ul>
</div>
</div>
</div>
</blockquote>
<div>For explicit time integration, Jacobian is not needed.</div>
<div><br class="">
</div>
<div>Hong</div>
<br class="">
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class="">
<div class=""><br class="">
</div>
</div>
<div class="">Thanks, any help is appreciated, if you see any errors or need more information please let me know,</div>
</div>
</div>
</blockquote>
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""><br class="">
</div>
<div class="">Regards,</div>
<div class=""><br class="">
</div>
<div class="">Manuel </div>
</div>
<br class="">
<div class="gmail_quote">
<div dir="ltr" class="gmail_attr">On Wed, Jun 26, 2019 at 9:54 PM Zhang, Hong <<a href="mailto:hongzhang@anl.gov" class="">hongzhang@anl.gov</a>> wrote:<br class="">
</div>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
<div style="overflow-wrap: break-word;" class=""><br class="">
<div class=""><br class="">
<blockquote type="cite" class="">
<div class="">On Jun 26, 2019, at 4:17 PM, Manuel Valera via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank" class="">petsc-users@mcs.anl.gov</a>> wrote:</div>
<br class="gmail-m_6853237684116060569Apple-interchange-newline">
<div class="">
<div dir="ltr" class="">Hi PETSc,<br class="">
<div class=""><br class="">
</div>
<div class="">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:</div>
<div class=""><br class="">
</div>
<div class="">
<ul class="">
<li class="">Interpolate u,v,w to the time advancing variable position.</li><li class="">Calculate nonlinear coefficients and advect velocities with a forward-backward shock capturing scheme.</li><li class="">Calculate the variable laplacian </li><li class="">Sum terms to create RHS (nonlinear advected velocities + laplacian)</li><li class="">Finally, the runge kutta integration is done in a typical way that looks like:</li></ul>
<div class="">                   newtemp(t+1) = prevtemp(t) + dt*RHS</div>
</div>
<div class=""><br class="">
</div>
<div class=""><br class="">
</div>
<div class="">So my questions are:</div>
<div class="">
<ul class="">
<li class="">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.</li></ul>
</div>
</div>
</div>
</blockquote>
<span class="">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. </span><span class=""><br class="">
</span>
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class="">
<ul class="">
<li class="">How do i know what are the appropriate routines i should be using here? so far i think i should use the following:</li></ul>
<div class="">call TSCreate(PETSC_COMM_WORLD,ts,ierr)<br class="">
call TSSetProblemType(ts,TS_LINEAR,ierr)<br class="">
call TSSetTimeStep(ts,dt,ierr)<br class="">
<br class="">
call TSSetFromOptions(ts,ierr)<br class="">
</div>
<div class=""><br class="">
</div>
<div class="">call TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,NULL,ierr)<br class="">
call TSSolve(ts,loctemperature,ierr)<br class="">
<br class="">
</div>
</div>
<div class="">Should i use call TSSetRHSJacobian for the temperature jacobian in this case?</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">I would suggest to write your own RHSFunction and set it to TS with TSSetRHSFunction().</div>
<div class=""><br class="">
</div>
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""> </div>
<div class=""><br class="">
</div>
<div class="">I am using <a href="https://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex4.c.html" target="_blank" class="">https://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex4.c.html</a> as a  general guide, is there a
 more appropriate example?</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
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.</div>
<div class=""><br class="">
</div>
<div class="">Hong (Mr.)</div>
<div class=""><br class="">
</div>
<div class="">
<blockquote type="cite" class="">
<div class="">
<div dir="ltr" class="">
<div class=""><br class="">
</div>
<div class="">The dt value and total timesteps are controlled by the model,</div>
<div class=""><br class="">
</div>
<div class="">Thanks for your help,</div>
<div class=""><br class="">
</div>
<div class=""><br class="">
</div>
</div>
</div>
</blockquote>
</div>
<br class="">
</div>
</blockquote>
</div>
</div>
</blockquote>
</div>
<br class="">
</body>
</html>