source: trunk/LMDZ.MARS/libf/dynlonlat_phylonlat/phymars/iniphysiq.F90 @ 1433

Last change on this file since 1433 was 1433, checked in by aslmd, 9 years ago

commit 1431 removed important initializations. this is corrected. only modifications in moldiff_red were kept from 1431.

File size: 5.2 KB
Line 
1subroutine iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep,           &
2                     rlatu,rlonv,aire,cu,cv,                             &
3                     prad,pg,pr,pcpp,iflag_phys)
4
5use dimphy, only : klev ! number of atmospheric levels
6use mod_grid_phy_lmdz, only : klon_glo ! number of atmospheric columns
7                                       ! (on full grid)
8use mod_phys_lmdz_para, only : klon_omp, & ! number of columns (on local omp grid)
9                               klon_omp_begin, & ! start index of local omp subgrid
10                               klon_omp_end, & ! end index of local omp subgrid
11                               klon_mpi_begin ! start indes of columns (on local mpi grid)
12
13use comgeomphy, only : initcomgeomphy, &
14                       airephy, & ! physics grid area (m2)
15                       cuphy, & ! cu coeff. (u_covariant = cu * u)
16                       cvphy, & ! cv coeff. (v_covariant = cv * v)
17                       rlond, & ! longitudes
18                       rlatd ! latitudes
19use infotrac, only : nqtot ! number of advected tracers
20use comgeomfi_h, only: ini_fillgeom
21
22implicit none
23
24include "iniprint.h"
25
26real,intent(in) :: prad ! radius of the planet (m)
27real,intent(in) :: pg ! gravitational acceleration (m/s2)
28real,intent(in) :: pr ! ! reduced gas constant R/mu
29real,intent(in) :: pcpp ! specific heat Cp
30real,intent(in) :: punjours ! length (in s) of a standard day
31!integer,intent(in) :: ngrid ! number of horizontal grid points in the physics (full grid)
32integer,intent(in) :: nlayer ! number of atmospheric layers
33integer,intent(in) :: ii ! number of atmospheric coulumns along longitudes
34integer,intent(in) :: jj  ! number of atompsheric columns along latitudes
35real,intent(in) :: rlatu(jj+1) ! latitudes of the dynamics U grid
36real,intent(in) :: rlonv(ii+1) ! longitudes of the dynamics V grid
37real,intent(in) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2)
38real,intent(in) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u)
39real,intent(in) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v)
40integer,intent(in) :: pdayref ! reference day of for the simulation
41real,intent(in) :: ptimestep !physics time step (s)
42integer,intent(in) :: iflag_phys ! type of physics to be called
43
44integer :: ibegin,iend,offset
45integer :: i,j
46character(len=20) :: modname='iniphysiq'
47character(len=80) :: abort_message
48real :: total_area_phy, total_area_dyn
49
50
51! global array, on full physics grid:
52real,allocatable :: latfi(:)
53real,allocatable :: lonfi(:)
54real,allocatable :: cufi(:)
55real,allocatable :: cvfi(:)
56real,allocatable :: airefi(:)
57
58IF (nlayer.NE.klev) THEN
59  write(*,*) 'STOP in ',trim(modname)
60  write(*,*) 'Problem with dimensions :'
61  write(*,*) 'nlayer     = ',nlayer
62  write(*,*) 'klev   = ',klev
63  abort_message = ''
64  CALL abort_gcm (modname,abort_message,1)
65ENDIF
66
67!IF (ngrid.NE.klon_glo) THEN
68!  write(*,*) 'STOP in ',trim(modname)
69!  write(*,*) 'Problem with dimensions :'
70!  write(*,*) 'ngrid     = ',ngrid
71!  write(*,*) 'klon   = ',klon_glo
72!  abort_message = ''
73!  CALL abort_gcm (modname,abort_message,1)
74!ENDIF
75
76! Generate global arrays on full physics grid
77allocate(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
78latfi(1)=rlatu(1)
79lonfi(1)=0.
80cufi(1) = cu(1)
81cvfi(1) = cv(1)
82DO j=2,jj
83  DO i=1,ii
84    latfi((j-2)*ii+1+i)= rlatu(j)
85    lonfi((j-2)*ii+1+i)= rlonv(i)
86    cufi((j-2)*ii+1+i) = cu((j-1)*(ii+1)+i)
87    cvfi((j-2)*ii+1+i) = cv((j-1)*(ii+1)+i)
88  ENDDO
89ENDDO
90latfi(klon_glo)= rlatu(jj+1)
91lonfi(klon_glo)= 0.
92cufi(klon_glo) = cu((ii+1)*jj+1)
93cvfi(klon_glo) = cv((ii+1)*jj-ii)
94
95! build airefi(), mesh area on physics grid
96allocate(airefi(klon_glo))
97CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi)
98! Poles are single points on physics grid
99airefi(1)=sum(aire(1:ii,1))
100airefi(klon_glo)=sum(aire(1:ii,jj+1))
101
102! Sanity check: do total planet area match between physics and dynamics?
103total_area_dyn=sum(aire(1:ii,1:jj+1))
104total_area_phy=sum(airefi(1:klon_glo))
105IF (total_area_dyn/=total_area_phy) THEN
106  WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
107  WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
108  WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
109  IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
110    ! stop here if the relative difference is more than 0.001%
111    abort_message = 'planet total surface discrepancy'
112    CALL abort_gcm(modname, abort_message, 1)
113  ENDIF
114ENDIF
115
116
117
118!$OMP PARALLEL
119! Now generate local lon/lat/cu/cv/area arrays
120call initcomgeomphy
121
122!!!!$OMP PARALLEL PRIVATE(ibegin,iend)
123!!!$OMP+         SHARED(parea,pcu,pcv,plon,plat)
124     
125offset=klon_mpi_begin-1
126airephy(1:klon_omp)=airefi(offset+klon_omp_begin:offset+klon_omp_end)
127cuphy(1:klon_omp)=cufi(offset+klon_omp_begin:offset+klon_omp_end)
128cvphy(1:klon_omp)=cvfi(offset+klon_omp_begin:offset+klon_omp_end)
129rlond(1:klon_omp)=lonfi(offset+klon_omp_begin:offset+klon_omp_end)
130rlatd(1:klon_omp)=latfi(offset+klon_omp_begin:offset+klon_omp_end)
131
132! copy some fundamental parameters to physics
133! and do some initializations
134call phys_state_var_init(klon_omp,nlayer,nqtot, &
135                         punjours,ptimestep,prad,pg,pr,pcpp)
136call ini_fillgeom(klon_omp,rlatd,rlond,airephy)
137call conf_phys(klon_omp,nlayer,nqtot)
138
139!$OMP END PARALLEL
140
141
142end subroutine iniphysiq
Note: See TracBrowser for help on using the repository browser.