source: codes/icosagcm/devel/src/dynamics/caldyn_adv.f90 @ 1050

Last change on this file since 1050 was 1050, checked in by dubos, 4 years ago

devel : cleanup USE ...

File size: 5.7 KB
Line 
1MODULE caldyn_adv_mod
2  USE icosa
3  IMPLICIT NONE
4  PRIVATE
5  PUBLIC :: init_caldyn, caldyn
6
7CONTAINS
8
9  SUBROUTINE init_caldyn
10  END SUBROUTINE init_caldyn
11
12  SUBROUTINE check_mass_conservation(f_ps,f_dps)
13    USE icosa
14    TYPE(t_field),POINTER :: f_ps(:)
15    TYPE(t_field),POINTER :: f_dps(:)
16    REAL(rstd),POINTER :: ps(:)
17    REAL(rstd),POINTER :: dps(:)
18    REAL(rstd) :: mass_tot,dmass_tot
19    INTEGER :: ind,i,j,ij
20
21    mass_tot=0
22    dmass_tot=0
23
24    CALL transfert_request(f_dps,req_i1)
25    CALL transfert_request(f_ps,req_i1)
26
27    DO ind=1,ndomain
28       IF (.NOT. assigned_domain(ind)) CYCLE
29       CALL swap_dimensions(ind)
30       CALL swap_geometry(ind)
31
32       ps=f_ps(ind)
33       dps=f_dps(ind)
34
35       DO j=jj_begin,jj_end
36          DO i=ii_begin,ii_end
37             ij=(j-1)*iim+i
38             IF (domain(ind)%own(i,j)) THEN
39                mass_tot=mass_tot+ps(ij)*Ai(ij)/g
40                dmass_tot=dmass_tot+dps(ij)*Ai(ij)/g
41             ENDIF
42          ENDDO
43       ENDDO
44
45    ENDDO
46    PRINT*, "mass_tot ", mass_tot,"      dmass_tot ",dmass_tot       
47
48  END SUBROUTINE check_mass_conservation
49
50
51  SUBROUTINE caldyn(write_out,f_phis, f_ps, f_theta_rhodz, f_u, f_q, &
52       f_hflux, f_wflux, f_dps, f_dtheta_rhodz, f_du)
53    USE icosa
54    USE output_field_mod
55    USE kinetic_mod
56    USE theta2theta_rhodz_mod
57    IMPLICIT NONE
58    LOGICAL,INTENT(IN)    :: write_out
59    TYPE(t_field),POINTER :: f_phis(:)
60    TYPE(t_field),POINTER :: f_ps(:)
61    TYPE(t_field),POINTER :: f_theta_rhodz(:)
62    TYPE(t_field),POINTER :: f_u(:)
63    TYPE(t_field),POINTER :: f_q(:)
64    TYPE(t_field),POINTER :: f_hflux(:), f_wflux(:)
65    TYPE(t_field) :: f_dps(:)
66    TYPE(t_field) :: f_dtheta_rhodz(:)
67    TYPE(t_field) :: f_du(:)
68
69    REAL(rstd),POINTER :: ps(:)
70    REAL(rstd),POINTER :: u(:,:)
71    REAL(rstd),POINTER :: dps(:)
72    REAL(rstd),POINTER :: hflux(:,:), wflux(:,:)
73    REAL(rstd),POINTER :: dtheta_rhodz(:,:), du(:,:)  ! set to 0
74    INTEGER :: ind
75
76    CALL transfert_request(f_ps,req_i1) 
77    CALL transfert_request(f_u,req_e1_vect)
78
79    DO ind=1,ndomain
80       IF (.NOT. assigned_domain(ind)) CYCLE
81       CALL swap_dimensions(ind)
82       CALL swap_geometry(ind)
83       ps=f_ps(ind)
84       u=f_u(ind)
85       dps=f_dps(ind)
86       hflux=f_hflux(ind)
87       wflux=f_wflux(ind)
88       dtheta_rhodz=f_dtheta_rhodz(ind)
89       du=f_du(ind)
90
91!       !$OMP PARALLEL DEFAULT(SHARED)
92       CALL compute_caldyn(ps,u,hflux, wflux, dps, dtheta_rhodz, du)
93!       !$OMP END PARALLEL
94    ENDDO
95
96    IF (write_out) THEN
97       CALL output_field("wflux",f_wflux)
98       CALL output_field("ps",f_ps)
99    ENDIF
100    !    CALL check_mass_conservation(f_ps,f_dps)
101
102  END SUBROUTINE caldyn
103
104  SUBROUTINE compute_caldyn(ps,u, hflux,wflux,dps, dtheta_rhodz,du)
105    USE icosa
106    USE disvert_mod
107    IMPLICIT NONE
108    REAL(rstd),INTENT(IN)  :: ps(iim*jjm)
109    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
110    REAL(rstd),INTENT(OUT) :: du(iim*3*jjm,llm), hflux(iim*3*jjm,llm) ! hflux in kg/s
111    REAL(rstd),INTENT(OUT) :: dtheta_rhodz(iim*jjm,llm)
112    REAL(rstd),INTENT(OUT) :: dps(iim*jjm)
113    REAL(rstd),INTENT(OUT) :: wflux(iim*jjm,llm+1) ! vertical mass flux (kg/m2/s)
114
115    REAL(rstd),ALLOCATABLE :: rhodz(:,:)
116    REAL(rstd),ALLOCATABLE :: divm(:,:)  ! mass flux divergence
117
118    INTEGER :: i,j,ij,l   
119    LOGICAL,SAVE :: first=.TRUE.
120
121    ALLOCATE(rhodz(iim*jjm,llm))
122    ALLOCATE(divm(iim*jjm,llm))    ! mass flux divergence
123
124    dtheta_rhodz(:,:)=0.
125    du(:,:)=0.
126
127!!! Compute mass
128    DO l = 1, llm
129       DO j=jj_begin-1,jj_end+1
130          DO i=ii_begin-1,ii_end+1
131             ij=(j-1)*iim+i
132             rhodz(ij,l) = (ap(l)-ap(l+1) + ps(ij)*(bp(l)-bp(l+1)) )/g
133          ENDDO
134       ENDDO
135    ENDDO
136
137    DO l = 1, llm
138!!!  Mass fluxes
139       DO j=jj_begin-1,jj_end+1
140          DO i=ii_begin-1,ii_end+1
141             ij=(j-1)*iim+i
142             hflux(ij+u_right,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_right,l))*u(ij+u_right,l)*le(ij+u_right)
143             hflux(ij+u_lup,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_lup,l))*u(ij+u_lup,l)*le(ij+u_lup)
144             hflux(ij+u_ldown,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_ldown,l))*u(ij+u_ldown,l)*le(ij+u_ldown)
145          ENDDO
146       ENDDO
147!!! Horizontal divergence of fluxes
148       DO j=jj_begin,jj_end
149          DO i=ii_begin,ii_end
150             ij=(j-1)*iim+i 
151             ! divm = +div(mass flux), sign convention as in Ringler et al. 2012, eq. 21
152             divm(ij,l)= 1./Ai(ij)*(ne(ij,right)*hflux(ij+u_right,l)   +  &
153                  ne(ij,rup)*hflux(ij+u_rup,l)       +  & 
154                  ne(ij,lup)*hflux(ij+u_lup,l)       +  & 
155                  ne(ij,left)*hflux(ij+u_left,l)     +  & 
156                  ne(ij,ldown)*hflux(ij+u_ldown,l)   +  & 
157                  ne(ij,rdown)*hflux(ij+u_rdown,l))
158          ENDDO
159       ENDDO
160    ENDDO
161
162!!! cumulate mass flux divergence from top to bottom
163    DO  l = llm-1, 1, -1
164       !$OMP DO
165       DO j=jj_begin,jj_end
166          DO i=ii_begin,ii_end
167             ij=(j-1)*iim+i
168             divm(ij,l) = divm(ij,l) + divm(ij,l+1)
169          ENDDO
170       ENDDO
171    ENDDO
172
173!!! Compute vertical mass flux
174    DO l = 1,llm-1
175       DO j=jj_begin,jj_end
176          DO i=ii_begin,ii_end
177             ij=(j-1)*iim+i
178             ! w = int(z,ztop,div(flux)dz) + B(eta)dps/dt
179             ! => w>0 for upward transport
180             wflux( ij, l+1 ) = divm( ij, l+1 ) - bp(l+1) * divm( ij, 1 )
181          ENDDO
182       ENDDO
183    ENDDO
184
185    ! compute dps, set vertical mass flux at the surface to 0
186    DO j=jj_begin,jj_end
187       DO i=ii_begin,ii_end
188          ij=(j-1)*iim+i
189          wflux(ij,1)  = 0.
190          ! dps/dt = -int(div flux)dz
191          dps(ij)=-divm(ij,1) * g 
192       ENDDO
193    ENDDO
194
195    DEALLOCATE(rhodz)
196    DEALLOCATE(divm)
197   
198  END SUBROUTINE compute_caldyn
199
200END MODULE caldyn_adv_mod
Note: See TracBrowser for help on using the repository browser.