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 !
5 
6 CONTAINS
7 
8 SUBROUTINE iniphysiq(iim,jjm,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, type_trac
23 ! USE comcstphy, ONLY: rradius, & ! planet radius (m)
24 ! rr, & ! recuced gas constant: R/molar mass of atm
25 ! rg, & ! gravity
26 ! rcpp ! specific heat of the atmosphere
27  USE inifis_mod, ONLY: inifis
28  USE phyaqua_mod, ONLY: iniaqua
31  east, west, north, south, &
35  USE nrtype, ONLY: pi
36  IMPLICIT NONE
37  !
38  !=======================================================================
39  ! Initialisation of the physical constants and some positional and
40  ! geometrical arrays for the physics
41  !=======================================================================
42 
43 
44  include "iniprint.h"
45 
46  REAL,INTENT(IN) :: prad ! radius of the planet (m)
47  REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
48  REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
49  REAL,INTENT(IN) :: pcpp ! specific heat Cp
50  REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
51  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
52  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
53  INTEGER, INTENT (IN) :: jjm ! number of atompsheric columns along latitudes
54  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
55  INTEGER, INTENT(IN) :: communicator ! MPI communicator
56  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
57  REAL, INTENT (IN) :: rlatv(jjm) ! latitude boundaries of the physics grid
58  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
59  REAL, INTENT (IN) :: rlonu(iim+1) ! longitude boundaries of the physics grid
60  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
61  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
62  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
63  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
64  REAL,INTENT(IN) :: ptimestep !physics time step (s)
65  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
66 
67  INTEGER :: ibegin,iend,offset
68  INTEGER :: i,j,k
69  CHARACTER (LEN=20) :: modname='iniphysiq'
70  CHARACTER (LEN=80) :: abort_message
71  REAL :: total_area_phy, total_area_dyn
72 
73  ! boundaries, on global grid
74  REAL,ALLOCATABLE :: boundslon_reg(:,:)
75  REAL,ALLOCATABLE :: boundslat_reg(:,:)
76 
77  ! global array, on full physics grid:
78  REAL,ALLOCATABLE :: latfi_glo(:)
79  REAL,ALLOCATABLE :: lonfi_glo(:)
80  REAL,ALLOCATABLE :: cufi_glo(:)
81  REAL,ALLOCATABLE :: cvfi_glo(:)
82  REAL,ALLOCATABLE :: airefi_glo(:)
83  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
84  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
85 
86  ! local arrays, on given MPI/OpenMP domain:
87  REAL,ALLOCATABLE,SAVE :: latfi(:)
88  REAL,ALLOCATABLE,SAVE :: lonfi(:)
89  REAL,ALLOCATABLE,SAVE :: cufi(:)
90  REAL,ALLOCATABLE,SAVE :: cvfi(:)
91  REAL,ALLOCATABLE,SAVE :: airefi(:)
92  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
93  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
94 !$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
95 
96  ! Initialize Physics distibution and parameters and interface with dynamics
97  CALL init_physics_distribution(regular_lonlat,4, &
98  nbp,iim,jjm+1,nlayer,communicator)
100 
101  ! init regular global longitude-latitude grid points and boundaries
102  ALLOCATE(boundslon_reg(iim,2))
103  ALLOCATE(boundslat_reg(jjm+1,2))
104 
105  DO i=1,iim
106  boundslon_reg(i,east)=rlonu(i)
107  boundslon_reg(i,west)=rlonu(i+1)
108  ENDDO
109 
110  boundslat_reg(1,north)= pi/2
111  boundslat_reg(1,south)= rlatv(1)
112  DO j=2,jjm
113  boundslat_reg(j,north)=rlatv(j-1)
114  boundslat_reg(j,south)=rlatv(j)
115  ENDDO
116  boundslat_reg(jjm+1,north)= rlatv(jjm)
117  boundslat_reg(jjm+1,south)= -pi/2
118 
119  ! Write values in module regular_lonlat_mod
120  CALL init_regular_lonlat(iim,jjm+1, rlonv(1:iim), rlatu, &
121  boundslon_reg, boundslat_reg)
122 
123  ! Generate global arrays on full physics grid
124  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
125  ALLOCATE(cufi_glo(klon_glo),cvfi_glo(klon_glo))
126  ALLOCATE(airefi_glo(klon_glo))
127  ALLOCATE(boundslonfi_glo(klon_glo,4))
128  ALLOCATE(boundslatfi_glo(klon_glo,4))
129 
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,jjm
144  DO i=1,iim
145  k=(j-2)*iim+1+i
146  latfi_glo(k)= rlatu(j)
147  lonfi_glo(k)= rlonv(i)
148  cufi_glo(k) = cu((j-1)*(iim+1)+i)
149  cvfi_glo(k) = cv((j-1)*(iim+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(jjm+1)
162  lonfi_glo(klon_glo)= 0.
163  cufi_glo(klon_glo) = cu((iim+1)*jjm+1)
164  cvfi_glo(klon_glo) = cv((iim+1)*jjm-iim)
165  boundslonfi_glo(klon_glo,north_east)= 0
166  boundslatfi_glo(klon_glo,north_east)= rlatv(jjm)
167  boundslonfi_glo(klon_glo,north_west)= 2*pi
168  boundslatfi_glo(klon_glo,north_west)= rlatv(jjm)
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,iim+1,jjm+1,klon_glo,aire,airefi_glo)
176  ! Poles are single points on physics grid
177  airefi_glo(1)=sum(aire(1:iim,1))
178  airefi_glo(klon_glo)=sum(aire(1:iim,jjm+1))
179 
180  ! Sanity check: do total planet area match between physics and dynamics?
181  total_area_dyn=sum(aire(1:iim,1:jjm+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 
194 !$OMP PARALLEL
195  ! Now generate local lon/lat/cu/cv/area/bounds arrays
196  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
197  ALLOCATE(airefi(klon_omp))
198  ALLOCATE(boundslonfi(klon_omp,4))
199  ALLOCATE(boundslatfi(klon_omp,4))
200 
201 
202  offset = klon_mpi_begin - 1
203  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
204  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
205  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
206  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
207  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
208  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
209  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
210 
211  ! copy over local grid longitudes and latitudes
212  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
213  airefi,cufi,cvfi)
214 
215  ! Initialize physical constants in physics:
216  CALL inifis(prad,pg,pr,pcpp)
217 
218  ! Initialize dimphy module
219  CALL init_dimphy(klon_omp,nlayer)
220 
221  ! Initialize tracer names, numbers, etc. for physics
223 
224  ! Additional initializations for aquaplanets
225  IF (iflag_phys>=100) THEN
226  CALL iniaqua(klon_omp,iflag_phys)
227  ENDIF
228 !$OMP END PARALLEL
229 
230 END SUBROUTINE iniphysiq
231 
232 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
!$Header!CDK comgeom COMMON comgeom rlatu
Definition: comgeom.h:25
subroutine inifis(prad, pg, pr, pcpp)
Definition: inifis_mod.F90:7
!$Id ysinus ok_gradsfile hybrid COMMON logici iflag_phys
Definition: logic.h:10
integer, save nqtot
Definition: infotrac.F90:6
subroutine iniaqua(nlon, iflag_phys)
Definition: phyaqua_mod.F90:11
integer, parameter east
subroutine init_dimphy(klon0, klev0)
Definition: dimphy.F90:20
integer, parameter west
integer, parameter north_west
Definition: nrtype.F90:1
!$Header!CDK comgeom COMMON comgeom rlonu
Definition: comgeom.h:25
!$Header!CDK comgeom COMMON comgeom rlatv
Definition: comgeom.h:25
character(len=4), save type_trac
Definition: infotrac.F90:40
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
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)
!$Header!CDK comgeom COMMON comgeom rlonv
Definition: comgeom.h:25