slim
index
/builds/dg/dg/build/slim.py

Second-generation Louvain-la-Neuve Ice-ocean model - User interface
 
Contact: jonathan.lambrechts@uclouvain.be
webpage: www.climate.be/slim
 
SLIM solves the governing equations on adaptative unstructured meshes using the discontinuous Galerkin finite element method.
 
To simulate the model, run the script my_run_script.py:
    - On a single core: rundgpy my_run_script.py
    - On multiple cores: mpirun -np N rundgpy my_run_script.py (replace N by the number of cores)
 
Equations: 
    - 2D:
        * Shallow water equations (see :class:`~ShallowWater2d`)
        * Tracer equation (see :class:`~slim.ShallowWaterTracer2d`)
        * Particle Transport (see class TODO)
    - 3D:
        * Shallow water equations (see class TODO)
        * Tracer equation (see class TODO)
        
Domain
    The domain must be represented by a mesh using the .msh format (generated by gmsh or unref)
    Various forms of physical domain can be used:
    
    - 2D domain:
        * domain which are projected by mean of cartographic transformation (such as utm projection, stereographic projection)
        * sphere. In this case, the local poblems are written in a local orthonormal curvilinear system. All the local problems are assembled in the global discrete algebraic system (for further details, see Comblen et al., 2009, A finite element method for solving the shallow water equations on the sphere).
    - 3D domain:
        * domain which are projected by mean of cartographic transformation (such as utm projection, stereographic projection)
 
Data: 
    Various file format can be used: 
        * netcdf file (data.nc) 
            Those files have a special structure: data are given to each nodes at one or several time(s). The file 'slimPre' helps generating those files.
        * msh file (data_COMP_0.msh)
             
        * idx file (data.idx or data_COMP_0.idx) 
             
    Plane domain: In this case, data such as initial condition, wind stress, forcing ... must be expressed in the local basis (x,y) which is defined by the projection.
    
    Sphere: In this case, data such as initial condition, wind stress, forcing ... must be expressed in the global orthonormal basis (x,y,z) with z pointing upwards.
 
partition_id:
    This function returns the number of the current partition as a string
 
partition_nb:
    This function returns the number of partitions used as a string

 
Modules
       
dgpy
dgpy.scripts.slim_private

 
Classes
       
builtins.object
Domain
Loop
ShallowWater2d
ShallowWater2dWD
ShallowWaterTracer2d

 
class Domain(builtins.object)
    Create the numerical domain
 
  Methods defined here:
__init__(self, mesh, bathy, g=None, density=1025, solve_on_sphere=False, order=1)
keyword arguments:
 
* mesh            
    path to the mesh file (.msh format). The partitioned mesh will automatically be loaded in multiprocessing
* bathy           
    bathymetric file [in meters, positive] (.msh, .idx or .nc format)
* g               
    map of the mean gravitational acceleration (.msh, .idx or .nc format) (default: 9.81 m/s^2)
* density         
    density of the liquid (default: 1025 kg/m^3)
* solve_on_sphere 
    boolean stating if you want to solve the shallow water equation on a sphere or on any 2d manifold instead of solving it on a plane (default: False)
* order
    Polynomial order of the discretization

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class Loop(builtins.object)
    Temporal solver
 
  Methods defined here:
__init__(self, maximum_time_step=3600, export_time=3600, path='./output', evaluateEveryTimeStep=False)
keyword arguments:
 
* maximum_time_step 
    maximum value for the time step allowed by the user in order to represent the physical processes. [in seconds] (default: 3600 s) 
 
    The real time step is computed from numerical restriction of the temporal scheme (if implicit, this value is used as time step).
* export_time       
    time step for results exportation [in seconds] (default: 3600 s)
* path 
    path to the output directory [string] (default: ./output)
* evaluateEveryTimeStep
    flag to define if quantities (export_value_at_points and compute_mass) are evaluated every time step (True) or only at export time (False) (default: False)
add_equation(self, equation)
Add equation to the temporal solver
 
keywords arguments:
 
* equation 
    equation to add to the temporal solver
compute_dt(self)
compute the optimal time step
evaluatePoints(self)
export_on_structured_grid(self, equation_list, x_min, x_max, y_min, y_max, d_x, d_y, file_name='output.nc')
Export the solution of equation_name on a structured grid in (longitude, latitude) in netcdf format
    
