source: codes/icosagcm/trunk/src/etat0.f90 @ 377

Last change on this file since 377 was 377, checked in by dubos, 8 years ago

New : positive advection option for theta

File size: 11.3 KB
Line 
1MODULE etat0_mod
2  USE icosa
3  IMPLICIT NONE         
4  PRIVATE
5
6    CHARACTER(len=255),SAVE :: etat0_type
7!$OMP THREADPRIVATE(etat0_type)
8
9    REAL(rstd) :: etat0_temp
10
11    PUBLIC :: etat0, etat0_type
12
13CONTAINS
14 
15  SUBROUTINE etat0(f_ps,f_mass,f_phis,f_theta_rhodz,f_u, f_geopot,f_w, f_q)
16    USE mpipara, ONLY : is_mpi_root
17    USE disvert_mod
18    ! Generic interface
19    USE etat0_dcmip1_mod, ONLY : getin_etat0_dcmip1=>getin_etat0
20    USE etat0_dcmip2_mod, ONLY : getin_etat0_dcmip2=>getin_etat0
21    USE etat0_dcmip4_mod, ONLY : getin_etat0_dcmip4=>getin_etat0
22    USE etat0_dcmip5_mod, ONLY : getin_etat0_dcmip5=>getin_etat0
23    USE etat0_bubble_mod, ONLY : getin_etat0_bubble=>getin_etat0
24    USE etat0_williamson_mod, ONLY : getin_etat0_williamson=>getin_etat0
25    USE etat0_temperature_mod, ONLY: getin_etat0_temperature=>getin_etat0
26    ! Ad hoc interfaces
27    USE etat0_academic_mod, ONLY : etat0_academic=>etat0 
28    USE etat0_heldsz_mod, ONLY : etat0_heldsz=>etat0 
29    USE etat0_venus_mod,  ONLY : etat0_venus=>etat0 
30    USE etat0_start_file_mod, ONLY : etat0_start_file=>etat0 
31
32    IMPLICIT NONE
33    TYPE(t_field),POINTER :: f_ps(:)
34    TYPE(t_field),POINTER :: f_mass(:)
35    TYPE(t_field),POINTER :: f_phis(:)
36    TYPE(t_field),POINTER :: f_theta_rhodz(:)
37    TYPE(t_field),POINTER :: f_u(:)
38    TYPE(t_field),POINTER :: f_geopot(:)
39    TYPE(t_field),POINTER :: f_w(:)
40    TYPE(t_field),POINTER :: f_q(:)
41   
42    REAL(rstd),POINTER :: ps(:), mass(:,:)
43    LOGICAL :: autoinit_mass, autoinit_geopot, collocated
44    INTEGER :: ind,i,j,ij,l
45
46    ! most etat0 routines set ps and not mass
47    ! in that case and if caldyn_eta == eta_lag
48    ! the initial distribution of mass is taken to be the same
49    ! as what the mass coordinate would dictate
50    ! however if etat0_XXX defines mass then the flag autoinit_mass must be set to .FALSE.
51    ! otherwise mass will be overwritten
52    autoinit_mass = (caldyn_eta == eta_lag)
53
54    etat0_type='jablonowsky06'
55    CALL getin("etat0",etat0_type)
56   
57    !------------------- Generic interface ---------------------
58    collocated=.TRUE.
59    SELECT CASE (TRIM(etat0_type))
60    CASE ('isothermal')
61       CALL getin_etat0_isothermal
62    CASE ('temperature_profile')
63       CALL getin_etat0_temperature
64    CASE ('jablonowsky06')
65    CASE ('dcmip1')
66        CALL getin_etat0_dcmip1
67    CASE ('dcmip2_mountain','dcmip2_schaer_noshear','dcmip2_schaer_shear')
68       CALL getin_etat0_dcmip2
69    CASE ('dcmip3')
70    CASE ('dcmip4')
71        CALL getin_etat0_dcmip4
72    CASE ('dcmip5')
73        CALL getin_etat0_dcmip5
74    CASE ('bubble')
75        CALL getin_etat0_bubble
76    CASE ('williamson91.6')
77       autoinit_mass=.FALSE.
78       CALL getin_etat0_williamson
79    CASE DEFAULT
80       collocated=.FALSE.
81       autoinit_mass = .FALSE.
82    END SELECT
83
84    !------------------- Ad hoc interfaces --------------------
85    SELECT CASE (TRIM(etat0_type))
86    CASE ('start_file')
87       CALL etat0_start_file(f_ps,f_phis,f_theta_rhodz,f_u, f_q)
88    CASE ('academic')
89       CALL etat0_academic(f_ps,f_phis,f_theta_rhodz,f_u, f_q)
90    CASE ('held_suarez')
91       PRINT *,"Held & Suarez (1994) test case"
92       CALL etat0_heldsz(f_ps,f_phis,f_theta_rhodz,f_u, f_q)
93    CASE ('venus')
94       CALL etat0_venus(f_ps, f_phis, f_theta_rhodz, f_u, f_q)
95       PRINT *, "Venus (Lebonnois et al., 2012) test case"
96   CASE DEFAULT
97      IF(collocated) THEN
98         CALL etat0_collocated(f_phis,f_ps,f_mass,f_theta_rhodz,f_u, f_geopot,f_W, f_q)
99      ELSE
100         PRINT*, 'Bad selector for variable etat0 <',etat0_type, &
101              '> options are <isothermal>, <temperature_profile>, <held_suarez>, &
102              &<bubble>, <venus>, <jablonowsky06>, <academic>, <dcmip[1-4]> '
103         STOP
104      END IF
105    END SELECT
106
107!       !$OMP BARRIER
108    IF(autoinit_mass) THEN
109       DO ind=1,ndomain
110          IF (.NOT. assigned_domain(ind)) CYCLE
111          CALL swap_dimensions(ind)
112          CALL swap_geometry(ind)
113          mass=f_mass(ind); ps=f_ps(ind)
114          CALL compute_rhodz(.TRUE., ps, mass) ! initialize mass distribution using ps
115       END DO
116    END IF
117 
118  END SUBROUTINE etat0
119
120  SUBROUTINE etat0_collocated(f_phis,f_ps,f_mass,f_theta_rhodz,f_u, f_geopot,f_W, f_q)
121    USE theta2theta_rhodz_mod
122    IMPLICIT NONE
123    TYPE(t_field),POINTER :: f_ps(:)
124    TYPE(t_field),POINTER :: f_mass(:)
125    TYPE(t_field),POINTER :: f_phis(:)
126    TYPE(t_field),POINTER :: f_theta_rhodz(:)
127    TYPE(t_field),POINTER :: f_u(:)
128    TYPE(t_field),POINTER :: f_geopot(:)
129    TYPE(t_field),POINTER :: f_W(:)
130    TYPE(t_field),POINTER :: f_q(:)
131 
132    TYPE(t_field),POINTER,SAVE :: f_temp(:)
133    REAL(rstd),POINTER :: ps(:)
134    REAL(rstd),POINTER :: mass(:,:)
135    REAL(rstd),POINTER :: phis(:)
136    REAL(rstd),POINTER :: theta_rhodz(:,:)
137    REAL(rstd),POINTER :: temp(:,:)
138    REAL(rstd),POINTER :: u(:,:)
139    REAL(rstd),POINTER :: geopot(:,:)
140    REAL(rstd),POINTER :: W(:,:)
141    REAL(rstd),POINTER :: q(:,:,:)
142    INTEGER :: ind
143
144    CALL allocate_field(f_temp,field_t,type_real,llm,name='temp')
145
146    DO ind=1,ndomain
147      IF (.NOT. assigned_domain(ind)) CYCLE
148      CALL swap_dimensions(ind)
149      CALL swap_geometry(ind)
150      ps=f_ps(ind)
151      mass=f_mass(ind)
152      phis=f_phis(ind)
153      theta_rhodz=f_theta_rhodz(ind)
154      temp=f_temp(ind)
155      u=f_u(ind)
156      geopot=f_geopot(ind)
157      w=f_w(ind)
158      q=f_q(ind)
159
160      IF( TRIM(etat0_type)=='williamson91.6' ) THEN
161       CALL compute_etat0_collocated(ps,mass, phis, theta_rhodz, u, geopot, W, q)
162      ELSE
163       CALL compute_etat0_collocated(ps,mass, phis, temp, u, geopot, W, q)
164      ENDIF
165    ENDDO
166   
167    IF( TRIM(etat0_type)/='williamson91.6' ) CALL temperature2theta_rhodz(f_ps,f_temp,f_theta_rhodz)
168   
169    CALL deallocate_field(f_temp)
170   
171  END SUBROUTINE etat0_collocated
172
173  SUBROUTINE compute_etat0_collocated(ps,mass,phis,temp_i,u, geopot,W, q)
174    USE wind_mod
175    USE disvert_mod
176    USE etat0_jablonowsky06_mod, ONLY : compute_jablonowsky06 => compute_etat0
177    USE etat0_dcmip1_mod, ONLY : compute_dcmip1 => compute_etat0
178    USE etat0_dcmip2_mod, ONLY : compute_dcmip2 => compute_etat0
179    USE etat0_dcmip3_mod, ONLY : compute_dcmip3 => compute_etat0
180    USE etat0_dcmip4_mod, ONLY : compute_dcmip4 => compute_etat0
181    USE etat0_dcmip5_mod, ONLY : compute_dcmip5 => compute_etat0
182    USE etat0_bubble_mod, ONLY : compute_bubble => compute_etat0 
183    USE etat0_williamson_mod, ONLY : compute_w91_6 => compute_etat0
184    USE etat0_temperature_mod, ONLY: compute_etat0_temperature => compute_etat0
185    IMPLICIT NONE
186    REAL(rstd),INTENT(INOUT) :: ps(iim*jjm)
187    REAL(rstd),INTENT(INOUT) :: mass(iim*jjm,llm)
188    REAL(rstd),INTENT(OUT) :: phis(iim*jjm)
189    REAL(rstd),INTENT(OUT) :: temp_i(iim*jjm,llm)
190    REAL(rstd),INTENT(OUT) :: u(3*iim*jjm,llm)
191    REAL(rstd),INTENT(OUT) :: W(iim*jjm,llm+1)
192    REAL(rstd),INTENT(OUT) :: geopot(iim*jjm,llm+1)
193    REAL(rstd),INTENT(OUT) :: q(iim*jjm,llm,nqtot)
194
195    REAL(rstd) :: ulon_i(iim*jjm,llm)
196    REAL(rstd) :: ulat_i(iim*jjm,llm)
197
198    REAL(rstd) :: ps_e(3*iim*jjm)
199    REAL(rstd) :: mass_e(3*iim*jjm,llm)
200    REAL(rstd) :: phis_e(3*iim*jjm)
201    REAL(rstd) :: temp_e(3*iim*jjm,llm)
202    REAL(rstd) :: geopot_e(3*iim*jjm,llm+1)
203    REAL(rstd) :: ulon_e(3*iim*jjm,llm)
204    REAL(rstd) :: ulat_e(3*iim*jjm,llm)
205    REAL(rstd) :: q_e(3*iim*jjm,llm,nqtot)
206
207    INTEGER :: l,i,j,ij
208    REAL :: p_ik, v_ik, mass_ik
209    LOGICAL :: autoinit_mass, autoinit_NH
210
211    ! For NH geopotential and vertical momentum must be initialized.
212    ! Unless autoinit_NH is set to .FALSE. , they will be initialized
213    ! to hydrostatic geopotential and zero
214    autoinit_mass = .TRUE.
215    autoinit_NH = .NOT. hydrostatic
216    w(:,:) = 0
217
218    !$OMP BARRIER
219
220    SELECT CASE (TRIM(etat0_type))
221    CASE ('isothermal')
222       CALL compute_etat0_isothermal(iim*jjm, phis, ps, temp_i, ulon_i, ulat_i, q)
223       CALL compute_etat0_isothermal(3*iim*jjm, phis_e, ps_e, temp_e, ulon_e, ulat_e, q_e)
224    CASE ('temperature_profile')
225       CALL compute_etat0_temperature(iim*jjm, phis, ps, temp_i, ulon_i, ulat_i, q)
226       CALL compute_etat0_temperature(3*iim*jjm, phis_e, ps_e, temp_e, ulon_e, ulat_e, q_e)
227    CASE('jablonowsky06')
228       CALL compute_jablonowsky06(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i)
229       CALL compute_jablonowsky06(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e)
230    CASE('dcmip1')
231       CALL compute_dcmip1(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i, q)
232       CALL compute_dcmip1(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e, q_e)
233    CASE ('dcmip2_mountain','dcmip2_schaer_noshear','dcmip2_schaer_shear')
234       CALL compute_dcmip2(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i)
235       CALL compute_dcmip2(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e)     
236    CASE('dcmip3')
237       CALL compute_dcmip3(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i, geopot, q)
238       CALL compute_dcmip3(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e, geopot_e, q_e)
239       autoinit_NH = .FALSE. ! compute_dcmip3 initializes geopot
240    CASE('dcmip4')
241       CALL compute_dcmip4(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i, q)
242       CALL compute_dcmip4(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e, q_e)
243    CASE('dcmip5')
244       CALL compute_dcmip5(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i, q)
245       CALL compute_dcmip5(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e, q_e)
246    CASE('bubble')
247       CALL compute_bubble(iim*jjm,lon_i,lat_i, phis, ps, temp_i, ulon_i, ulat_i, geopot, q)
248       CALL compute_bubble(3*iim*jjm,lon_e,lat_e, phis_e, ps_e, temp_e, ulon_e, ulat_e, geopot_e, q_e)
249!       autoinit_NH = .FALSE. ! compute_bubble initializes geopot
250    CASE('williamson91.6')
251       CALL compute_w91_6(iim*jjm,lon_i,lat_i, phis, mass(:,1), temp_i(:,1), ulon_i(:,1), ulat_i(:,1))
252       CALL compute_w91_6(3*iim*jjm,lon_e,lat_e, phis_e, mass_e(:,1), temp_e(:,1), ulon_e(:,1), ulat_e(:,1))
253       autoinit_mass = .FALSE. ! do not overwrite mass
254    END SELECT
255
256    IF(autoinit_mass) CALL compute_rhodz(.TRUE., ps, mass) ! initialize mass distribution using ps
257    IF(autoinit_NH) THEN
258       geopot(:,1) = phis(:) ! surface geopotential
259       DO l = 1, llm
260          DO ij=1,iim*jjm
261             ! hybrid pressure coordinate
262             p_ik = ptop + mass_ak(l) + mass_bk(l)*ps(ij)
263             mass_ik = (mass_dak(l) + mass_dbk(l)*ps(ij))/g
264             ! v=R.T/p, R=kappa*cpp
265             v_ik = kappa*cpp*temp_i(ij,l)/p_ik
266             geopot(ij,l+1) = geopot(ij,l) + mass_ik*v_ik*g
267          END DO
268       END DO
269    END IF
270
271    !$OMP BARRIER
272
273    CALL compute_wind_perp_from_lonlat_compound(ulon_e, ulat_e, u)
274
275  END SUBROUTINE compute_etat0_collocated
276
277!----------------------------- Resting isothermal state --------------------------------
278
279  SUBROUTINE getin_etat0_isothermal
280    etat0_temp=300
281    CALL getin("etat0_isothermal_temp",etat0_temp)
282  END SUBROUTINE getin_etat0_isothermal
283
284  SUBROUTINE compute_etat0_isothermal(ngrid, phis, ps, temp, ulon, ulat, q)
285    IMPLICIT NONE 
286    INTEGER, INTENT(IN)    :: ngrid
287    REAL(rstd),INTENT(OUT) :: phis(ngrid)
288    REAL(rstd),INTENT(OUT) :: ps(ngrid)
289    REAL(rstd),INTENT(OUT) :: temp(ngrid,llm)
290    REAL(rstd),INTENT(OUT) :: ulon(ngrid,llm)
291    REAL(rstd),INTENT(OUT) :: ulat(ngrid,llm)
292    REAL(rstd),INTENT(OUT) :: q(ngrid,llm,nqtot)
293    phis(:)=0
294    ps(:)=preff
295    temp(:,:)=etat0_temp
296    ulon(:,:)=0
297    ulat(:,:)=0
298    q(:,:,:)=0
299  END SUBROUTINE compute_etat0_isothermal
300
301END MODULE etat0_mod
Note: See TracBrowser for help on using the repository browser.