source: codes/icosagcm/trunk/src/euler_scheme.f90 @ 360

Last change on this file since 360 was 360, checked in by dubos, 9 years ago

Towards HEVI time stepping

File size: 4.7 KB
Line 
1MODULE euler_scheme_mod
2  USE field_mod
3  IMPLICIT NONE
4  PRIVATE
5
6  TYPE(t_field),POINTER,SAVE,PUBLIC :: f_geopot(:), f_q(:), &
7       f_rhodz(:), f_mass(:), f_massm1(:), f_massm2(:), f_dmass(:), &
8       f_phis(:), f_ps(:),f_psm1(:), f_psm2(:), f_dps(:), &
9       f_u(:),f_um1(:),f_um2(:), f_du(:), & ! previous state
10       f_theta_rhodz(:),f_theta_rhodzm1(:),f_theta_rhodzm2(:), f_dtheta_rhodz(:), & ! trends
11       f_hflux(:), f_wflux(:), f_hfluxt(:), f_wfluxt(:), & ! accumulated mass fluxes
12                                ! slow/fast trend arrays for HEVI time scheme
13       f_dmass_slow(:,:), f_dps_slow(:,:), f_dtheta_rhodz_slow(:,:), &
14       f_du_slow(:,:), f_du_fast(:,:)
15
16  INTEGER, PARAMETER, PUBLIC :: explicit=1, hevi=2, euler=1, rk4=2, mlf=3, rk25=4, ark23=6, ark33=7
17
18  INTEGER,SAVE, PUBLIC :: nb_stage, matsuno_period, scheme, scheme_family
19  !$OMP THREADPRIVATE(nb_stage, matsuno_period, scheme, scheme_family)
20
21
22  PUBLIC :: euler_scheme, accumulate_fluxes
23CONTAINS
24
25  SUBROUTINE Euler_scheme(with_dps,fluxt_zero)
26    USE icosa
27    USE disvert_mod
28    USE omp_para
29    USE trace
30    LOGICAL :: with_dps
31    LOGICAL, OPTIONAL :: fluxt_zero(ndomain) ! set to .TRUE. to start accumulating fluxes in time
32    REAL(rstd),POINTER :: ps(:), dps(:)
33    REAL(rstd),POINTER :: u(:,:) , du(:,:)
34    REAL(rstd),POINTER :: mass(:,:), dmass(:,:)
35    REAL(rstd),POINTER :: theta_rhodz(:,:), dtheta_rhodz(:,:)
36    REAL(rstd),POINTER :: hflux(:,:),wflux(:,:),hfluxt(:,:),wfluxt(:,:)
37    INTEGER :: ind
38    INTEGER :: i,j,ij,l
39    CALL trace_start("Euler_scheme") 
40
41    DO ind=1,ndomain
42       IF (.NOT. assigned_domain(ind)) CYCLE
43       CALL swap_dimensions(ind)
44       CALL swap_geometry(ind)
45
46       IF(with_dps) THEN ! update ps/mass
47          IF(caldyn_eta==eta_mass) THEN ! update ps
48             ps=f_ps(ind) ; dps=f_dps(ind) ;             
49             IF (is_omp_first_level) THEN
50                !$SIMD
51                DO ij=ij_begin,ij_end
52                   ps(ij)=ps(ij)+dt*dps(ij)
53                ENDDO
54             ENDIF
55          ELSE ! update mass
56             mass=f_mass(ind) ; dmass=f_dmass(ind) ;             
57             DO l=ll_begin,ll_end
58                !$SIMD
59                DO ij=ij_begin,ij_end
60                   mass(ij,l)=mass(ij,l)+dt*dmass(ij,l)
61                ENDDO
62             END DO
63          END IF
64
65          hflux=f_hflux(ind);     hfluxt=f_hfluxt(ind)
66          wflux=f_wflux(ind);     wfluxt=f_wfluxt(ind)
67          CALL accumulate_fluxes(hflux,wflux,hfluxt,wfluxt,dt,fluxt_zero(ind))
68       END IF ! update ps/mass
69
70       u=f_u(ind) ; theta_rhodz=f_theta_rhodz(ind)
71       du=f_du(ind) ; dtheta_rhodz=f_dtheta_rhodz(ind)
72
73       DO l=ll_begin,ll_end
74          !$SIMD
75          DO ij=ij_begin,ij_end
76             u(ij+u_right,l)=u(ij+u_right,l)+dt*du(ij+u_right,l)
77             u(ij+u_lup,l)=u(ij+u_lup,l)+dt*du(ij+u_lup,l)
78             u(ij+u_ldown,l)=u(ij+u_ldown,l)+dt*du(ij+u_ldown,l)
79             theta_rhodz(ij,l)=theta_rhodz(ij,l)+dt*dtheta_rhodz(ij,l)
80          ENDDO
81       ENDDO
82    ENDDO
83
84    CALL trace_end("Euler_scheme") 
85
86  END SUBROUTINE Euler_scheme
87
88  SUBROUTINE accumulate_fluxes(hflux,wflux, hfluxt,wfluxt, tau,fluxt_zero)
89    USE dimensions
90    USE omp_para
91    USE disvert_mod
92    IMPLICIT NONE
93    REAL(rstd), INTENT(IN)    :: hflux(3*iim*jjm,llm), wflux(iim*jjm,llm+1)
94    REAL(rstd), INTENT(INOUT) :: hfluxt(3*iim*jjm,llm), wfluxt(iim*jjm,llm+1)
95    REAL(rstd), INTENT(IN) :: tau
96    LOGICAL, INTENT(INOUT) :: fluxt_zero
97    INTEGER :: l,i,j,ij
98
99    IF(fluxt_zero) THEN
100
101       fluxt_zero=.FALSE.
102
103       DO l=ll_begin,ll_end
104          !$SIMD
105          DO ij=ij_begin_ext,ij_end_ext
106             hfluxt(ij+u_right,l) = tau*hflux(ij+u_right,l)
107             hfluxt(ij+u_lup,l) = tau*hflux(ij+u_lup,l)
108             hfluxt(ij+u_ldown,l) = tau*hflux(ij+u_ldown,l)
109          ENDDO
110       ENDDO
111
112       IF(caldyn_eta==eta_mass) THEN ! no need for vertical fluxes if eta_lag
113          DO l=ll_begin,ll_endp1
114             !$SIMD
115             DO ij=ij_begin,ij_end
116                wfluxt(ij,l) = tau*wflux(ij,l)
117             ENDDO
118          ENDDO
119       END IF
120
121    ELSE
122
123       DO l=ll_begin,ll_end
124          !$SIMD
125          DO ij=ij_begin_ext,ij_end_ext
126             hfluxt(ij+u_right,l) = hfluxt(ij+u_right,l)+tau*hflux(ij+u_right,l)
127             hfluxt(ij+u_lup,l) = hfluxt(ij+u_lup,l)+tau*hflux(ij+u_lup,l)
128             hfluxt(ij+u_ldown,l) = hfluxt(ij+u_ldown,l)+tau*hflux(ij+u_ldown,l)
129          ENDDO
130       ENDDO
131
132       IF(caldyn_eta==eta_mass) THEN ! no need for vertical fluxes if eta_lag
133          DO l=ll_begin,ll_endp1
134             !$SIMD
135             DO ij=ij_begin,ij_end
136                wfluxt(ij,l) = wfluxt(ij,l)+tau*wflux(ij,l)
137             ENDDO
138          ENDDO
139       END IF
140
141    END IF
142
143  END SUBROUTINE accumulate_fluxes
144
145END MODULE euler_scheme_mod
Note: See TracBrowser for help on using the repository browser.