source: LMDZ5/trunk/libf/dynlonlat_phylonlat/phymar/iniphysiq_mod.F90 @ 2351

Last change on this file since 2351 was 2351, checked in by Ehouarn Millour, 9 years ago

More on physics/dynamics separation and cleanup:

  • Set things up so that all physics-related initializations are done via iniphysiq.
  • Created a "geometry_mod.F90" module in phy_common to store information on the loacl grid (i.e. replaces comgeomphy) and moreover give these variables more obvious names (e.g.: rlond => longitude, rlatd => latitude, airephy => cell_area).
  • removed obsolete comgeomphy.h and comgeomphy.F90

EM

File size: 9.8 KB
Line 
1!
2! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $
3!
4MODULE iniphysiq_mod
5
6CONTAINS
7
8SUBROUTINE iniphysiq(iim,jjm,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  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 comcstphy, ONLY: rradius, & ! planet radius (m)
22                       rr, & ! recuced gas constant: R/molar mass of atm
23                       rg, & ! gravity
24                       rcpp  ! specific heat of the atmosphere
25!  USE phyaqua_mod, ONLY: iniaqua
26  USE physics_distribution_mod, ONLY : init_physics_distribution
27  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
28                                 east, west, north, south, &
29                                 north_east, north_west, &
30                                 south_west, south_east
31  USE mod_interface_dyn_phys, ONLY :  init_interface_dyn_phys
32  USE nrtype, ONLY: pi
33  IMPLICIT NONE
34  !
35  !=======================================================================
36  !   Initialisation of the physical constants and some positional and
37  !   geometrical arrays for the physics
38  !=======================================================================
39 
40 
41  include "iniprint.h"
42
43  REAL,INTENT(IN) :: prad ! radius of the planet (m)
44  REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
45  REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
46  REAL,INTENT(IN) :: pcpp ! specific heat Cp
47  REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
48  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
49  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
50  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
51  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
52  INTEGER, INTENT(IN) :: communicator ! MPI communicator
53  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
54  REAL, INTENT (IN) :: rlatv(jjm) ! latitude boundaries of the physics grid
55  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
56  REAL, INTENT (IN) :: rlonu(iim+1) ! longitude boundaries of the physics grid
57  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
58  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
59  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
60  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
61  REAL,INTENT(IN) :: ptimestep !physics time step (s)
62  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
63
64  INTEGER :: ibegin,iend,offset
65  INTEGER :: i,j,k
66  CHARACTER (LEN=20) :: modname='iniphysiq'
67  CHARACTER (LEN=80) :: abort_message
68  REAL :: total_area_phy, total_area_dyn
69
70  ! boundaries, on global grid
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
83  ! local arrays, on given MPI/OpenMP domain:
84  REAL,ALLOCATABLE,SAVE :: latfi(:)
85  REAL,ALLOCATABLE,SAVE :: lonfi(:)
86  REAL,ALLOCATABLE,SAVE :: cufi(:)
87  REAL,ALLOCATABLE,SAVE :: cvfi(:)
88  REAL,ALLOCATABLE,SAVE :: airefi(:)
89  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
90  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
91!$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
92
93  ! Initialize Physics distibution and parameters and interface with dynamics
94  CALL init_physics_distribution(regular_lonlat,4, &
95                                 nbp,ii,jj+1,nlayer,communicator)
96  CALL init_interface_dyn_phys
97 
98  ! init regular global longitude-latitude grid points and boundaries
99  ALLOCATE(boundslon_reg(iim,2))
100  ALLOCATE(boundslat_reg(jjm+1,2))
101 
102  DO i=1,iim
103   boundslon_reg(i,east)=rlonu(i)
104   boundslon_reg(i,west)=rlonu(i+1)
105  ENDDO
106
107  boundslat_reg(1,north)= PI/2
108  boundslat_reg(1,south)= rlatv(1)
109  DO j=2,jjm
110   boundslat_reg(j,north)=rlatv(j-1)
111   boundslat_reg(j,south)=rlatv(j)
112  ENDDO
113  boundslat_reg(jjm+1,north)= rlatv(jjm)
114  boundslat_reg(jjm+1,south)= -PI/2
115
116  ! Write values in module regular_lonlat_mod
117  CALL init_regular_lonlat(iim,jjm+1, rlonv(1:iim), rlatu, &
118                           boundslon_reg, boundslat_reg)
119
120  ! Generate global arrays on full physics grid
121  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
122  ALLOCATE(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,jjm
142      DO i=1,iim
143        k=(j-2)*iim+1+i
144        latfi_glo((j-2)*iim+1+i)= rlatu(j)
145        lonfi_glo((j-2)*iim+1+i)= rlonv(i)
146        cufi_glo((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
147        cvfi_glo((j-2)*iim+1+i) = cv((j-1)*iim+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(jjm+1)
160    lonfi_glo(klon_glo)= 0.
161    cufi_glo(klon_glo) = cu((iim+1)*jjm+1)
162    cvfi_glo(klon_glo) = cv((iim+1)*jjm-iim)
163    boundslonfi_glo(klon_glo,north_east)= 0
164    boundslatfi_glo(klon_glo,north_east)= rlatv(jjm)
165    boundslonfi_glo(klon_glo,north_west)= 2*PI
166    boundslatfi_glo(klon_glo,north_west)= rlatv(jjm)
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,iim+1,jjm+1,klon_glo,aire,airefi_glo)
174    ! Poles are single points on physics grid
175    airefi_glo(1)=sum(aire(1:iim,1))
176    airefi_glo(klon_glo)=sum(aire(1:iim,jjm+1))
177
178    ! Sanity check: do total planet area match between physics and dynamics?
179    total_area_dyn=sum(aire(1:iim,1:jjm+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
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 local grid longitudes and latitudes
226  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
227                     airefi,cufi,cvfi)
228
229
230! copy some fundamental parameters to physics
231  rradius=prad
232  rg=pg
233  rr=pr
234  rcpp=pcpp
235
236!$OMP END PARALLEL
237
238!      print*,'ATTENTION !!! TRAVAILLER SUR INIPHYSIQ'
239!      print*,'CONTROLE DES LATITUDES, LONGITUDES, PARAMETRES ...'
240
241! Additional initializations for aquaplanets
242!!$OMP PARALLEL
243!      if (iflag_phys>=100) then
244!        call iniaqua(klon_omp,rlatd,rlond,iflag_phys)
245!      endif
246!!$OMP END PARALLEL
247
248END SUBROUTINE iniphysiq
249
250END MODULE iniphysiq_mod
Note: See TracBrowser for help on using the repository browser.