source: LMDZ5/branches/AI-cosp/libf/dynphy_lonlat/phylmd/iniphysiq_mod.F90 @ 2932

Last change on this file since 2932 was 2396, checked in by Laurent Fairhead, 9 years ago

Modifications from 2274 that were somehow lost
LF

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 12.4 KB
Line 
1!
2! $Id: iniphysiq_mod.F90 2396 2015-11-18 14:12:46Z fairhead $
3!
4MODULE iniphysiq_mod
5
6CONTAINS
7
8SUBROUTINE iniphysiq(ii,jj,nlayer, &
9                     nbp, communicator, &
10                     punjours, pdayref,ptimestep, &
11                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv,       &
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                               nbp_lon, nbp_lat, nbp_lev
17  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
18                                klon_omp_begin, & ! start index of local omp subgrid
19                                klon_omp_end, & ! end index of local omp subgrid
20                                klon_mpi_begin ! start indes of columns (on local mpi grid)
21  USE geometry_mod, ONLY : init_geometry
22  USE vertical_layers_mod, ONLY : init_vertical_layers
23  USE infotrac, ONLY: nqtot,nqo,nbtr,tname,ttext,type_trac,&
24                      niadv,conv_flg,pbl_flg,solsym,&
25                      nqfils,nqdesc,nqdesc_tot,iqfils,iqpere,&
26                      ok_isotopes,ok_iso_verif,ok_isotrac,&
27                      ok_init_iso,niso_possibles,tnat,&
28                      alpha_ideal,use_iso,iqiso,iso_num,&
29                      iso_indnum,zone_num,phase_num,&
30                      indnum_fn_num,index_trac,&
31                      niso,ntraceurs_zone,ntraciso
32#ifdef REPROBUS
33  USE CHEM_REP, ONLY : Init_chem_rep_phys
34#endif
35  USE control_mod, ONLY: dayref,anneeref,day_step,nday,offline, iphysiq
36  USE inifis_mod, ONLY: inifis
37  USE time_phylmdz_mod, ONLY: init_time
38  USE infotrac_phy, ONLY: init_infotrac_phy
39  USE phystokenc_mod, ONLY: init_phystokenc
40  USE phyaqua_mod, ONLY: iniaqua
41  USE physics_distribution_mod, ONLY : init_physics_distribution
42  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
43                                 east, west, north, south, &
44                                 north_east, north_west, &
45                                 south_west, south_east
46  USE mod_interface_dyn_phys, ONLY :  init_interface_dyn_phys
47#ifdef INCA
48  USE indice_sol_mod, ONLY: nbsrf, is_oce, is_sic, is_ter, is_lic
49  USE parallel_lmdz, ONLY : mpi_size
50  USE mod_const_mpi, ONLY : COMM_LMDZ
51  USE bands, ONLY : distrib_phys
52  USE mod_phys_lmdz_omp_data, ONLY: klon_omp
53#endif
54  IMPLICIT NONE
55
56  ! =======================================================================
57  ! Initialisation of the physical constants and some positional and
58  ! geometrical arrays for the physics
59  ! =======================================================================
60
61  include "dimensions.h"
62  include "comvert.h"
63  include "comconst.h"
64  include "iniprint.h"
65  include "temps.h"
66  include "tracstoke.h"
67
68  REAL, INTENT (IN) :: prad ! radius of the planet (m)
69  REAL, INTENT (IN) :: pg ! gravitational acceleration (m/s2)
70  REAL, INTENT (IN) :: pr ! ! reduced gas constant R/mu
71  REAL, INTENT (IN) :: pcpp ! specific heat Cp
72  REAL, INTENT (IN) :: punjours ! length (in s) of a standard day
73  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
74  INTEGER, INTENT (IN) :: ii ! number of atmospheric columns along longitudes
75  INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes
76  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
77  INTEGER, INTENT(IN) :: communicator ! MPI communicator
78  REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid
79  REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid
80  REAL, INTENT (IN) :: rlonv(ii+1) ! longitudes of the physics grid
81  REAL, INTENT (IN) :: rlonu(ii+1) ! longitude boundaries of the physics grid
82  REAL, INTENT (IN) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2)
83  REAL, INTENT (IN) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
84  REAL, INTENT (IN) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
85  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
86  REAL, INTENT (IN) :: ptimestep !physics time step (s)
87  INTEGER, INTENT (IN) :: iflag_phys ! type of physics to be called
88
89  INTEGER :: ibegin, iend, offset
90  INTEGER :: i,j,k
91  CHARACTER (LEN=20) :: modname = 'iniphysiq'
92  CHARACTER (LEN=80) :: abort_message
93  REAL :: total_area_phy, total_area_dyn
94
95  ! boundaries, on global grid
96  REAL,ALLOCATABLE :: boundslon_reg(:,:)
97  REAL,ALLOCATABLE :: boundslat_reg(:,:)
98
99  ! global array, on full physics grid:
100  REAL,ALLOCATABLE :: latfi_glo(:)
101  REAL,ALLOCATABLE :: lonfi_glo(:)
102  REAL,ALLOCATABLE :: cufi_glo(:)
103  REAL,ALLOCATABLE :: cvfi_glo(:)
104  REAL,ALLOCATABLE :: airefi_glo(:)
105  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
106  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
107
108  ! local arrays, on given MPI/OpenMP domain:
109  REAL,ALLOCATABLE,SAVE :: latfi(:)
110  REAL,ALLOCATABLE,SAVE :: lonfi(:)
111  REAL,ALLOCATABLE,SAVE :: cufi(:)
112  REAL,ALLOCATABLE,SAVE :: cvfi(:)
113  REAL,ALLOCATABLE,SAVE :: airefi(:)
114  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
115  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
116!$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
117
118  ! Initialize Physics distibution and parameters and interface with dynamics
119  IF (ii*jj>1) THEN ! general 3D case
120    CALL init_physics_distribution(regular_lonlat,4, &
121                                 nbp,ii,jj+1,nlayer,communicator)
122  ELSE ! For 1D model
123    CALL init_physics_distribution(regular_lonlat,4, &
124                                 1,1,1,nlayer,communicator) 
125  ENDIF
126  CALL init_interface_dyn_phys
127 
128  ! init regular global longitude-latitude grid points and boundaries
129  ALLOCATE(boundslon_reg(ii,2))
130  ALLOCATE(boundslat_reg(jj+1,2))
131 
132  DO i=1,ii
133   boundslon_reg(i,east)=rlonu(i)
134   boundslon_reg(i,west)=rlonu(i+1)
135  ENDDO
136
137  boundslat_reg(1,north)= PI/2
138  boundslat_reg(1,south)= rlatv(1)
139  DO j=2,jj
140   boundslat_reg(j,north)=rlatv(j-1)
141   boundslat_reg(j,south)=rlatv(j)
142  ENDDO
143  boundslat_reg(jj+1,north)= rlatv(jj)
144  boundslat_reg(jj+1,south)= -PI/2
145
146  ! Write values in module regular_lonlat_mod
147  CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, &
148                           boundslon_reg, boundslat_reg)
149
150  ! Generate global arrays on full physics grid
151  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
152  ALLOCATE(cufi_glo(klon_glo),cvfi_glo(klon_glo))
153  ALLOCATE(airefi_glo(klon_glo))
154  ALLOCATE(boundslonfi_glo(klon_glo,4))
155  ALLOCATE(boundslatfi_glo(klon_glo,4))
156 
157
158  IF (klon_glo>1) THEN ! general case
159    ! North pole
160    latfi_glo(1)=rlatu(1)
161    lonfi_glo(1)=0.
162    cufi_glo(1) = cu(1)
163    cvfi_glo(1) = cv(1)
164    boundslonfi_glo(1,north_east)=0
165    boundslatfi_glo(1,north_east)=PI/2
166    boundslonfi_glo(1,north_west)=2*PI
167    boundslatfi_glo(1,north_west)=PI/2
168    boundslonfi_glo(1,south_west)=2*PI
169    boundslatfi_glo(1,south_west)=rlatv(1)
170    boundslonfi_glo(1,south_east)=0
171    boundslatfi_glo(1,south_east)=rlatv(1)
172    DO j=2,jj
173      DO i=1,ii
174        k=(j-2)*ii+1+i
175        latfi_glo(k)= rlatu(j)
176        lonfi_glo(k)= rlonv(i)
177        cufi_glo(k) = cu((j-1)*(ii+1)+i)
178        cvfi_glo(k) = cv((j-1)*(ii+1)+i)
179        boundslonfi_glo(k,north_east)=rlonu(i)
180        boundslatfi_glo(k,north_east)=rlatv(j-1)
181        boundslonfi_glo(k,north_west)=rlonu(i+1)
182        boundslatfi_glo(k,north_west)=rlatv(j-1)
183        boundslonfi_glo(k,south_west)=rlonu(i+1)
184        boundslatfi_glo(k,south_west)=rlatv(j)
185        boundslonfi_glo(k,south_east)=rlonu(i)
186        boundslatfi_glo(k,south_east)=rlatv(j)
187      ENDDO
188    ENDDO
189    ! South pole
190    latfi_glo(klon_glo)= rlatu(jj+1)
191    lonfi_glo(klon_glo)= 0.
192    cufi_glo(klon_glo) = cu((ii+1)*jj+1)
193    cvfi_glo(klon_glo) = cv((ii+1)*jj-ii)
194    boundslonfi_glo(klon_glo,north_east)= 0
195    boundslatfi_glo(klon_glo,north_east)= rlatv(jj)
196    boundslonfi_glo(klon_glo,north_west)= 2*PI
197    boundslatfi_glo(klon_glo,north_west)= rlatv(jj)
198    boundslonfi_glo(klon_glo,south_west)= 2*PI
199    boundslatfi_glo(klon_glo,south_west)= -PI/2
200    boundslonfi_glo(klon_glo,south_east)= 0
201    boundslatfi_glo(klon_glo,south_east)= -Pi/2
202
203    ! build airefi_glo(), mesh area on physics grid
204    CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi_glo)
205    ! Poles are single points on physics grid
206    airefi_glo(1)=sum(aire(1:ii,1))
207    airefi_glo(klon_glo)=sum(aire(1:ii,jj+1))
208
209    ! Sanity check: do total planet area match between physics and dynamics?
210    total_area_dyn=sum(aire(1:ii,1:jj+1))
211    total_area_phy=sum(airefi_glo(1:klon_glo))
212    IF (total_area_dyn/=total_area_phy) THEN
213      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
214      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
215      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
216      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
217        ! stop here if the relative difference is more than 0.001%
218        abort_message = 'planet total surface discrepancy'
219        CALL abort_gcm(modname, abort_message, 1)
220      ENDIF
221    ENDIF
222  ELSE ! klon_glo==1, running the 1D model
223    ! just copy over input values
224    latfi_glo(1)=rlatu(1)
225    lonfi_glo(1)=rlonv(1)
226    cufi_glo(1)=cu(1)
227    cvfi_glo(1)=cv(1)
228    airefi_glo(1)=aire(1,1)
229    boundslonfi_glo(1,north_east)=rlonu(1)
230    boundslatfi_glo(1,north_east)=PI/2
231    boundslonfi_glo(1,north_west)=rlonu(2)
232    boundslatfi_glo(1,north_west)=PI/2
233    boundslonfi_glo(1,south_west)=rlonu(2)
234    boundslatfi_glo(1,south_west)=rlatv(1)
235    boundslonfi_glo(1,south_east)=rlonu(1)
236    boundslatfi_glo(1,south_east)=rlatv(1)
237  ENDIF ! of IF (klon_glo>1)
238
239!$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/)
240  ! Now generate local lon/lat/cu/cv/area/bounds arrays
241  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
242  ALLOCATE(airefi(klon_omp))
243  ALLOCATE(boundslonfi(klon_omp,4))
244  ALLOCATE(boundslatfi(klon_omp,4))
245
246
247  offset = klon_mpi_begin - 1
248  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
249  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
250  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
251  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
252  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
253  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
254  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
255
256  ! copy over local grid longitudes and latitudes
257  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
258                     airefi,cufi,cvfi)
259
260  ! copy over preff , ap(), bp(), etc
261  CALL init_vertical_layers(nlayer,preff,scaleheight, &
262                            ap,bp,presnivs,pseudoalt)
263
264  ! Initialize physical constants in physics:
265  CALL inifis(punjours,prad,pg,pr,pcpp)
266
267  CALL init_time(annee_ref,day_ref,day_ini,start_time,nday,ptimestep)
268
269  ! Initialize dimphy module (unless in 1D where it has already been done)
270  IF (klon_glo>1) CALL Init_dimphy(klon_omp,nlayer)
271
272  ! Copy over "offline" settings
273  CALL init_phystokenc(offline,istphy)
274
275  ! Initialize tracer names, numbers, etc. for physics
276  CALL init_infotrac_phy(nqtot,nqo,nbtr,tname,ttext,type_trac,&
277                         niadv,conv_flg,pbl_flg,solsym,&
278                         nqfils,nqdesc,nqdesc_tot,iqfils,iqpere,&
279                         ok_isotopes,ok_iso_verif,ok_isotrac,&
280                         ok_init_iso,niso_possibles,tnat,&
281                         alpha_ideal,use_iso,iqiso,iso_num,&
282                         iso_indnum,zone_num,phase_num,&
283                         indnum_fn_num,index_trac,&
284                         niso,ntraceurs_zone,ntraciso)
285
286  ! Initializations for Reprobus
287  IF (type_trac == 'repr') THEN
288#ifdef REPROBUS
289    CALL Init_chem_rep_phys(klon_omp,nlayer)
290#endif
291  ENDIF
292!$OMP END PARALLEL
293
294  IF (type_trac == 'inca') THEN
295#ifdef INCA
296     call init_const_lmdz( &
297          anneeref,dayref, &
298          iphysiq,day_step,nday,  &
299          nbsrf, is_oce,is_sic, &
300          is_ter,is_lic, calend)
301     call init_inca_para( &
302          nbp_lon,nbp_lat,nbp_lev,klon_glo,mpi_size, &
303          distrib_phys,COMM_LMDZ)
304#endif
305  END IF
306!$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/)
307
308  ! Additional initializations for aquaplanets
309  IF (iflag_phys>=100) THEN
310    CALL iniaqua(klon_omp,iflag_phys)
311  END IF
312
313  IF (type_trac == 'inca') THEN
314#ifdef INCA
315     CALL init_inca_dim(klon_omp,nbp_lev,nbp_lon,nbp_lat - 1, &
316          rlonu,rlatu,rlonv,rlatv)
317#endif
318  END IF
319
320!$OMP END PARALLEL
321
322END SUBROUTINE iniphysiq
323
324END MODULE iniphysiq_mod
Note: See TracBrowser for help on using the repository browser.