source: trunk/LMDZ.GENERIC/libf/dynlonlat_phylonlat/phystd/iniphysiq_mod.F90 @ 1529

Last change on this file since 1529 was 1529, checked in by emillour, 9 years ago

Generic GCM: Towards a cleaner separation between physics and dynamics

  • Got rid of references to "dimensions.h" from physics packages: use nbp_lon (=iim), nbp_lat (=jjp1) and nbp_lev from module mod_grid_phy_lmdz (in phy_common) instead.
  • Removed module "comhdiff_mod.F90", as it is only used by module surf_heat_transp_mod.F90, moved module variables there.
  • Addedin "surf_heat_transp_mod" local versions of some arrays and routines (from dyn3d) required to compute gradient, divergence, etc. on the global dynamics grid. As before, the slab ocean only works in serial.

EM

File size: 7.1 KB
Line 
1MODULE iniphysiq_mod
2
3CONTAINS
4
5subroutine iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep,           &
6               rlatudyn,rlatvdyn,rlonudyn,rlonvdyn,airedyn,cudyn,cvdyn,  &
7                     prad,pg,pr,pcpp,iflag_phys)
8
9use dimphy, only : klev ! number of atmospheric levels
10use mod_grid_phy_lmdz, only : klon_glo ! number of atmospheric columns
11                                       ! (on full grid)
12use mod_phys_lmdz_para, only : klon_omp, & ! number of columns (on local omp grid)
13                               klon_omp_begin, & ! start index of local omp subgrid
14                               klon_omp_end, & ! end index of local omp subgrid
15                               klon_mpi_begin ! start indes of columns (on local mpi grid)
16use control_mod, only: nday
17use comgeomphy, only : initcomgeomphy, &
18                       airephy, & ! physics grid area (m2)
19                       cuphy, & ! cu coeff. (u_covariant = cu * u)
20                       cvphy, & ! cv coeff. (v_covariant = cv * v)
21                       rlond, & ! longitudes
22                       rlatd ! latitudes
23use surf_heat_transp_mod, only: ini_surf_heat_transp
24use infotrac, only : nqtot ! number of advected tracers
25use planete_mod, only: ini_planete_mod
26USE comvert_mod, ONLY: ap,bp,preff
27use inifis_mod, only: inifis
28use regular_lonlat_mod, only: init_regular_lonlat, &
29                              east, west, north, south, &
30                              north_east, north_west, &
31                              south_west, south_east
32use ioipsl_getin_p_mod, only: getin_p
33
34implicit none
35include "dimensions.h"
36include "paramet.h"
37include "comgeom.h"
38include "iniprint.h"
39
40real,intent(in) :: prad ! radius of the planet (m)
41real,intent(in) :: pg ! gravitational acceleration (m/s2)
42real,intent(in) :: pr ! ! reduced gas constant R/mu
43real,intent(in) :: pcpp ! specific heat Cp
44real,intent(in) :: punjours ! length (in s) of a standard day
45!integer,intent(in) :: ngrid ! number of horizontal grid points in the physics (full grid)
46integer,intent(in) :: nlayer ! number of atmospheric layers
47integer,intent(in) :: ii ! number of atmospheric coulumns along longitudes
48integer,intent(in) :: jj  ! number of atompsheric columns along latitudes
49real,intent(in) :: rlatudyn(jj+1) ! latitudes of the physics grid
50real,intent(in) :: rlatvdyn(jj) ! latitude boundaries of the physics grid
51real,intent(in) :: rlonvdyn(ii+1) ! longitudes of the physics grid
52real,intent(in) :: rlonudyn(ii+1) ! longitude boundaries of the physics grid
53real,intent(in) :: airedyn(ii+1,jj+1) ! area of the dynamics grid (m2)
54real,intent(in) :: cudyn((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
55real,intent(in) :: cvdyn((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
56integer,intent(in) :: pdayref ! reference day of for the simulation
57real,intent(in) :: ptimestep !physics time step (s)
58integer,intent(in) :: iflag_phys ! type of physics to be called
59
60integer :: ibegin,iend,offset
61integer :: i,j
62character(len=20) :: modname='iniphysiq'
63character(len=80) :: abort_message
64real :: total_area_phy, total_area_dyn
65real :: pi
66logical :: ok_slab_ocean
67
68! boundaries, on global grid
69real,allocatable :: boundslon_reg(:,:)
70real,allocatable :: boundslat_reg(:,:)
71
72! global array, on full physics grid:
73real,allocatable :: latfi(:)
74real,allocatable :: lonfi(:)
75real,allocatable :: cufi(:)
76real,allocatable :: cvfi(:)
77real,allocatable :: airefi(:)
78
79pi=2.*asin(1.0)
80
81IF (nlayer.NE.klev) THEN
82  write(*,*) 'STOP in ',trim(modname)
83  write(*,*) 'Problem with dimensions :'
84  write(*,*) 'nlayer     = ',nlayer
85  write(*,*) 'klev   = ',klev
86  abort_message = ''
87  CALL abort_gcm (modname,abort_message,1)
88ENDIF
89
90!IF (ngrid.NE.klon_glo) THEN
91!  write(*,*) 'STOP in ',trim(modname)
92!  write(*,*) 'Problem with dimensions :'
93!  write(*,*) 'ngrid     = ',ngrid
94!  write(*,*) 'klon   = ',klon_glo
95!  abort_message = ''
96!  CALL abort_gcm (modname,abort_message,1)
97!ENDIF
98
99!call init_phys_lmdz(iim,jjm+1,llm,1,(/(jjm-1)*iim+2/))
100
101! init regular global longitude-latitude grid points and boundaries
102ALLOCATE(boundslon_reg(ii,2))
103ALLOCATE(boundslat_reg(jj+1,2))
104 
105DO i=1,ii
106   boundslon_reg(i,east)=rlonudyn(i)
107   boundslon_reg(i,west)=rlonudyn(i+1)
108ENDDO
109
110boundslat_reg(1,north)= PI/2
111boundslat_reg(1,south)= rlatvdyn(1)
112DO j=2,jj
113   boundslat_reg(j,north)=rlatvdyn(j-1)
114   boundslat_reg(j,south)=rlatvdyn(j)
115ENDDO
116boundslat_reg(jj+1,north)= rlatvdyn(jj)
117boundslat_reg(jj+1,south)= -PI/2
118
119! Write values in module regular_lonlat_mod
120CALL init_regular_lonlat(ii,jj+1, rlonvdyn(1:ii), rlatudyn, &
121                         boundslon_reg, boundslat_reg)
122
123! Generate global arrays on full physics grid
124allocate(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
125latfi(1)=rlatudyn(1)
126lonfi(1)=0.
127cufi(1) = cudyn(1)
128cvfi(1) = cvdyn(1)
129DO j=2,jj
130  DO i=1,ii
131    latfi((j-2)*ii+1+i)= rlatudyn(j)
132    lonfi((j-2)*ii+1+i)= rlonvdyn(i)
133    cufi((j-2)*ii+1+i) = cudyn((j-1)*(ii+1)+i)
134    cvfi((j-2)*ii+1+i) = cvdyn((j-1)*(ii+1)+i)
135  ENDDO
136ENDDO
137latfi(klon_glo)= rlatudyn(jj+1)
138lonfi(klon_glo)= 0.
139cufi(klon_glo) = cudyn((ii+1)*jj+1)
140cvfi(klon_glo) = cvdyn((ii+1)*jj-ii)
141
142! build airefi(), mesh area on physics grid
143allocate(airefi(klon_glo))
144CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,airedyn,airefi)
145! Poles are single points on physics grid
146airefi(1)=sum(airedyn(1:ii,1))
147airefi(klon_glo)=sum(airedyn(1:ii,jj+1))
148
149! Sanity check: do total planet area match between physics and dynamics?
150total_area_dyn=sum(airedyn(1:ii,1:jj+1))
151total_area_phy=sum(airefi(1:klon_glo))
152IF (total_area_dyn/=total_area_phy) THEN
153  WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
154  WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
155  WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
156  IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
157    ! stop here if the relative difference is more than 0.001%
158    abort_message = 'planet total surface discrepancy'
159    CALL abort_gcm(modname, abort_message, 1)
160  ENDIF
161ENDIF
162
163
164!$OMP PARALLEL
165! Now generate local lon/lat/cu/cv/area arrays
166call initcomgeomphy
167
168!!!!$OMP PARALLEL PRIVATE(ibegin,iend) &
169!!!     !$OMP SHARED(airefi,cufi,cvfi,lonfi,latfi)
170
171offset=klon_mpi_begin-1
172airephy(1:klon_omp)=airefi(offset+klon_omp_begin:offset+klon_omp_end)
173cuphy(1:klon_omp)=cufi(offset+klon_omp_begin:offset+klon_omp_end)
174cvphy(1:klon_omp)=cvfi(offset+klon_omp_begin:offset+klon_omp_end)
175rlond(1:klon_omp)=lonfi(offset+klon_omp_begin:offset+klon_omp_end)
176rlatd(1:klon_omp)=latfi(offset+klon_omp_begin:offset+klon_omp_end)
177
178! copy over preff , ap() and bp()
179call ini_planete_mod(nlayer,preff,ap,bp)
180
181! for slab ocean, copy over some arrays
182ok_slab_ocean=.false. ! default value
183call getin_p("ok_slab_ocean",ok_slab_ocean)
184if (ok_slab_ocean) then
185  call ini_surf_heat_transp(ip1jm,ip1jmp1,unsairez,fext,unsaire, &
186                            cu,cuvsurcv,cv,cvusurcu,aire,apoln,apols, &
187                            aireu,airev)
188endif
189
190! copy some fundamental parameters to physics
191! and do some initializations
192call inifis(klon_omp,nlayer,nqtot,pdayref,punjours,nday,ptimestep, &
193            rlatd,rlond,airephy,prad,pg,pr,pcpp)
194
195!$OMP END PARALLEL
196
197
198end subroutine iniphysiq
199
200
201END MODULE iniphysiq_mod
Note: See TracBrowser for help on using the repository browser.