!> \file bmelt-ant-regions_mod.f90 !! Module qui calcule la fusion basale (grounded ou ice shelves). !< !> \namespace BMELT_ANT_REGIONS !! Module qui calcule la fusion basale (grounded ou ice shelves) !! \author CatRitz !! \date juillet !! @note pour les ice shelves Antarctique, tient compte des différentes régions !! A choisir dans le module_choix !! @note Used module !! @note - use module3D_phy !< MODULE BMELT_ANT_REGIONS ! cat juillet 2005 USE module3D_phy implicit none REAL,dimension(nx,ny) :: bmgrz !< fusion basale a la grounding zone real,dimension(nx,ny) :: bmshelf !< fusion basale sous shelf REAL,dimension(nx,ny) :: dist_talu !< distance du point au talu continental REAL,dimension(nx,ny) :: typeshelf !< Type de shelf Ronne->1 Ross ->2 .... !< utilises pour moduler la fusion sous le shelf integer, dimension(10) :: region !< pour écrire dans le fichier param character(len=30),dimension(10) :: regname !< nom des régions real :: bsupshelf integer :: grd !< pour une sortie real :: bmelt_Ross,bmgrz_Ross !< fusion basale Ross real :: bmelt_FRis,bmgrz_FRis !< fusion basale filchner-Ronne ice shelf real :: bmelt_Amery,bmgrz_Amery !< Amery ice shelves real :: bmelt_PIG,bmgrz_PIG !< Pine Island glacier real :: bmelt_Pen,bmgrz_Pen !< Petits ice shelves Peninsule real :: bmelt_other,bmgrz_other !< Autre ice shelves real :: bmelt_talus,bmgrz_talus !< Au dela du talus continental CONTAINS !------------------------------------------------------------------------------- subroutine init_bmelt ! Cette routine fait l'initialisation pour la fusion basale. ! Elle est appelée par inputfile-vec-0.5.f90 namelist/bmelt_ant_reg/bmelt_Ross,bmgrz_Ross,bmelt_FRis,bmgrz_FRis, & bmelt_Amery,bmgrz_Amery,bmelt_PIG,bmgrz_PIG, & bmelt_Pen,bmgrz_Pen,bmelt_other,bmgrz_other, & bmelt_talus,bmgrz_talus rewind(num_param) ! pour revenir au debut du fichier param_list.dat read(num_param,bmelt_ant_reg) ! ecriture dans le fichier parametres write(num_rep_42,*)'fusion basale sous les ice shelves : bmelt-ant-regions_mod' write(num_rep_42,*)'----------------------------------------------------------' ! lecture du fichier contenant les distances au talu continental (m) open(88,file=TRIM(DIRNAMEINP)//'distance_talu-40km.xy') do j=1,ny do i=1,nx read(88,'(i3,1x,i3,1x,f10.2)') k,k,dist_talu(i,j) enddo enddo close(88) ! les bords de la grille sont mis a dist_talus=0 de force dist_talu(1,:)=0. !a gauche dist_talu(2,:)=0. dist_talu(nx-1,:)=0. ! a droite dist_talu(nx,:)=0. dist_talu(:,1)=0. ! en bas dist_talu(:,2)=0. dist_talu(:,ny-1)=0. ! en haut dist_talu(:,ny)=0. debug_3D(:,:,32)=dist_talu(:,:) ! lecture du fichier contenant les types de shelves ! Ronne-Flichner ->1, Ross -> 2 , Amery -> 3, ! PIG-> 4, les petits shelves au dessus de PIG -> 5 open(88,file=TRIM(DIRNAMEINP)//'numer-ice-shelves-juil07.dat',status='OLD') typeshelf(:,:)=100 do k=1,nx*ny read(88,*,end=36) i,j,typeshelf(i,j) end do 36 close(88) region(:)=0 regname(1)='Ronne-Fichner' regname(2)='Ross' regname(3)='Amery' regname(4)='Pig' regname(5)='Petits ice shelves peninsule' regname(6)='autres ice shelves' regname(7)='Au dela du talus continental' bms_init: do j=1,ny do i=1,nx talus: if (dist_talu(i,j).GT.5000.) then if (nint(typeshelf(i,j)).eq.1) then ! Ronne-Filchner FRis bmshelf(i,j)=bmelt_FRis bmgrz(i,j)=bmgrz_FRis if (region(1).eq.0) then ! pour l'impression dans 42 region(1)=1 write(num_rep_42,80)regname(1),bmshelf(i,j),bmgrz(i,j) endif 80 format(a32,' bmshelf=',f8.4,' bmgrz=',f8.4) else if (nint(typeshelf(i,j)).eq.2) then ! Ross bmshelf(i,j)=bmelt_Ross bmgrz(i,j)=bmgrz_Ross if (region(2).eq.0) then region(2)=1 write(num_rep_42,80)regname(2),bmshelf(i,j),bmgrz(i,j) endif else if (nint(typeshelf(i,j)).eq.3) then ! Amery bmshelf(i,j)=bmelt_Amery bmgrz(i,j)=bmgrz_Amery if (region(3).eq.0) then region(3)=1 write(num_rep_42,80)regname(3),bmshelf(i,j),bmgrz(i,j) endif else if (nint(typeshelf(i,j)).eq.4) then ! Pig bmshelf(i,j)=bmelt_PIG bmgrz(i,j)=bmgrz_PIG if (region(4).eq.0) then region(4)=1 write(num_rep_42,80)regname(4),bmshelf(i,j),bmgrz(i,j) endif else if (nint(typeshelf(i,j)).eq.5) then ! petits shelves peninsule bmshelf(i,j)=bmelt_Pen bmgrz(i,j)=bmgrz_Pen if (region(5).eq.0) then region(5)=1 write(num_rep_42,80)regname(5),bmshelf(i,j),bmgrz(i,j) endif else ! typeshelf(i,j)=6 bmshelf(i,j)=bmelt_other bmgrz(i,j)=bmgrz_other if (region(6).eq.0) then region(6)=1 write(num_rep_42,80)regname(6),bmshelf(i,j),bmgrz(i,j) endif endif else ! au dela du talus continental bmshelf(i,j)=bmelt_talus bmgrz(i,j)=bmgrz_talus if (region(7).eq.0) then region(7)=1 write(num_rep_42,80)regname(7),bmshelf(i,j),bmgrz(i,j) endif endif talus ! fin du test sur dist_talus enddo enddo bms_init debug_3D(:,:,33)=typeshelf(:,:) debug_3D(:,:,34)=bmshelf(:,:) return end subroutine init_bmelt !________________________________________________________________________________ subroutine bmeltshelf ! cette routine calcule la fusion basale proprement dite integer :: ngr ! nombre de voisins flottants real :: coef_talus ! pour ne pas changer la fusion au dessus de l'ocean profond do j=2,ny-1 do i=2,nx-1 talus_nochange : if (typeshelf(i,j).lt.10) then coef_talus = coefbmshelf else coef_talus = 1 ! pas de changement au dela du talus continental endif talus_nochange shelf: if (flot(i,j)) then ! partie flottante bmelt(i,j)=coef_talus*bmshelf(i,j) ! bloc bsupshelf à enlever ! if (time.gt.5000.) then ! bmelt(i,j)=bmelt(i,j)+bsupshelf ! endif if (fbm(i,j)) then bmelt(i,j)=coef_talus*bmgrz(i,j) endif ! if (time.gt.5000.) then ! bmelt(i,j)=bmelt(i,j)+bsupshelf ! endif ! ATTENTION le bloc suivant sert a determiner la fusion basale d'equilibre ! pour les shelves stationnaires if (igrdline.eq.1) then corrbmelt(i,j)=hdot(i,j)+bmelt(i,j) ! le bmelt d'equilibre debug_3D(i,j,28)=corrbmelt(i,j) endif else ! point posé, on compte le nombre de voisins flottants ngr=0 if (flot(i+1,j)) ngr=ngr+1 if (flot(i-1,j)) ngr=ngr+1 if (flot(i,j+1)) ngr=ngr+1 if (flot(i,j-1)) ngr=ngr+1 ! la fusion des points limites est une combinaison entre valeur posée et valeur flottante ! en fonction du nombre de points flottants bmelt(i,j)= ngr/4.*bmgrz(i,j)*coef_talus+(1.-ngr/4.)*bmelt(i,j) endif shelf end do end do ! les bords de la grille sont mis a dist_talus=0 de force ! attention, pourrait poser des problemes avec AGRIF bmelt(1,:)=bmelt_talus !a gauche bmelt(2,:)=bmelt_talus bmelt(nx-1,:)=bmelt_talus ! a droite bmelt(nx,:)=bmelt_talus bmelt(:,1)=bmelt_talus ! en bas bmelt(:,2)=bmelt_talus bmelt(:,ny-1)=bmelt_talus ! en haut bmelt(:,ny)=bmelt_talus if (igrdline.eq.1) then bmelt(:,:)=0. ! hdot donne alors le -bmelt endif ! ecriture des bmelt pour en faire des statistiques !-------------------------------------------------- !open(144,file='bmelt-test') !do j=1,ny ! do i=1,nx ! if (fbm(i,j)) then ! grd=1 ! else ! grd=0 ! endif ! if ((flot(i,j)).and.(H(i,j).gt.1.)) then ! write(144,999) i,j,grd,nint(typeshelf(i,j)), & ! bmelt(i,j),bmgrz(i,j),bmshelf(i,j),hdot(i,j),dist_talu(i,j) ! endif ! end do !end do !999 format(4(i0,2x),5(e12.3,2x)) !close(144) return end subroutine bmeltshelf END MODULE BMELT_ANT_REGIONS