source: trunk/LMDZ.GENERIC/libf/dynlonlat_phylonlat/phystd/iniphysiq_mod.F90 @ 1523

Last change on this file since 1523 was 1523, checked in by emillour, 9 years ago

All models: More updates to make planetary codes (+Earth) setups converge.

  • in dyn3d_common:
  • convmas.F => convmas.F90
  • enercin.F => enercin.F90
  • flumass.F => flumass.F90
  • massbar.F => massbar.F90
  • tourpot.F => tourpot.F90
  • vitvert.F => vitvert.F90
  • in misc:
  • move "q_sat" from "dyn3d_common" to "misc" (in Earth model, it is also called by the physics)
  • move "write_field" from "dyn3d_common" to "misc"(may be called from physics or dynamics and depends on neither).
  • in phy_common:
  • move "write_field_phy" here since it may be called from any physics package)
  • add module "regular_lonlat_mod" to store global information on lon-lat grid
  • in dynlonlat_phylonlat/phy*:
  • turn "iniphysiq.F90" into module "iniphysiq_mod.F90" (and of course adapt gcm.F[90] and 1D models accordingly)

EM

File size: 6.5 KB
Line 
1MODULE iniphysiq_mod
2
3CONTAINS
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
9use dimphy, only : klev ! number of atmospheric levels
10use mod_grid_phy_lmdz, only : klon_glo ! number of atmospheric columns
11                                       ! (on full grid)
12use 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
17use comgeomphy, only : initcomgeomphy, &
18                       airephy, & ! physics grid area (m2)
19                       cuphy, & ! cu coeff. (u_covariant = cu * u)
20                       cvphy, & ! cv coeff. (v_covariant = cv * v)
21                       rlond, & ! longitudes
22                       rlatd ! latitudes
23use infotrac, only : nqtot ! number of advected tracers
24use planete_mod, only: ini_planete_mod
25USE comvert_mod, ONLY: ap,bp,preff
26use regular_lonlat_mod, only: init_regular_lonlat, &
27                              east, west, north, south, &
28                              north_east, north_west, &
29                              south_west, south_east
30
31implicit none
32include "dimensions.h"
33include "iniprint.h"
34
35real,intent(in) :: prad ! radius of the planet (m)
36real,intent(in) :: pg ! gravitational acceleration (m/s2)
37real,intent(in) :: pr ! ! reduced gas constant R/mu
38real,intent(in) :: pcpp ! specific heat Cp
39real,intent(in) :: punjours ! length (in s) of a standard day
40!integer,intent(in) :: ngrid ! number of horizontal grid points in the physics (full grid)
41integer,intent(in) :: nlayer ! number of atmospheric layers
42integer,intent(in) :: ii ! number of atmospheric coulumns along longitudes
43integer,intent(in) :: jj  ! number of atompsheric columns along latitudes
44real,intent(in) :: rlatu(jj+1) ! latitudes of the physics grid
45real,intent(in) :: rlatv(jj) ! latitude boundaries of the physics grid
46real,intent(in) :: rlonv(ii+1) ! longitudes of the physics grid
47real,intent(in) :: rlonu(ii+1) ! longitude boundaries of the physics grid
48real,intent(in) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2)
49real,intent(in) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
50real,intent(in) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
51integer,intent(in) :: pdayref ! reference day of for the simulation
52real,intent(in) :: ptimestep !physics time step (s)
53integer,intent(in) :: iflag_phys ! type of physics to be called
54
55integer :: ibegin,iend,offset
56integer :: i,j
57character(len=20) :: modname='iniphysiq'
58character(len=80) :: abort_message
59real :: total_area_phy, total_area_dyn
60real :: pi
61
62! boundaries, on global grid
63real,allocatable :: boundslon_reg(:,:)
64real,allocatable :: boundslat_reg(:,:)
65
66! global array, on full physics grid:
67real,allocatable :: latfi(:)
68real,allocatable :: lonfi(:)
69real,allocatable :: cufi(:)
70real,allocatable :: cvfi(:)
71real,allocatable :: airefi(:)
72
73pi=2.*asin(1.0)
74
75IF (nlayer.NE.klev) THEN
76  write(*,*) 'STOP in ',trim(modname)
77  write(*,*) 'Problem with dimensions :'
78  write(*,*) 'nlayer     = ',nlayer
79  write(*,*) 'klev   = ',klev
80  abort_message = ''
81  CALL abort_gcm (modname,abort_message,1)
82ENDIF
83
84!IF (ngrid.NE.klon_glo) THEN
85!  write(*,*) 'STOP in ',trim(modname)
86!  write(*,*) 'Problem with dimensions :'
87!  write(*,*) 'ngrid     = ',ngrid
88!  write(*,*) 'klon   = ',klon_glo
89!  abort_message = ''
90!  CALL abort_gcm (modname,abort_message,1)
91!ENDIF
92
93!call init_phys_lmdz(iim,jjm+1,llm,1,(/(jjm-1)*iim+2/))
94
95! init regular global longitude-latitude grid points and boundaries
96ALLOCATE(boundslon_reg(ii,2))
97ALLOCATE(boundslat_reg(jj+1,2))
98 
99DO i=1,ii
100   boundslon_reg(i,east)=rlonu(i)
101   boundslon_reg(i,west)=rlonu(i+1)
102ENDDO
103
104boundslat_reg(1,north)= PI/2
105boundslat_reg(1,south)= rlatv(1)
106DO j=2,jj
107   boundslat_reg(j,north)=rlatv(j-1)
108   boundslat_reg(j,south)=rlatv(j)
109ENDDO
110boundslat_reg(jj+1,north)= rlatv(jj)
111boundslat_reg(jj+1,south)= -PI/2
112
113! Write values in module regular_lonlat_mod
114CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, &
115                         boundslon_reg, boundslat_reg)
116
117! Generate global arrays on full physics grid
118allocate(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
119latfi(1)=rlatu(1)
120lonfi(1)=0.
121cufi(1) = cu(1)
122cvfi(1) = cv(1)
123DO j=2,jj
124  DO i=1,ii
125    latfi((j-2)*ii+1+i)= rlatu(j)
126    lonfi((j-2)*ii+1+i)= rlonv(i)
127    cufi((j-2)*ii+1+i) = cu((j-1)*(ii+1)+i)
128    cvfi((j-2)*ii+1+i) = cv((j-1)*(ii+1)+i)
129  ENDDO
130ENDDO
131latfi(klon_glo)= rlatu(jj+1)
132lonfi(klon_glo)= 0.
133cufi(klon_glo) = cu((ii+1)*jj+1)
134cvfi(klon_glo) = cv((ii+1)*jj-ii)
135
136! build airefi(), mesh area on physics grid
137allocate(airefi(klon_glo))
138CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi)
139! Poles are single points on physics grid
140airefi(1)=sum(aire(1:ii,1))
141airefi(klon_glo)=sum(aire(1:ii,jj+1))
142
143! Sanity check: do total planet area match between physics and dynamics?
144total_area_dyn=sum(aire(1:ii,1:jj+1))
145total_area_phy=sum(airefi(1:klon_glo))
146IF (total_area_dyn/=total_area_phy) THEN
147  WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
148  WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
149  WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
150  IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
151    ! stop here if the relative difference is more than 0.001%
152    abort_message = 'planet total surface discrepancy'
153    CALL abort_gcm(modname, abort_message, 1)
154  ENDIF
155ENDIF
156
157
158!$OMP PARALLEL
159! Now generate local lon/lat/cu/cv/area arrays
160call initcomgeomphy
161
162!!!!$OMP PARALLEL PRIVATE(ibegin,iend) &
163!!!     !$OMP SHARED(airefi,cufi,cvfi,lonfi,latfi)
164
165offset=klon_mpi_begin-1
166airephy(1:klon_omp)=airefi(offset+klon_omp_begin:offset+klon_omp_end)
167cuphy(1:klon_omp)=cufi(offset+klon_omp_begin:offset+klon_omp_end)
168cvphy(1:klon_omp)=cvfi(offset+klon_omp_begin:offset+klon_omp_end)
169rlond(1:klon_omp)=lonfi(offset+klon_omp_begin:offset+klon_omp_end)
170rlatd(1:klon_omp)=latfi(offset+klon_omp_begin:offset+klon_omp_end)
171
172! copy over preff , ap() and bp()
173call ini_planete_mod(nlayer,preff,ap,bp)
174
175! copy some fundamental parameters to physics
176! and do some initializations
177call inifis(klon_omp,nlayer,nqtot,pdayref,punjours,ptimestep, &
178            rlatd,rlond,airephy,prad,pg,pr,pcpp)
179
180!$OMP END PARALLEL
181
182
183end subroutine iniphysiq
184
185
186END MODULE iniphysiq_mod
Note: See TracBrowser for help on using the repository browser.