source: codes/icosagcm/trunk/src/etat0_dcmip3.f90 @ 295

Last change on this file since 295 was 295, checked in by ymipsl, 10 years ago

Merging OpenMP parallisme mode : by subdomain and on vertical level.
This feature is actually experimental but may be retro-compatible with the last method based only on subdomain

YM

File size: 6.0 KB
Line 
1MODULE etat0_dcmip3_mod
2
3! test cases DCMIP 2012, category 3 : Non-hydrostatic gravity waves
4
5! Questions
6! Replace ps0 by preff ??
7
8  USE genmod
9  USE dcmip_initial_conditions_test_1_2_3
10
11  PRIVATE
12
13  PUBLIC  etat0
14
15CONTAINS
16
17
18  SUBROUTINE etat0(f_ps,f_phis,f_theta_rhodz,f_u, f_q)
19  USE icosa
20  USE theta2theta_rhodz_mod
21  IMPLICIT NONE
22    TYPE(t_field),POINTER :: f_ps(:)
23    TYPE(t_field),POINTER :: f_phis(:)
24    TYPE(t_field),POINTER :: f_theta_rhodz(:)
25    TYPE(t_field),POINTER :: f_u(:)
26    TYPE(t_field),POINTER :: f_q(:)
27    TYPE(t_field),POINTER,SAVE :: f_temp(:)
28   
29    REAL(rstd),POINTER :: ps(:)
30    REAL(rstd),POINTER :: phis(:)
31    REAL(rstd),POINTER :: u(:,:)
32    REAL(rstd),POINTER :: Temp(:,:)
33    REAL(rstd),POINTER :: q(:,:,:)
34
35    INTEGER :: ind
36
37    CALL allocate_field(f_temp,field_t,type_real,llm,name='temp') 
38
39    DO ind=1,ndomain
40      IF (.NOT. assigned_domain(ind)) CYCLE
41      CALL swap_dimensions(ind)
42      CALL swap_geometry(ind)
43
44      ps=f_ps(ind)
45      phis=f_phis(ind)
46      u=f_u(ind)
47      q=f_q(ind)
48      temp=f_temp(ind)
49      CALL compute_etat0_DCMIP3(ps,phis,u,Temp,q)
50    ENDDO
51
52    CALL temperature2theta_rhodz(f_ps,f_temp,f_theta_rhodz)
53    CALL deallocate_field(f_temp)
54           
55  END SUBROUTINE etat0
56 
57
58  SUBROUTINE compute_etat0_DCMIP3(ps, phis, u, temp,q)
59  USE icosa
60  USE pression_mod
61  USE theta2theta_rhodz_mod
62  USE wind_mod
63  IMPLICIT NONE
64  REAL(rstd),PARAMETER :: u0=20.         ! Maximum amplitude of the zonal wind (m.s-1)
65  REAL(rstd),PARAMETER :: N=0.01         ! Brunt-Vaisala frequency (s-1)
66  REAL(rstd),PARAMETER :: Teq=300.       ! Surface temperature at the equator (K)
67  REAL(rstd),PARAMETER :: Peq=1e5        ! Reference surface pressure at the equator (hPa)
68  REAL(rstd),PARAMETER :: d=5000.        ! Witdth parameter for theta
69  REAL(rstd),PARAMETER :: lonc=2*pi/3    ! Longitudinal centerpoint of theta
70  REAL(rstd),PARAMETER :: latc=0         ! Longitudinal centerpoint of theta
71  REAL(rstd),PARAMETER :: dtheta=1.      ! Maximum amplitude of theta (K)
72  REAL(rstd),PARAMETER :: Lz=20000.      ! Vertical wave lenght of the theta perturbation
73
74  REAL(rstd), INTENT(OUT) :: ps(iim*jjm)
75  REAL(rstd), INTENT(OUT) :: phis(iim*jjm)
76  REAL(rstd), INTENT(OUT) :: u(3*iim*jjm,llm)
77  REAL(rstd), INTENT(OUT) :: Temp(iim*jjm,llm)
78  REAL(rstd), INTENT(OUT) :: q(iim*jjm,llm,nqtot)
79 
80  REAL(rstd) :: Ts(iim*jjm)
81  REAL(rstd) :: s(iim*jjm)
82  REAL(rstd) :: p(iim*jjm,llm+1)
83  REAL(rstd) :: theta(iim*jjm,llm)
84  REAL(rstd) :: ulon(3*iim*jjm,llm)
85  REAL(rstd) :: ulat(3*iim*jjm,llm)
86 
87 
88  INTEGER :: i,j,l,ij
89  REAL(rstd) :: Rd        ! gas constant of dry air, P=rho.Rd.T
90  REAL(rstd) :: alpha, rm
91  REAL(rstd) :: C0, C1, GG
92  REAL(rstd) :: p0psk, pspsk,r,zz, thetab, thetap
93  REAL(rstd) :: dummy, pp
94  LOGICAL    :: use_dcmip_routine
95 
96  Rd=cpp*kappa
97 
98  GG=(g/N)**2/cpp
99  C0=0.25*u0*(u0+2.*Omega*radius)
100 
101  q(:,:,:)=0
102 
103!  use_dcmip_routine=.TRUE.
104  use_dcmip_routine=.FALSE.
105  dummy=0.
106
107  pp=peq
108  DO j=jj_begin,jj_end
109     DO i=ii_begin,ii_end
110        ij=(j-1)*iim+i
111
112        IF(use_dcmip_routine) THEN
113           CALL test3_gravity_wave(lon_i(ij),lat_i(ij),pp,dummy,0, dummy,dummy,dummy,dummy,phis(ij),ps(ij),dummy,dummy)
114        ELSE
115           C1=C0*(cos(2*lat_i(ij))-1)
116           
117           !--- GROUND GEOPOTENTIAL
118           phis(ij)=0.
119           
120           !--- GROUND TEMPERATURE
121           Ts(ij) = GG+(Teq-GG)*EXP(-C1*(N/g)**2)
122           
123           !--- GROUND PRESSURE
124           Ps(ij) = peq*EXP(C1/GG/Rd)*(Ts(ij)/Teq)**(1/kappa)
125           
126           
127           r=radius*acos(sin(latc)*sin(lat_i(ij))+cos(latc)*cos(lat_i(ij))*cos(lon_i(ij)-lonc))
128           s(ij)= d**2/(d**2+r**2)
129        END IF
130     END DO
131  END DO
132 
133!$OMP BARRIER
134  CALL compute_pression(ps,p,0)
135!$OMP BARRIER
136 
137  DO l=1,llm
138     DO j=jj_begin,jj_end
139        DO i=ii_begin,ii_end
140           ij=(j-1)*iim+i
141           pp=0.5*(p(ij,l+1)+p(ij,l))  ! full-layer pressure
142           IF(use_dcmip_routine) THEN
143              CALL test3_gravity_wave(lon_i(ij),lat_i(ij),pp,dummy,0, &
144                   dummy,dummy,dummy,Temp(ij,l),dummy,dummy,dummy,dummy)
145           ELSE
146              pspsk=(pp/ps(ij))**kappa
147              p0psk=(Peq/ps(ij))**kappa
148              thetab = Ts(ij)*p0psk / ( Ts(ij) / GG * ( pspsk-1) +1)  ! background pot. temp.
149              zz     = -g/N**2*log( Ts(ij)/GG * (pspsk -1)+1)         ! altitude
150              thetap = dtheta *sin(2*Pi*zz/Lz) * s(ij)                ! perturbation pot. temp.
151              theta(ij,l) = thetab + thetap
152              Temp(ij,l) = theta(ij,l)* ((pp/peq)**kappa)
153              ! T(ij,l) = Ts(ij)*pspsk / ( Ts(ij) / GG * ( pspsk-1) +1)  ! background temp.
154           END IF
155        ENDDO
156     ENDDO
157  ENDDO
158 
159!  IF(use_dcmip_routine) THEN
160!     CALL compute_temperature2theta_rhodz(ps,T,theta_rhodz,0)
161!  ELSE
162!     CALL compute_temperature2theta_rhodz(ps,T,theta_rhodz,0)
163!  END IF
164 
165  pp=peq
166  DO l=1,llm
167     DO j=jj_begin-1,jj_end+1
168        DO i=ii_begin-1,ii_end+1
169           ij=(j-1)*iim+i
170           IF(use_dcmip_routine) THEN
171              CALL test3_gravity_wave(lon_e(ij+u_right),lat_e(ij+u_right), &
172                   pp,0.,0, ulon(ij+u_right,l),ulat(ij+u_right,l),&
173                   dummy,dummy,dummy,dummy,dummy,dummy)
174              CALL test3_gravity_wave(lon_e(ij+u_lup),lat_e(ij+u_lup), &
175                   pp,0.,0, ulon(ij+u_lup,l),ulat(ij+u_lup,l),&
176                   dummy,dummy,dummy,dummy,dummy,dummy)
177              CALL test3_gravity_wave(lon_e(ij+u_ldown),lat_e(ij+u_ldown), &
178                   pp,0.,0, ulon(ij+u_ldown,l),ulat(ij+u_ldown,l),&
179                   dummy,dummy,dummy,dummy,dummy,dummy)
180           ELSE
181              ulon(ij+u_right,l) = u0*cos(lat_e(ij+u_right))
182              ulat(ij+u_right,l) = 0
183              ulon(ij+u_lup,l)   = u0*cos(lat_e(ij+u_lup))
184              ulat(ij+u_lup,l)   = 0
185              ulon(ij+u_ldown,l) = u0*cos(lat_e(ij+u_ldown))
186              ulat(ij+u_ldown,l) = 0
187           END IF
188        ENDDO
189     ENDDO
190  ENDDO
191 
192  CALL compute_wind_perp_from_lonlat_compound(ulon,ulat,u)   
193 
194END SUBROUTINE compute_etat0_DCMIP3
195
196
197END MODULE etat0_DCMIP3_mod
Note: See TracBrowser for help on using the repository browser.