class SharpClawSolver3D(SharpClawSolver): # ======================================================================== """ Three Dimensional SharpClawSolver """ __doc__ += add_parent_doc(SharpClawSolver) def __init__(self,riemann_solver=None,claw_package=None): r""" Create 3D SharpClaw solver See :class:`SharpClawSolver3D` for more info. """ self.num_dim = 3 self.reflect_index = [1,2,3] super(SharpClawSolver3D,self).__init__(riemann_solver,claw_package) def teardown(self): r""" Deallocate F90 module arrays. Also delete Fortran objects, which otherwise tend to persist in Python sessions. """ if self.kernel_language=='Fortran': self.fmod.workspace.dealloc_workspace(self.char_decomp) self.fmod.reconstruct.dealloc_recon_workspace(self.fmod.clawparams.lim_type,self.fmod.clawparams.char_decomp) self.fmod.clawparams.dealloc_clawparams() del self.fmod def dq_hyperbolic(self,state): """Compute dq/dt * (delta t) for the hyperbolic hyperbolic system Note that the capa array, if present, should be located in the aux variable. Indexing works like this (here num_ghost=2 as an example):: 0 1 2 3 4 mx+num_ghost-2 mx+num_ghost mx+num_ghost+2 | mx+num_ghost-1 | mx+num_ghost+1 | | | | | ... | | | | | 0 1 | 2 3 mx+num_ghost-2 |mx+num_ghost mx+num_ghost-1 mx+num_ghost+1 The top indices represent the values that are located on the grid cell boundaries such as waves, s and other Riemann problem values, the bottom for the cell centered values such as q. In particular the ith grid cell boundary has the following related information:: i-1 i i+1 | | | | i-1 | i | | | | Again, grid cell boundary quantities are at the top, cell centered values are in the cell. """ self._apply_bcs(state) q = self.qbc grid = state.grid num_ghost=self.num_ghost mx=grid.num_cells[0] my=grid.num_cells[1] mz=grid.num_cells[2] maxm = max(mx,my,mz) if self.kernel_language=='Fortran': rpn3 = self.rp.rpn3._cpointer if self.tfluct_solver: tfluct3 = self.tfluct.tfluct3._cpointer else: tfluct3 = self.tfluct dq,cfl=self.fmod.flux3(q,self.auxbc,self.dt,state.t,num_ghost,maxm,mx,my,mz,rpn3,tfluct3) else: raise Exception('Only Fortran kernels are supported in 3D.') self.cfl.update_global_max(cfl) return dq[:,num_ghost:-num_ghost,num_ghost:-num_ghost,num_ghost:-num_ghost]
class SharpClawSolver2D(SharpClawSolver): # ======================================================================== """ Two Dimensional SharpClawSolver """ __doc__ += add_parent_doc(SharpClawSolver) def __init__(self,riemann_solver=None,claw_package=None): r""" Create 2D SharpClaw solver See :class:`SharpClawSolver2D` for more info. """ self.num_dim = 2 self.reflect_index = [1,2] super(SharpClawSolver2D,self).__init__(riemann_solver,claw_package) def dq_hyperbolic(self,state): """Compute dq/dt * (delta t) for the hyperbolic hyperbolic system Note that the capa array, if present, should be located in the aux variable. Indexing works like this (here num_ghost=2 as an example):: 0 1 2 3 4 mx+num_ghost-2 mx+num_ghost mx+num_ghost+2 | mx+num_ghost-1 | mx+num_ghost+1 | | | | | ... | | | | | 0 1 | 2 3 mx+num_ghost-2 |mx+num_ghost mx+num_ghost-1 mx+num_ghost+1 The top indices represent the values that are located on the grid cell boundaries such as waves, s and other Riemann problem values, the bottom for the cell centered values such as q. In particular the ith grid cell boundary has the following related information:: i-1 i i+1 | | | | i-1 | i | | | | Again, grid cell boundary quantities are at the top, cell centered values are in the cell. """ self._apply_bcs(state) q = self.qbc grid = state.grid num_ghost=self.num_ghost mx=grid.num_cells[0] my=grid.num_cells[1] maxm = max(mx,my) if self.kernel_language=='Fortran': rpn2 = self.rp.rpn2._cpointer if self.tfluct_solver: tfluct2 = self.tfluct.tfluct2._cpointer else: tfluct2 = self.tfluct dq,cfl=self.fmod.flux2(q,self.auxbc,self.dt,state.t,num_ghost,maxm,mx,my,rpn2,tfluct2) else: raise Exception('Only Fortran kernels are supported in 2D.') self.cfl.update_global_max(cfl) return dq[:,num_ghost:-num_ghost,num_ghost:-num_ghost]
class SharpClawSolver1D(SharpClawSolver): # ======================================================================== """ SharpClaw solver for one-dimensional problems. Used to solve 1D hyperbolic systems using the SharpClaw algorithms, which are based on WENO reconstruction and Runge-Kutta time stepping. """ __doc__ += add_parent_doc(SharpClawSolver) def __init__(self,riemann_solver=None,claw_package=None): r""" See :class:`SharpClawSolver1D` for more info. """ self.num_dim = 1 self.reflect_index = [1] super(SharpClawSolver1D,self).__init__(riemann_solver,claw_package) def dq_hyperbolic(self,state): r""" Compute dq/dt * (delta t) for the hyperbolic hyperbolic system. Note that the capa array, if present, should be located in the aux variable. Indexing works like this (here num_ghost=2 as an example):: 0 1 2 3 4 mx+num_ghost-2 mx+num_ghost mx+num_ghost+2 | mx+num_ghost-1 | mx+num_ghost+1 | | | | | ... | | | | | 0 1 | 2 3 mx+num_ghost-2 |mx+num_ghost mx+num_ghost-1 mx+num_ghost+1 The top indices represent the values that are located on the grid cell boundaries such as waves, s and other Riemann problem values, the bottom for the cell centered values such as q. In particular the ith grid cell boundary has the following related information:: i-1 i i+1 | | | | i-1 | i | | | | Again, grid cell boundary quantities are at the top, cell centered values are in the cell. """ import numpy as np self._apply_bcs(state) q = self.qbc aux = self.auxbc grid = state.grid mx = grid.num_cells[0] ixy=1 if self.kernel_language=='Fortran': rp1 = self.rp.rp1._cpointer if self.tfluct_solver: tfluct1 = self.tfluct.tfluct1._cpointer else: tfluct1 = self.tfluct dq,cfl=self.fmod.flux1(q,self.auxbc,self.dt,state.t,ixy,mx,self.num_ghost,mx,rp1,tfluct1) elif self.kernel_language=='Python': dtdx = np.zeros( (mx+2*self.num_ghost) ,order='F') dq = np.zeros( (state.num_eqn,mx+2*self.num_ghost) ,order='F') # Find local value for dt/dx if state.index_capa>=0: dtdx = self.dt / (grid.delta[0] * aux[state.index_capa,:]) else: dtdx += self.dt/grid.delta[0] if aux.shape[0]>0: aux_l=aux[:,:-1] aux_r=aux[:,1: ] else: aux_l = None aux_r = None #Reconstruct (wave reconstruction uses a Riemann solve) if self.lim_type==-1: #1st-order Godunov ql=q; qr=q elif self.lim_type==0: #Unlimited reconstruction raise NotImplementedError('Unlimited reconstruction not implemented') elif self.lim_type==1: #TVD Reconstruction raise NotImplementedError('TVD reconstruction not implemented') elif self.lim_type==2: #WENO Reconstruction if self.char_decomp==0: #No characteristic decomposition ql,qr=recon.weno(5,q) elif self.char_decomp==1: #Wave-based reconstruction q_l=q[:,:-1] q_r=q[:,1: ] wave,s,amdq,apdq = self.rp(q_l,q_r,aux_l,aux_r,state.problem_data) ql,qr=recon.weno5_wave(q,wave,s) elif self.char_decomp==2: #Characteristic-wise reconstruction raise NotImplementedError # Solve Riemann problem at each interface q_l=qr[:,:-1] q_r=ql[:,1: ] wave,s,amdq,apdq = self.rp(q_l,q_r,aux_l,aux_r,state.problem_data) # Loop limits for local portion of grid # THIS WON'T WORK IN PARALLEL! LL = self.num_ghost - 1 UL = grid.num_cells[0] + self.num_ghost + 1 # Compute maximum wave speed cfl = 0.0 for mw in range(self.num_waves): smax1 = np.max( dtdx[LL :UL] *s[mw,LL-1:UL-1]) smax2 = np.max(-dtdx[LL-1:UL-1]*s[mw,LL-1:UL-1]) cfl = max(cfl,smax1,smax2) #Find total fluctuation within each cell wave,s,amdq2,apdq2 = self.rp(ql,qr,aux,aux,state.problem_data) # Compute dq for m in range(state.num_eqn): dq[m,LL:UL] = -dtdx[LL:UL]*(amdq[m,LL:UL] + apdq[m,LL-1:UL-1] \ + apdq2[m,LL:UL] + amdq2[m,LL:UL]) else: raise Exception('Unrecognized value of solver.kernel_language.') self.cfl.update_global_max(cfl) return dq[:,self.num_ghost:-self.num_ghost]
class ClawSolver3D(ClawSolver): r""" 3D Classic (Clawpack) solver. Solve using the wave propagation algorithms of Randy LeVeque's Clawpack code (www.clawpack.org). In addition to the attributes of ClawSolver, ClawSolver3D also has the following options: .. attribute:: dimensional_split If True, use dimensional splitting (Godunov splitting). Dimensional splitting with Strang splitting is not supported at present but could easily be enabled if necessary. If False, use unsplit Clawpack algorithms, possibly including transverse Riemann solves. .. attribute:: transverse_waves If dimensional_split is True, this option has no effect. If dim_plit is False, then transverse_waves should be one of the following values: ClawSolver3D.no_trans: Transverse Riemann solver not used. The stable CFL for this algorithm is 0.5. Not recommended. ClawSolver3D.trans_inc: Transverse increment waves are computed and propagated. ClawSolver3D.trans_cor: Transverse increment waves and transverse correction waves are computed and propagated. Note that only Fortran routines are supported for now in 3D -- there is no pure-python version. """ __doc__ += add_parent_doc(ClawSolver) no_trans = 0 trans_inc = 11 trans_cor = 22 def __init__(self, riemann_solver=None, claw_package=None): r""" Create 3d Clawpack solver See :class:`ClawSolver3D` for more info. """ # Add the functions as required attributes self.dimensional_split = True self.transverse_waves = self.trans_cor self.num_dim = 3 self.aux1 = None self.aux2 = None self.aux3 = None self.work = None super(ClawSolver3D, self).__init__(riemann_solver, claw_package) # ========== Setup routine ============================= def _allocate_workspace(self, solution): r""" Allocate auxN and work arrays for use in Fortran subroutines. """ import numpy as np state = solution.states[0] num_eqn, num_aux, num_waves, num_ghost, aux = state.num_eqn, state.num_aux, self.num_waves, self.num_ghost, state.aux #The following is a hack to work around an issue #with f2py. It involves wastefully allocating three arrays. #f2py seems not able to handle multiple zero-size arrays being passed. # it appears the bug is related to f2py/src/fortranobject.c line 841. if (aux == None): num_aux = 1 grid = state.grid maxmx, maxmy, maxmz = grid.num_cells[0], grid.num_cells[ 1], grid.num_cells[2] maxm = max(maxmx, maxmy, maxmz) # These work arrays really ought to live inside a fortran module # as is done for sharpclaw self.aux1 = np.empty((num_aux, maxm + 2 * num_ghost, 3), order='F') self.aux2 = np.empty((num_aux, maxm + 2 * num_ghost, 3), order='F') self.aux3 = np.empty((num_aux, maxm + 2 * num_ghost, 3), order='F') mwork = (maxm + 2 * num_ghost) * (31 * num_eqn + num_waves + num_eqn * num_waves) self.work = np.empty((mwork), order='F') # ========== Hyperbolic Step ===================================== def step_hyperbolic(self, solution): r""" Take a step on the homogeneous hyperbolic system using the Clawpack algorithm. Clawpack is based on the Lax-Wendroff method, combined with Riemann solvers and TVD limiters applied to waves. """ if (self.kernel_language == 'Fortran'): state = solution.states[0] grid = state.grid dx, dy, dz = grid.delta mx, my, mz = grid.num_cells maxm = max(mx, my, mz) self._apply_q_bcs(state) if state.num_aux > 0: self._apply_aux_bcs(state) qnew = self.qbc qold = qnew.copy('F') rpn3 = self.rp.rpn3._cpointer rpt3 = self.rp.rpt3._cpointer rptt3 = self.rp.rptt3._cpointer if self.dimensional_split: #Right now only Godunov-dimensional-splitting is implemented. #Strang-dimensional-splitting could be added following dimsp2.f in Clawpack. q, cfl_x = self.fmod.step3ds(maxm,self.num_ghost,mx,my,mz, \ qold,qnew,self.auxbc,dx,dy,dz,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,1,rpn3,rpt3,rptt3) q, cfl_y = self.fmod.step3ds(maxm,self.num_ghost,mx,my,mz, \ q,q,self.auxbc,dx,dy,dz,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,2,rpn3,rpt3,rptt3) q, cfl_z = self.fmod.step3ds(maxm,self.num_ghost,mx,my,mz, \ q,q,self.auxbc,dx,dy,dz,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,3,rpn3,rpt3,rptt3) cfl = max(cfl_x, cfl_y, cfl_z) else: q, cfl = self.fmod.step3(maxm,self.num_ghost,mx,my,mz, \ qold,qnew,self.auxbc,dx,dy,dz,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,rpn3,rpt3,rptt3) self.cfl.update_global_max(cfl) state.set_q_from_qbc(self.num_ghost, self.qbc) if state.num_aux > 0: state.set_aux_from_auxbc(self.num_ghost, self.auxbc) else: raise NotImplementedError( "No python implementation for step_hyperbolic in 3D.")
class ClawSolver2D(ClawSolver): r""" 2D Classic (Clawpack) solver. Solve using the wave propagation algorithms of Randy LeVeque's Clawpack code (www.clawpack.org). In addition to the attributes of ClawSolver1D, ClawSolver2D also has the following options: .. attribute:: dimensional_split If True, use dimensional splitting (Godunov splitting). Dimensional splitting with Strang splitting is not supported at present but could easily be enabled if necessary. If False, use unsplit Clawpack algorithms, possibly including transverse Riemann solves. .. attribute:: transverse_waves If dimensional_split is True, this option has no effect. If dim_plit is False, then transverse_waves should be one of the following values: ClawSolver2D.no_trans: Transverse Riemann solver not used. The stable CFL for this algorithm is 0.5. Not recommended. ClawSolver2D.trans_inc: Transverse increment waves are computed and propagated. ClawSolver2D.trans_cor: Transverse increment waves and transverse correction waves are computed and propagated. Note that only the fortran routines are supported for now in 2D. """ __doc__ += add_parent_doc(ClawSolver) no_trans = 0 trans_inc = 1 trans_cor = 2 def __init__(self, riemann_solver=None, claw_package=None): r""" Create 2d Clawpack solver See :class:`ClawSolver2D` for more info. """ self.dimensional_split = True self.transverse_waves = self.trans_inc self.num_dim = 2 self.aux1 = None self.aux2 = None self.aux3 = None self.work = None super(ClawSolver2D, self).__init__(riemann_solver, claw_package) def _check_cfl_settings(self): if (not self.dimensional_split) and (self.transverse_waves == 0): cfl_recommended = 0.5 else: cfl_recommended = 1.0 if self.cfl_max > cfl_recommended: import warnings warnings.warn( 'cfl_max is set higher than the recommended value of %s' % cfl_recommended) warnings.warn(str(self.cfl_desired)) def _allocate_workspace(self, solution): r""" Pack parameters into format recognized by Clawpack (Fortran) code. Sets the method array and the cparam common block for the Riemann solver. """ import numpy as np state = solution.state num_eqn, num_aux, num_waves, num_ghost, aux = state.num_eqn, state.num_aux, self.num_waves, self.num_ghost, state.aux #The following is a hack to work around an issue #with f2py. It involves wastefully allocating three arrays. #f2py seems not able to handle multiple zero-size arrays being passed. # it appears the bug is related to f2py/src/fortranobject.c line 841. if (aux == None): num_aux = 1 grid = state.grid maxmx, maxmy = grid.num_cells[0], grid.num_cells[1] maxm = max(maxmx, maxmy) # These work arrays really ought to live inside a fortran module # as is done for sharpclaw self.aux1 = np.empty((num_aux, maxm + 2 * num_ghost), order='F') self.aux2 = np.empty((num_aux, maxm + 2 * num_ghost), order='F') self.aux3 = np.empty((num_aux, maxm + 2 * num_ghost), order='F') mwork = (maxm + 2 * num_ghost) * (5 * num_eqn + num_waves + num_eqn * num_waves) self.work = np.empty((mwork), order='F') # ========== Hyperbolic Step ===================================== def step_hyperbolic(self, solution): r""" Take a step on the homogeneous hyperbolic system using the Clawpack algorithm. Clawpack is based on the Lax-Wendroff method, combined with Riemann solvers and TVD limiters applied to waves. """ if (self.kernel_language == 'Fortran'): state = solution.states[0] grid = state.grid dx, dy = grid.delta mx, my = grid.num_cells maxm = max(mx, my) self._apply_q_bcs(state) if state.num_aux > 0: self._apply_aux_bcs(state) qold = self.qbc.copy('F') rpn2 = self.rp.rpn2._cpointer rpt2 = self.rp.rpt2._cpointer if self.dimensional_split: #Right now only Godunov-dimensional-splitting is implemented. #Strang-dimensional-splitting could be added following dimsp2.f in Clawpack. self.qbc, cfl_x = self.fmod.step2ds(maxm,self.num_ghost,mx,my, \ qold,self.qbc,self.auxbc,dx,dy,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,1,self.fwave,rpn2,rpt2) self.qbc, cfl_y = self.fmod.step2ds(maxm,self.num_ghost,mx,my, \ self.qbc,self.qbc,self.auxbc,dx,dy,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,2,self.fwave,rpn2,rpt2) cfl = max(cfl_x, cfl_y) else: self.qbc, cfl = self.fmod.step2(maxm,self.num_ghost,mx,my, \ qold,self.qbc,self.auxbc,dx,dy,self.dt,self._method,self._mthlim,\ self.aux1,self.aux2,self.aux3,self.work,self.fwave,rpn2,rpt2) self.cfl.update_global_max(cfl) state.set_q_from_qbc(self.num_ghost, self.qbc) if state.num_aux > 0: state.set_aux_from_auxbc(self.num_ghost, self.auxbc) else: raise NotImplementedError( "No python implementation for step_hyperbolic in 2D.")
class ClawSolver1D(ClawSolver): r""" Clawpack evolution routine in 1D This class represents the 1d clawpack solver on a single grid. Note that there are routines here for interfacing with the fortran time stepping routines and the Python time stepping routines. The ones used are dependent on the argument given to the initialization of the solver (defaults to python). """ __doc__ += add_parent_doc(ClawSolver) def __init__(self, riemann_solver=None, claw_package=None): r""" Create 1d Clawpack solver Output: - (:class:`ClawSolver1D`) - Initialized 1d clawpack solver See :class:`ClawSolver1D` for more info. """ self.num_dim = 1 super(ClawSolver1D, self).__init__(riemann_solver, claw_package) # ========== Homogeneous Step ===================================== def step_hyperbolic(self, solution): r""" Take one time step on the homogeneous hyperbolic system. :Input: - *solution* - (:class:`~pyclaw.solution.Solution`) Solution that will be evolved """ import numpy as np state = solution.states[0] grid = state.grid self._apply_q_bcs(state) if state.num_aux > 0: self._apply_aux_bcs(state) num_eqn, num_ghost = state.num_eqn, self.num_ghost if (self.kernel_language == 'Fortran'): mx = grid.num_cells[0] dx, dt = grid.delta[0], self.dt dtdx = np.zeros((mx + 2 * num_ghost)) + dt / dx rp1 = self.rp.rp1._cpointer self.qbc, cfl = self.fmod.step1(num_ghost, mx, self.qbc, self.auxbc, dx, dt, self._method, self._mthlim, self.fwave, rp1) elif (self.kernel_language == 'Python'): q = self.qbc aux = self.auxbc # Limiter to use in the pth family limiter = np.array(self._mthlim, ndmin=1) dtdx = np.zeros((2 * self.num_ghost + grid.num_cells[0])) # Find local value for dt/dx if state.index_capa >= 0: dtdx = self.dt / (grid.delta[0] * state.aux[state.index_capa, :]) else: dtdx += self.dt / grid.delta[0] # Solve Riemann problem at each interface q_l = q[:, :-1] q_r = q[:, 1:] if state.aux is not None: aux_l = aux[:, :-1] aux_r = aux[:, 1:] else: aux_l = None aux_r = None wave, s, amdq, apdq = self.rp(q_l, q_r, aux_l, aux_r, state.problem_data) # Update loop limits, these are the limits for the Riemann solver # locations, which then update a grid cell value # We include the Riemann problem just outside of the grid so we can # do proper limiting at the grid edges # LL | | UL # | LL | | | | ... | | | UL | | # | | LL = self.num_ghost - 1 UL = self.num_ghost + grid.num_cells[0] + 1 # Update q for Godunov update for m in xrange(num_eqn): q[m, LL:UL] -= dtdx[LL:UL] * apdq[m, LL - 1:UL - 1] q[m, LL - 1:UL - 1] -= dtdx[LL - 1:UL - 1] * amdq[m, LL - 1:UL - 1] # Compute maximum wave speed cfl = 0.0 for mw in xrange(wave.shape[1]): smax1 = np.max(dtdx[LL:UL] * s[mw, LL - 1:UL - 1]) smax2 = np.max(-dtdx[LL - 1:UL - 1] * s[mw, LL - 1:UL - 1]) cfl = max(cfl, smax1, smax2) # If we are doing slope limiting we have more work to do if self.order == 2: # Initialize flux corrections f = np.zeros((num_eqn, grid.num_cells[0] + 2 * self.num_ghost)) # Apply Limiters to waves if (limiter > 0).any(): wave = tvd.limit(state.num_eqn, wave, s, limiter, dtdx) # Compute correction fluxes for second order q_{xx} terms dtdxave = 0.5 * (dtdx[LL - 1:UL - 1] + dtdx[LL:UL]) if self.fwave: for mw in xrange(wave.shape[1]): sabs = np.abs(s[mw, LL - 1:UL - 1]) om = 1.0 - sabs * dtdxave[:UL - LL] ssign = np.sign(s[mw, LL - 1:UL - 1]) for m in xrange(num_eqn): f[m, LL:UL] += 0.5 * ssign * om * wave[m, mw, LL - 1:UL - 1] else: for mw in xrange(wave.shape[1]): sabs = np.abs(s[mw, LL - 1:UL - 1]) om = 1.0 - sabs * dtdxave[:UL - LL] for m in xrange(num_eqn): f[m, LL:UL] += 0.5 * sabs * om * wave[m, mw, LL - 1:UL - 1] # Update q by differencing correction fluxes for m in xrange(num_eqn): q[m, LL:UL - 1] -= dtdx[LL:UL - 1] * (f[m, LL + 1:UL] - f[m, LL:UL - 1]) else: raise Exception( "Unrecognized kernel_language; choose 'Fortran' or 'Python'") self.cfl.update_global_max(cfl) state.set_q_from_qbc(num_ghost, self.qbc) if state.num_aux > 0: state.set_aux_from_auxbc(num_ghost, self.auxbc)