keywords arguments:
 
* equation_list
    list of equation(s) whose the solution will be exported on a structured grid 
* x_min
    minimum x
* x_max
    maximum x
* y_min
    minimum y
* y_max
    maximum y
* d_x
    step between x
* d_y
    step between y
* file_name     
    name of the output file (default: output.nc)
export_solutions(self, temporal_scheme)
Export the solutions
 
keywords arguments:
 
* temporal_scheme 
    scheme for the temporal integration (default: None)
                       
    * "explicit" 
        explicit Euler scheme
    * "implicit" 
        implicit order 2 Runge-Kutta
get_time_step(self)
Get the time step of the equation
iterate(self, temporal_scheme, offline, dt=None)
Time integration loop with exports
 
keywords arguments:
 
* temporal_scheme 
    scheme for the temporal integration (default: None)
    
    * "explicit" 
        explicit Euler scheme
    * "implicit" 
        implicit order 2 Runge-Kutta
* offline         
    boolean stating if you want to solve the shallow water tracer equation offline
* dt              
    you can set a time step different from the automatically computed time step
restart(self, index, time)
Restart run from a time-step already resolved
keywords arguments:
 
* index 
    index used for restart
* time
    time corresponding to this index
run(self)
Time integration loop with exports
set_time_step_offline(self, dt, periodic=False, index_start=1, n_index_per_period=-1)
Set the time step for the tracer equation if the hydrodynamics has already been computed
 
keywords arguments:
 
* dt 
    time step for the tracer equation with offline hydrodynamics
* periodic 
    bool  stating if the offline solution is periodic (default: False)       
* index_start 
    index of the solution exported at each time step starting the period (defined if periodic = True)
* n_index_per_period
    number of file per period for offline method (defined if periodic = True)

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class ShallowWater2d(builtins.object)
    Create the shallow water equation with various options
 
  Methods defined here:
__init__(self, domain, temporal_scheme, linear_equation=False, wetting_drying=None, export_every_sub_time_step=False, initial_time=None, final_time=None)
keyword arguments:
 
* domain                     
    object of class Domain
* temporal_scheme 
    scheme for the temporal integration
    
    * "explicit" 
        explicit Euler scheme
    * "implicit" 
        implicit order 2 Runge-Kutta
    * "multirate"
        multirate scheme (see Seny et al. 2012, 2014)
        
        Warning: This temporal scheme is not yet implemented for tracers
* linear_equation            
    boolean stating if you want tos solve the linear version of the  shallow water equation (Default: False)
* wetting_drying             
    value for the wetting & drying algorithm (see Karna et al. 2011) (Default: None)
* export_every_sub_time_step 
    boolean stating if you want to export the hydrodynamics every sub time step for further use of offline equation (Default: False)
* initial_time               
    initial time of the simulation (specific to each equation) (Default: 0)
    
    Format: 
 
    - float [in seconds]
    - date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-02-15 10:05:00")
* final_time
    final time of the simulation (default: 0)
    
    Format:
        
    - float [in seconds]
    - date string in format "year-month-day hours:minutes:seconds" (e.g. "2007-06-25 22:42:17")
compute_mass(self, output_name)
compute the total mass of the domain
 
keyword arguments:
 
* output_name 
    name of the output file
export_value_at_point(self, output_name, x, y, z=None)
set the coordinates where a time serie of the solution is exported
 
keyword arguments:
 
* output_name 
    name of the output file
* coord_x     
    vector containing the x coordinates where the solution is evaluated
* coord_y     
    vector containing the y coordinates where the solution is evaluated
* coord_z     
    vector containing the z coordinates where the solution is evaluated (Default: None, replaced by vector of 0)
set_boundary_coast(self, physical_tag, slip=False)
Boundary condition: normal speed set to zero and tangential speed set to zero or free (no slip or free slip condition)
 
keyword arguments:
 
* physical_tag 
    tag of the part of the boundary on which this boundary condition is applied.
* slip         
    flag wheter applying slip (True) or no slip (False) condition at the boundary (Default: False)
set_boundary_open(self, physical_tag, discharge=None, sse=None, ux=None, uy=None, uz=None, transport_flux=False, ramp_period=0)
Open boundary condition: one can impose the discharge or the sea surface elevation and/or the velocity (ux and uy have to be defined at the same time but they can be zero).
 
