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

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

Some missing threadprivate...

YM

File size: 9.8 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 
92
93  CALL init_physics_distribution(regular_lonlat, 4, nbp, ii, jj+1, nlayer, communicator)
94  CALL init_interface_dyn_phys
95
96
97!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
98! init regular longitude-latitude grid
99
100  ALLOCATE(boundslon_reg(ii,2))
101  ALLOCATE(boundslat_reg(jj+1,2))
102 
103  DO i=1,ii
104   boundslon_reg(i,east)=rlonu(i)
105   boundslon_reg(i,west)=rlonu(i+1)
106  ENDDO
107
108  boundslat_reg(1,north)= PI/2
109  boundslat_reg(1,south)= rlatv(1)
110  DO j=2,jj
111   boundslat_reg(j,north)=rlatv(j-1)
112   boundslat_reg(j,south)=rlatv(j)
113  ENDDO
114  boundslat_reg(jj+1,north)= rlatv(jj)
115  boundslat_reg(jj+1,south)= -PI/2
116 
117  CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, boundslon_reg, boundslat_reg)
118 
119!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
120 
121  ! Generate global arrays on full physics grid
122  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo),cufi_glo(klon_glo),cvfi_glo(klon_glo))
123  ALLOCATE(airefi_glo(klon_glo))
124  ALLOCATE(boundslonfi_glo(klon_glo,4))
125  ALLOCATE(boundslatfi_glo(klon_glo,4))
126
127  IF (klon_glo>1) THEN ! general case
128    ! North pole
129    latfi_glo(1)=rlatu(1)
130    lonfi_glo(1)=0.
131    cufi_glo(1) = cu(1)
132    cvfi_glo(1) = cv(1)
133    boundslonfi_glo(1,north_east)=0
134    boundslatfi_glo(1,north_east)=PI/2
135    boundslonfi_glo(1,north_west)=2*PI
136    boundslatfi_glo(1,north_west)=PI/2
137    boundslonfi_glo(1,south_west)=2*PI
138    boundslatfi_glo(1,south_west)=rlatv(1)
139    boundslonfi_glo(1,south_east)=0
140    boundslatfi_glo(1,south_east)=rlatv(1)
141    DO j=2,jj
142      DO i=1,ii
143        k=(j-2)*ii+1+i
144        latfi_glo(k)= rlatu(j)
145        lonfi_glo(k)= rlonv(i)
146        cufi_glo(k) = cu((j-1)*ii+1+i)
147        cvfi_glo(k) = cv((j-1)*ii+1+i)
148        boundslonfi_glo(k,north_east)=rlonu(i)
149        boundslatfi_glo(k,north_east)=rlatv(j-1)
150        boundslonfi_glo(k,north_west)=rlonu(i+1)
151        boundslatfi_glo(k,north_west)=rlatv(j-1)
152        boundslonfi_glo(k,south_west)=rlonu(i+1)
153        boundslatfi_glo(k,south_west)=rlatv(j)
154        boundslonfi_glo(k,south_east)=rlonu(i)
155        boundslatfi_glo(k,south_east)=rlatv(j)
156      ENDDO
157    ENDDO
158    ! South pole
159    latfi_glo(klon_glo)= rlatu(jj+1)
160    lonfi_glo(klon_glo)= 0.
161    cufi_glo(klon_glo) = cu((ii+1)*jj+1)
162    cvfi_glo(klon_glo) = cv((ii+1)*jj-ii)
163    boundslonfi_glo(klon_glo,north_east)= 0
164    boundslatfi_glo(klon_glo,north_east)= rlatv(jj)
165    boundslonfi_glo(klon_glo,north_west)= 2*PI
166    boundslatfi_glo(klon_glo,north_west)= rlatv(jj)
167    boundslonfi_glo(klon_glo,south_west)= 2*PI
168    boundslatfi_glo(klon_glo,south_west)= -PI/2
169    boundslonfi_glo(klon_glo,south_east)= 0
170    boundslatfi_glo(klon_glo,south_east)= -Pi/2
171
172    ! build airefi(), mesh area on physics grid
173    CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi_glo)
174    ! Poles are single points on physics grid
175    airefi_glo(1)=sum(aire(1:ii,1))
176    airefi_glo(klon_glo)=sum(aire(1:ii,jj+1))
177
178    ! Sanity check: do total planet area match between physics and dynamics?
179    total_area_dyn=sum(aire(1:ii,1:jj+1))
180    total_area_phy=sum(airefi_glo(1:klon_glo))
181    IF (total_area_dyn/=total_area_phy) THEN
182      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
183      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
184      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
185      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
186        ! stop here if the relative difference is more than 0.001%
187        abort_message = 'planet total surface discrepancy'
188        CALL abort_gcm(modname, abort_message, 1)
189      ENDIF
190    ENDIF
191  ELSE ! klon_glo==1, running the 1D model
192    ! just copy over input values
193    latfi_glo(1)=rlatu(1)
194    lonfi_glo(1)=rlonv(1)
195    cufi_glo(1)=cu(1)
196    cvfi_glo(1)=cv(1)
197    airefi_glo(1)=aire(1,1)
198    boundslonfi_glo(1,north_east)=rlonu(1)
199    boundslatfi_glo(1,north_east)=PI/2
200    boundslonfi_glo(1,north_west)=rlonu(2)
201    boundslatfi_glo(1,north_west)=PI/2
202    boundslonfi_glo(1,south_west)=rlonu(2)
203    boundslatfi_glo(1,south_west)=rlatv(1)
204    boundslonfi_glo(1,south_east)=rlonu(1)
205    boundslatfi_glo(1,south_east)=rlatv(1)
206  ENDIF ! of IF (klon_glo>1)
207
208!$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/)
209  ! Now generate local lon/lat/cu/cv/area/bounds arrays
210  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
211  ALLOCATE(airefi(klon_omp))
212  ALLOCATE(boundslonfi(klon_omp,4))
213  ALLOCATE(boundslatfi(klon_omp,4))
214
215
216  offset = klon_mpi_begin - 1
217  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
218  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
219  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
220  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
221  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
222  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
223  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
224
225  ! copy over global grid longitudes and latitudes
226  CALL init_geometry(lonfi, latfi, boundslonfi, boundslatfi, airefi, cufi, cvfi)
227
228 
229  ! copy over preff , ap(), bp(), etc
230  CALL init_vertical_layers(nlayer,preff,ap,bp,presnivs,pseudoalt)
231
232  ! Initialize tracer names, numbers, etc. for physics
233  CALL init_infotrac_phy(nqtot,nqo,nbtr,tname,ttext,type_trac,&
234                         niadv,conv_flg,pbl_flg,solsym)
235
236  ! transfer some flags/infos from dynamics to physics
237 
238  CALL inifis(punjours,prad,pg,pr,pcpp)
239
240  CALL init_time(annee_ref, day_ref, day_ini, start_time, nday, ptimestep)
241 
242  ! Additional initializations for aquaplanets
243  IF (iflag_phys>=100) THEN
244    CALL iniaqua(klon_omp, iflag_phys)
245  END IF
246!$OMP END PARALLEL
247
248END SUBROUTINE iniphysiq
Note: See TracBrowser for help on using the repository browser.