mirror of
				https://github.com/smyalygames/FiniteVolumeGPU.git
				synced 2025-11-04 04:59:49 +01:00 
			
		
		
		
	
		
			
				
	
	
		
			160 lines
		
	
	
		
			5.6 KiB
		
	
	
	
		
			Python
		
	
	
	
	
	
			
		
		
	
	
			160 lines
		
	
	
		
			5.6 KiB
		
	
	
	
		
			Python
		
	
	
	
	
	
# -*- coding: utf-8 -*-
 | 
						|
 | 
						|
"""
 | 
						|
This python module implements the Kurganov-Petrova numerical scheme 
 | 
						|
for the shallow water equations, described in 
 | 
						|
A. Kurganov & Guergana Petrova
 | 
						|
A Second-Order Well-Balanced Positivity Preserving Central-Upwind
 | 
						|
Scheme for the Saint-Venant System Communications in Mathematical
 | 
						|
Sciences, 5 (2007), 133-160. 
 | 
						|
 | 
						|
Copyright (C) 2016  SINTEF ICT
 | 
						|
 | 
						|
This program is free software: you can redistribute it and/or modify
 | 
						|
it under the terms of the GNU General Public License as published by
 | 
						|
the Free Software Foundation, either version 3 of the License, or
 | 
						|
(at your option) any later version.
 | 
						|
 | 
						|
This program is distributed in the hope that it will be useful,
 | 
						|
but WITHOUT ANY WARRANTY; without even the implied warranty of
 | 
						|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 | 
						|
GNU General Public License for more details.
 | 
						|
 | 
						|
You should have received a copy of the GNU General Public License
 | 
						|
along with this program.  If not, see <http://www.gnu.org/licenses/>.
 | 
						|
"""
 | 
						|
 | 
						|
#Import packages we need
 | 
						|
import numpy as np
 | 
						|
 | 
						|
import pycuda.compiler as cuda_compiler
 | 
						|
import pycuda.gpuarray
 | 
						|
import pycuda.driver as cuda
 | 
						|
 | 
						|
from SWESimulators import Common
 | 
						|
 | 
						|
 | 
						|
 | 
						|
 | 
						|
 | 
						|
"""
 | 
						|
Class that solves the SW equations using the dimentionally split KP07 scheme
 | 
						|
"""
 | 
						|
