Ignore:
Timestamp:
01/09/19 17:09:26 (6 years ago)
Author:
aquiquet
Message:

Sealevel is now treated as a 2D variable (sealevel_2d while sealevel remains the eustatic sea level), results should remain identical as sealevel_2d is equal to sealevel in this revision.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/SOURCES/bilan_eau_mod.f90

    r192 r237  
    4545subroutine init_bilan_eau 
    4646 ! initialisation des variables 
    47         diff_H=0. 
    48         sum_H_old = sum(H(2:nx-1,2:ny-1),mask=ice(2:nx-1,2:ny-1)==1) 
    49         tot_water(:,:)=0. 
    50         bm_dt(:,:)=0. 
    51         bmelt_dt(:,:)=0. 
     47        diff_H=0. 
     48        sum_H_old = sum(H(2:nx-1,2:ny-1),mask=ice(2:nx-1,2:ny-1)==1) 
     49        tot_water(:,:)=0. 
     50        bm_dt(:,:)=0. 
     51        bmelt_dt(:,:)=0. 
    5252        alpha_flot=ro/row 
    53         archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel 
     53        archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel_2d(:,:) 
    5454        gr_line(:,:)=0 
    5555        do j=1,ny 
    5656           do i=1,nx 
    57               !afq    if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel)) then ! grounded with ice 
     57              !afq    if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel_2d(i,j))) then ! grounded with ice 
    5858              if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.)) then ! grounded with ice 
    5959                 if (archimtab(i-1,j).LT.0..and.Uxbar(i,j).LT.0..and..not.flot_marais(i-1,j)) gr_line(i,j)=1 
     
    6767end subroutine init_bilan_eau 
    6868 
    69          
     69 
    7070 
    7171 
     
    101101 
    102102where (ice(2:nx-1,2:ny-1).eq.1) 
    103         Bm_dtt(2:nx-1,2:ny-1) = Bm_dtt(2:nx-1,2:ny-1) + Bm_dt(2:nx-1,2:ny-1) !* dt ! somme Bm sur dt 
    104         bmelt_dtt(2:nx-1,2:ny-1) = bmelt_dtt(2:nx-1,2:ny-1) + bmelt_dt(2:nx-1,2:ny-1) ! * dt ! somme bmelt sur dt 
     103        Bm_dtt(2:nx-1,2:ny-1) = Bm_dtt(2:nx-1,2:ny-1) + Bm_dt(2:nx-1,2:ny-1) !* dt ! somme Bm sur dt 
     104        bmelt_dtt(2:nx-1,2:ny-1) = bmelt_dtt(2:nx-1,2:ny-1) + bmelt_dt(2:nx-1,2:ny-1) ! * dt ! somme bmelt sur dt 
    105105endwhere 
    106106 
    107 archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel 
     107archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel_2d(:,:) 
    108108gr_line(:,:)=0 
    109109do j=1,ny 
    110110  do i=1,nx 
    111      !afq    if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel)) then ! grounded with ice 
     111     !afq    if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel_2d(i,j))) then ! grounded with ice 
    112112     if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.)) then ! grounded with ice 
    113113      if (archimtab(i-1,j).LT.0..and.Uxbar(i,j).LT.0..and..not.flot_marais(i-1,j)) gr_line(i,j)=1 
     
    141141    
    142142! bilan d'eau sur la grille : 
    143         water_bilan=sum(tot_water(:,:)) 
    144         diff_H = diff_H/dtt 
     143        water_bilan=sum(tot_water(:,:)) 
     144        diff_H = diff_H/dtt 
    145145 
    146146!999 format(f0.2,1x,e15.8,1x,i10,8(1x,e15.8)) 
    147147!       write(6,999),time,sum_H,count(ice(:,:)==1),diff_H,water_bilan,sum(calv_dtt(:,:))/dtt,sum(ablbord_dtt(:,:))/dtt,sum(bmelt_dtt(:,:),mask=ice(:,:)==1)/dtt,sum(bm(:,:),mask=ice(:,:)==1),sum(Bm_dtt(:,:))/dtt,sum(bmelt_dtt(:,:))/dtt 
    148         diff_H_water_bilan(2:nx-1,2:ny-1)=tot_water(2:nx-1,2:ny-1)-diff_H_2D(2:nx-1,2:ny-1) 
     148        diff_H_water_bilan(2:nx-1,2:ny-1)=tot_water(2:nx-1,2:ny-1)-diff_H_2D(2:nx-1,2:ny-1) 
    149149 
    150150endif 
Note: See TracChangeset for help on using the changeset viewer.