|
1 | 1 | !=============================================================================================================================== |
2 | 2 | !**************** routines for horizontal tracer advection *********************** |
3 | | -module oce_adv_tra_hor_interfaces |
4 | | - interface |
| 3 | +module oce_adv_tra_hor_module |
| 4 | + use MOD_MESH |
| 5 | + use MOD_TRACER |
| 6 | + USE MOD_PARTIT |
| 7 | + USE MOD_PARSUP |
| 8 | + use g_comm_auto |
| 9 | + use O_PARAM, only: r_earth |
| 10 | + implicit none |
| 11 | + private |
| 12 | + public :: adv_tra_hor_upw1, adv_tra_hor_muscl, adv_tra_hor_mfct |
| 13 | +contains |
5 | 14 | ! (low order upwind) |
6 | 15 | ! returns flux given at edges which contributes with |
7 | 16 | ! plus sign into 1st. node and with the minus sign into the 2nd node |
8 | 17 | ! IF init_zero=.TRUE. : flux will be set to zero before computation |
9 | 18 | ! IF init_zero=.FALSE. : flux=flux-input flux |
10 | 19 | ! flux is not multiplied with dt |
11 | | - subroutine adv_tra_hor_upw1(vel, ttf, partit, mesh, flux, o_init_zero) |
12 | | - use MOD_MESH |
13 | | - use MOD_TRACER |
14 | | - USE MOD_PARTIT |
15 | | - USE MOD_PARSUP |
16 | | - type(t_partit),intent(in), target :: partit |
17 | | - type(t_mesh), intent(in), target :: mesh |
18 | | - real(kind=WP), intent(in) :: ttf( mesh%nl-1, partit%myDim_nod2D+partit%eDim_nod2D) |
19 | | - real(kind=WP), intent(in) :: vel(2, mesh%nl-1, partit%myDim_elem2D+partit%eDim_elem2D) |
20 | | - real(kind=WP), intent(inout) :: flux( mesh%nl-1, partit%myDim_edge2D) |
21 | | - logical, optional :: o_init_zero |
22 | | - end subroutine adv_tra_hor_upw1 |
23 | | -!=============================================================================== |
24 | | -! MUSCL |
25 | | -! returns flux given at edges which contributes with |
26 | | -! plus sign into 1st. node and with the minus sign into the 2nd node |
27 | | -! IF init_zero=.TRUE. : flux will be set to zero before computation |
28 | | -! IF init_zero=.FALSE. : flux=flux-input flux |
29 | | -! flux is not multiplied with dt |
30 | | - subroutine adv_tra_hor_muscl(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_grad, nboundary_lay, o_init_zero) |
31 | | - use MOD_MESH |
32 | | - USE MOD_PARTIT |
33 | | - USE MOD_PARSUP |
34 | | - type(t_partit),intent(in), target :: partit |
35 | | - type(t_mesh), intent(in), target :: mesh |
36 | | - real(kind=WP), intent(in) :: num_ord ! num_ord is the fraction of fourth-order contribution in the solution |
37 | | - real(kind=WP), intent(in) :: ttf( mesh%nl-1, partit%myDim_nod2D+partit%eDim_nod2D) |
38 | | - real(kind=WP), intent(in) :: vel(2, mesh%nl-1, partit%myDim_elem2D+partit%eDim_elem2D) |
39 | | - real(kind=WP), intent(inout) :: flux( mesh%nl-1, partit%myDim_edge2D) |
40 | | - integer, intent(in) :: nboundary_lay(partit%myDim_nod2D+partit%eDim_nod2D) |
41 | | - real(kind=WP), intent(in) :: edge_up_dn_grad(4, mesh%nl-1, partit%myDim_edge2D) |
42 | | - logical, optional :: o_init_zero |
43 | | - end subroutine adv_tra_hor_muscl |
44 | | -! a not stable version of MUSCL (reconstruction in the vicinity of bottom topography is not upwind) |
45 | | -! it runs with FCT option only |
46 | | - subroutine adv_tra_hor_mfct(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_grad, o_init_zero) |
47 | | - use MOD_MESH |
48 | | - USE MOD_PARTIT |
49 | | - USE MOD_PARSUP |
50 | | - type(t_partit),intent(inout), target :: partit |
51 | | - type(t_mesh), intent(in), target :: mesh |
52 | | - real(kind=WP), intent(in) :: num_ord ! num_ord is the fraction of fourth-order contribution in the solution |
53 | | - real(kind=WP), intent(in) :: ttf( mesh%nl-1, partit%myDim_nod2D+partit%eDim_nod2D) |
54 | | - real(kind=WP), intent(in) :: vel(2, mesh%nl-1, partit%myDim_elem2D+partit%eDim_elem2D) |
55 | | - real(kind=WP), intent(inout) :: flux( mesh%nl-1, partit%myDim_edge2D) |
56 | | - real(kind=WP), intent(in) :: edge_up_dn_grad(4, mesh%nl-1, partit%myDim_edge2D) |
57 | | - logical, optional :: o_init_zero |
58 | | - end subroutine adv_tra_hor_mfct |
59 | | - end interface |
60 | | -end module oce_adv_tra_hor_interfaces |
61 | | -! |
62 | | -! |
63 | | -!=============================================================================== |
64 | | -subroutine adv_tra_hor_upw1(vel, ttf, partit, mesh, flux, o_init_zero) |
65 | | - use MOD_MESH |
66 | | - use O_PARAM, only: r_earth |
67 | | - USE MOD_PARTIT |
68 | | - USE MOD_PARSUP |
69 | | - use g_comm_auto |
70 | | - implicit none |
| 20 | + subroutine adv_tra_hor_upw1(vel, ttf, partit, mesh, flux, o_init_zero) |
71 | 21 | type(t_partit),intent(in), target :: partit |
72 | 22 | type(t_mesh), intent(in), target :: mesh |
73 | 23 | real(kind=WP), intent(in) :: ttf( mesh%nl-1, partit%myDim_nod2D+partit%eDim_nod2D) |
@@ -255,17 +205,15 @@ subroutine adv_tra_hor_upw1(vel, ttf, partit, mesh, flux, o_init_zero) |
255 | 205 | #else |
256 | 206 | !$ACC END PARALLEL LOOP |
257 | 207 | #endif |
258 | | -end subroutine adv_tra_hor_upw1 |
259 | | -! |
260 | | -! |
| 208 | + end subroutine adv_tra_hor_upw1 |
261 | 209 | !=============================================================================== |
262 | | -subroutine adv_tra_hor_muscl(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_grad, nboundary_lay, o_init_zero) |
263 | | - use MOD_MESH |
264 | | - use MOD_TRACER |
265 | | - USE MOD_PARTIT |
266 | | - USE MOD_PARSUP |
267 | | - use g_comm_auto |
268 | | - implicit none |
| 210 | +! MUSCL |
| 211 | +! returns flux given at edges which contributes with |
| 212 | +! plus sign into 1st. node and with the minus sign into the 2nd node |
| 213 | +! IF init_zero=.TRUE. : flux will be set to zero before computation |
| 214 | +! IF init_zero=.FALSE. : flux=flux-input flux |
| 215 | +! flux is not multiplied with dt |
| 216 | + subroutine adv_tra_hor_muscl(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_grad, nboundary_lay, o_init_zero) |
269 | 217 | type(t_partit),intent(in), target :: partit |
270 | 218 | type(t_mesh), intent(in), target :: mesh |
271 | 219 | real(kind=WP), intent(in) :: num_ord ! num_ord is the fraction of fourth-order contribution in the solution |
@@ -540,17 +488,11 @@ subroutine adv_tra_hor_muscl(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_g |
540 | 488 | end do |
541 | 489 | !$OMP END DO |
542 | 490 | !$OMP END PARALLEL |
543 | | -end subroutine adv_tra_hor_muscl |
| 491 | + end subroutine adv_tra_hor_muscl |
544 | 492 | ! |
545 | 493 | ! |
546 | 494 | !=============================================================================== |
547 | | - subroutine adv_tra_hor_mfct(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_grad, o_init_zero) |
548 | | - use MOD_MESH |
549 | | - use MOD_TRACER |
550 | | - USE MOD_PARTIT |
551 | | - USE MOD_PARSUP |
552 | | - use g_comm_auto |
553 | | - implicit none |
| 495 | + subroutine adv_tra_hor_mfct(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_grad, o_init_zero) |
554 | 496 | type(t_partit),intent(inout), target :: partit |
555 | 497 | type(t_mesh), intent(in), target :: mesh |
556 | 498 | real(kind=WP), intent(in) :: num_ord ! num_ord is the fraction of fourth-order contribution in the solution |
@@ -832,4 +774,6 @@ subroutine adv_tra_hor_mfct(vel, ttf, partit, mesh, num_ord, flux, edge_up_dn_gr |
832 | 774 | #else |
833 | 775 | !$ACC END PARALLEL LOOP |
834 | 776 | #endif |
835 | | -end subroutine adv_tra_hor_mfct |
| 777 | + end subroutine adv_tra_hor_mfct |
| 778 | + |
| 779 | +end module oce_adv_tra_hor_module |
0 commit comments