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

Last change on this file since 3184 was 3184, checked in by afalco, 10 months ago

Pluto PCM:
Added LMDZ.PLUTO, a copy of the generic model,
cleaned from some unnecessary modules (water, ...)
AF

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