source: codes/icosagcm/trunk/src/physics.f90 @ 196

Last change on this file since 196 was 196, checked in by dubos, 10 years ago

First draft of generic dynamics-physics interface - works with DCMIP5.1

File size: 6.2 KB
RevLine 
[81]1MODULE physics_mod
2
[196]3  USE field_mod
4
5  PRIVATE
6
7  INTEGER, PARAMETER :: phys_none=0, phys_HS94=1, phys_DCMIP=2
8
9  INTEGER :: nb_extra_physics_2D, nb_extra_physics_3D, phys_type
10  TYPE(t_field),POINTER :: f_extra_physics_2D(:), f_extra_physics_3D(:)
11
[170]12  CHARACTER(LEN=255) :: physics_type="automatic"
[186]13!$OMP THREADPRIVATE(physics_type)
[81]14
[196]15  PUBLIC :: physics, init_physics
[81]16CONTAINS
17
[98]18  SUBROUTINE init_physics
[178]19    USE mpipara
20    USE etat0_mod
[170]21    USE icosa
[196]22    USE physics_interface_mod
23    USE physics_dcmip_mod, ONLY : init_physics_dcmip=>init_physics
[170]24    IMPLICIT NONE
25
[178]26    physics_type='automatic'
[81]27    CALL getin("physics",physics_type)
[170]28
[81]29    SELECT CASE(TRIM(physics_type))
[170]30    CASE ('automatic')
[178]31       etat0_type='jablonowsky06'
32       CALL getin("etat0",etat0_type)
33       SELECT CASE(TRIM(etat0_type))
34       CASE('held_suarez')
[196]35          phys_type = phys_HS94
[178]36       CASE DEFAULT
37          IF(is_mpi_root) PRINT*,"NO PHYSICAL PACKAGE USED"
[196]38          phys_type = phys_none
[178]39       END SELECT
[149]40
[170]41    CASE ('dcmip')
42       CALL init_physics_dcmip
[196]43       nb_extra_physics_2D=1 ! precl
44       nb_extra_physics_3D=0
45       phys_type = phys_DCMIP
[170]46    CASE DEFAULT
[196]47       IF(is_mpi_root) PRINT*, 'init_physics : Bad selector for variable physics <',&
48            TRIM(physics_type), '> options are <automatic>, <dcmip>, <dry>'
[170]49       STOP
[81]50    END SELECT
[170]51
[196]52    IF(is_mpi_root) THEN
53       PRINT *, 'phys_type = ',phys_type
54       PRINT *, 'nb_extra_physics_2D = ', nb_extra_physics_2D
55       PRINT *, 'nb_extra_physics_3D = ', nb_extra_physics_3D
56    END IF
57    IF(nb_extra_physics_2D>0) CALL allocate_field(f_extra_physics_2D,field_t,type_real,nb_extra_physics_2D)
58    IF(nb_extra_physics_3D>0) CALL allocate_field(f_extra_physics_3D,field_t,type_real,llm,nb_extra_physics_3D)
59
[81]60  END SUBROUTINE init_physics
61
[149]62  SUBROUTINE physics(it,jD_cur,jH_cur,f_phis, f_ps, f_theta_rhodz, f_ue, f_q)
[170]63    USE icosa
[196]64    USE physics_interface_mod
65 !   USE etat0_mod
[170]66    USE etat0_heldsz_mod
67    IMPLICIT NONE
[99]68    INTEGER, INTENT(IN)   :: it
[149]69    REAL(rstd),INTENT(IN)::jD_cur,jH_cur
[81]70    TYPE(t_field),POINTER :: f_phis(:)
71    TYPE(t_field),POINTER :: f_ps(:)
72    TYPE(t_field),POINTER :: f_theta_rhodz(:)
73    TYPE(t_field),POINTER :: f_ue(:)
74    TYPE(t_field),POINTER :: f_q(:)
[196]75    REAL(rstd),POINTER :: phis(:)
76    REAL(rstd),POINTER :: ps(:)
77    REAL(rstd),POINTER :: theta_rhodz(:,:)
78    REAL(rstd),POINTER :: ue(:,:)
79    REAL(rstd),POINTER :: q(:,:,:)
80
[149]81    LOGICAL:: firstcall,lastcall
[196]82    INTEGER :: ind
83    TYPE(physics_inout) :: args
[170]84
[196]85    SELECT CASE(phys_type)
86    CASE (phys_none)
87       ! No physics, do nothing
88    CASE(phys_HS94)
89       CALL held_suarez(f_ps,f_theta_rhodz,f_ue) 
90    CASE DEFAULT
91       CALL transfert_request(f_ps,req_i1)
92       CALL transfert_request(f_theta_rhodz,req_i1)
93       CALL transfert_request(f_ue,req_e1_vect)
94       CALL transfert_request(f_q,req_i1)
95       
96       args%dt_phys = dt
97       DO ind=1,ndomain
98          IF (.NOT. assigned_domain(ind)) CYCLE
99          CALL swap_dimensions(ind)
100          CALL swap_geometry(ind)
101     
102          phis=f_phis(ind)
103          ps=f_ps(ind)
104          theta_rhodz=f_theta_rhodz(ind)
105          ue=f_ue(ind)
106          q=f_q(ind)
[149]107
[196]108          IF(nb_extra_physics_2D>0) args%extra_2D=f_extra_physics_2D(ind)
109          IF(nb_extra_physics_3D>0) args%extra_3D=f_extra_physics_3D(ind)
110          CALL physics_column(args, phis, ps, theta_rhodz, ue, q)
111       ENDDO
[149]112
[196]113       IF (mod(it,itau_out)==0 ) THEN
114          IF(nb_extra_physics_2D>0) CALL writefield("extra_physics_2D",f_extra_physics_2D)
115          IF(nb_extra_physics_3D>0) CALL writefield("extra_physics_3D",f_extra_physics_3D)
116       ENDIF
117    END SELECT
[170]118
[196]119  END SUBROUTINE physics
[170]120
[196]121  SUBROUTINE physics_column(args, phis, ps, theta_rhodz, ue, q)
122    USE icosa
123    USE wind_mod
124    USE pression_mod
125    USE theta2theta_rhodz_mod
126    USE physics_interface_mod
127    USE physics_dcmip_mod
128    IMPLICIT NONE
129    TYPE(physics_inout) :: args   
130    REAL(rstd) :: phis(iim*jjm)
131    REAL(rstd) :: ps(iim*jjm)
132    REAL(rstd) :: theta_rhodz(iim*jjm,llm)
133    REAL(rstd) :: ue(3*iim*jjm,llm)
134    REAL(rstd), TARGET :: q(iim*jjm,llm,nqtot)
135    ! local arrays
136    REAL(rstd), TARGET :: lat(iim*jjm)
137    REAL(rstd), TARGET :: lon(iim*jjm)
138    REAL(rstd), TARGET :: p(iim*jjm,llm+1)
139    REAL(rstd), TARGET :: Temp(iim*jjm,llm)
140    REAL(rstd), TARGET :: ulon(iim*jjm,llm)
141    REAL(rstd), TARGET :: ulat(iim*jjm,llm)
142    REAL(rstd), TARGET :: dTemp(iim*jjm,llm)
143    REAL(rstd), TARGET :: dulon(iim*jjm,llm)
144    REAL(rstd), TARGET :: dulat(iim*jjm,llm)
145    REAL(rstd), TARGET :: dq(iim*jjm,llm,nqtot)
146    REAL(rstd) :: uc(iim*jjm,3,llm)  ! 3D velocity at cell centers
147
148    INTEGER :: i,j,ij,l
149    REAL(rstd) :: due, dt2
150
151    DO j=jj_begin,jj_end
152      DO i=ii_begin,ii_end
153        ij=(j-1)*iim+i
154        CALL xyz2lonlat(xyz_i(ij,:),lon(ij),lat(ij)) 
155      ENDDO
156    ENDDO
157
158    ! Reconstruct wind vector at hexagons
159    CALL compute_pression(ps,p,0)
160    CALL compute_theta_rhodz2temperature(ps,theta_rhodz,Temp,0)
161    CALL compute_wind_centered(ue,uc)
162    CALL compute_wind_centered_lonlat_compound(uc, ulon, ulat)
163    args%ngrid = iim*jjm
164    args%lon => lon
165    args%lat => lat
166    args%p => p
167    args%Temp => Temp
168    args%ulon => ulon
169    args%ulat => ulat
170    args%q => q
171    args%dTemp => dTemp
172    args%dulon => dulon
173    args%dulat => dulat
174    args%dq => dq
175
176    SELECT CASE(phys_type)
177    CASE (phys_DCMIP)
178       CALL compute_phys_wrap(args)
[81]179    END SELECT
[170]180
[196]181    q = q + args%dt_phys * dq
182    Temp = Temp + args%dt_phys * dTemp
183    CALL compute_temperature2theta_rhodz(ps,Temp,theta_rhodz,0)
184   
185    ! Reconstruct wind tendencies at edges and add
186    CALL compute_wind_centered_from_lonlat_compound(dulon,dulat,uc)
187    dt2=.5*args%dt_phys
188    DO l=1,llm
189      DO j=jj_begin,jj_end
190        DO i=ii_begin,ii_end
191          ij=(j-1)*iim+i
192          due = sum((uc(ij,:,l) + uc(ij+t_right,:,l))*ep_e(ij+u_right,:) )
193          ue(ij+u_right,l) = ue(ij+u_right,l) + dt2*due
[81]194
[196]195          due = sum( (uc(ij,:,l) + uc(ij+t_lup,:,l))*ep_e(ij+u_lup,:) )
196          ue(ij+u_lup,l)=ue(ij+u_lup,l) + dt2*due
197
198          due = sum( (uc(ij,:,l) + uc(ij+t_ldown,:,l))*ep_e(ij+u_ldown,:) )
199          ue(ij+u_ldown,l)=ue(ij+u_ldown,l) + dt2*due
200        ENDDO
201      ENDDO
202    ENDDO
203
204  END SUBROUTINE physics_column
205
[81]206END MODULE physics_mod
Note: See TracBrowser for help on using the repository browser.