source: trunk/LMDZ.GENERIC/libf/phystd/phyetat0_mod.F90 @ 3809

Last change on this file since 3809 was 3773, checked in by afalco, 7 weeks ago

Generic: allow to write controle in XIOS output.
AF

File size: 15.1 KB
RevLine 
[1669]1module phyetat0_mod
2
3implicit none
4
[3522]5real, save :: tab_cntrl_mod(100)
6
[3552]7integer,save :: nid_start ! NetCDF file identifier for startfi.nc file
8
[3522]9!$OMP THREADPRIVATE(tab_cntrl_mod)
10
[1669]11contains
12
13subroutine phyetat0 (startphy_file, &
14                     ngrid,nlayer,fichnom,tab0,Lmodif,nsoil,nq, &
[3335]15                     day_ini,time,tsurf,tsoil,emis,albedo, &
16                     q2,qsurf,cloudfrac,totcloudfrac,hice, &
[3397]17                     rnat,pctsrf_sic,tslab,tsea_ice,tice,sea_ice)
[787]18
[1709]19! to use  'getin_p'
20      use ioipsl_getin_p_mod, only: getin_p
[3100]21!!
22  use write_field_phy, only: Writefield_phy
23!!
[3773]24  use tabfi_mod, only: tabfi,ini_tab_controle_dyn_xios
[3100]25  USE tracer_h, ONLY: noms, igcm_h2o_vap
[3335]26  USE radinc_h, ONLY: L_NSPECTV
[1216]27  USE surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam, zthe
28  use iostart, only: nid_start, open_startphy, close_startphy, &
29                     get_field, get_var, inquire_field, &
30                     inquire_dimension, inquire_dimension_length
[3100]31!  use slab_ice_h, only: noceanmx
32  USE ocean_slab_mod, ONLY: nslay
[2557]33  use callkeys_mod, only: CLFvarying,surfalbedo,surfemis, callsoil
[2299]34  use nonoro_gwd_ran_mod, only: du_nonoro_gwd, dv_nonoro_gwd, &
35                                east_gwstress, west_gwstress
[1216]36  implicit none
[787]37
[1216]38!======================================================================
[135]39!  Arguments:
40!  ---------
41!  inputs:
[1669]42  logical,intent(in) :: startphy_file ! .true. if reading start file
[1216]43  integer,intent(in) :: ngrid
[1308]44  integer,intent(in) :: nlayer
[1216]45  character*(*),intent(in) :: fichnom ! "startfi.nc" file
46  integer,intent(in) :: tab0
47  integer,intent(in) :: Lmodif
48  integer,intent(in) :: nsoil ! # of soil layers
49  integer,intent(in) :: nq
[1669]50  integer,intent(out) :: day_ini
51  real,intent(out) :: time
[135]52
53!  outputs:
[1216]54  real,intent(out) :: tsurf(ngrid) ! surface temperature
55  real,intent(out) :: tsoil(ngrid,nsoil) ! soil temperature
56  real,intent(out) :: emis(ngrid) ! surface emissivity
[3335]57  real,intent(out) :: albedo(ngrid,L_NSPECTV) ! albedo of the surface
[3522]58  real,intent(out) :: q2(ngrid,nlayer+1) !
[1216]59  real,intent(out) :: qsurf(ngrid,nq) ! tracers on surface
60! real co2ice(ngrid) ! co2 ice cover
[1308]61  real,intent(out) :: cloudfrac(ngrid,nlayer)
[1216]62  real,intent(out) :: hice(ngrid), totcloudfrac(ngrid)
[3522]63  real,intent(out) :: pctsrf_sic(ngrid),tslab(ngrid,nslay)
[1297]64  real,intent(out) :: tsea_ice(ngrid),sea_ice(ngrid)
65  real,intent(out) :: rnat(ngrid)
[3397]66  real,intent(out) :: tice(ngrid)
[135]67
68!======================================================================
69!  Local variables:
70
71!      INTEGER radpas
72!      REAL co2_ppm
73!      REAL solaire
74
75      real xmin,xmax ! to display min and max of a field
[1216]76!
[135]77      INTEGER ig,iq,lmax
78      INTEGER nid, nvarid
79      INTEGER ierr, i, nsrf
[3522]80!      integer isoil
[135]81!      INTEGER length
82!      PARAMETER (length=100)
83      CHARACTER*7 str7
84      CHARACTER*2 str2
85      CHARACTER*1 yes
[1216]86!
[135]87      REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec
88      INTEGER nqold
89
90! flag which identifies if 'startfi.nc' file is using old names (qsurf01,...)
[1216]91!      logical :: oldtracernames=.false.
[135]92      integer :: count
93      character(len=30) :: txt ! to store some text
[3522]94
[1216]95      INTEGER :: indextime=1 ! index of selected time, default value=1
96      logical :: found
[3522]97
[1669]98      character(len=8) :: modname="phyetat0"
[135]99
[1216]100!
101! ALLOCATE ARRAYS IN surfdat_h
102!
103IF (.not. ALLOCATED(albedodat)) ALLOCATE(albedodat(ngrid))
104IF (.not. ALLOCATED(phisfi)) ALLOCATE(phisfi(ngrid))
105IF (.not. ALLOCATED(zmea)) ALLOCATE(zmea(ngrid))
106IF (.not. ALLOCATED(zstd)) ALLOCATE(zstd(ngrid))
107IF (.not. ALLOCATED(zsig)) ALLOCATE(zsig(ngrid))
108IF (.not. ALLOCATED(zgam)) ALLOCATE(zgam(ngrid))
109IF (.not. ALLOCATED(zthe)) ALLOCATE(zthe(ngrid))
[787]110
[1669]111if (startphy_file) then
112  ! open physics initial state file:
[3552]113  call open_startphy(fichnom,nid_start)
[1297]114
[1669]115  ! possibility to modify tab_cntrl in tabfi
116  write(*,*)
117  write(*,*) 'TABFI in phyeta0: Lmodif=',Lmodif," tab0=",tab0
118  call tabfi (ngrid,nid_start,Lmodif,tab0,day_ini,lmax,p_rad, &
119                   p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
[135]120
[1669]121else ! "academic" initialization of planetary parameters via tabfi
122  call tabfi (ngrid,0,0,0,day_ini,lmax,p_rad, &
[1216]123                   p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
[1669]124endif ! of if (startphy_file)
[135]125
[1669]126if (startphy_file) then
127  ! Load surface geopotential:
[3552]128  call get_field(nid_start,"phisfi",phisfi,found)
[1669]129  if (.not.found) then
130    call abort_physic(modname,"Failed loading <phisfi>",1)
131  endif
[1216]132else
[1709]133  phisfi(:)=0.
[1669]134endif ! of if (startphy_file)
135write(*,*) "phyetat0: surface geopotential <phisfi> range:", &
136               minval(phisfi), maxval(phisfi)
[253]137
[1669]138if (startphy_file) then
[3335]139  ! Load bare ground albedo: (will be stored in surfdat_h)
[3552]140  call get_field(nid_start,"albedodat",albedodat,found)
[1669]141  if (.not.found) then
142    call abort_physic(modname,"Failed loading <albedodat>",1)
143  endif
[1216]144else
[1709]145  ! If no startfi file, use parameter surfalbedo in def file
146  surfalbedo=0.5
147  call getin_p("surfalbedo",surfalbedo)
148  print*,"surfalbedo",surfalbedo
149  albedodat(:)=surfalbedo
[1669]150endif ! of if (startphy_file)
151write(*,*) "phyetat0: Bare ground albedo <albedodat> range:", &
[1216]152             minval(albedodat), maxval(albedodat)
[253]153
[3335]154if (startphy_file) then
155  ! Load surface albedo (for now assume it is spectrally homogeneous)
[3552]156  call get_field(nid_start,"albedo",albedo(:,1),found)
[3335]157  if (.not.found) then
158    write(*,*) modname//": Failed loading <albedo>"
159    write(*,*) " setting it to bare ground albedo"
160    albedo(1:ngrid,1)=albedodat(1:ngrid)
161  endif
162  ! copy value to all spectral bands
163  do i=2,L_NSPECTV
164    albedo(1:ngrid,i)=albedo(1:ngrid,1)
165  enddo
166else
167  ! If no startfi file, use bare ground value
168  do i=1,L_NSPECTV
169    albedo(1:ngrid,i)=albedodat(1:ngrid)
170  enddo
171endif ! of if (startphy_file)
172write(*,*) "phyetat0: Surface albedo <albedo> range:", &
173             minval(albedo), maxval(albedo)
174
[1216]175! ZMEA
[1669]176if (startphy_file) then
[3552]177  call get_field(nid_start,"ZMEA",zmea,found)
[1669]178  if (.not.found) then
179    call abort_physic(modname,"Failed loading <ZMEA>",1)
180  endif
[1216]181else
[1709]182  zmea(:)=0.
[1669]183endif ! of if (startphy_file)
184write(*,*) "phyetat0: <ZMEA> range:", &
[1216]185             minval(zmea), maxval(zmea)
[253]186
[1216]187! ZSTD
[1669]188if (startphy_file) then
[3552]189  call get_field(nid_start,"ZSTD",zstd,found)
[1669]190  if (.not.found) then
191    call abort_physic(modname,"Failed loading <ZSTD>",1)
192  endif
[1216]193else
[1709]194  zstd(:)=0.
[1669]195endif ! of if (startphy_file)
196write(*,*) "phyetat0: <ZSTD> range:", &
[1216]197             minval(zstd), maxval(zstd)
[253]198
[1216]199! ZSIG
[1669]200if (startphy_file) then
[3552]201  call get_field(nid_start,"ZSIG",zsig,found)
[1669]202  if (.not.found) then
203    call abort_physic(modname,"Failed loading <ZSIG>",1)
204  endif
[1216]205else
[1709]206  zsig(:)=0.
[1669]207endif ! of if (startphy_file)
208write(*,*) "phyetat0: <ZSIG> range:", &
[1216]209             minval(zsig), maxval(zsig)
[253]210
[1216]211! ZGAM
[1669]212if (startphy_file) then
[3552]213  call get_field(nid_start,"ZGAM",zgam,found)
[1669]214  if (.not.found) then
215    call abort_physic(modname,"Failed loading <ZGAM>",1)
216  endif
[1216]217else
[1709]218  zgam(:)=0.
[1669]219endif ! of if (startphy_file)
220write(*,*) "phyetat0: <ZGAM> range:", &
[1216]221             minval(zgam), maxval(zgam)
[253]222
[1216]223! ZTHE
[1669]224if (startphy_file) then
[3552]225  call get_field(nid_start,"ZTHE",zthe,found)
[1669]226  if (.not.found) then
227    call abort_physic(modname,"Failed loading <ZTHE>",1)
228  endif
[1216]229else
[1709]230  zthe(:)=0.
[1669]231endif ! of if (startphy_file)
232write(*,*) "phyetat0: <ZTHE> range:", &
[1216]233             minval(zthe), maxval(zthe)
[253]234
[1216]235! Surface temperature :
[1669]236if (startphy_file) then
[3552]237  call get_field(nid_start,"tsurf",tsurf,found,indextime)
[1669]238  if (.not.found) then
239    call abort_physic(modname,"Failed loading <tsurf>",1)
240  endif
[1216]241else
[1709]242  tsurf(:)=0. ! will be updated afterwards in physiq !
[1669]243endif ! of if (startphy_file)
244write(*,*) "phyetat0: Surface temperature <tsurf> range:", &
[1216]245             minval(tsurf), maxval(tsurf)
[135]246
[1216]247! Surface emissivity
[1669]248if (startphy_file) then
[3552]249  call get_field(nid_start,"emis",emis,found,indextime)
[1669]250  if (.not.found) then
251    call abort_physic(modname,"Failed loading <emis>",1)
252  endif
[1216]253else
[1709]254  ! If no startfi file, use parameter surfemis in def file
255  surfemis=1.0
256  call getin_p("surfemis",surfemis)
257  print*,"surfemis",surfemis
258  emis(:)=surfemis
[1669]259endif ! of if (startphy_file)
260write(*,*) "phyetat0: Surface emissivity <emis> range:", &
[1216]261             minval(emis), maxval(emis)
[135]262
[1216]263! Cloud fraction (added by BC 2010)
[1421]264if (CLFvarying) then
[1669]265  if (startphy_file) then
[3552]266    call get_field(nid_start,"cloudfrac",cloudfrac,found,indextime)
[1669]267    if (.not.found) then
268      call abort_physic(modname,"Failed loading <cloudfrac>",1)
269    endif
270  else
[1709]271    cloudfrac(:,:)=0.0
[1669]272  endif ! of if (startphy_file)
[1216]273  write(*,*) "phyetat0: Cloud fraction <cloudfrac> range:", &
[1669]274               minval(cloudfrac), maxval(cloudfrac)
[1421]275else
[1669]276  cloudfrac(:,:)=0.0
277endif ! of if (CLFvarying)
[1216]278
279! Total cloud fraction (added by BC 2010)
[1421]280if (CLFvarying) then
[1669]281  if (startphy_file) then
[3552]282    call get_field(nid_start,"totcloudfrac",totcloudfrac,found,indextime)
[1669]283    if (.not.found) then
284      call abort_physic(modname,"Failed loading <totcloudfrac>",1)
285    endif
286  else
[1709]287    totcloudfrac(:)=0.0
[1669]288  endif ! of if (startphy_file)
[1216]289  write(*,*) "phyetat0: Total cloud fraction <totcloudfrac> range:", &
290             minval(totcloudfrac), maxval(totcloudfrac)
[1421]291else
[1669]292  totcloudfrac(:)=0.0
293endif ! of if (CLFvarying)
[1216]294
295! Height of oceanic ice (added by BC 2010)
[1669]296if (startphy_file) then
[3552]297  call get_field(nid_start,"hice",hice,found,indextime)
[1669]298  if (.not.found) then
299    write(*,*) "phyetat0: Failed loading <hice>"
300    !  call abort
[1709]301    hice(:)=0.
[1669]302  endif
[1216]303else
[1709]304  hice(:)=0.
[1669]305endif ! of if (startphy_file)
306write(*,*) "phyetat0: Height of oceanic ice <hice> range:", &
[1216]307             minval(hice), maxval(hice)
308
[1297]309! SLAB OCEAN (added by BC 2014)
[1669]310if (startphy_file) then
311  ! nature of the surface
[3552]312  call get_field(nid_start,"rnat",rnat,found,indextime)
[1669]313  if (.not.found) then
314    write(*,*) "phyetat0: Failed loading <rnat>"
315        rnat(1:ngrid)=1.
316  else
[1297]317      do ig=1,ngrid
318        if((nint(rnat(ig)).eq.2).or.(nint(rnat(ig)).eq.0))then
319          rnat(ig)=0.
320        else
321          rnat(ig)=1.
[3522]322        endif
[1297]323      enddo
[1669]324  endif ! of if (.not.found)
325else
[1709]326  rnat(:)=1.
[1669]327endif ! of if (startphy_file)
328write(*,*) "phyetat0: Nature of surface <rnat> range:", &
329             minval(rnat), maxval(rnat)
[1297]330
[1669]331if (startphy_file) then
332  ! Pourcentage of sea ice cover
[3552]333  call get_field(nid_start,"pctsrf_sic",pctsrf_sic,found,indextime)
[1669]334  if (.not.found) then
335    write(*,*) "phyetat0: Failed loading <pctsrf_sic>"
336    pctsrf_sic(1:ngrid)=0.
337  endif
[1297]338else
[1709]339  pctsrf_sic(:)=0.
[1669]340endif ! of if (startphy_file)
341write(*,*) "phyetat0: Pourcentage of sea ice cover <pctsrf_sic> range:", &
[1297]342             minval(pctsrf_sic), maxval(pctsrf_sic)
[1669]343
344if (startphy_file) then
345  ! Slab ocean temperature (2 layers)
[3552]346  call get_field(nid_start,"tslab",tslab,found,indextime)
[1669]347  if (.not.found) then
348    write(*,*) "phyetat0: Failed loading <tslab>"
[3100]349    do iq=1,nslay
[1669]350      tslab(1:ngrid,iq)=tsurf(1:ngrid)
351    enddo
352  endif
[1297]353else
[3100]354  do iq=1,nslay
[1669]355    tslab(1:ngrid,iq)=tsurf(1:ngrid)
356  enddo
357endif ! of if (startphy_file)
358write(*,*) "phyetat0: Slab ocean temperature <tslab> range:", &
[1297]359             minval(tslab), maxval(tslab)
[1669]360
361if (startphy_file) then
362  ! Oceanic ice temperature
[3552]363  call get_field(nid_start,"tsea_ice",tsea_ice,found,indextime)
[1669]364  if (.not.found) then
365    write(*,*) "phyetat0: Failed loading <tsea_ice>"
366    tsea_ice(1:ngrid)=273.15-1.8
367  endif
[1297]368else
[1669]369  tsea_ice(1:ngrid)=273.15-1.8
370endif ! of if (startphy_file)
[3397]371write(*,*) "phyetat0: Oceanic ice/snow temperature <tsea_ice> range:", &
[1297]372             minval(tsea_ice), maxval(tsea_ice)
[1669]373
[3397]374if (startphy_file) then
[3100]375  ! Oceanic ice temperature
[3552]376  call get_field(nid_start,"tice",tice,found,indextime)
[3397]377  if (.not.found) then
378    write(*,*) "phyetat0: Failed loading <tice>"
379    tice(1:ngrid)=273.15-1.8
380  endif
381else
382  tice(1:ngrid)=273.15-1.8
383endif ! of if (startphy_file)
384write(*,*) "phyetat0: Oceanic ice surface temperature <tice> range:", &
385             minval(tice), maxval(tice)
[3100]386
[1669]387if (startphy_file) then
388  !  Oceanic ice quantity (kg/m^2)
[3552]389  call get_field(nid_start,"sea_ice",sea_ice,found,indextime)
[1669]390  if (.not.found) then
391    write(*,*) "phyetat0: Failed loading <sea_ice>"
[1763]392    sea_ice(1:ngrid)=0.
[1669]393  endif
[1297]394else
[1763]395  sea_ice(1:ngrid)=0.
[1669]396endif ! of if (startphy_file)
397write(*,*) "phyetat0: Oceanic ice quantity <sea_ice> range:", &
[1297]398             minval(sea_ice), maxval(sea_ice)
399
400
[1216]401! pbl wind variance
[1669]402if (startphy_file) then
[3552]403  call get_field(nid_start,"q2",q2,found,indextime)
[1669]404  if (.not.found) then
405    call abort_physic(modname,"Failed loading <q2>",1)
406  endif
[1216]407else
[1709]408  q2(:,:)=0.
[1669]409endif ! of if (startphy_file)
410write(*,*) "phyetat0: PBL wind variance <q2> range:", &
[1216]411             minval(q2), maxval(q2)
412
413! tracer on surface
414if (nq.ge.1) then
415  do iq=1,nq
[1621]416    txt=noms(iq)
[1669]417    if (startphy_file) then
[3552]418      call get_field(nid_start,txt,qsurf(:,iq),found,indextime)
[1669]419      if (.not.found) then
420        write(*,*) "phyetat0: Failed loading <",trim(txt),">"
421        write(*,*) "         ",trim(txt)," is set to zero"
422        qsurf(:,iq) = 0.
423      endif
[1216]424    else
[1709]425      qsurf(:,iq)=0.
[1669]426    endif ! of if (startphy_file)
427    write(*,*) "phyetat0: Surface tracer <",trim(txt),"> range:", &
[1216]428                 minval(qsurf(:,iq)), maxval(qsurf(:,iq))
[1669]429  enddo! of do iq=1,nq
[1216]430endif ! of if (nq.ge.1)
431
[3100]432!!call WriteField_phy("in_phyetat0_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1)
[1297]433
[1669]434if (startphy_file) then
435  ! Call to soil_settings, in order to read soil temperatures,
436  ! as well as thermal inertia and volumetric heat capacity
[2557]437  if (callsoil) then
438    call soil_settings(nid_start,ngrid,nsoil,tsurf,tsoil,indextime)
439  endif
[1669]440endif ! of if (startphy_file)
[2299]441
[3522]442! Non-orographic gravity waves
[2299]443if (startphy_file) then
[3552]444   call get_field(nid_start,"du_nonoro_gwd",du_nonoro_gwd,found,indextime)
[2299]445   if (.not.found) then
446      write(*,*) "phyetat0: <du_nonoro_gwd> not in file"
447      du_nonoro_gwd(:,:)=0.
448   else
449      write(*,*) "phyetat0: Memory of zonal wind tendency due to non-orographic GW"
450      write(*,*) " <du_nonoro_gwd> range:", &
451             minval(du_nonoro_gwd), maxval(du_nonoro_gwd)
452   endif
453endif ! of if (startphy_file)
454if (startphy_file) then
[3552]455   call get_field(nid_start,"dv_nonoro_gwd",dv_nonoro_gwd,found,indextime)
[2299]456   if (.not.found) then
457      write(*,*) "phyetat0: <dv_nonoro_gwd> not in file"
458      dv_nonoro_gwd(:,:)=0.
459   else
460      write(*,*) "phyetat0: Memory of meridional wind tendency due to non-orographic GW"
461      write(*,*) " <dv_nonoro_gwd> range:", &
462             minval(dv_nonoro_gwd), maxval(dv_nonoro_gwd)
463   endif
464endif ! of if (startphy_file)
465if (startphy_file) then
[3552]466   call get_field(nid_start,"east_gwstress",east_gwstress,found,indextime)
[2299]467   if (.not.found) then
468      write(*,*) "phyetat0: <east_gwstress> not in file"
469      east_gwstress(:,:)=0.
470   else
471      write(*,*) "phyetat0: Memory of eastward stress profile due to non-orographic GW"
472      write(*,*) " <east_gwstress> range:", &
473             minval(east_gwstress), maxval(east_gwstress)
474   endif
475endif ! of if (startphy_file)
476if (startphy_file) then
[3552]477   call get_field(nid_start,"west_gwstress",west_gwstress,found,indextime)
[2299]478   if (.not.found) then
479      write(*,*) "phyetat0: <west_gwstress> not in file"
480      west_gwstress(:,:)=0.
481   else
482      write(*,*) "phyetat0: Memory of westward stress profile due to non-orographic GW"
483      write(*,*) " <west_gwstress> range:", &
484             minval(west_gwstress), maxval(west_gwstress)
485   endif
486endif ! of if (startphy_file)
[3773]487
488call ini_tab_controle_dyn_xios(day_ini)
489
490
[1216]491!
492! close file:
493!
[3552]494if (startphy_file) call close_startphy(nid_start)
[1216]495
[1669]496end subroutine phyetat0
497
498end module phyetat0_mod
Note: See TracBrowser for help on using the repository browser.