source: dynamico_lmdz/aquaplanet/LMDZ5/libf/dynlonlat_phylonlat/phylmd/iniphysiq.F90 @ 3831

Last change on this file since 3831 was 3831, checked in by ymipsl, 9 years ago

module reorganisation for a cleaner dyn-phys interface
YM

File size: 9.9 KB
Line 
1
2! $Id: iniphysiq.F90 2242 2015-03-24 08:08:43Z emillour $
3
4
5SUBROUTINE iniphysiq(ii, jj, nbp, communicator, nlayer,punjours, pdayref,ptimestep,         &
6                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv,               &
7                     prad,pg,pr,pcpp,iflag_phys)
8 
9  USE dimphy, ONLY: klev ! number of atmospheric levels
10  USE mod_grid_phy_lmdz, ONLY: klon_glo,  & ! number of atmospheric columns (on full grid)
11                               regular_lonlat  ! regular longitude-latitude grid type
12  USE 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)
16  USE geometry_mod, ONLY : init_geometry
17  USE vertical_layers_mod, ONLY : init_vertical_layers
18  USE misc_mod, ONLY: debug
19  USE infotrac, ONLY: nqtot,nqo,nbtr,tname,ttext,type_trac,&
20                      niadv,conv_flg,pbl_flg,solsym
21  USE control_mod, ONLY: dayref,anneeref,day_step,iphysiq,nday,&
22                         raz_date,offline
23  USE inifis_mod, ONLY: inifis
24  USE time_phylmdz_mod, ONLY: init_time
25  USE infotrac_phy, ONLY: init_infotrac_phy
26  USE phyaqua_mod, ONLY: iniaqua
27  USE physics_distribution_mod, ONLY : init_physics_distribution
28  USE regular_lonlat_mod, ONLY : init_regular_lonlat, east, west, north, south, north_east, north_west, south_west, south_east
29  USE mod_interface_dyn_phys, ONLY :  init_interface_dyn_phys
30  IMPLICIT NONE
31
32  ! =======================================================================
33  ! Initialisation of the physical constants and some positional and
34  ! geometrical arrays for the physics
35  ! =======================================================================
36
37!  include "YOMCST.h"
38  include "dimensions.h"
39  include "comvert.h"
40  include "comconst.h"
41  include "iniprint.h"
42  include "temps.h"
43
44  REAL, INTENT (IN) :: prad ! radius of the planet (m)
45  REAL, INTENT (IN) :: pg ! gravitational acceleration (m/s2)
46  REAL, INTENT (IN) :: pr ! ! reduced gas constant R/mu
47  REAL, INTENT (IN) :: pcpp ! specific heat Cp
48  REAL, INTENT (IN) :: punjours ! length (in s) of a standard day
49  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
50  INTEGER, INTENT (IN) :: ii ! number of atmospheric columns along longitudes
51  INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes
52  INTEGER, INTENT (IN) :: nbp ! number of physics points (local)
53  INTEGER, INTENT (IN) :: communicator ! mpi communicator
54  REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid
55  REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid
56  REAL, INTENT (IN) :: rlonv(ii+1) ! longitudes of the physics grid
57  REAL, INTENT (IN) :: rlonu(ii+1) ! longitude boundaries of the physics grid
58  REAL, INTENT (IN) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2)
59  REAL, INTENT (IN) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
60  REAL, INTENT (IN) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
61  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
62  REAL, INTENT (IN) :: ptimestep !physics time step (s)
63  INTEGER, INTENT (IN) :: iflag_phys ! type of physics to be called
64
65  INTEGER :: ibegin, iend, offset
66  INTEGER :: i,j,k
67  CHARACTER (LEN=20) :: modname = 'iniphysiq'
68  CHARACTER (LEN=80) :: abort_message
69  REAL :: total_area_phy, total_area_dyn
70
71  REAL,ALLOCATABLE :: boundslon_reg(:,:)
72  REAL,ALLOCATABLE :: boundslat_reg(:,:)
73
74  ! global array, on full physics grid:
75  REAL,ALLOCATABLE :: latfi_glo(:)
76  REAL,ALLOCATABLE :: lonfi_glo(:)
77  REAL,ALLOCATABLE :: cufi_glo(:)
78  REAL,ALLOCATABLE :: cvfi_glo(:)
79  REAL,ALLOCATABLE :: airefi_glo(:)
80  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
81  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
82  ! local arrays, on given MPI/OpenMP domain:
83  REAL,ALLOCATABLE,SAVE :: latfi(:)
84  REAL,ALLOCATABLE,SAVE :: lonfi(:)
85  REAL,ALLOCATABLE,SAVE :: cufi(:)
86  REAL,ALLOCATABLE,SAVE :: cvfi(:)
87  REAL,ALLOCATABLE,SAVE :: airefi(:)
88  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
89  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
90!$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
91  INTEGER :: itaufin_phy
92 
93
94  CALL init_physics_distribution(regular_lonlat, 4, nbp, ii, jj+1, nlayer, communicator)
95  CALL init_interface_dyn_phys
96
97
98!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99! init regular longitude-latitude grid
100
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  CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, boundslon_reg, boundslat_reg)
119 
120!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
121 
122  ! Generate global arrays on full physics grid
123  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo),cufi_glo(klon_glo),cvfi_glo(klon_glo))
124  ALLOCATE(airefi_glo(klon_glo))
125  ALLOCATE(boundslonfi_glo(klon_glo,4))
126  ALLOCATE(boundslatfi_glo(klon_glo,4))
127
128  IF (klon_glo>1) THEN ! general case
129    ! North pole
130    latfi_glo(1)=rlatu(1)
131    lonfi_glo(1)=0.
132    cufi_glo(1) = cu(1)
133    cvfi_glo(1) = cv(1)
134    boundslonfi_glo(1,north_east)=0
135    boundslatfi_glo(1,north_east)=PI/2
136    boundslonfi_glo(1,north_west)=2*PI
137    boundslatfi_glo(1,north_west)=PI/2
138    boundslonfi_glo(1,south_west)=2*PI
139    boundslatfi_glo(1,south_west)=rlatv(1)
140    boundslonfi_glo(1,south_east)=0
141    boundslatfi_glo(1,south_east)=rlatv(1)
142    DO j=2,jj
143      DO i=1,ii
144        k=(j-2)*ii+1+i
145        latfi_glo(k)= rlatu(j)
146        lonfi_glo(k)= rlonv(i)
147        cufi_glo(k) = cu((j-1)*ii+1+i)
148        cvfi_glo(k) = cv((j-1)*ii+1+i)
149        boundslonfi_glo(k,north_east)=rlonu(i)
150        boundslatfi_glo(k,north_east)=rlatv(j-1)
151        boundslonfi_glo(k,north_west)=rlonu(i+1)
152        boundslatfi_glo(k,north_west)=rlatv(j-1)
153        boundslonfi_glo(k,south_west)=rlonu(i+1)
154        boundslatfi_glo(k,south_west)=rlatv(j)
155        boundslonfi_glo(k,south_east)=rlonu(i)
156        boundslatfi_glo(k,south_east)=rlatv(j)
157      ENDDO
158    ENDDO
159    ! South pole
160    latfi_glo(klon_glo)= rlatu(jj+1)
161    lonfi_glo(klon_glo)= 0.
162    cufi_glo(klon_glo) = cu((ii+1)*jj+1)
163    cvfi_glo(klon_glo) = cv((ii+1)*jj-ii)
164    boundslonfi_glo(klon_glo,north_east)= 0
165    boundslatfi_glo(klon_glo,north_east)= rlatv(jj)
166    boundslonfi_glo(klon_glo,north_west)= 2*PI
167    boundslatfi_glo(klon_glo,north_west)= rlatv(jj)
168    boundslonfi_glo(klon_glo,south_west)= 2*PI
169    boundslatfi_glo(klon_glo,south_west)= -PI/2
170    boundslonfi_glo(klon_glo,south_east)= 0
171    boundslatfi_glo(klon_glo,south_east)= -Pi/2
172
173    ! build airefi(), mesh area on physics grid
174    CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi_glo)
175    ! Poles are single points on physics grid
176    airefi_glo(1)=sum(aire(1:ii,1))
177    airefi_glo(klon_glo)=sum(aire(1:ii,jj+1))
178
179    ! Sanity check: do total planet area match between physics and dynamics?
180    total_area_dyn=sum(aire(1:ii,1:jj+1))
181    total_area_phy=sum(airefi_glo(1:klon_glo))
182    IF (total_area_dyn/=total_area_phy) THEN
183      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
184      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
185      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
186      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
187        ! stop here if the relative difference is more than 0.001%
188        abort_message = 'planet total surface discrepancy'
189        CALL abort_gcm(modname, abort_message, 1)
190      ENDIF
191    ENDIF
192  ELSE ! klon_glo==1, running the 1D model
193    ! just copy over input values
194    latfi_glo(1)=rlatu(1)
195    lonfi_glo(1)=rlonv(1)
196    cufi_glo(1)=cu(1)
197    cvfi_glo(1)=cv(1)
198    airefi_glo(1)=aire(1,1)
199    boundslonfi_glo(1,north_east)=rlonu(1)
200    boundslatfi_glo(1,north_east)=PI/2
201    boundslonfi_glo(1,north_west)=rlonu(2)
202    boundslatfi_glo(1,north_west)=PI/2
203    boundslonfi_glo(1,south_west)=rlonu(2)
204    boundslatfi_glo(1,south_west)=rlatv(1)
205    boundslonfi_glo(1,south_east)=rlonu(1)
206    boundslatfi_glo(1,south_east)=rlatv(1)
207  ENDIF ! of IF (klon_glo>1)
208
209!$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/)
210  ! Now generate local lon/lat/cu/cv/area/bounds arrays
211  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
212  ALLOCATE(airefi(klon_omp))
213  ALLOCATE(boundslonfi(klon_omp,4))
214  ALLOCATE(boundslatfi(klon_omp,4))
215
216
217  offset = klon_mpi_begin - 1
218  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
219  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
220  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
221  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
222  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
223  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
224  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
225
226  ! copy over global grid longitudes and latitudes
227  CALL init_geometry(lonfi, latfi, boundslonfi, boundslatfi, airefi, cufi, cvfi)
228
229 
230  ! copy over preff , ap(), bp(), etc
231  CALL init_vertical_layers(nlayer,preff,ap,bp,presnivs,pseudoalt)
232
233  ! Initialize tracer names, numbers, etc. for physics
234  CALL init_infotrac_phy(nqtot,nqo,nbtr,tname,ttext,type_trac,&
235                         niadv,conv_flg,pbl_flg,solsym)
236
237  ! transfer some flags/infos from dynamics to physics
238 
239  itaufin_phy=itaufin/iphysiq
240
241  CALL inifis(punjours,prad,pg,pr,pcpp)
242
243  CALL init_time(annee_ref, day_ref, day_ini, start_time, nday, ptimestep)
244 
245  ! Additional initializations for aquaplanets
246  IF (iflag_phys>=100) THEN
247    CALL iniaqua(klon_omp, iflag_phys)
248  END IF
249!$OMP END PARALLEL
250
251END SUBROUTINE iniphysiq
Note: See TracBrowser for help on using the repository browser.