LMDZ
iniphysiq_mod.F90
Go to the documentation of this file.
1 !
2 ! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $
3 !
4 MODULE iniphysiq_mod
5 
6 CONTAINS
7 
8 SUBROUTINE iniphysiq(ii,jj,nlayer, &
9  nbp, communicator, &
10  punjours, pdayref,ptimestep, &
12  prad,pg,pr,pcpp,iflag_phys)
13  USE dimphy, ONLY: init_dimphy
14  USE mod_grid_phy_lmdz, ONLY: klon_glo, & ! number of atmospheric columns (on full grid)
15  regular_lonlat ! regular longitude-latitude grid type
16  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
17  klon_omp_begin, & ! start index of local omp subgrid
18  klon_omp_end, & ! end index of local omp subgrid
19  klon_mpi_begin ! start indes of columns (on local mpi grid)
20  USE geometry_mod, ONLY : init_geometry
21  USE infotrac, ONLY: nqtot
22  USE comcstphy, ONLY: rradius, & ! planet radius (m)
23  rr, & ! recuced gas constant: R/molar mass of atm
24  rg, & ! gravity
25  rcpp ! specific heat of the atmosphere
26 ! USE phyaqua_mod, ONLY: iniaqua
29  east, west, north, south, &
34  USE nrtype, ONLY: pi
35  IMPLICIT NONE
36  !
37  !=======================================================================
38  ! Initialisation of the physical constants and some positional and
39  ! geometrical arrays for the physics
40  !=======================================================================
41 
42 
43  include "iniprint.h"
44 
45  REAL,INTENT(IN) :: prad ! radius of the planet (m)
46  REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
47  REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
48  REAL,INTENT(IN) :: pcpp ! specific heat Cp
49  REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
50  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
51  INTEGER, INTENT (IN) :: ii ! number of atmospheric coulumns along longitudes
52  INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes
53  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
54  INTEGER, INTENT(IN) :: communicator ! MPI communicator
55  REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid
56  REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid
57  REAL, INTENT (IN) :: rlonv(ii+1) ! longitudes of the physics grid
58  REAL, INTENT (IN) :: rlonu(ii+1) ! longitude boundaries of the physics grid
59  REAL, INTENT (IN) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2)
60  REAL, INTENT (IN) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
61  REAL, INTENT (IN) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
62  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
63  REAL,INTENT(IN) :: ptimestep !physics time step (s)
64  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
65 
66  INTEGER :: ibegin,iend,offset
67  INTEGER :: i,j,k
68  CHARACTER (LEN=20) :: modname='iniphysiq'
69  CHARACTER (LEN=80) :: abort_message
70  REAL :: total_area_phy, total_area_dyn
71 
72  ! boundaries, on global grid
73  REAL,ALLOCATABLE :: boundslon_reg(:,:)
74  REAL,ALLOCATABLE :: boundslat_reg(:,:)
75 
76  ! global array, on full physics grid:
77  REAL,ALLOCATABLE :: latfi_glo(:)
78  REAL,ALLOCATABLE :: lonfi_glo(:)
79  REAL,ALLOCATABLE :: cufi_glo(:)
80  REAL,ALLOCATABLE :: cvfi_glo(:)
81  REAL,ALLOCATABLE :: airefi_glo(:)
82  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
83  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
84 
85  ! local arrays, on given MPI/OpenMP domain:
86  REAL,ALLOCATABLE,SAVE :: latfi(:)
87  REAL,ALLOCATABLE,SAVE :: lonfi(:)
88  REAL,ALLOCATABLE,SAVE :: cufi(:)
89  REAL,ALLOCATABLE,SAVE :: cvfi(:)
90  REAL,ALLOCATABLE,SAVE :: airefi(:)
91  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
92  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
93 !$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
94 
95  ! Initialize Physics distibution and parameters and interface with dynamics
96  CALL init_physics_distribution(regular_lonlat,4, &
97  nbp,ii,jj+1,nlayer,communicator)
99 
100  ! init regular global longitude-latitude grid points and boundaries
101  ALLOCATE(boundslon_reg(ii,2))
102  ALLOCATE(boundslat_reg(jj+1,2))
103 
104  DO i=1,ii
105  boundslon_reg(i,east)=rlonu(i)
106  boundslon_reg(i,west)=rlonu(i+1)
107  ENDDO
108 
109  boundslat_reg(1,north)= pi/2
110  boundslat_reg(1,south)= rlatv(1)
111  DO j=2,jj
112  boundslat_reg(j,north)=rlatv(j-1)
113  boundslat_reg(j,south)=rlatv(j)
114  ENDDO
115  boundslat_reg(jj+1,north)= rlatv(jj)
116  boundslat_reg(jj+1,south)= -pi/2
117 
118  ! Write values in module regular_lonlat_mod
119  CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, &
120  boundslon_reg, boundslat_reg)
121 
122  ! Generate global arrays on full physics grid
123  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
124  ALLOCATE(cufi_glo(klon_glo),cvfi_glo(klon_glo))
125  ALLOCATE(airefi_glo(klon_glo))
126  ALLOCATE(boundslonfi_glo(klon_glo,4))
127  ALLOCATE(boundslatfi_glo(klon_glo,4))
128 
129  IF (klon_glo>1) THEN ! general case
130  ! North pole
131  latfi_glo(1)=rlatu(1)
132  lonfi_glo(1)=0.
133  cufi_glo(1) = cu(1)
134  cvfi_glo(1) = cv(1)
135  boundslonfi_glo(1,north_east)=0
136  boundslatfi_glo(1,north_east)=pi/2
137  boundslonfi_glo(1,north_west)=2*pi
138  boundslatfi_glo(1,north_west)=pi/2
139  boundslonfi_glo(1,south_west)=2*pi
140  boundslatfi_glo(1,south_west)=rlatv(1)
141  boundslonfi_glo(1,south_east)=0
142  boundslatfi_glo(1,south_east)=rlatv(1)
143  DO j=2,jj
144  DO i=1,ii
145  k=(j-2)*ii+1+i
146  latfi_glo((j-2)*ii+1+i)= rlatu(j)
147  lonfi_glo((j-2)*ii+1+i)= rlonv(i)
148  cufi_glo((j-2)*ii+1+i) = cu((j-1)*(ii+1)+i)
149  cvfi_glo((j-2)*ii+1+i) = cv((j-1)*(ii+1)+i)
150  boundslonfi_glo(k,north_east)=rlonu(i)
151  boundslatfi_glo(k,north_east)=rlatv(j-1)
152  boundslonfi_glo(k,north_west)=rlonu(i+1)
153  boundslatfi_glo(k,north_west)=rlatv(j-1)
154  boundslonfi_glo(k,south_west)=rlonu(i+1)
155  boundslatfi_glo(k,south_west)=rlatv(j)
156  boundslonfi_glo(k,south_east)=rlonu(i)
157  boundslatfi_glo(k,south_east)=rlatv(j)
158  ENDDO
159  ENDDO
160  ! South pole
161  latfi_glo(klon_glo)= rlatu(jj+1)
162  lonfi_glo(klon_glo)= 0.
163  cufi_glo(klon_glo) = cu((ii+1)*jj+1)
164  cvfi_glo(klon_glo) = cv((ii+1)*jj-ii)
165  boundslonfi_glo(klon_glo,north_east)= 0
166  boundslatfi_glo(klon_glo,north_east)= rlatv(jj)
167  boundslonfi_glo(klon_glo,north_west)= 2*pi
168  boundslatfi_glo(klon_glo,north_west)= rlatv(jj)
169  boundslonfi_glo(klon_glo,south_west)= 2*pi
170  boundslatfi_glo(klon_glo,south_west)= -pi/2
171  boundslonfi_glo(klon_glo,south_east)= 0
172  boundslatfi_glo(klon_glo,south_east)= -pi/2
173 
174  ! build airefi(), mesh area on physics grid
175  CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi_glo)
176  ! Poles are single points on physics grid
177  airefi_glo(1)=sum(aire(1:ii,1))
178  airefi_glo(klon_glo)=sum(aire(1:ii,jj+1))
179 
180  ! Sanity check: do total planet area match between physics and dynamics?
181  total_area_dyn=sum(aire(1:ii,1:jj+1))
182  total_area_phy=sum(airefi_glo(1:klon_glo))
183  IF (total_area_dyn/=total_area_phy) THEN
184  WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
185  WRITE (lunout, *) ' in the dynamics total_area_dyn=', total_area_dyn
186  WRITE (lunout, *) ' but in the physics total_area_phy=', total_area_phy
187  IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
188  ! stop here if the relative difference is more than 0.001%
189  abort_message = 'planet total surface discrepancy'
190  CALL abort_gcm(modname, abort_message, 1)
191  ENDIF
192  ENDIF
193  ELSE ! klon_glo==1, running the 1D model
194  ! just copy over input values
195  latfi_glo(1)=rlatu(1)
196  lonfi_glo(1)=rlonv(1)
197  cufi_glo(1)=cu(1)
198  cvfi_glo(1)=cv(1)
199  airefi_glo(1)=aire(1,1)
200  boundslonfi_glo(1,north_east)=rlonu(1)
201  boundslatfi_glo(1,north_east)=pi/2
202  boundslonfi_glo(1,north_west)=rlonu(2)
203  boundslatfi_glo(1,north_west)=pi/2
204  boundslonfi_glo(1,south_west)=rlonu(2)
205  boundslatfi_glo(1,south_west)=rlatv(1)
206  boundslonfi_glo(1,south_east)=rlonu(1)
207  boundslatfi_glo(1,south_east)=rlatv(1)
208  ENDIF ! of IF (klon_glo>1)
209 
210 !$OMP PARALLEL
211  ! Now generate local lon/lat/cu/cv/area/bounds arrays
212  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
213  ALLOCATE(airefi(klon_omp))
214  ALLOCATE(boundslonfi(klon_omp,4))
215  ALLOCATE(boundslatfi(klon_omp,4))
216 
217 
218  offset = klon_mpi_begin - 1
219  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
220  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
221  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
222  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
223  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
224  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
225  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
226 
227  ! copy over local grid longitudes and latitudes
228  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
229  airefi,cufi,cvfi)
230 
231 
232  ! Initialize tracer names, numbers, etc. for physics
234 
235 ! copy some fundamental parameters to physics
236  rradius=prad
237  rg=pg
238  rr=pr
239  rcpp=pcpp
240 
241 !$OMP END PARALLEL
242 
243 ! print*,'ATTENTION !!! TRAVAILLER SUR INIPHYSIQ'
244 ! print*,'CONTROLE DES LATITUDES, LONGITUDES, PARAMETRES ...'
245 
246 ! Additional initializations for aquaplanets
247 !!$OMP PARALLEL
248 ! if (iflag_phys>=100) then
249 ! call iniaqua(klon_omp,rlatd,rlond,iflag_phys)
250 ! endif
251 !!$OMP END PARALLEL
252 
253 END SUBROUTINE iniphysiq
254 
255 END MODULE iniphysiq_mod
subroutine iniphysiq(iim, jjm, nlayer, nbp, communicator, punjours, pdayref, ptimestep, rlatu, rlatv, rlonu, rlonv, aire, cu, cv, prad, pg, pr, pcpp, iflag_phys)
integer, parameter north
integer, parameter south
integer, save klon_glo
subroutine abort_gcm(modname, message, ierr)
Definition: abort_gcm.F:7
!$Id mode_top_bound COMMON comconstr && pi
Definition: comconst.h:7
!$Header!CDK comgeom COMMON comgeom aire
Definition: comgeom.h:25
c c $Id c c calculs statistiques distribution nuage ftion du regime dynamique c c Ce calcul doit etre fait a partir de valeurs mensuelles CALL nbregdyn DO k
Definition: calcul_REGDYN.h:12
!$Header!CDK comgeom COMMON comgeom rlatu
Definition: comgeom.h:25
!$Id ysinus ok_gradsfile hybrid COMMON logici iflag_phys
Definition: logic.h:10
!$Id klon initialisation mois suivants day_rain itap ENDIF!Calcul fin de nday_rain calcul nday_rain itap DO i
Definition: calcul_divers.h:24
integer, save nqtot
Definition: infotrac.F90:6
integer, parameter east
subroutine init_dimphy(klon0, klev0)
Definition: dimphy.F90:20
real rradius
Definition: comcstphy.F90:3
integer, parameter west
integer, parameter north_west
!$Id Turb_fcg!implicit none!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!cette routine permet d obtenir hq et ainsi de!pouvoir calculer la convergence et le cisaillement dans la physiq!cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc INTEGER klev REAL j
Definition: 1Dconv.h:27
Definition: nrtype.F90:1
!$Header!CDK comgeom COMMON comgeom rlonu
Definition: comgeom.h:25
!$Header!CDK comgeom COMMON comgeom rlatv
Definition: comgeom.h:25
subroutine gr_dyn_fi(nfield, im, jm, ngrid, pdyn, pfi)
Definition: gr_dyn_fi.F:5
!$Id!Parameters for parameters that control the rate of approach!to quasi equilibrium noff nlm real tlcrit real entp real sigd real coeffs real dtmax real cu real betad real damp real delta COMMON cvparam nlm tlcrit sigd coeffs cu
Definition: cvparam.h:12
subroutine init_geometry(klon, longitude_, latitude_, boundslon_, boundslat_, cell_area_, dx_, dy_)
subroutine init_infotrac_phy(nqtot_, type_trac_)
!$Header!CDK comgeom COMMON comgeom cv
Definition: comgeom.h:25
subroutine init_regular_lonlat(nbp_lon, nbp_lat, lon_reg_, lat_reg_, boundslon_reg_, boundslat_reg_)
integer, parameter north_east
real rr
Definition: comcstphy.h:1
integer, parameter south_east
Definition: dimphy.F90:1
integer, parameter south_west
!$Header!gestion des impressions de sorties et de débogage la sortie standard prt_level COMMON comprint lunout
Definition: iniprint.h:7
subroutine init_physics_distribution(grid_type, nvertex, nbp, nbp_lon, nbp_lat, nbp_lev, communicator)
real rg
Definition: comcstphy.h:1
!$Header!CDK comgeom COMMON comgeom rlonv
Definition: comgeom.h:25