source: LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/iniphysiq.F90 @ 2346

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

Physics/dynamics separation:

  • remove all references to dimensions.h from physics. nbp_lon (==iim) , nbp_lat (==jjm+1) and nbp_lev (==llm) from mod_grid_phy_lmdz should be used instead.
  • added module regular_lonlat_mod in phy_common to store information about the global (lon-lat) grid cell boundaries and centers.

EM

  • 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: 8.3 KB
Line 
1
2! $Id: iniphysiq.F90 2346 2015-08-21 15:13:46Z emillour $
3
4
5SUBROUTINE iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep, &
6                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv,       &
7                     prad,pg,pr,pcpp,iflag_phys)
8  USE dimphy, ONLY: klev ! number of atmospheric levels
9  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
10                                        ! (on full grid)
11  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
12                                klon_omp_begin, & ! start index of local omp subgrid
13                                klon_omp_end, & ! end index of local omp subgrid
14                                klon_mpi_begin ! start indes of columns (on local mpi grid)
15  USE vertical_layers_mod, ONLY : init_vertical_layers
16  USE infotrac, ONLY: nqtot,nqo,nbtr,tname,ttext,type_trac,&
17                      niadv,conv_flg,pbl_flg,solsym,&
18                      nqfils,nqdesc,nqdesc_tot,iqfils,iqpere,&
19                      ok_isotopes,ok_iso_verif,ok_isotrac,&
20                      ok_init_iso,niso_possibles,tnat,&
21                      alpha_ideal,use_iso,iqiso,iso_num,&
22                      iso_indnum,zone_num,phase_num,&
23                      indnum_fn_num,index_trac,&
24                      niso,ntraceurs_zone,ntraciso
25  USE control_mod, ONLY: dayref,anneeref,day_step,nday,offline
26  USE comgeomphy, ONLY: initcomgeomphy, &
27                        airephy, & ! physics grid area (m2)
28                        cuphy, & ! cu coeff. (u_covariant = cu * u)
29                        cvphy, & ! cv coeff. (v_covariant = cv * v)
30                        rlond, & ! longitudes
31                        rlatd ! latitudes
32  USE inifis_mod, ONLY: inifis
33  USE time_phylmdz_mod, ONLY: init_time
34  USE infotrac_phy, ONLY: init_infotrac_phy
35  USE phystokenc_mod, ONLY: init_phystokenc
36  USE phyaqua_mod, ONLY: iniaqua
37  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
38                                 east, west, north, south, &
39                                 north_east, north_west, &
40                                 south_west, south_east
41  IMPLICIT NONE
42
43  ! =======================================================================
44  ! Initialisation of the physical constants and some positional and
45  ! geometrical arrays for the physics
46  ! =======================================================================
47
48  include "dimensions.h"
49  include "comvert.h"
50  include "comconst.h"
51  include "iniprint.h"
52  include "temps.h"
53  include "tracstoke.h"
54
55  REAL, INTENT (IN) :: prad ! radius of the planet (m)
56  REAL, INTENT (IN) :: pg ! gravitational acceleration (m/s2)
57  REAL, INTENT (IN) :: pr ! ! reduced gas constant R/mu
58  REAL, INTENT (IN) :: pcpp ! specific heat Cp
59  REAL, INTENT (IN) :: punjours ! length (in s) of a standard day
60  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
61  INTEGER, INTENT (IN) :: ii ! number of atmospheric columns along longitudes
62  INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes
63  REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid
64  REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid
65  REAL, INTENT (IN) :: rlonv(ii+1) ! longitudes of the physics grid
66  REAL, INTENT (IN) :: rlonu(ii+1) ! longitude boundaries of the physics grid
67  REAL, INTENT (IN) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2)
68  REAL, INTENT (IN) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
69  REAL, INTENT (IN) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
70  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
71  REAL, INTENT (IN) :: ptimestep !physics time step (s)
72  INTEGER, INTENT (IN) :: iflag_phys ! type of physics to be called
73
74  INTEGER :: ibegin, iend, offset
75  INTEGER :: i,j
76  CHARACTER (LEN=20) :: modname = 'iniphysiq'
77  CHARACTER (LEN=80) :: abort_message
78  REAL :: total_area_phy, total_area_dyn
79
80  ! boundaries, on global grid
81  REAL,ALLOCATABLE :: boundslon_reg(:,:)
82  REAL,ALLOCATABLE :: boundslat_reg(:,:)
83
84  ! global array, on full physics grid:
85  REAL,ALLOCATABLE :: latfi(:)
86  REAL,ALLOCATABLE :: lonfi(:)
87  REAL,ALLOCATABLE :: cufi(:)
88  REAL,ALLOCATABLE :: cvfi(:)
89  REAL,ALLOCATABLE :: airefi(:)
90
91  IF (nlayer/=klev) THEN
92    WRITE (lunout, *) 'nlayer     = ', nlayer
93    WRITE (lunout, *) 'klev   = ', klev
94    CALL abort_gcm(modname, 'Problem with dimensions', 1)
95  END IF
96
97  !call init_phys_lmdz(ii,jj+1,llm,1,(/(jj-1)*ii+2/))
98 
99  ! init regular global longitude-latitude grid points and boundaries
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  ! Write values in module regular_lonlat_mod
118  CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, &
119                           boundslon_reg, boundslat_reg)
120
121  ! Generate global arrays on full physics grid
122  ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
123  ALLOCATE(airefi(klon_glo))
124
125  IF (klon_glo>1) THEN ! general case
126    ! North pole
127    latfi(1)=rlatu(1)
128    lonfi(1)=0.
129    cufi(1) = cu(1)
130    cvfi(1) = cv(1)
131    DO j=2,jj
132      DO i=1,ii
133        latfi((j-2)*ii+1+i)= rlatu(j)
134        lonfi((j-2)*ii+1+i)= rlonv(i)
135        cufi((j-2)*ii+1+i) = cu((j-1)*(ii+1)+i)
136        cvfi((j-2)*ii+1+i) = cv((j-1)*(ii+1)+i)
137      ENDDO
138    ENDDO
139    ! South pole
140    latfi(klon_glo)= rlatu(jj+1)
141    lonfi(klon_glo)= 0.
142    cufi(klon_glo) = cu((ii+1)*jj+1)
143    cvfi(klon_glo) = cv((ii+1)*jj-ii)
144
145    ! build airefi(), mesh area on physics grid
146    CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi)
147    ! Poles are single points on physics grid
148    airefi(1)=sum(aire(1:ii,1))
149    airefi(klon_glo)=sum(aire(1:ii,jj+1))
150
151    ! Sanity check: do total planet area match between physics and dynamics?
152    total_area_dyn=sum(aire(1:ii,1:jj+1))
153    total_area_phy=sum(airefi(1:klon_glo))
154    IF (total_area_dyn/=total_area_phy) THEN
155      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
156      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
157      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
158      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
159        ! stop here if the relative difference is more than 0.001%
160        abort_message = 'planet total surface discrepancy'
161        CALL abort_gcm(modname, abort_message, 1)
162      ENDIF
163    ENDIF
164  ELSE ! klon_glo==1, running the 1D model
165    ! just copy over input values
166    latfi(1)=rlatu(1)
167    lonfi(1)=rlonv(1)
168    cufi(1)=cu(1)
169    cvfi(1)=cv(1)
170    airefi(1)=aire(1,1)
171  ENDIF ! of IF (klon_glo>1)
172
173!$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/)
174  ! Initialize physical constants in physics:
175  CALL inifis(punjours,prad,pg,pr,pcpp)
176  CALL init_time(annee_ref,day_ref,day_ini,start_time,nday,ptimestep)
177
178  ! Copy over "offline" settings
179  CALL init_phystokenc(offline,istphy)
180
181  ! copy over preff , ap(), bp(), etc
182  CALL init_vertical_layers(nlayer,preff,scaleheight, &
183                            ap,bp,presnivs,pseudoalt)
184
185  ! Initialize tracer names, numbers, etc. for physics
186  CALL init_infotrac_phy(nqtot,nqo,nbtr,tname,ttext,type_trac,&
187                         niadv,conv_flg,pbl_flg,solsym,&
188                         nqfils,nqdesc,nqdesc_tot,iqfils,iqpere,&
189                         ok_isotopes,ok_iso_verif,ok_isotrac,&
190                         ok_init_iso,niso_possibles,tnat,&
191                         alpha_ideal,use_iso,iqiso,iso_num,&
192                         iso_indnum,zone_num,phase_num,&
193                         indnum_fn_num,index_trac,&
194                         niso,ntraceurs_zone,ntraciso)
195
196  ! Now generate local lon/lat/cu/cv/area arrays
197  CALL initcomgeomphy
198
199  offset = klon_mpi_begin - 1
200  airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
201  cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
202  cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
203  rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
204  rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
205
206  ! Additional initializations for aquaplanets
207  IF (iflag_phys>=100) THEN
208    CALL iniaqua(klon_omp, rlatd, rlond, iflag_phys)
209  END IF
210!$OMP END PARALLEL
211
212END SUBROUTINE iniphysiq
Note: See TracBrowser for help on using the repository browser.