Ignore:
Timestamp:
Jun 15, 2024, 5:17:08 PM (11 days ago)
Author:
crisi
Message:

suppress isotope_params.def + update physiq_mod + proof of concept of 3rd dimension with reevap routine

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/trunk/libf/phylmdiso/isotopes_routines_mod.F90

    r4491 r4982  
    15251525#endif       
    15261526       pxtfra=max(min(pxtfra,alpha_max),0.0)
    1527 
    15281527
    15291528      end subroutine fractcalk_liq
     
    1592215921
    1592315922      ! verif
    15924 !      text="phyisoetat0 67"
    15925 !      write(*,*) 'snow(8,1)=',snow(8,1)
    15926 !      write(*,*) 'xtsnow(4,8,1)=',xtsnow(4,8,1)
    1592715923#ifdef ISOVERIF
    1592815924      do i=1,klon
     
    1593415930         enddo !do ixt=1,niso
    1593515931      enddo !do i=1,klon
    15936 #endif     
    15937 #ifdef ISOVERIF
    1593815932      do i=1,klon
    1593915933         if (iso_eau.gt.0) then
     
    1602116015         endif
    1602216016       enddo !do i=1,klon
    16023 
    1602416017#endif
    1602516018      !end verif
     
    1612816121          deltaD_run_off_lic_0(ixt)=deltaD_sol(ixt)
    1612916122          deltaD_land_ice(ixt)=deltaD_snow(ixt)
    16130           call fractcalk_liq(ixt, 283.0, alpha(ixt))           
     16123          call fractcalk_liq(ixt, 283.0, alpha(ixt))   
    1613116124        enddo !do ixt=1,niso
    1613216125        call calcul_kcin(2.0,kcin)
     
    1883018823        if ((iso_HDO.gt.0).and.(ixt.eq.iso_HDO)) then
    1883118824            if (q.gt.ridicule) then
     18825                    write(*,*) 'xt,q=',xt,q
     18826                    write(*,*) 'alpha=',alpha
     18827                    write(*,*) 'toce,kcin,h0=',toce,kcin,h0
     18828                    write(*,*) 'RMerlivat=',RMerlivat
    1883218829                call iso_verif_aberrant_encadre( xt/q, 'isotopes_routines_mod 18930b: iso_init_ideal')
    1883318830            endif
     
    1890218899end subroutine appel_stewart_debug
    1890318900
     18901
     18902subroutine dispatch(klon,klev,qx,q_seri,xt_seri,ql_seri,xtl_seri,qs_seri,xts_seri)
     18903
     18904use infotrac_phy, ONLY: nqtot,nqo,ivap,iliq,isol,iqIsoPha,ntraciso=>ntiso
     18905implicit none
     18906
     18907! inputs
     18908integer, intent(in) :: klon,klev
     18909real,dimension(klon,klev,nqtot), intent(in) ::qx
     18910
     18911! outputs
     18912real,dimension(klon,klev), intent(out) ::q_seri,ql_seri,qs_seri
     18913real,dimension(ntraciso,klon,klev), intent(out) :: xt_seri,xtl_seri,xts_seri
     18914
     18915! locals
     18916integer :: i,k,ixt
     18917
     18918do k=1,klev
     18919do i=1,klon
     18920    q_seri(i,k)  = qx(i,k,ivap)
     18921    ql_seri(i,k) = qx(i,k,iliq)
     18922    IF (nqo.EQ.2) THEN             !--vapour and liquid only
     18923             qs_seri(i,k) = 0.
     18924    ELSE IF (nqo.ge.3) THEN        !--vapour, liquid and ice
     18925             qs_seri(i,k) = qx(i,k,isol)
     18926    ENDIF
     18927    do ixt=1,ntraciso
     18928          xt_seri(ixt,i,k)  = qx(i,k,iqIsoPha(ixt,ivap))
     18929          xtl_seri(ixt,i,k) = qx(i,k,iqIsoPha(ixt,iliq))
     18930          if (nqo.eq.2) then
     18931             xts_seri(ixt,i,k) = 0.
     18932          else if (nqo.eq.3) then
     18933             xts_seri(ixt,i,k) = qx(i,k,iqIsoPha(ixt,isol))
     18934          endif
     18935    enddo !do ixt=1,niso
     18936
     18937enddo
     18938enddo
     18939
     18940end subroutine dispatch
     18941
     18942subroutine together(klon,klev,qx,q_seri,xt_seri,ql_seri,xtl_seri,qs_seri,xts_seri)
     18943
     18944use infotrac_phy, ONLY: nqtot,nqo,ivap,iliq,isol,iqIsoPha,ntraciso=>ntiso
     18945implicit none
     18946
     18947! inputs
     18948integer, intent(in) :: klon,klev
     18949real,dimension(klon,klev), intent(in) ::q_seri,ql_seri,qs_seri
     18950real,dimension(ntraciso,klon,klev), intent(in) :: xt_seri,xtl_seri,xts_seri
     18951
     18952! inputs
     18953real,dimension(klon,klev,nqtot), intent(out) ::qx
     18954
     18955! locals
     18956integer :: i,k,ixt
     18957
     18958do k=1,klev
     18959do i=1,klon 
     18960    qx(i,k,ivap)  = q_seri(i,k)
     18961    qx(i,k,iliq) = ql_seri(i,k)
     18962    IF (nqo.ge.3) THEN        !--vapour, liquid and ice
     18963             qx(i,k,isol) = qs_seri(i,k)
     18964    ENDIF
     18965    do ixt=1,ntraciso
     18966          qx(i,k,iqIsoPha(ixt,ivap)) = xt_seri(ixt,i,k) 
     18967          qx(i,k,iqIsoPha(ixt,iliq)) = xtl_seri(ixt,i,k)
     18968          if (nqo.ge.3) then
     18969             qx(i,k,iqIsoPha(ixt,isol)) = xts_seri(ixt,i,k)
     18970          endif
     18971    enddo !do ixt=1,niso
     18972
     18973enddo
     18974enddo
     18975
     18976end subroutine together
     18977
     18978
    1890418979END MODULE isotopes_routines_mod
    1890518980#endif
Note: See TracChangeset for help on using the changeset viewer.