keyword arguments:
 
* physical_tag  
    tag of the part of the boundary on which this boundary condition is applied.
* discharge    
    netcdf or .msh file containing the discharge which will be applied at the boundary [in m^3s^-1] (default: None)
* sse            
    netcdf or .msh file containing the sea surface elevation which will be applied at the boundary [in m] (default: None, i.e. the sea surface elevation computed by the model is set to this boundary)
* ux             
    netcdf or .msh file containing the velocity along the x axis of the local basis  or the velocity along the x axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uy             
    netcdf or .msh file containing the velocity along the y axis of the local basis  or the velocity along the y axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uz             
    netcdf or .msh file containing the velocity along the z axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* transport_flux 
    set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False) 
* ramp_period  
    period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
set_coriolis(self, coriolis)
Add a coriolis term (2*Omega*cos(latitude)) in the shallow water equation 
 
 keyword argument: 
 
* coriolis 
     netcdf or .msh file containing the coriolis term for the whole domain.
set_dissipation(self, mode, coefficient=None)
Dissipation due to bottom friction: linear, basic quadratic scheme or Chezy-Manning scheme
 
keyword arguments:
 
* mode
    scheme used for the bottom friction
                              
    * "linear"    
        linear dissipation tau_b/(rho*H) = gamma*u (default: gamma = 1e-6 s^-1)
    * "quadratic" 
        basic quadratic dissipation tau_b/(rho*H) = (C_d/H)*|u|*u (default: C_d = 2.5e-3)
    * "manning"   
        Chezy-Manning-Strickler formulation tau_b/(rho*H) = (n*n*g/(H^4/3))*|u|*u (default: n = 0.03 sm^(-1/3))                
    * "slope"   
        Bottom slope formulation tau_b/(rho*H) = C_s*||u.grad(h)||*u/H (default: C_s = 10)
* coefficient 
    netcdf or .msh file containing the coefficient value (gamma, C_d, n or C_s depending on the mode) over the whole domain (default: None, the default value of the mode will be used)
set_forcing(self, forcing_x, forcing_y, forcing_z=None)
Add a momentum source term derived from a gravitational potential in the shallow water equation.
 
keyword arguments: 
 
* forcing_x 
    netcdf or .msh file containing the forcing term along the x-axis in the local basis 
    
    (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].  
* forcing_y 
    netcdf or .msh file containing the forcing term along the y-axis in the local basis 
    
    (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_z 
    netcdf or .msh file containing the forcing term along the z-axis in the local basis 
    
    (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
set_initial_condition(self, sse, ux, uy, uz=None)
Initial conditions
 
keyword arguments:
 
* sse 
    sea surface elevation [in m] (.msh or .nc format) 
* ux  
    velocity along the x-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) 
* uy  
    velocity along the y-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) 
* uz  
    velocity along the z-axis in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) (default: None)
set_nudging(self, coefficient, external_sse, external_ux, external_uy, external_uz=None, ramp_period=0, transport_flux=False)
Add a nudging term as a source term in the shallow water equation.
 
keyword arguments: 
 
* coefficient    
    netcdf or .msh file containing the nudging coefficient for the whole domain [in s^-1]
* external_sse   
    netcdf or .msh file containing the external sea surface elevation for the whole domain [in m]
* external_ux    
    netcdf or .msh file containing the external velocity along the x-axis expressed in the local basis (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in m/s]
* external_uy    
    netcdf or .msh file containing the external velocity along the y-axis expressed in the local basis (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in m/s]
* external_uz    
    netcdf or .msh file containing the external velocity along the z-axis expressed in the global basis (x,y,z) for the whole domain [in m/s]
* ramp_period    
    period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp) 
* transport_flux 
    set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False)
set_viscosity(self, mode, smagorinsky_coefficient=0.1, constant_viscosity=1e-06)
Turbulent viscosity : Smagorinsky scheme or constant value
 
keyword arguments:
 
* mode       
    type of turbulent viscosity
    
    * "smagorinsky" 
        Smagorinsky scheme (see Smagorinsky 1962)
    * "constant"    
        constant value 
* smagorinsky_coefficient 
    coefficient value for Smagorinsky scheme (default: 0.1)
* constant_viscosity      
    constant value [in m^2/s] (default: 1e-6 m^2/s)
