namespace mom_tracer_advect¶
Overview¶
This module contains the subroutines that advect tracers along coordinate surfaces. More…
namespace mom_tracer_advect { // global variables integer id_clock_advect; // global functions subroutine, public advect_tracer( h_end h_end, uhtr uhtr, vhtr vhtr, OBC OBC, dt dt, G G, GV GV, US US, CS CS, Reg Reg, h_prev_opt h_prev_opt, max_iter_in max_iter_in, x_first_in x_first_in, uhr_out uhr_out, vhr_out vhr_out, h_out h_out ); subroutine, public tracer_advect_init(Time Time, G G, param_file param_file, diag diag, CS CS); subroutine, public tracer_advect_end(CS CS); } // namespace mom_tracer_advect
Detailed Documentation¶
This module contains the subroutines that advect tracers along coordinate surfaces.
This program contains the subroutines that advect tracers horizontally (i.e. along layers).
section_mom_advect_intro¶
advect_tracer advects tracer concentrations using a combination of the modified flux advection scheme from Easter (Mon. Wea. Rev., 1993) with tracer distributions given by the monotonic modified van Leer scheme proposed by Lin et al. (Mon. Wea. Rev., 1994). This scheme conserves the total amount of tracer while avoiding spurious maxima and minima of the tracer concentration. If a higher order accuracy scheme is needed, suggest monotonic piecewise parabolic method, as described in Carpenter et al. (MWR, 1990).
advect_tracer has 4 arguments, described below. This subroutine determines the volume of a layer in a grid cell at the previous instance when the tracer concentration was changed, so it is essential that the volume fluxes should be correct. It is also important that the tracer advection occurs before each calculation of the diabatic forcing.
Global Functions¶
subroutine, public advect_tracer( h_end h_end, uhtr uhtr, vhtr vhtr, OBC OBC, dt dt, G G, GV GV, US US, CS CS, Reg Reg, h_prev_opt h_prev_opt, max_iter_in max_iter_in, x_first_in x_first_in, uhr_out uhr_out, vhr_out vhr_out, h_out h_out )
This routine time steps the tracer concentration using a monotonic, conservative, weakly diffusive scheme.
Parameters:
g |
ocean grid structure |
gv |
ocean vertical grid structure |
h_end |
layer thickness after advection [H ~> m or kg m-2] |
uhtr |
accumulated volume/mass flux through zonal face [H L2 ~> m3 or kg] |
vhtr |
accumulated volume/mass flux through merid face [H L2 ~> m3 or kg] |
obc |
specifies whether, where, and what OBCs are used |
dt |
time increment [s] |
us |
A dimensional unit scaling type |
cs |
control structure for module |
reg |
pointer to tracer registry |
h_prev_opt |
layer thickness before advection [H ~> m or kg m-2] |
max_iter_in |
The maximum number of iterations |
x_first_in |
If present, indicate whether to update first in the x- or y-direction. |
uhr_out |
accumulated volume/mass flux through zonal face |
vhr_out |
accumulated volume/mass flux through merid face |
h_out |
layer thickness before advection [H ~> m or kg m-2] |
subroutine, public tracer_advect_init( Time Time, G G, param_file param_file, diag diag, CS CS )
Initialize lateral tracer advection module.
Parameters:
time |
current model time |
g |
ocean grid structure |
param_file |
open file to parse for model parameters |
diag |
regulates diagnostic output |
cs |
module control structure |
subroutine, public tracer_advect_end(CS CS)
Close the tracer advection module.
Parameters:
cs |
module control structure |