source: trunk/LMDZ.PLUTO/libf/phypluto/phyetat0_mod.F90 @ 3502

Last change on this file since 3502 was 3502, checked in by afalco, 13 days ago

Pluto: better thermal inertia name when no startfi given.
AF

File size: 8.1 KB
Line 
1module phyetat0_mod
2
3implicit none
4
5contains
6
7subroutine phyetat0 (startphy_file, &
8                     ngrid,nlayer,fichnom,tab0,Lmodif,nsoil,nq, &
9                     day_ini,time,tsurf,tsoil, &
10                     emis,q2,qsurf,therm_inertia)
11                    !  ,cloudfrac,totcloudfrac,hice, &
12                    !  rnat,pctsrf_sic)
13                    !  ,tslab,tsea_ice,sea_ice) !AF24
14
15! to use  'getin_p'
16      use ioipsl_getin_p_mod, only: getin_p
17!!
18  use write_field_phy, only: Writefield_phy
19!!
20  use comsoil_h, only: nsoilmx
21  use tabfi_mod, only: tabfi
22  USE tracer_h, ONLY: noms
23  USE surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam, zthe
24  use iostart, only: nid_start, open_startphy, close_startphy, &
25                     get_field, get_var, inquire_field, &
26                     inquire_dimension, inquire_dimension_length
27  use callkeys_mod, only: surfalbedo,surfemis, callsoil
28  implicit none
29
30!======================================================================
31!  Arguments:
32!  ---------
33!  inputs:
34  logical,intent(in) :: startphy_file ! .true. if reading start file
35  integer,intent(in) :: ngrid
36  integer,intent(in) :: nlayer
37  character*(*),intent(in) :: fichnom ! "startfi.nc" file
38  integer,intent(in) :: tab0
39  integer,intent(in) :: Lmodif
40  integer,intent(in) :: nsoil ! # of soil layers
41  integer,intent(in) :: nq
42  integer,intent(out) :: day_ini
43  real,intent(out) :: time
44
45!  outputs:
46  real,intent(out) :: tsurf(ngrid) ! surface temperature
47  real,intent(out) :: tsoil(ngrid,nsoil) ! soil temperature
48  real,intent(out) :: emis(ngrid) ! surface emissivity
49  real,intent(out) :: q2(ngrid,nlayer+1) !
50  real,intent(out) :: qsurf(ngrid,nq) ! tracers on surface
51  real,intent(out) :: therm_inertia(ngrid,nsoilmx) ! thermal inertia
52! real n2ice(ngrid) ! n2 ice cover
53
54!======================================================================
55!  Local variables:
56
57!      INTEGER radpas
58!      REAL n2_ppm
59!      REAL solaire
60
61      real xmin,xmax ! to display min and max of a field
62!
63      INTEGER ig,iq,lmax
64      INTEGER nid, nvarid
65      INTEGER ierr, i, nsrf
66!      integer isoil
67!      INTEGER length
68!      PARAMETER (length=100)
69      CHARACTER*7 str7
70      CHARACTER*2 str2
71      CHARACTER*1 yes
72!
73      REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec,thermal_inertia_if_no_startfi
74      INTEGER nqold
75
76! flag which identifies if 'startfi.nc' file is using old names (qsurf01,...)
77!      logical :: oldtracernames=.false.
78      integer :: count
79      character(len=30) :: txt ! to store some text
80
81      INTEGER :: indextime=1 ! index of selected time, default value=1
82      logical :: found
83
84      character(len=8) :: modname="phyetat0"
85
86!
87! ALLOCATE ARRAYS IN surfdat_h
88!
89IF (.not. ALLOCATED(albedodat)) ALLOCATE(albedodat(ngrid))
90IF (.not. ALLOCATED(phisfi)) ALLOCATE(phisfi(ngrid))
91IF (.not. ALLOCATED(zmea)) ALLOCATE(zmea(ngrid))
92IF (.not. ALLOCATED(zstd)) ALLOCATE(zstd(ngrid))
93IF (.not. ALLOCATED(zsig)) ALLOCATE(zsig(ngrid))
94IF (.not. ALLOCATED(zgam)) ALLOCATE(zgam(ngrid))
95IF (.not. ALLOCATED(zthe)) ALLOCATE(zthe(ngrid))
96
97if (startphy_file) then
98  ! open physics initial state file:
99  call open_startphy(fichnom)
100
101  ! possibility to modify tab_cntrl in tabfi
102  write(*,*)
103  write(*,*) 'TABFI in phyeta0: Lmodif=',Lmodif," tab0=",tab0
104  call tabfi (ngrid,nid_start,Lmodif,tab0,day_ini,lmax,p_rad, &
105                   p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
106
107else ! "academic" initialization of planetary parameters via tabfi
108  call tabfi (ngrid,0,0,0,day_ini,lmax,p_rad, &
109                   p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
110endif ! of if (startphy_file)
111
112if (startphy_file) then
113  ! Load surface geopotential:
114  call get_field("phisfi",phisfi,found)
115  if (.not.found) then
116    call abort_physic(modname,"Failed loading <phisfi>",1)
117  endif
118else
119  phisfi(:)=0.
120endif ! of if (startphy_file)
121write(*,*) "phyetat0: surface geopotential <phisfi> range:", &
122               minval(phisfi), maxval(phisfi)
123
124if (startphy_file) then
125  ! Load bare ground albedo:
126  call get_field("albedodat",albedodat,found)
127  if (.not.found) then
128    call abort_physic(modname,"Failed loading <albedodat>",1)
129  endif
130else
131  ! If no startfi file, use parameter surfalbedo in def file
132  surfalbedo=0.5
133  call getin_p("surfalbedo",surfalbedo)
134  print*,"surfalbedo",surfalbedo
135  albedodat(:)=surfalbedo
136endif ! of if (startphy_file)
137write(*,*) "phyetat0: Bare ground albedo <albedodat> range:", &
138             minval(albedodat), maxval(albedodat)
139
140! ZMEA
141if (startphy_file) then
142  call get_field("ZMEA",zmea,found)
143  if (.not.found) then
144    call abort_physic(modname,"Failed loading <ZMEA>",1)
145  endif
146else
147  zmea(:)=0.
148endif ! of if (startphy_file)
149write(*,*) "phyetat0: <ZMEA> range:", &
150             minval(zmea), maxval(zmea)
151
152! ZSTD
153if (startphy_file) then
154  call get_field("ZSTD",zstd,found)
155  if (.not.found) then
156    call abort_physic(modname,"Failed loading <ZSTD>",1)
157  endif
158else
159  zstd(:)=0.
160endif ! of if (startphy_file)
161write(*,*) "phyetat0: <ZSTD> range:", &
162             minval(zstd), maxval(zstd)
163
164! ZSIG
165if (startphy_file) then
166  call get_field("ZSIG",zsig,found)
167  if (.not.found) then
168    call abort_physic(modname,"Failed loading <ZSIG>",1)
169  endif
170else
171  zsig(:)=0.
172endif ! of if (startphy_file)
173write(*,*) "phyetat0: <ZSIG> range:", &
174             minval(zsig), maxval(zsig)
175
176! ZGAM
177if (startphy_file) then
178  call get_field("ZGAM",zgam,found)
179  if (.not.found) then
180    call abort_physic(modname,"Failed loading <ZGAM>",1)
181  endif
182else
183  zgam(:)=0.
184endif ! of if (startphy_file)
185write(*,*) "phyetat0: <ZGAM> range:", &
186             minval(zgam), maxval(zgam)
187
188! ZTHE
189if (startphy_file) then
190  call get_field("ZTHE",zthe,found)
191  if (.not.found) then
192    call abort_physic(modname,"Failed loading <ZTHE>",1)
193  endif
194else
195  zthe(:)=0.
196endif ! of if (startphy_file)
197write(*,*) "phyetat0: <ZTHE> range:", &
198             minval(zthe), maxval(zthe)
199
200! Surface temperature :
201if (startphy_file) then
202  call get_field("tsurf",tsurf,found,indextime)
203  if (.not.found) then
204    call abort_physic(modname,"Failed loading <tsurf>",1)
205  endif
206else
207  tsurf(:)=0. ! will be updated afterwards in physiq !
208endif ! of if (startphy_file)
209write(*,*) "phyetat0: Surface temperature <tsurf> range:", &
210             minval(tsurf), maxval(tsurf)
211
212! Surface emissivity
213if (startphy_file) then
214  call get_field("emis",emis,found,indextime)
215  if (.not.found) then
216    call abort_physic(modname,"Failed loading <emis>",1)
217  endif
218else
219  ! If no startfi file, use parameter surfemis in def file
220  surfemis=1.0
221  call getin_p("surfemis",surfemis)
222  print*,"surfemis",surfemis
223  emis(:)=surfemis
224endif ! of if (startphy_file)
225write(*,*) "phyetat0: Surface emissivity <emis> range:", &
226             minval(emis), maxval(emis)
227
228! !AF24 removed ocean stuff
229
230! pbl wind variance
231if (startphy_file) then
232  call get_field("q2",q2,found,indextime)
233  if (.not.found) then
234    call abort_physic(modname,"Failed loading <q2>",1)
235  endif
236else
237  q2(:,:)=0.
238endif ! of if (startphy_file)
239write(*,*) "phyetat0: PBL wind variance <q2> range:", &
240             minval(q2), maxval(q2)
241
242! tracer on surface
243if (nq.ge.1) then
244  do iq=1,nq
245    txt=noms(iq)
246    if (startphy_file) then
247      call get_field(txt,qsurf(:,iq),found,indextime)
248      if (.not.found) then
249        write(*,*) "phyetat0: Failed loading <",trim(txt),">"
250        write(*,*) "         ",trim(txt)," is set to zero"
251        qsurf(:,iq) = 0.
252      endif
253    else
254      qsurf(:,iq)=0.
255    endif ! of if (startphy_file)
256    write(*,*) "phyetat0: Surface tracer <",trim(txt),"> range:", &
257                 minval(qsurf(:,iq)), maxval(qsurf(:,iq))
258  enddo! of do iq=1,nq
259endif ! of if (nq.ge.1)
260
261!!call WriteField_phy("in_phyetat0_qsurf",qsurf(1:ngrid,igcm_h2o_gas),1)
262
263if (startphy_file) then
264  ! Call to soil_settings, in order to read soil temperatures,
265  ! as well as thermal inertia and volumetric heat capacity
266  if (callsoil) then
267    call soil_settings(nid_start,ngrid,nsoil,tsurf,tsoil,indextime)
268  endif
269else
270  thermal_inertia_if_no_startfi=400 ! default value
271  call getin_p("thermal_inertia_if_no_startfi",thermal_inertia_if_no_startfi)
272  therm_inertia(:,:) = thermal_inertia_if_no_startfi
273    !AF24
274endif ! of if (startphy_file)
275
276!
277! close file:
278!
279if (startphy_file) call close_startphy
280
281end subroutine phyetat0
282
283end module phyetat0_mod
Note: See TracBrowser for help on using the repository browser.