Changeset 2020 for trunk


Ignore:
Timestamp:
Oct 19, 2018, 11:44:55 AM (6 years ago)
Author:
mlefevre
Message:

Inclusion of large scale forcing (temperature and water) initialization for LES generic model

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/MESOSCALE/LMD_MM_MARS/SRC/LES/modif_mars/module_initialize_les.F

    r1768 r2020  
    118118 REAL, DIMENSION(nl_max) :: prescribed_sw,prescribed_lw,prescribed_dyn
    119119 REAL, DIMENSION(nl_max) :: hrsw,hrlw,hrdyn
     120 REAL, DIMENSION(nl_max) :: lsf_dt,lsf_dq,lsfdt,lsfdq
    120121 REAL, DIMENSION(nl_max) :: venus_hrdyn
    121122 REAL, DIMENSION(nl_max) :: altitude
     
    955956      ENDDO
    956957      close(unit=20)
     958    ENDIF
     959
     960    IF (planet.eq."generic") THEN
     961      call read_lsf(lsfdt,lsfdq,nl_in)
     962      open(unit=17,file="lsf.txt",action="write")
     963      DO k=1,kte!-1
     964        p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
     965        lsf_dt = interp_0_log( lsfdt, pd_in, p_level, nl_in )
     966        lsf_dq = interp_0_log( lsfdq, pd_in, p_level, nl_in )
     967        write (17,*) lsf_dt(k),lsf_dq(k)
     968      ENDDO
     969    ENDIF
     970
    957971
    958972    open(unit=21,file="altitude.txt",action="write")
     
    963977    ENDDO
    964978    close(unit=21)
    965     ENDIF
    966979
    967980!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     
    13171330      end subroutine read_hr
    13181331
     1332      subroutine read_lsf(dt,dq,n)
     1333      implicit none
     1334      integer n
     1335      real dt(n+1),dq(n+1)
     1336      logical end_of_file
     1337
     1338      integer k
     1339
     1340! first element is the surface
     1341
     1342      open(unit=12,file='input_lsf',form='formatted',status='old')
     1343      rewind(12)
     1344      end_of_file = .false.
     1345      k = 0
     1346      do while (.not. end_of_file)
     1347
     1348        read(12,*,end=103) dt(k+1),dq(k+1)
     1349        write(*,*) k,dt(k+1),dq(k+1)
     1350        k = k+1
     1351        go to 114
     1352 103    end_of_file = .true.
     1353 114    continue
     1354      enddo
     1355
     1356      close(unit=12,status = 'keep')
     1357
     1358      end subroutine read_lsf
     1359
    13191360END MODULE module_initialize_ideal
Note: See TracChangeset for help on using the changeset viewer.