set_wind_stress(self, mode, wind_x, wind_y, density_air=-1, wind_z=None)
Add a wind stress term in the shallow water equation.
 
keyword arguments:
 
* mode                  
    type of wind forcing
    
    * "speed"   
        wind speed given in m s^-1. wind_stress will be computed with Smith-Banke formula.
          stress = C_D(speed)*density_air*speed
        Setting the air density is mandatory 
    * "stress" 
        wind stress given in kg m^-1 s^-1. (air density will not be used)
* wind_x 
    netcdf or .msh file containing the wind term along the x-axis in the local basis 
    
    (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.  
* wind_y 
    netcdf or .msh file containing the wind term along the y-axis in the local basis 
    
    (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.  
* density_air  
    density of the ambiant air (only necessary for "speed" mode) (default: -1)
* wind_z 
    netcdf or .msh file containing the wind term along the z-axis in the global basis (x,y,z) for the whole domain (only for equation on a sphere). (default: None)

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class ShallowWater2dWD(builtins.object)
    Create the shallow water equation with various options
 
  Methods defined here:
__init__(self, domain, temporal_scheme, export_every_sub_time_step=False, initial_time=None, final_time=None, hlimt=1, hlimd=0.1, linDrag=0.1, manning=0.02)
keyword arguments:
 
* domain                     
    object of class Domain
* temporal_scheme 
    scheme for the temporal integration
    
    * "explicit" 
        explicit Euler scheme
    * "implicit" 
        implicit order 2 Runge-Kutta
    * "multirate"
        multirate scheme (see Seny et al. 2012, 2014)
        
        Warning: This temporal scheme is not yet implemented for tracers
* export_every_sub_time_step 
    boolean stating if you want to export the hydrodynamics every sub time step for further use of offline equation (Default: False)
* initial_time               
    initial time of the simulation (specific to each equation) (Default: 0)
    
    Format: 
 
    - float [in seconds]
    - date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-02-15 10:05:00")
* final_time
    final time of the simulation (default: 0)
    
    Format:
        
    - float [in seconds]
    - date string in format "year-month-day hours:minutes:seconds" (e.g. "2007-06-25 22:42:17")
compute_mass(self, output_name)
compute the total mass of the domain
 
keyword arguments:
 
* output_name 
    name of the output file
export_value_at_point(self, output_name, x, y, z=None)
set the coordinates where a time serie of the solution is exported
 
keyword arguments:
 
* output_name 
    name of the output file
* coord_x     
    vector containing the x coordinates where the solution is evaluated
* coord_y     
    vector containing the y coordinates where the solution is evaluated
* coord_z     
    vector containing the z coordinates where the solution is evaluated (Default: None, replaced by vector of 0)
set_boundary_coast(self, physical_tag)
Boundary condition: normal speed set to zero and tangential speed set to zero or free (no slip or free slip condition)
 
keyword arguments:
 
* physical_tag 
    tag of the part of the boundary on which this boundary condition is applied.
set_boundary_open(self, physical_tag, discharge=None, sse=None, ux=None, uy=None, uz=None, transport_flux=False, ramp_period=0)
Open boundary condition: one can impose the discharge or the sea surface elevation and/or the velocity (ux and uy have to be defined at the same time but they can be zero).
 
keyword arguments:
 
* physical_tag  
    tag of the part of the boundary on which this boundary condition is applied.
* discharge    
    netcdf or .msh file containing the discharge which will be applied at the boundary [in m^3s^-1] (default: None)
* sse            
    netcdf or .msh file containing the sea surface elevation which will be applied at the boundary [in m] (default: None, i.e. the sea surface elevation computed by the model is set to this boundary)
* ux             
    netcdf or .msh file containing the velocity along the x axis of the local basis  or the velocity along the x axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uy             
    netcdf or .msh file containing the velocity along the y axis of the local basis  or the velocity along the y axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uz             
    netcdf or .msh file containing the velocity along the z axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* transport_flux 
    set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False) 
* ramp_period  
    period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
set_coriolis(self, coriolis)
Add a coriolis term (2*Omega*cos(latitude)) in the shallow water equation 
 
 keyword argument: 
 
* coriolis 
     netcdf or .msh file containing the coriolis term for the whole domain.
set_initial_condition(self, sse, ux, uy, uz=None)
Initial conditions
 
keyword arguments:
 
* sse 
    sea surface elevation [in m] (.msh or .nc format) 
* ux  
    velocity along the x-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) 