class KP07_dimsplit:
 | 
						|
 | 
						|
    """
 | 
						|
    Initialization routine
 | 
						|
    h0: Water depth incl ghost cells, (nx+1)*(ny+1) cells
 | 
						|
    hu0: Initial momentum along x-axis incl ghost cells, (nx+1)*(ny+1) cells
 | 
						|
    hv0: Initial momentum along y-axis incl ghost cells, (nx+1)*(ny+1) cells
 | 
						|
    nx: Number of cells along x-axis
 | 
						|
    ny: Number of cells along y-axis
 | 
						|
    dx: Grid cell spacing along x-axis (20 000 m)
 | 
						|
    dy: Grid cell spacing along y-axis (20 000 m)
 | 
						|
    dt: Size of each timestep (90 s)
 | 
						|
    g: Gravitational accelleration (9.81 m/s^2)
 | 
						|
    """
 | 
						|
    def __init__(self, \
 | 
						|
                 context, \
 | 
						|
                 h0, hu0, hv0, \
 | 
						|
                 nx, ny, \
 | 
						|
                 dx, dy, dt, \
 | 
						|
                 g, \
 | 
						|
                 theta=1.3, \
 | 
						|
                 block_width=16, block_height=16):
 | 
						|
        #Create a CUDA stream
 | 
						|
        self.stream = cuda.Stream()
 | 
						|
 | 
						|
        #Get kernels
 | 
						|
        self.kp07_dimsplit_module = context.get_kernel("KP07_dimsplit_kernel.cu", block_width, block_height)
 | 
						|
        self.kp07_dimsplit_kernel = self.kp07_dimsplit_module.get_function("KP07DimsplitKernel")
 | 
						|
        self.kp07_dimsplit_kernel.prepare("iifffffiPiPiPiPiPiPi")
 | 
						|
        
 | 
						|
        #Create data by uploading to device
 | 
						|
        ghost_cells_x = 2
 | 
						|
        ghost_cells_y = 2
 | 
						|
        self.data = Common.SWEDataArakawaA(self.stream, nx, ny, ghost_cells_x, ghost_cells_y, h0, hu0, hv0)
 | 
						|
        
 | 
						|
        #Save input parameters
 | 
						|
        #Notice that we need to specify them in the correct dataformat for the
 | 
						|
        #OpenCL kernel
 | 
						|
        self.nx = np.int32(nx)
 | 
						|
        self.ny = np.int32(ny)
 | 
						|
        self.dx = np.float32(dx)
 | 
						|
        self.dy = np.float32(dy)
 | 
						|
        self.dt = np.float32(dt)
 | 
						|
        self.g = np.float32(g)
 | 
						|
        self.theta = np.float32(theta)
 | 
						|
        
 | 
						|
        #Initialize time
 | 
						|
        self.t = np.float32(0.0)
 | 
						|
        
 | 
						|
        #Compute kernel launch parameters
 | 
						|
        self.local_size = (block_width, block_height, 1) 
 | 
						|
        self.global_size = ( \
 | 
						|
                       int(np.ceil(self.nx / float(self.local_size[0]))), \
 | 
						|
                       int(np.ceil(self.ny / float(self.local_size[1]))) \
 | 
						|
                      ) 
 | 
						|
    
 | 
						|
    
 | 
						|
    def __str__(self):
 | 
						|
        return "Kurganov-Petrova 2007 dimensionally split"
 | 
						|
    
 | 
						|
    
 | 
						|
    """
 | 
						|
    Function which steps n timesteps
 | 
						|
    """
 | 
						|
    def step(self, t_end=0.0):
 | 
						|
        n = int(t_end / (2.0*self.dt) + 1)
 | 
						|
        
 | 
						|
        for i in range(0, n): 
 | 
						|
            #Dimensional splitting: second order accurate for every other timestep,
 | 
						|
            #thus run two timesteps in a go
 | 
						|
            
 | 
						|
            #Compute timestep
 | 
						|
            local_dt = np.float32(0.5*min(2*self.dt, t_end-2*i*self.dt))
 | 
						|
            if (local_dt <= 0.0):
 | 
						|
                break
 | 
						|
                
 | 
						|
            #Along X, then Y
 | 
						|
            self.kp07_dimsplit_kernel.prepared_async_call(self.global_size, self.local_size, self.stream, \
 | 
						|
                    self.nx, self.ny, \
 | 
						|
                    self.dx, self.dy, local_dt, \
 | 
						|
                    self.g, \
 | 
						|
                    self.theta, \
 | 
						|
                    np.int32(0), \
 | 
						|
                    self.data.h0.data.gpudata,  self.data.h0.pitch, \
 | 
						|
                    self.data.hu0.data.gpudata, self.data.hu0.pitch, \
 | 
						|
                    self.data.hv0.data.gpudata, self.data.hv0.pitch, \
 | 
						|
                    self.data.h1.data.gpudata,  self.data.h1.pitch, \
 | 
						|
                    self.data.hu1.data.gpudata, self.data.hu1.pitch, \
 | 
						|
                    self.data.hv1.data.gpudata, self.data.hv1.pitch)
 | 
						|
            self.data.swap()
 | 
						|
            
 | 
						|
            #Along Y, then X
 | 
						|
            self.kp07_dimsplit_kernel.prepared_async_call(self.global_size, self.local_size, self.stream, \
 | 
						|
                    self.nx, self.ny, \
 | 
						|
                    self.dx, self.dy, local_dt, \
 | 
						|
                    self.g, \
 | 
						|
                    self.theta, \
 | 
						|
                    np.int32(1), \
 | 
						|
                    self.data.h0.data.gpudata,  self.data.h0.pitch, \
 | 
						|
                    self.data.hu0.data.gpudata, self.data.hu0.pitch, \
 | 
						|
                    self.data.hv0.data.gpudata, self.data.hv0.pitch, \
 | 
						|
                    self.data.h1.data.gpudata,  self.data.h1.pitch, \
 | 
						|
                    self.data.hu1.data.gpudata, self.data.hu1.pitch, \
 | 
						|
                    self.data.hv1.data.gpudata, self.data.hv1.pitch)
 | 
						|
            self.data.swap()
 | 
						|
            
 | 
						|
            self.t += 2.0*local_dt
 | 
						|
            
 | 
						|
        
 | 
						|
        return self.t
 | 
						|
    
 | 
						|
    
 | 
						|
    
 | 
						|
    
 | 
						|
    def download(self):
 | 
						|
        return self.data.download(self.stream)
 | 
						|
 |