Changeset 146 for branches/iLoveclim/SOURCES/bilan_eau_mod.f90
- Timestamp:
- 10/20/17 09:31:39 (7 years ago)
- Location:
- branches/iLoveclim
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/iLoveclim
- Property svn:mergeinfo changed
/trunk merged: 124,127-145
- Property svn:mergeinfo changed
-
branches/iLoveclim/SOURCES/bilan_eau_mod.f90
r126 r146 31 31 real,dimension(nx,ny) :: Bmelt_dtt !< basal melting on ice points accumulated during dtt 32 32 real,dimension(nx,ny) :: calv_dtt !< calving sur dtt (pour calcul bilan d'eau) 33 real,dimension(nx,ny) :: archimtab ! point pose si > 0 34 real,dimension(nx,ny) :: grline_dtt ! grounding line flux during dtt 33 35 34 36 real :: sum_H_old 35 37 real :: diff_H 38 real :: alpha_flot ! ro/row 36 39 37 40 real,dimension(nx,ny) :: bm_dt,bmelt_dt … … 42 45 subroutine init_bilan_eau 43 46 ! initialisation des variables 44 diff_H=0. 45 sum_H_old = sum(H(2:nx-1,2:ny-1),mask=ice(2:nx-1,2:ny-1)==1) 46 tot_water(:,:)=0. 47 bm_dt(:,:)=0. 48 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. 52 alpha_flot=ro/row 49 53 ! iLOVECLIM initialisation of water conservation related variables 50 51 52 53 54 trendWAC=0. 55 smbWAC(:,:)=0. 56 bmeltWAC(:,:)=0. 57 calvingWAC(:,:)=0. 54 58 end subroutine init_bilan_eau 55 59 56 60 57 61 58 62 … … 88 92 89 93 where (ice(2:nx-1,2:ny-1).eq.1) 90 91 94 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 95 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 92 96 endwhere 97 98 archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel 99 gr_line(:,:)=0 100 do j=1,ny 101 do i=1,nx 102 !afq if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel)) then ! grounded with ice 103 if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.)) then ! grounded with ice 104 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 105 if (archimtab(i+1,j).LT.0..and.Uxbar(i+1,j).GT.0..and..not.flot_marais(i+1,j)) gr_line(i,j)=1 106 if (archimtab(i,j-1).LT.0..and.Uybar(i,j).LT.0..and..not.flot_marais(i,j-1)) gr_line(i,j)=1 107 if (archimtab(i,j+1).LT.0..and.Uybar(i,j+1).GT.0..and..not.flot_marais(i,j+1)) gr_line(i,j)=1 108 endif 109 enddo 110 enddo 111 112 where (gr_line(:,:).eq.1) 113 !~ grline_dtt(:,:)= (((uxbar(:,:)+eoshift(uxbar(:,:),shift=1,boundary=0.,dim=1))**2+ & 114 !~ (uybar(:,:)+eoshift(uybar(:,:),shift=1,boundary=0.,dim=2))**2)**0.5)*0.5 & 115 !~ *H(:,:) + grline_dtt(:,:) 116 grline_dtt(:,:)= - sqrt( & 117 ( (uxbar(:,:)+eoshift(uxbar(:,:),shift=1,boundary=0.,dim=1))*dy/2. )**2 + & 118 ( (uybar(:,:)+eoshift(uybar(:,:),shift=1,boundary=0.,dim=2))*dx/2. )**2 ) & 119 * H(:,:) * dt + grline_dtt(:,:) 120 endwhere 121 93 122 94 123 if (isynchro.eq.1) then … … 103 132 104 133 ! bilan d'eau sur la grille : 105 106 134 water_bilan=sum(tot_water(:,:)) 135 diff_H = diff_H/dtt 107 136 108 137 !999 format(f0.2,1x,e15.8,1x,i10,8(1x,e15.8)) 109 138 ! 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 110 139 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) 111 140 112 141 endif
Note: See TracChangeset
for help on using the changeset viewer.