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