* uy  
    velocity along the y-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) 
* uz  
    velocity along the z-axis in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) (default: None)

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class ShallowWaterTracer2d(builtins.object)
    Create the shallow water equation for a tracer
 
  Methods defined here:
__init__(self, domain, temporal_scheme, hydro_sol, linear_equation=False, name='Tracer', offline=False, wetting_drying=None, initial_time=None, final_time=None)
keyword arguments:
 
* domain                     
    object of class Domain
* temporal_scheme 
    scheme for the temporal integration
                       
    * "explicit" 
        explicit Euler scheme
    * "implicit" 
        implicit order 2 Runge-Kutta
* hydro_sol
    hydrodynamics solution
    
    Format:
    
    - equation : ShallowWater2d object (online mode)
    - string : path to the folder containing the hydrodynamics folder 'offline_sw2d' (offline mode).
* linear_equation            
    boolean stating linear shallow water tracer equation (default: False)
* name                       
    name of the Tracer (default: "Tracer")
* offline                    
    resolve shallow water equation (False) or use file to provide hydrodynamics (True) (default: False)
* wetting_drying             
    value for the wetting & drying algorithm. It has to be the same as for the shallow water equation (see Karna et al. 2011) (default: None)
* initial_time               
    initial time of the simulation (default: 0)
  
    Format:
    
    - float [in seconds]
    - date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-06-20 10:05:00")
* final_time                 
    final time of the simulation (default: 0)
    
    Format:
    
    - float [in seconds]
    - date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-06-20 22:42:17")
* is_sediment                    
    use  shallow water tracer law to predict the sediment transport (default: False)
compute_mass(self, output_name)
Compute the integral of the tracer over the whole domain
 
keyword arguments:
 
* output_name 
    name of the output file
compute_sediment(self, dt)
Compute the sediment flux
 
keyword arguments:
 
* dt 
    time step
set_boundary_coast(self, physical_tag)
Boundary condition: no flux condition at the boundary
 
keyword argument:
 
* physical_tag 
    tag of the part of the boundary on which this boundary condition is applied.
set_boundary_flux(self, physical_tag, flux)
Boundary condition for a tracer flux: one can impose the flux throw the boundary condition, in [c s^-1], (with the tracer units [c]).
 
keyword arguments:
 
* physical_tag             
    tag of the part of the boundary on which this boundary condition is applied.
* flux            
    netcdf or .msh file containing the flux which will be applied at the boundary in [c s^-1]
set_boundary_open(self, physical_tag, concentration)
Open boundary condition: one can impose the external concentration
 
keyword arguments:
 
* physical_tag             
    tag of the part of the boundary on which this boundary condition is applied.
* concentration            
    netcdf or .msh file containing the concentration which will be applied at the boundary [in c]
set_diffusivity(self, mode, okubo_coefficient=0.03, constant_diffusivity=1e-06)
Diffusivity
 
keyword arguments:
 
* mode                  
    type of diffusivity
    
    * "okubo"   
        Okubo scheme 
    * "constant" 
        constant value 
* okubo_coefficient   
    coefficient for Okubo scheme [in m^0.85/s] (default: 0.03 m^0.85/s)
* constant_diffusivity 
    constant value [in m^2/s] (default: 1e-6 m^2/s)
set_initial_condition(self, initial_condition)
Initial conditions
 
keyword argument:
 
* tracer 
    tracer initial condition [-] (.msh or .nc format) (default: 0 )
set_sediment(self, initial_bottom_concentration, windU, windV)
Define wind velocity for sediment transport.
 
keyword argument: 
 
*initial_bottom_concentration
    netcdf or .msh file containing the initial sediment concentration at the bottom of the domain for the whole domain [in c/m^2].
* windU
    netcdf or .msh file containing the surface wind velocity along the x-axis in the local basis [in m*s^-1].
* windV
    netcdf or .msh file containing the surface wind velocity along the y-axis in the local basis [in m*s^-1].
set_source(self, tracer_source)
Add a tracer source term in the tracer equation.
 
keyword argument: 
 
* tracer_source 
    netcdf or .msh file containing the source term for the whole domain [in c*s^-1].

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
Functions
       
partition_id()
Return the number of the current partition as a string
partition_nb()
Return the number of partitions used as a string