Ignore:
Timestamp:
Aug 5, 2025, 1:29:32 PM (5 weeks ago)
Author:
tbertrand
Message:

LMDZ.PLUTO:
For newstart : reading in the intial albedo map and the new topography map
TB

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.PLUTO/libf/dynphy_lonlat/phypluto/lect_start_archive.F

    r3732 r3878  
    11      SUBROUTINE lect_start_archive(ngrid,nlayer,
    2      &     date,tsurf,tsoil,emis,q2,
     2     &     date,tsurf,tsoil,emis,alb,q2,
    33     &     t,ucov,vcov,ps,h,phisold_newgrid,
    44     &     q,qsurf,surfith,nid)
     
    6868      REAL n2ice(ngrid) ! N2 ice layer
    6969      REAL,INTENT(OUT) :: emis(ngrid)
     70      REAL,INTENT(OUT) :: alb(ngrid)
    7071      REAL,INTENT(OUT) :: q2(ngrid,llm+1),qsurf(ngrid,nqtot)
    7172      ! REAL,INTENT(OUT) :: tslab(ngrid,nslay)
     
    100101      real n2iceS(iip1,jjp1)
    101102      real emisS(iip1,jjp1)
     103      real albS(iip1,jjp1)
    102104      REAL q2S(iip1,jjp1,llm+1),qsurfS(iip1,jjp1,nqtot)
    103105      ! real tslabS(iip1,jjp1,nslay)
     
    135137      real, dimension(:,:), allocatable :: tsurfold
    136138      real, dimension(:,:), allocatable :: emisold
     139      real, dimension(:,:), allocatable :: albold
    137140      real, dimension(:,:,:,:), allocatable :: qold
    138141      ! real, dimension(:,:,:), allocatable :: tslabold
     
    284287      allocate(tsurfold(imold+1,jmold+1))
    285288      allocate(emisold(imold+1,jmold+1))
     289      allocate(albold(imold+1,jmold+1))
    286290      allocate(q2old(imold+1,jmold+1,lmold+1))
    287291!      allocate(tsoilold(imold+1,jmold+1,nsoilmx))
     
    645649      IF (ierr .NE. NF_NOERR) THEN
    646650         PRINT*, "lect_start_archive: Failed loading <emis>"
     651         CALL abort
     652      ENDIF
     653c
     654      ierr = NF_INQ_VARID (nid, "albedodat", nvarid)
     655      IF (ierr .NE. NF_NOERR) THEN
     656         PRINT*, "lect_start_archive: Field <albedodat> not found"
     657         CALL abort
     658      ENDIF
     659#ifdef NC_DOUBLE
     660      ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,albold)
     661#else
     662      ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,albold)
     663#endif
     664      IF (ierr .NE. NF_NOERR) THEN
     665         PRINT*, "lect_start_archive: Failed loading <albedodat>"
    647666         CALL abort
    648667      ENDIF
     
    10861105      call gr_dyn_fi (1,iim+1,jjm+1,ngrid,emiss,emis)
    10871106c     write(46,*) 'emis',emis
     1107c Albedo de la surface
     1108      call interp_horiz (albold,albs,imold,jmold,iim,jjm,1,
     1109     &                   rlonuold,rlatvold,rlonu,rlatv)
     1110      call gr_dyn_fi (1,iim+1,jjm+1,ngrid,albs,alb)
     1111c     write(46,*) 'alb',alb
    10881112
    10891113
     
    15031527      deallocate(tsurfold)
    15041528      deallocate(emisold)
     1529      deallocate(albold)
    15051530      deallocate(q2old)
    15061531      deallocate(tsoilold)
Note: See TracChangeset for help on using the changeset viewer.