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

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

Make iniphysiq a module.
Fix call to iniphysiq in lmdz1d (missing arguments and arrays of wrong sizes).
EM

File size: 7.2 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,punjours, pdayref,ptimestep, &
9                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv,         &
10                     prad,pg,pr,pcpp,iflag_phys)
11  USE dimphy, ONLY: klev ! number of atmospheric levels
12  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
13                                        ! (on full grid)
14  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
15                                klon_omp_begin, & ! start index of local omp subgrid
16                                klon_omp_end, & ! end index of local omp subgrid
17                                klon_mpi_begin ! start indes of columns (on local mpi grid)
18  USE comgeomphy, ONLY: initcomgeomphy, &
19                        airephy, & ! physics grid area (m2)
20                        cuphy, & ! cu coeff. (u_covariant = cu * u)
21                        cvphy, & ! cv coeff. (v_covariant = cv * v)
22                        rlond, & ! longitudes
23                        rlatd ! latitudes
24  USE comcstphy, ONLY: rradius, & ! planet radius (m)
25                       rr, & ! recuced gas constant: R/molar mass of atm
26                       rg, & ! gravity
27                       rcpp  ! specific heat of the atmosphere
28!  USE phyaqua_mod, ONLY: iniaqua
29  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
30                                 east, west, north, south, &
31                                 north_east, north_west, &
32                                 south_west, south_east
33  USE nrtype, ONLY: pi
34  IMPLICIT NONE
35  !
36  !=======================================================================
37  !   Initialisation of the physical constants and some positional and
38  !   geometrical arrays for the physics
39  !=======================================================================
40 
41 
42  include "iniprint.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) :: iim ! number of atmospheric coulumns along longitudes
51  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
52  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
53  REAL, INTENT (IN) :: rlatv(jjm) ! latitude boundaries of the physics grid
54  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
55  REAL, INTENT (IN) :: rlonu(iim+1) ! longitude boundaries of the physics grid
56  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
57  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
58  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
59  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
60  REAL,INTENT(IN) :: ptimestep !physics time step (s)
61  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
62
63  INTEGER :: ibegin,iend,offset
64  INTEGER :: i,j
65  CHARACTER (LEN=20) :: modname='iniphysiq'
66  CHARACTER (LEN=80) :: abort_message
67  REAL :: total_area_phy, total_area_dyn
68
69  ! boundaries, on global grid
70  REAL,ALLOCATABLE :: boundslon_reg(:,:)
71  REAL,ALLOCATABLE :: boundslat_reg(:,:)
72
73  ! global array, on full physics grid:
74  REAL,ALLOCATABLE :: latfi(:)
75  REAL,ALLOCATABLE :: lonfi(:)
76  REAL,ALLOCATABLE :: cufi(:)
77  REAL,ALLOCATABLE :: cvfi(:)
78  REAL,ALLOCATABLE :: airefi(:)
79 
80  IF (nlayer.NE.klev) THEN
81    WRITE(lunout,*) 'STOP in ',trim(modname)
82    WRITE(lunout,*) 'Problem with dimensions :'
83    WRITE(lunout,*) 'nlayer     = ',nlayer
84    WRITE(lunout,*) 'klev   = ',klev
85    abort_message = ''
86    CALL abort_gcm (modname,abort_message,1)
87  ENDIF
88
89  ! init regular global longitude-latitude grid points and boundaries
90  ALLOCATE(boundslon_reg(iim,2))
91  ALLOCATE(boundslat_reg(jjm+1,2))
92 
93  DO i=1,iim
94   boundslon_reg(i,east)=rlonu(i)
95   boundslon_reg(i,west)=rlonu(i+1)
96  ENDDO
97
98  boundslat_reg(1,north)= PI/2
99  boundslat_reg(1,south)= rlatv(1)
100  DO j=2,jjm
101   boundslat_reg(j,north)=rlatv(j-1)
102   boundslat_reg(j,south)=rlatv(j)
103  ENDDO
104  boundslat_reg(jjm+1,north)= rlatv(jjm)
105  boundslat_reg(jjm+1,south)= -PI/2
106
107  ! Write values in module regular_lonlat_mod
108  CALL init_regular_lonlat(iim,jjm+1, rlonv(1:iim), rlatu, &
109                           boundslon_reg, boundslat_reg)
110
111  ! Generate global arrays on full physics grid
112  ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
113  ALLOCATE(airefi(klon_glo))
114
115  IF (klon_glo>1) THEN ! general case
116    ! North pole
117    latfi(1)=rlatu(1)
118    lonfi(1)=0.
119    cufi(1) = cu(1)
120    cvfi(1) = cv(1)
121    DO j=2,jjm
122      DO i=1,iim
123        latfi((j-2)*iim+1+i)= rlatu(j)
124        lonfi((j-2)*iim+1+i)= rlonv(i)
125        cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
126        cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
127      ENDDO
128    ENDDO
129    ! South pole
130    latfi(klon_glo)= rlatu(jjm+1)
131    lonfi(klon_glo)= 0.
132    cufi(klon_glo) = cu((iim+1)*jjm+1)
133    cvfi(klon_glo) = cv((iim+1)*jjm-iim)
134
135    ! build airefi(), mesh area on physics grid
136    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
137    ! Poles are single points on physics grid
138    airefi(1)=sum(aire(1:iim,1))
139    airefi(klon_glo)=sum(aire(1:iim,jjm+1))
140
141    ! Sanity check: do total planet area match between physics and dynamics?
142    total_area_dyn=sum(aire(1:iim,1:jjm+1))
143    total_area_phy=sum(airefi(1:klon_glo))
144    IF (total_area_dyn/=total_area_phy) THEN
145      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
146      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
147      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
148      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
149        ! stop here if the relative difference is more than 0.001%
150        abort_message = 'planet total surface discrepancy'
151        CALL abort_gcm(modname, abort_message, 1)
152      ENDIF
153    ENDIF
154  ELSE ! klon_glo==1, running the 1D model
155    ! just copy over input values
156    latfi(1)=rlatu(1)
157    lonfi(1)=rlonv(1)
158    cufi(1)=cu(1)
159    cvfi(1)=cv(1)
160    airefi(1)=aire(1,1)
161  ENDIF ! of IF (klon_glo>1)
162
163!$OMP PARALLEL
164  ! Now generate local lon/lat/cu/cv/area arrays
165  CALL initcomgeomphy
166     
167  offset = klon_mpi_begin - 1
168  airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
169  cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
170  cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
171  rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
172  rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
173
174! copy some fundamental parameters to physics
175  rradius=prad
176  rg=pg
177  rr=pr
178  rcpp=pcpp
179
180!$OMP END PARALLEL
181
182!      print*,'ATTENTION !!! TRAVAILLER SUR INIPHYSIQ'
183!      print*,'CONTROLE DES LATITUDES, LONGITUDES, PARAMETRES ...'
184
185! Additional initializations for aquaplanets
186!!$OMP PARALLEL
187!      if (iflag_phys>=100) then
188!        call iniaqua(klon_omp,rlatd,rlond,iflag_phys)
189!      endif
190!!$OMP END PARALLEL
191
192END SUBROUTINE iniphysiq
193
194END MODULE iniphysiq_mod
Note: See TracBrowser for help on using the repository browser.