! Construit les cartes qui servent dans le recul de grounding line ! et teste les methodes de recul. ! la diminution de H est lineaire avec le recul module toy_retreat_mod use io_netcdf_grisli use declar_toy_retreat implicit none contains !---------------------------------------------------------------------------------- ! subroutine time_step_recu : ce qui se fait a chaque pas de temps ! subroutine init_proto_recul : lit les donnees et initialise ! subroutine init_masques : initialise les masques ! subroutine pos_gr_line_init : Calcule les positions initiales de Grounding line ! subroutine update_retreat_rates : en fonction du temps ! subroutine calc_delHdt_maj : calcul le dHdt et le new H sur chaque maille ! subroutine calc_H_float : calcule la hauteur de flottaison ! subroutine update_flottants : points qui deviennent flottant ! subroutine update_a_traiter : points qui deviennent a traiter ! subroutine prescribe_Hp : renvoie Hp et I_Hp vers GRISLI !---------------------------------------------------------------------------------- !---------------------------------------------------------------------------------- ! time_step_recul : ce qui se fait a chaque pas de temps !---------------------------------------------------------------------------------- subroutine time_step_recul implicit none call update_retreat_rates ! update les vitesses de retrait (autorise que si time > time_dep) call calc_eps_max ! pour le sanity check call calc_delHdt_maj ! calcule delHdt et nouveau H pour chaque point call calc_H_float ! calcule la valeur de H_flot call Schoof_flux ! calcul coef_schoof, bf_x_schoof, bf_y_schoof call update_flottants ! update les nouveaux points flottants et xgl, delHdx associes ! call update_ramollo ! ramolli eventuellement les ice shelves apres le temps de demarrage call prescribe_Hp ! renvoie Hp et I_Hp vers GRISLI call update_a_traiter ! update le masque des points a traiter ! quelques sorties debug_3D(:,:,94) = time_float(:,:) debug_3D(:,:,95) = ux_schoof(:,:) debug_3D(:,:,96) = uy_schoof(:,:) end subroutine time_step_recul !---------------------------------------------------------------------------------- ! init_proto_recul : lit les donnees et initialise !---------------------------------------------------------------------------------- subroutine init_proto_recul use declar_toy_retreat implicit none namelist/retreat_ice2sea/region_file,file_exp rewind(num_param) ! pour revenir au debut du fichier param_list.dat 428 format(A) read(num_param,retreat_ice2sea) write(num_rep_42,428) '!___________________________________________________________' write(num_rep_42,428) '! read parameters of retreat' write(num_rep_42,retreat_ice2sea) write(num_rep_42,428) '! file_exp est le fichier qui contient toutes les experiences' write(num_rep_42,428) '! le numero de l experience correspond au numero du job' write(num_rep_42,428) '!___________________________________________________________' call lect_retreat_files ! call lect_retreat_nolinfiles ! pour pouvoir lire aussi l'exposant de la loi de frottement call init_masques ! initialise mk_traiter, retreat0_x, retreat0_y, Deltatot_H0, call pos_gr_line_init ! calcule x_gl_mx, y_gl_my, delHdx_mx, delHdy_my a l'initialisation ! ensuite cette position est updated dans la boucle call calc_eps_max delHdt_sanity(:,:) = 5000. call init_Schoof_flux end subroutine init_proto_recul !------------------------------------------------------------------------- !< init_masques : initialise les masques !------------------------------------------------------------------------- subroutine init_masques !-------------------------------------------------------------------------- ! Pour chaque noeud pose dont le socle est sous le niveau des mers ! ! definition des noeuds a traiter !------------------------------------ ! ! mk_traiter = 2 : pose loin de gl ! mk_traiter = 1 : en cours de recul ! mk_traiter = -1: flottant implique dans le recul ! mk_traiter = -2: flottant loin de GL ! mk_traiter = 0 : pas a traiter use declar_toy_retreat implicit none integer :: som_voisins !< pour des tests ! definition du masque a traiter where ((Mk_gr(:,:).eq.1).and.(Bsoc(:,:).lt.0.)) ! points à traiter mk_traiter(:,:) = 1 H_float(:,:) = Bsoc(:,:) / Coef_bflot elsewhere ((Mk_gr(:,:).eq.1).and.(Bsoc(:,:).ge.0.)) ! pose et stable mk_traiter(:,:) = 0 H_float(:,:) = -10000. elsewhere (Mk_gr(:,:).eq.0) mk_traiter(:,:) = -1 H_float (:,:) = Bsoc(:,:) / Coef_bflot end where ! calcul de la hauteur de flottaison pour le socle test where ((Mk_gr(:,:).eq.1).and.(B_test(:,:).lt.0.)) ! points à traiter H_float_test(:,:) = B_test(:,:) / Coef_bflot elsewhere ((Mk_gr(:,:).eq.1).and.(B_test(:,:).ge.0.)) ! pose et stable H_float_test(:,:) = -10000. elsewhere (Mk_gr(:,:).eq.0) H_float_test(:,:) = B_test(:,:) / Coef_bflot end where do j=2,ny-1 do i =2,nx-1 if (mk_traiter(i,j).eq.1) then ! pose et a traiter if (any(mk_gr(i-1:i+1,j-1:j+1).eq.0)) then mk_traiter(i,j) = 1 ! tout de suite else mk_traiter(i,j) = 2 ! apres propagation end if end if if (mk_traiter(i,j).eq.-1) then if (any(mk_gr(i-1:i+1,j-1:j+1).eq.1)) then mk_traiter(i,j) = -1 ! proche de la grounding line else mk_traiter(i,j) = -2 ! loin de la grounding line end if end if end do end do call init_retreat0_alpha call init_time_depart call update_retreat_rates ! retreat rate depend du temps ramollo(:,:) = 0. ! les ice shelves ont la viscosite standard mk_traiter0(:,:) = mk_traiter(:,:) end subroutine init_masques !--------------------------------------------------------------------------- !< pos_gr_line : Calcule les positions initiales de Grounding line !--------------------------------------------------------------------------- subroutine pos_gr_line_init ! Calcul des positions de Grounding line !__________________________________________ ! ! les points sont sur les mailles mineures. ! la coordonnee est positive quand le point grounded est à l'Est ou au Nord ! ! ~--------------x------------o ! i-1 i x_gl_mx >0 ! float grounded ! <- West -> East ! ! si deux points flottants: x_gl_mx = -10 dx ! si deux points poses : x_gl_mx = 10 dx use declar_toy_retreat implicit none ! variables locales real :: H_0 ! epaisseur en amont gl real :: H_1 ! epaisseur en aval gl real :: B_0 ! socle en amont gl real :: B_1 ! socle en aval gl real :: dyy ! variable de travail longueur de la maille en diagonale dyy = dx*(2.**0.5) !longueur de la diagonale do j=2,ny-1 do i=2,nx-1 ! recherche position sous maille selon x ------------------------------------------- x_gl: if ((mk_gr(i,j).eq.1).and.(mk_gr(i-1,j).eq.0)) then ! grounded a l'Est H_1 = H(i-1,j) H_0 = H(i,j) B_1 = Bsoc(i-1,j) B_0 = Bsoc(i,j) call calc_xgl(dx,coef_Bflot,H_0,B_0,H_1,B_1, x_gl_mx(i,j),delHdx_mx(i,j)) else if ((mk_gr(i,j).eq.0).and.(mk_gr(i-1,j).eq.1)) then ! grounded a l'West H_1 = H(i,j) H_0 = H(i-1,j) B_1 = Bsoc(i,j) B_0 = Bsoc(i-1,j) call calc_xgl(dx,coef_Bflot,H_0,B_0,H_1,B_1, x_gl_mx(i,j),delHdx_mx(i,j)) x_gl_mx(i,j) = - x_gl_mx(i,j) ! on change le signe de x_gl else if ((mk_gr(i,j).eq.0).and.(mk_gr(i-1,j).eq.0)) then ! tout le monde flottant x_gl_mx(i,j) = -10.* dx delHdx_mx(i,j) = 0. else ! tout le monde pose x_gl_mx(i,j) = 10.*dx delHdx_mx(i,j) = 0. end if x_gl ! recherche position sous maille selon y ------------------------------------------- y_gl: if ((mk_gr(i,j).eq.1).and.(mk_gr(i,j-1).eq.0)) then ! grounded au Nord H_1 = H(i,j-1) H_0 = H(i,j) B_1 = Bsoc(i,j-1) B_0 = Bsoc(i,j) call calc_xgl(dx,coef_Bflot,H_0,B_0,H_1,B_1, y_gl_my(i,j),delHdy_my(i,j)) else if ((mk_gr(i,j).eq.0).and.(mk_gr(i,j-1).eq.1)) then ! grounded au Sud H_1 = H(i,j) H_0 = H(i,j-1) B_1 = Bsoc(i,j) B_0 = Bsoc(i,j-1) call calc_xgl(dx,coef_Bflot,H_0,B_0,H_1,B_1, y_gl_my(i,j),delHdy_my(i,j)) y_gl_my(i,j) = - y_gl_my(i,j) ! on change le signe de y_gl else if ((mk_gr(i,j).eq.0).and.(mk_gr(i,j-1).eq.0)) then ! tout le monde flottant y_gl_my(i,j) = -10.*dx delHdy_my(i,j) = 0. else ! tout le monde pose y_gl_my(i,j) = 10.*dx delHdy_my(i,j) = 0. end if y_gl ! recherche position sous maille selon diagonale SW-NE ------------------------------------------- SW_gl: if ((mk_gr(i,j).eq.1).and.(mk_gr(i-1,j-1).eq.0)) then ! grounded au nord-est H_1 = H(i-1,j-1) H_0 = H(i,j) B_1 = Bsoc(i-1,j-1) B_0 = Bsoc(i,j) call calc_xgl(dyy,coef_Bflot,H_0,B_0,H_1,B_1, SW_gl_m(i,j),delHdx_SW(i,j)) else if ((mk_gr(i,j).eq.0).and.(mk_gr(i-1,j-1).eq.1)) then ! grounded au Sud-west H_1 = H(i,j) H_0 = H(i-1,j-1) B_1 = Bsoc(i,j) B_0 = Bsoc(i-1,j-1) call calc_xgl(dyy,coef_Bflot,H_0,B_0,H_1,B_1, SW_gl_m(i,j),delHdx_SW(i,j)) SW_gl_m(i,j) = - SW_gl_m(i,j) ! on change le signe de SW_gl else if ((mk_gr(i,j).eq.0).and.(mk_gr(i-1,j-1).eq.0)) then ! tout le monde flottant SW_gl_m(i,j) = -10.*dx delHdx_SW(i,j) = 0. else ! tout le monde pose SW_gl_m(i,j) = 10.*dx delHdx_SW(i,j) = 0. end if SW_gl ! recherche position sous maille selon diagonale SE-NW ------------------------------------------- SE_gl: if ((mk_gr(i,j).eq.1).and.(mk_gr(i+1,j-1).eq.0)) then ! grounded au nord-west H_1 = H(i+1,j-1) H_0 = H(i,j) B_1 = Bsoc(i+1,j-1) B_0 = Bsoc(i,j) call calc_xgl(dyy,coef_Bflot,H_0,B_0,H_1,B_1,SE_gl_m(i,j),delHdx_SE(i,j)) else if ((mk_gr(i,j).eq.0).and.(mk_gr(i+1,j-1).eq.1)) then ! grounded au Sud-Est H_1 = H(i,j) H_0 = H(i+1,j-1) B_1 = Bsoc(i,j) B_0 = Bsoc(i+1,j-1) call calc_xgl(dyy,coef_Bflot,H_0,B_0,H_1,B_1,SE_gl_m(i,j),delHdx_SE(i,j)) SE_gl_m(i,j) = - SE_gl_m(i,j) ! on change le signe de SE_gl else if ((mk_gr(i,j).eq.0).and.(mk_gr(i+1,j-1).eq.0)) then ! tout le monde flottant SE_gl_m(i,j) = -10.*dx delHdx_SE(i,j) = 0. else ! tout le monde pose SE_gl_m(i,j) = 10.*dx delHdx_SE(i,j) = 0. end if SE_gl end do end do return end subroutine pos_gr_line_init !------------------------------------------------------------------------- !< update_retreat_rates : en fonction du temps !------------------------------------------------------------------------- subroutine update_retreat_rates use declar_toy_retreat implicit none call init_retreat0_alpha ! reinitialise retreat0 pour tenir compte des variations beta where (time_dep_mx(:,:).lt.time) retreat_x(:,:) = retreat0_x(:,:) elsewhere retreat_x(:,:) = 0. end where where (time_dep_my(:,:).lt.time) retreat_y(:,:) = retreat0_y(:,:) elsewhere retreat_y(:,:) = 0. end where where (time_dep_SW(:,:).le.time) retreat_SW(:,:) = retreat0_SW(:,:) elsewhere retreat_SW(:,:) = 0. end where where (time_dep_SE(:,:).le.time) retreat_SE(:,:) = retreat0_SE(:,:) elsewhere retreat_SE(:,:) = 0. end where return end subroutine update_retreat_rates !------------------------------------------------------------------------ !< pour ramollir eventuellement les ice shelves !------------------------------------------------------------------------ subroutine update_ramollo use declar_toy_retreat implicit none where ((time_dep(:,:).lt.time).and.(flot(:,:))) ramollo(:,:) = 1. elsewhere ramollo(:,:) = 0. end where return end subroutine update_ramollo !------------------------------------------------------------------------- !< calc_delHdt_maj : calcul le dHdt et le new H sur chaque maille !------------------------------------------------------------------------- subroutine calc_delHdt_maj use declar_toy_retreat implicit none real,dimension(2) :: dH_x_voisin ! quand on balaye les voisins real,dimension(2) :: dH_y_voisin ! quand on balaye les voisins real,dimension(2) :: dH_SW_voisin ! quand on balaye les voisins real,dimension(2) :: dH_SE_voisin ! quand on balaye les voisins real :: max_x ! valeur max de dH_x_voisin real :: max_y ! valeur max de dH_y_voisin real :: max_SW ! valeur max de dH_SW_voisin real :: max_SE ! valeur max de dH_SE_voisin integer :: som_voisins !< pour des tests ! conditions grounding line du point de vue du noeud majeur i ! ~--------------x--------------o---------------x-----------------~ ! i-1 i i+1 ! 0 < x_gl(i) < 10 dx -10 dx < x_gl(i+1) < 0 do j=2,ny-1 do i=2,nx-1 point_a_traiter: if (mk_traiter(i,j).eq.1) then ! Point à traiter som_voisins = 0 ! balaye les voisins en croix dH_x_voisin(:) = 0. dH_y_voisin(:) = 0. dH_SW_voisin(:) = 0. ! et en diagonale dH_SE_voisin(:) = 0. ! if ((mk_traiter(i-1,j).eq.-1).and. & ! voisin west--------- (x_gl_mx(i,j).lt.10.*dx).and.(x_gl_mx(i,j).ge.0.)) then dH_x_voisin(1) = delHdx_mx(i,j)*retreat_x(i,j) som_voisins = som_voisins + 1 end if if ((mk_traiter(i+1,j).eq.-1).and. & ! voisin East--------- (x_gl_mx(i+1,j).gt.-10.*dx).and.(x_gl_mx(i+1,j).le.0.)) then dH_x_voisin(2) = delHdx_mx(i+1,j)*retreat_x(i+1,j) som_voisins = som_voisins + 1 end if if ((mk_traiter(i,j-1).eq.-1).and. & ! voisin Sud--------- (y_gl_my(i,j).lt.10.*dx).and.(y_gl_my(i,j).ge.0.)) then dH_y_voisin(1) = delHdy_my(i,j)*retreat_y(i,j) som_voisins = som_voisins + 1 end if if ((mk_traiter(i,j+1).eq.-1).and. & ! voisin Nord--------- (y_gl_my(i,j+1).gt.-10.*dx).and.(y_gl_my(i,j+1).le.0.)) then dH_y_voisin(2) = delHdy_my(i,j+1)*retreat_y(i,j+1) som_voisins = som_voisins + 1 endif ! en diagonale if ((mk_traiter(i-1,j-1).eq.-1).and. & ! voisin SW--------- (SW_gl_m(i,j).lt.10.*dx).and.(SW_gl_m(i,j).ge.0.)) then dH_SW_voisin(1) = delHdx_SW(i,j)*retreat_SW(i,j) som_voisins = som_voisins + 1 end if if ((mk_traiter(i+1,j+1).eq.-1).and. & ! voisin NE--------- (SW_gl_m(i+1,j+1).gt.-10.*dx).and.(SW_gl_m(i+1,j+1).le.0.)) then dH_SW_voisin(2) = delHdx_SW(i+1,j+1)*retreat_SW(i+1,j+1) som_voisins = som_voisins + 1 end if if ((mk_traiter(i+1,j-1).le.-1).and. & ! voisin SE--------- (SE_gl_m(i,j).lt.10.*dx).and.(SE_gl_m(i,j).ge.0.)) then dH_SE_voisin(1) = delHdx_SE(i,j)*retreat_SE(i,j) som_voisins = som_voisins + 1 end if if ((mk_traiter(i-1,j+1).le.-1).and. & ! voisin NW--------- (SE_gl_m(i-1,j+1).gt.-10.*dx).and.(SE_gl_m(i-1,j+1).le.0.)) then dH_SE_voisin(2) = delHdx_SE(i-1,j+1)*retreat_SE(i-1,j+1) som_voisins = som_voisins + 1 endif if (som_voisins.eq.0) then write(6,'("probleme pour i,j=",11(i3,x))'), i,j,mk_traiter(i,j) write(6,*)'j-1', mk_gr(i-1:i+1,j-1) write(6,*)'j ', mk_gr(i-1:i+1,j) write(6,*)'j+1', mk_gr(i-1:i+1,j+1) end if ! seuls les voisins au dessus peuvent venir ! pour les differentes version VERIFER LA LECTURE de B_TEST et B_VOISIN ! utilise un autre socle que celui du modele pour faire le test ! call teste_socle_Btest(i,j,dH_x_voisin,dH_y_voisin,dH_SW_voisin,dH_SE_voisin) call teste_Btest_schoof(i,j,dH_x_voisin,dH_y_voisin,dH_SW_voisin,dH_SE_voisin) ! utilise les socles min et max pour faire les tests ! call teste_socle_Bminmax(i,j,dH_x_voisin,dH_y_voisin,dH_SW_voisin,dH_SE_voisin) ! utilise le socle du model pour le test ! call teste_socle_Bsoc(i,j,dH_x_voisin,dH_y_voisin,dH_SW_voisin,dH_SE_voisin) ! only for tests ! calcul de la variation d'epaisseur pour ce point max_x = maxval(dH_x_voisin) max_y = maxval(dH_y_voisin) max_SW = maxval(dH_SW_voisin) max_SE = maxval(dH_SE_voisin) delHdt(i,j) = max(max_x,max_y,max_SE,max_SW) ! call sanity_check call sanity_check(i,j) ! call no_sanity_check(i,j) delHdt(i,j) = min(delHdt(i,j),delHdt_sanity(i,j)) ! pour tester un recul maxi : le sanity check est aussi utilis comme minimum ! delHdt(i,j) = max(delHdt(i,j),delHdt_sanity(i,j)*0.01) ! delHdt(i,j) = max(delHdt(i,j),0.) H(i,j) = H(i,j)-delHdt(i,j)*dt ! attention peut etre sous-flottaison ! sanity check end if point_a_traiter end do end do end subroutine calc_delHdt_maj !------------------------------------------------------------------------- !< calc_H_float : calcule la hauteur de flottaison !------------------------------------------------------------------------- subroutine calc_H_float use declar_toy_retreat implicit none H_float(:,:) = 0. where ((Mk_gr(:,:).gt.0).and.(Bsoc(:,:).lt.0.)) H_float(:,:) = Bsoc(:,:) / Coef_bflot elsewhere ((Mk_gr(:,:).gt.0).and.(Bsoc(:,:).ge.0.)) H_float(:,:) = -10000. elsewhere (Mk_gr(:,:).eq.0) H_float(:,:) = Bsoc(:,:) / Coef_bflot end where end subroutine calc_H_float !------------------------------------------------------------------------- !< update_flottants : points qui deviennent flottant !------------------------------------------------------------------------- subroutine update_flottants use declar_toy_retreat implicit none integer :: som_voisins !< pour des tests real :: Hf !< Hf = 0 si socle > 0 sinon H_float real :: dyy ! variable de travail longueur de la maille en diagonale dyy = dx*(2.**0.5) !longueur de la diagonale ! update des masques : points qui deviennent flottant do j=2,ny-1 do i=2,nx-1 new_float: if ((mk_traiter(i,j).eq.1).and.(H(i,j).le.H_float(i,j))) then ! le noeud est passe flottant mk_traiter(i,j) = -1 time_float(i,j) = time H(i,j) = H_float(i,j)-20. mk_gr(i,j) = 0 ! update vers les noeuds voisins if (mk_gr(i-1,j).eq.1) then ! noeud west est grounded x_gl_mx(i,j) = -(dx-1.) ! position gl proche de i,j Hf = max(0.,H_float(i-1,j)) delHdx_mx(i,j) = ( H(i-1,j) - Hf )/dx ! taux d'amincissement de (i-1,j) else x_gl_mx(i,j) = -10.*dx delHdx_mx(i,j) = 0. end if if (mk_gr(i+1,j).eq.1) then ! noeud Est est grounded x_gl_mx(i+1,j) = dx-1. ! position gl proche de i,j Hf = max(0.,H_float(i+1,j)) delHdx_mx(i+1,j) = ( H(i+1,j) - Hf )/dx ! taux d'amincissement de (i+1,j) else x_gl_mx(i+1,j) = -10.*dx delHdx_mx(i+1,j) = 0. end if if (mk_gr(i,j-1).eq.1) then ! noeud Sud est grounded y_gl_my(i,j) = -(dx-1.) ! position gl proche de i,j Hf = max(0.,H_float(i,j-1)) delHdy_my(i,j) = ( H(i,j-1) - Hf )/dx ! taux d'amincissement de (i-1,j) else y_gl_my(i,j) = -10.*dx delHdy_my(i,j) = 0. end if if (mk_gr(i,j+1).eq.1) then ! noeud Nord est grounded y_gl_my(i,j+1) = dx-1. ! position gl proche de i,j Hf = max(0.,H_float(i,j+1)) delHdy_my(i,j+1) = ( H(i,j+1) - Hf )/dx ! taux d'amincissement de (i+1,j) else y_gl_my(i,j+1) = -10.*dx delHdy_my(i,j+1) = 0. end if ! en diagonale if (mk_gr(i-1,j-1).eq.1) then ! noeud SW grounded majeur (i-1,j-1) SW_gl_m(i,j) = -(dyy-1.) ! position gl proche de i,j <0 Hf = max(0.,H_float(i-1,j-1)) ! noeud mineur i,j delHdx_SW(i,j) = ( H(i-1,j-1) - Hf )/dyy ! taux d'amincissement de (i-1,j-1) else SW_gl_m(i,j) = -10.*dx delHdx_SW(i,j) = 0. end if if (mk_gr(i+1,j+1).eq.1) then ! noeud NE grounded majeur (i+1,j+1) SW_gl_m(i+1,j+1) = dyy-1. ! position gl proche de i,j Hf = max(0.,H_float(i+1,j+1)) ! noeud mineur (i+1,j+1) delHdx_SW(i+1,j+1) = ( H(i+1,j+1) - Hf )/dyy ! taux d'amincissement de (i+1,j+1) else SW_gl_m(i+1,j+1) = -10.*dx delHdx_SW(i+1,j+1) = 0. end if if (mk_gr(i+1,j-1).eq.1) then ! noeud SE grounded majeur (i+1,j-1) SE_gl_m(i,j) = -(dyy-1.) ! position gl proche de i,j <0 Hf = max(0.,H_float(i+1,j-1)) ! noeud mineur i,j delHdx_SE(i,j) = ( H(i+1,j-1) - Hf )/dyy ! taux d'amincissement de (i+1,j-1) else SE_gl_m(i,j) = -10.*dx delHdx_SE(i,j) = 0. end if if (mk_gr(i-1,j+1).eq.1) then ! noeud NW grounded majeur (i-1,j+1) SE_gl_m(i-1,j+1) = dyy-1. ! position gl proche de i,j Hf = max(0.,H_float(i-1,j+1)) ! noeud mineur (i-1,j+1) delHdx_SE(i-1,j+1) = ( H(i-1,j+1) - Hf )/dyy ! taux d'amincissement de (i-1,j+1) else SE_gl_m(i-1,j+1) = -10.*dx delHdx_SE(i-1,j+1) = 0. end if end if new_float end do end do end subroutine update_flottants !------------------------------------------------------------------------- !< update_a_traiter : points qui deviennent a traiter !------------------------------------------------------------------------- subroutine update_a_traiter use declar_toy_retreat implicit none integer :: som_voisins !< pour des tests ! les nouveaux points a traiter do j=2,ny-1 do i =2,nx-1 ! formulation en croix ! som_voisins = mk_gr(i-1,j)+mk_gr(i+1,j)+mk_gr(i,j-1)+mk_gr(i,j+1) ! if ((mk_traiter(i,j).eq.2).and.(som_voisins.gt.0).and.(som_voisins.lt.4)) then ! mk_traiter(i,j) = 1 ! end if ! Formulation avec diagonales if ((mk_traiter(i,j).eq.2).and.(any(mk_gr(i-1:i+1,j-1:j+1).eq.0))) then mk_traiter(i,j) = 1 end if end do end do end subroutine update_a_traiter !------------------------------------------------------------------------- !< prescribe_Hp : renvoie Hp et I_Hp vers GRISLI !------------------------------------------------------------------------- subroutine prescribe_Hp use declar_toy_retreat implicit none ! remet les masques Hp et i_Hp aux valeurs fixes i_HP(:,:) = i_Hp0(:,:) Hp(:,:) = Hp0(:,:) i=208 j=176 where (mk_gr(:,:).eq.0) ! points flottants i_Hp(:,:) = 1 Hp (:,:) = H(:,:) ! les points flottants ne sont pas modifies end where where(mk_traiter(:,:).eq.1) ! points imposes Epaisseur en cours de decroissance i_Hp(:,:) = 1 Hp (:,:) = H(:,:) end where end subroutine prescribe_Hp !------------------------------------------------------------------------- !< calc_eps_max : calcule le epsilon max (sanity check) !------------------------------------------------------------------------- subroutine calc_eps_max use declar_toy_retreat implicit none real :: gamma !< coefficient de flottaison gamma = ro*g*(1.+coef_Bflot) epsmax(:,:) = Abar(:,:)*H(:,:)**3*(gamma/4.)**3 end subroutine calc_eps_max !------------------------------------------------------------------------- !< calc_xgl : calcule la position sous maille de la grounding line !< en supposant que l'epaisseur varie lineairement !------------------------------------------------------------------------- subroutine calc_xgl(dy,alpha,H_0,B_0,H_1,B_1,xpos,delHdx) implicit none ! dummy real,intent(in) :: dy !< longueur de la maille real,intent(in) :: alpha !< coefficient de flottaison = coef_Bflot real,intent(in) :: H_0 !< epaisseur au point pose real,intent(in) :: B_0 !< altitude socle au point pose real,intent(in) :: H_1 !< epaisseur au point flottant real,intent(in) :: B_1 !< altitude socle au point flottant real,intent(out) :: xpos !< position de la ligne (en distance depuis le point pose real,intent(out) :: delHdx !< variation d'epaisseur au noeud pose en fonction d'une variation de xpos real :: Cgl ! variable de travail Cgl = (alpha * (H_1-H_0) - (B_1-B_0)) if (abs(Cgl).gt.1.e-5) then xpos = dy * (B_0 - alpha * H_0) / Cgl else xpos = dy-1. ! verifier end if ! la variation d'epaisseur est proportionnelle au recul if (xpos.gt.1.) then delHdx = (H_0 - B_0 /alpha) / xpos else delHdx = (H_0 - B_0 /alpha) end if end subroutine calc_xgl !---------------------------------------------------------------------- ! lect_retreat_files : lit les fichiers specifiques au retrait !---------------------------------------------------------------------- subroutine lect_retreat_files use netcdf use io_netcdf_grisli use declar_toy_retreat implicit none integer :: num ! pour la lecture integer :: num_job ! numero du job integer :: lenrun ! longueur du run name integer :: lenb ! longueur du nom de fichier integer :: pos ! position d'un sous_string real*8, dimension(:,:), pointer :: tab !< tableau 2d real ecrit dans le fichier character(len=4) :: char_num ! pour la lecture character(len=80) :: file_B_test !< file name of the bedrock map to test topo instability character(len=80) :: file_voisin !< file name of the bedrock map min or max to test topo instability ! recupere le numero du job dans le runname lenrun = len(trim(runname)) char_num = runname(lenrun-3:lenrun) read(char_num,*) num_job write(6,*) num_job ! Nom du fichier d'experience file_exp = trim(dirnameinp)//trim(file_exp) ! Read the experiment list and keep the line corresponding to the job number open(99,file=file_exp) !read(99,*) ! lit la ligne de titre do i=1,2000 ! each line is a parameter set read(99,*,end=200) num, alpha_lim_1,alpha_lim_2,retreat_1,retreat_2,(time_region(j),j=1,nb_regions),file_B_test if (num.eq.num_job) exit ! stop at the job number end do close(99) ! read the bedrock map on which we test topo instability file_B_test = trim(dirnameinp)//'Socles_Tony/'//file_B_test call Read_Ncdf_var('z',trim(file_B_test),tab) ! lit la variable 'z' B_test(:,:) = tab(:,:) ! lit le fichier region (en general le même que celui des sorties) region_file = trim(dirnameinp)//region_file call Read_Ncdf_var('z',region_file,tab) ! lit la variable 'z' map_region(:,:) = nint(tab(:,:)) write(6,*) 'runname ',runname,' num',num write(6,*) 'alpha_lim_1,alpha_lim_2,retreat1,retreat2', alpha_lim_1,alpha_lim_2,retreat_1,retreat_2 write(6,*) 'time',time_region(:) write(6,*) 'B_test file :', file_B_test !write(6,*) 'B_voisin file :', file_voisin write(6,*) '----------------------------------------------------------' return 200 write(6,*) ' This job number',num,' runname',runname,' is not in the file',file_exp close(99) end subroutine lect_retreat_files !---------------------------------------------------------------------- ! init_retreat_rates : calcule les taux de retrait en fonction de divers criteres ! ici avec une limite une variation lieaire en alpha !---------------------------------------------------------------------- subroutine init_retreat0_alpha use declar_toy_retreat implicit none real :: A_alpha !pour calcul retreat = A_alpha * alpha + B_alpha real :: B_alpha ! A_alpha = (retreat_2 - retreat_1)/(alpha_lim_2 - alpha_lim_1) B_alpha = retreat_1 - A_alpha * alpha_lim_1 ! critere sur le alpha (log10 de beta) !------------------------ travail_centre(:,:) = log10(max(beta_centre(:,:),1.e-5)) call moyennes_demi_mailles where (travail_mx(:,:).gt.alpha_lim_2) !---- en x retreat0_x(:,:) = retreat_2*1000. elsewhere (travail_mx(:,:).lt.alpha_lim_1) retreat0_x(:,:) = retreat_1*1000. elsewhere retreat0_x(:,:) = A_alpha * travail_mx(:,:) + B_alpha end where where (travail_my(:,:).gt.alpha_lim_2) !---- en y retreat0_y(:,:) = retreat_2*1000. elsewhere (travail_my(:,:).lt.alpha_lim_1) retreat0_y(:,:) = retreat_1*1000. elsewhere retreat0_y(:,:) = A_alpha * travail_my(:,:) + B_alpha end where where (travail_SW(:,:).gt.alpha_lim_2) !---- en SW retreat0_SW(:,:) = retreat_2*1000. elsewhere (travail_SW(:,:).lt.alpha_lim_1) retreat0_SW(:,:) = retreat_1*1000. elsewhere retreat0_SW(:,:) = A_alpha * travail_SW(:,:) + B_alpha end where where (travail_SE(:,:).gt.alpha_lim_2) !---- en SE retreat0_SE(:,:) = retreat_2*1000. elsewhere (travail_SE(:,:).lt.alpha_lim_1) retreat0_SE(:,:) = retreat_1*1000. elsewhere retreat0_SE(:,:) = A_alpha * travail_SE(:,:) + B_alpha end where ! enleve les regions qui ne sont pas a traiter do j=2,ny-1 do i =2,nx-1 if (mk_traiter(i,j)*mk_traiter(i-1,j).eq.0) then !---- en x retreat0_x(i,j) = 0. end if if (mk_traiter(i,j)*mk_traiter(i,j-1).eq.0) then !---- en y retreat0_y(i,j) = 0. end if if (mk_traiter(i,j)*mk_traiter(i-1,j-1).eq.0) then !---- en SW retreat0_SW(i,j) = 0. end if if (mk_traiter(i,j)*mk_traiter(i+1,j-1).eq.0) then !---- en SE retreat0_SE(i,j) = 0. end if end do end do end subroutine init_retreat0_alpha !---------------------------------------------------------------------- ! init_retreat_rates : calcule les taux de retrait en fonction de divers criteres ! ici avec une limite abrupte en beta !---------------------------------------------------------------------- subroutine init_retreat0_beta use declar_toy_retreat implicit none ! critere sur le beta !------------------------ travail_centre(:,:) = beta_centre(:,:) call moyennes_demi_mailles where (travail_mx(:,:).gt.beta_lim_ret) !---- en x retreat0_x(:,:) = retreat_2*1000. elsewhere retreat0_x(:,:) = retreat_1*1000. end where where (travail_my(:,:).gt.beta_lim_ret) !---- en y retreat0_y(:,:) = retreat_2*1000. elsewhere retreat0_y(:,:) = retreat_1*1000. end where where (travail_SW(:,:).gt.beta_lim_ret) !---- en SW retreat0_SW(:,:) = retreat_2*1000. elsewhere retreat0_SW(:,:) = retreat_1*1000. end where where (travail_SE(:,:).gt.beta_lim_ret) !---- en SE retreat0_SE(:,:) = retreat_2*1000. elsewhere retreat0_SE(:,:) = retreat_1*1000. end where ! enleve les regions qui ne sont pas a traiter do j=2,ny-1 do i =2,nx-1 if (mk_traiter(i,j)*mk_traiter(i-1,j).eq.0) then !---- en x retreat0_x(i,j) = 0. end if if (mk_traiter(i,j)*mk_traiter(i,j-1).eq.0) then !---- en y retreat0_y(i,j) = 0. end if if (mk_traiter(i,j)*mk_traiter(i-1,j-1).eq.0) then !---- en SW retreat0_SW(i,j) = 0. end if if (mk_traiter(i,j)*mk_traiter(i+1,j-1).eq.0) then !---- en SE retreat0_SE(i,j) = 0. end if end do end do end subroutine init_retreat0_beta !---------------------------------------------------------------------- ! moyennes_demi_mailles : calcule les moyennes d'un tableau ! travail !---------------------------------------------------------------------- subroutine moyennes_demi_mailles use declar_toy_retreat implicit none do j=2,ny-1 do i =2,nx-1 travail_mx(i,j) = 0.5* ( travail_centre(i,j) + travail_centre(i-1,j) ) travail_my(i,j) = 0.5* ( travail_centre(i,j) + travail_centre(i,j-1) ) travail_SW(i,j) = 0.5* ( travail_centre(i,j) +travail_centre(i-1,j-1)) travail_SE(i,j) = 0.5* ( travail_centre(i,j) +travail_centre(i+1,j-1)) end do end do end subroutine moyennes_demi_mailles !---------------------------------------------------------------------- ! max_demi_mailles : calcule les max d'un tableau ! travail !---------------------------------------------------------------------- subroutine max_demi_mailles use declar_toy_retreat implicit none do j=2,ny-1 do i =2,nx-1 travail_mx(i,j) = max ( travail_centre(i,j) , travail_centre(i-1,j) ) travail_my(i,j) = max ( travail_centre(i,j) , travail_centre(i,j-1) ) travail_SW(i,j) = max ( travail_centre(i,j) , travail_centre(i-1,j-1)) travail_SE(i,j) = max ( travail_centre(i,j) , travail_centre(i+1,j-1)) end do end do end subroutine max_demi_mailles !---------------------------------------------------------------------- ! init_time_depart : donne les temps de depart en fonction des regions !---------------------------------------------------------------------- ! init_time_depart : donne les temps de depart en fonction des regions !---------------------------------------------------------------------- subroutine init_time_depart use declar_toy_retreat implicit none real :: time_inacc = 10000. ! temps donne aux noeuds inaccessibles real :: t_already_float = -10. ! temps donne aux noeuds qui flottent deja real :: time_run_begin = 0. ! debut du run real :: time_acc = 1000. ! time for the regions below sea level (to see in graphic output) ! region par region time_region(0) = t_already_float-time_run_begin do j=2,ny-1 do i=2,nx-1 time_dep(i,j) = (time_region(map_region(i,j)) - time_run_begin) end do end do where (mk_traiter(:,:).eq.0) time_dep(:,:) = time_inacc end where ! calcul sur les demi-mailles travail_centre(:,:) = time_dep(:,:) ! call moyennes_demi_mailles call max_demi_mailles time_dep_mx(:,:) = travail_mx(:,:) time_dep_my(:,:) = travail_my(:,:) time_dep_SW(:,:) = travail_SW(:,:) time_dep_SE(:,:) = travail_SE(:,:) ! le temps est tres eleve pour les zones inaccessibles do j=2,ny-1 do i =2,nx-1 if (mk_traiter(i,j)*mk_traiter(i-1,j).eq.0) then !---- en x time_dep_mx(i,j) = time_inacc end if if (mk_traiter(i,j)*mk_traiter(i,j-1).eq.0) then !---- en y time_dep_my(i,j) = time_inacc end if if (mk_traiter(i,j)*mk_traiter(i-1,j-1).eq.0) then !---- en SW time_dep_SW(i,j) = time_inacc end if if (mk_traiter(i,j)*mk_traiter(i+1,j-1).eq.0) then !---- en SE time_dep_SE(i,j) = time_inacc end if end do end do ! pour bien voir les régions dans time_float where(mk_traiter(:,:).eq.2) time_float(:,:) = time_acc elsewhere (mk_traiter(:,:).eq.0) time_float(:,:) = time_inacc elsewhere (mk_traiter(:,:).lt.0) time_float(:,:) = t_already_float end where end subroutine init_time_depart !---------------------------------------------------------------------- !< teste_socle_Btest : pour un point i,j : teste si les voisins sont plus haut !---------------------------------------------------------------------- subroutine teste_socle_Btest(ii,jj,A_x,A_y,A_SW,A_SE) use declar_toy_retreat implicit none integer,intent(in) :: ii,jj ! indice du point considere real,dimension(2),intent(inout) :: A_x ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_y ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SW ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SE ! quand on balaye les voisins ! Btest est plus haut = Bsoc + Bsoc_sigma ! Mais pas suffisamment if (B_test(ii-1,jj ).le.B_test(ii,jj)) A_x(1) = 0. !---- en x if (B_test(ii+1,jj ).le.B_test(ii,jj)) A_x(2) = 0. if (B_test(ii ,jj-1).le.B_test(ii,jj)) A_y(1) = 0. !---- en y if (B_test(ii ,jj+1).le.B_test(ii,jj)) A_y(2) = 0. if (B_test(ii-1,jj-1).le.B_test(ii,jj)) A_SW(1) = 0. !---- en SW if (B_test(ii+1,jj+1).le.B_test(ii,jj)) A_SW(2) = 0. if (B_test(ii+1,jj-1).le.B_test(ii,jj)) A_SE(1) = 0. !---- en SE if (B_test(ii-1,jj+1).le.B_test(ii,jj)) A_SE(2) = 0. end subroutine teste_socle_Btest !---------------------------------------------------------------------- !< teste_socle_Btest : pour un point i,j : teste si les voisins sont plus haut !---------------------------------------------------------------------- subroutine teste_Btest_schoof(ii,jj,A_x,A_y,A_SW,A_SE) use declar_toy_retreat implicit none integer,intent(in) :: ii,jj ! indice du point considere real,dimension(2),intent(inout) :: A_x ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_y ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SW ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SE ! quand on balaye les voisins real :: U2_schoof ! variable de calcul ! test associe sur la pente du socle et le flux de schoof ! le fait de faire le test de schoof sur les vitesses sous-entends ! qu'on fait le flux avec l'epaisseur H_float_test(i,j) ce qui est ! la condiition la plus dure ! si la pente remonte vers le noeud et le flux de Schoof < flux dynamique. Arret. if ((B_test(ii-1,jj ).le.B_test(ii,jj)) .and. & !---- en x (Ux_schoof(ii,jj).lt.abs(Uxbar(ii,jj)))) & ! test schoof A_x(1) = 0. ! noeud mineur_x ii,jj if ((B_test(ii+1,jj ).le.B_test(ii,jj)) .and. & (Ux_schoof(ii,jj).lt.abs(Uxbar(ii+1,jj)))) & ! test schoof A_x(2) = 0. ! noeud mineur_x ii+1 if ((B_test(ii ,jj-1).le.B_test(ii,jj)).and. & !---- en y (Uy_schoof(ii,jj).lt.abs(Uybar(ii,jj)) )) & A_y(1) = 0. ! noeud mineur_y ii,jj if ((B_test(ii ,jj+1).le.B_test(ii,jj)).and. & (Uy_schoof(ii,jj).lt.abs(Uybar(ii,jj+1)) )) & ! test schoof A_y(2) = 0. ! noeud mineur_y ii,jj+1 U2_schoof= Ux_schoof(ii,jj)**2+Uy_schoof(ii,jj)**2 if ((B_test(ii-1,jj-1).le.B_test(ii,jj)) .and. & !---- en SW (U2_schoof.lt.(Uxbar(ii,jj)**2+Uybar(ii,jj)**2))) & ! test Schoof A_SW(1) = 0. ! noeuds mineur_SW (ii,jj) if ((B_test(ii+1,jj+1).le.B_test(ii,jj)).and. & !---- en NE (U2_schoof.lt.(Uxbar(ii+1,jj)**2+Uybar(ii,jj+1)**2))) & ! test Schoof A_SW(2) = 0. ! noeuds mineurs_SW ii+1,j et ii,jj+1 if ((B_test(ii+1,jj-1).le.B_test(ii,jj)).and. & !---- en SE (U2_schoof.lt.(Uxbar(ii+1,jj)**2+Uybar(ii,jj)**2))) & ! test Schoof A_SE(1) = 0. ! noeuds mineurs_SE ii+1,j et ii,jj if ((B_test(ii-1,jj+1).le.B_test(ii,jj)).and. & !--- en NW (U2_schoof.lt.(Uxbar(ii,jj)**2+Uybar(ii,jj+1)**2))) & ! test Schoof A_SE(2) = 0. ! noeuds mineurs_SE ii,j et ii,jj+1 end subroutine teste_Btest_schoof !---------------------------------------------------------------------- !---------------------------------------------------------------------- !< teste_socle_Btest : pour un point i,j : teste si les voisins sont plus haut !---------------------------------------------------------------------- subroutine teste_socle_Bminmax(ii,jj,A_x,A_y,A_SW,A_SE) use declar_toy_retreat implicit none integer,intent(in) :: ii,jj ! indice du point considere real,dimension(2),intent(inout) :: A_x ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_y ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SW ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SE ! quand on balaye les voisins ! B_voisin est le socle test du point voisin (par ex. le max de la grille 15) ! B_test est le socle test du point considere (par ex. le min de la grille 15) ! Btest est plus haut = Bsoc + Bsoc_sigma ! Mais pas suffisamment if (B_voisin(ii-1,jj ).le.B_test(ii,jj)) A_x(1) = 0. !---- en x if (B_voisin(ii+1,jj ).le.B_test(ii,jj)) A_x(2) = 0. if (B_voisin(ii ,jj-1).le.B_test(ii,jj)) A_y(1) = 0. !---- en y if (B_voisin(ii ,jj+1).le.B_test(ii,jj)) A_y(2) = 0. if (B_voisin(ii-1,jj-1).le.B_test(ii,jj)) A_SW(1) = 0. !---- en SW if (B_voisin(ii+1,jj+1).le.B_test(ii,jj)) A_SW(2) = 0. if (B_voisin(ii+1,jj-1).le.B_test(ii,jj)) A_SE(1) = 0. !---- en SE if (B_voisin(ii-1,jj+1).le.B_test(ii,jj)) A_SE(2) = 0. end subroutine teste_socle_Bminmax !---------------------------------------------------------------------- ! idem previous subroutine but the test is done on Bsoc ! only for tests !< teste_socle_Btest : pour un point i,j : teste si les voisins sont plus haut !---------------------------------------------------------------------- subroutine teste_socle_Bsoc(ii,jj,A_x,A_y,A_SW,A_SE) use declar_toy_retreat implicit none integer,intent(in) :: ii,jj ! indice du point considere real,dimension(2),intent(inout) :: A_x ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_y ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SW ! quand on balaye les voisins real,dimension(2),intent(inout) :: A_SE ! quand on balaye les voisins ! Btest est plus haut = Bsoc + Bsoc_sigma ! Mais pas suffisamment if (B_test(ii-1,jj ).le.Bsoc(ii,jj)) A_x(1) = 0. !---- en x if (B_test(ii+1,jj ).le.Bsoc(ii,jj)) A_x(2) = 0. if (B_test(ii ,jj-1).le.Bsoc(ii,jj)) A_y(1) = 0. !---- en y if (B_test(ii ,jj+1).le.Bsoc(ii,jj)) A_y(2) = 0. if (B_test(ii-1,jj-1).le.Bsoc(ii,jj)) A_SW(1) = 0. !---- en SW if (B_test(ii+1,jj+1).le.Bsoc(ii,jj)) A_SW(2) = 0. if (B_test(ii+1,jj-1).le.Bsoc(ii,jj)) A_SE(1) = 0. !---- en SE if (B_test(ii-1,jj+1).le.Bsoc(ii,jj)) A_SE(2) = 0. end subroutine teste_socle_Bsoc !---------------------------------------------------------------------- !< sanity_check : pour un point i,j limite le dhdt !---------------------------------------------------------------------- subroutine sanity_check(ii,jj) use declar_toy_retreat implicit none integer,intent(in) :: ii,jj ! indice du point considere real :: dH_x ! variable de travail real :: dH_y ! variable de travail DelHdt_sanity(ii,jj) = 0. ! en x ----------------------------------------------------------------------------------- if ((Uxbar(ii,jj).ge.0.).and.(Uxbar(ii+1,jj).ge.0.)) then ! ecoulement ---> dH_x = Uxbar(ii,jj)*(H(ii,jj)-H(ii-1,jj)) / dx dH_x = Dh_x + epsmax(ii,jj) * H(ii,jj) else if ((Uxbar(ii,jj).le.0.).and.(Uxbar(ii+1,jj).le.0.)) then ! ecoulement <--- dH_x = Uxbar(ii+1,jj)*(H(ii+1,jj)-H(ii,jj)) / dx dH_x = Dh_x + epsmax(ii,jj) * H(ii,jj) else if ((Uxbar(ii,jj).le.0.).and.(Uxbar(ii+1,jj).ge.0.)) then ! ecoulement <------> dH_x = Dh_x + epsmax(ii,jj) * H(ii,jj) else ! ecoulement ---><----- dH_x = Uxbar(ii,jj)*(H(ii,jj)-H(ii-1,jj)) / dx + Uxbar(ii+1,jj)*(H(ii+1,jj)-H(ii,jj)) / dx ! mais ici pas de terme en epsmax car confine endif ! en y ----------------------------------------------------------------------------------- if ((Uybar(ii,jj).ge.0.).and.(Uybar(ii,jj+1).ge.0.)) then ! ecoulement ^ dH_y = Uybar(ii,jj)*(H(ii,jj)-H(ii,jj-1)) / dx dH_y = Dh_y + epsmax(ii,jj) * H(ii,jj) else if ((Uybar(ii,jj).le.0.).and.(Uybar(ii,jj+1).le.0.)) then ! ecoulement v dH_y = Uybar(ii,jj+1)*(H(ii,jj+1)-H(ii,jj)) / dx dH_y = Dh_y + epsmax(ii,jj) * H(ii,jj) else if ((Uybar(ii,jj).le.0.).and.(Uybar(ii,jj+1).ge.0.)) then ! ecoulement divergent dH_y = Dh_y + epsmax(ii,jj) * H(ii,jj) else ! ecoulement convergent dH_y = Uybar(ii,jj)*(H(ii,jj)-H(ii,jj-1)) / dx + Uybar(ii,jj+1)*(H(ii,jj+1)-H(ii,jj)) / dx ! mais ici pas de terme en epsmax car confine endif DelHdt_sanity(ii,jj) = max(dH_x + dH_y,0.) end subroutine sanity_check !----------------------------------------------------------------------------------- ! no_sanity_check pour tester l'influence du sanity check !----------------------------------------------------------------------------------- subroutine no_sanity_check(ii,jj) use declar_toy_retreat implicit none integer,intent(in) :: ii,jj ! indice du point considere real :: dH_x ! variable de travail real :: dH_y ! variable de travail DelHdt_sanity(ii,jj) = 5000. end subroutine no_sanity_check !-------------------------------------------------------------------- ! Schoof_flux : calcule le coefficient et le flux de Schoof !-------------------------------------------------------------------- subroutine init_Schoof_flux use declar_toy_retreat implicit none nschoof = 3 ! loi de Glen mschoof = 1./alpha_drag ! exposant loi de friction = alpha_drag m_scoof = 1/alpha_drag m1_schoof = 1/(mschoof+1) exp_schoof = (nschoof+3+mschoof)* m1_schoof ! exposant de l'epaisseur cst_schoof = ( (ro*g)**(nschoof+1) * (1+coef_Bflot)**nschoof / 4**nschoof ) **(m1_schoof) ! formulation no lin end subroutine Init_Schoof_flux !-------------------------------------------------------------------- ! Schoof_flux : calcule le coefficient et le flux de Schoof !-------------------------------------------------------------------- subroutine Schoof_flux use declar_toy_retreat implicit none where (H_float_test(:,:).gt.0.) ! epaisseur de flottaison avec le socle test ! coef_schoof(:,:) = (Abar(:,:)/max(1.,beta_centre(:,:)))**m1_schoof ! coef_schoof(:,:) = min(1000.,beta_centre(:,:)) ! coef_schoof(:,:) = beta_centre(:,:) ! seulement vrai pour dragging lineaire coef_schoof(:,:) = coef_drag(:,:) ! coef_drag calcule dans init_dragging coef_schoof(:,:) = (Abar(:,:)/max(1.,coef_schoof(:,:)))**m1_schoof coef_schoof(:,:) = coef_schoof(:,:) * cst_schoof * H_float_test(:,:)**exp_schoof tauf_schoof(:,:) = 0.5 * ro * g * H_float_test(:,:) * (1+coef_Bflot) elsewhere coef_schoof(:,:) = 0. tauf_schoof(:,:) = 0. end where debug_3D(:,:,105) = coef_drag(:,:) ! calcul du terme de confinement. A revoir !!$do j=2,ny-1 !!$ do i=2,nx-1 !!$ if (H_float_test(i,j).gt.0.) then !!$ !!$ ! 2 lignes en dessous, H et pas H_float_test parce que c'est pour tirer la viscosite moyenne !!$ if (H(i,j).gt.0.) then !!$ bf_x_schoof(i,j) = 2.*(Uxbar(i+1,j)-Uxbar(i,j)) /dx *pvi(i,j) /H(i,j) ! tau'xx, !!$ bf_y_schoof(i,j) = 2.*(Uybar(i,j+1)-Uybar(i,j)) /dx *pvi(i,j) /H(i,j) ! tau'yy !!$ end if !!$ !!$ if (tauf_schoof(i,j).gt.0.) then !!$ bf_x_schoof(i,j) = abs(bf_x_schoof(i,j))/tauf_schoof(i,j) ! rapport tauxx/tauf !!$ bf_y_schoof(i,j) = abs(bf_y_schoof(i,j))/tauf_schoof(i,j) ! rapport tauyy/tauf !!$ else !!$ bf_x_schoof(i,j) = 1. !!$ bf_y_schoof(i,j) = 1. !!$ end if !!$ !!$ if (bf_x_schoof(i,j).gt.0) then ! a la puissance nschoof*m1_schoof !!$ bf_x_schoof(i,j) = bf_x_schoof(i,j) **(nschoof*m1_schoof) !!$ else !!$ bf_x_schoof(i,j) = 0. !!$ end if !!$ !!$ if (bf_y_schoof(i,j).gt.0) then ! a la puissance nschoof*m1_schoof !!$ bf_y_schoof(i,j) = bf_y_schoof(i,j) **(nschoof*m1_schoof) !!$ else !!$ bf_y_schoof(i,j) = 0. !!$ end if !!$ else ! pas de flottaison !!$ bf_x_schoof(i,j) = 0. !!$ bf_y_schoof(i,j) = 0. !!$ end if !!$ !!$ end do !!$end do bf_x_schoof(:,:) = 1. bf_y_schoof(:,:) = 1. where ((H_float_test(:,:).gt.0.).and.(H(:,:).gt.2.)) Ux_schoof(:,:) = bf_x_schoof(:,:)*coef_schoof(:,:)/H_float_test(:,:) Uy_schoof(:,:) = bf_y_schoof(:,:)*coef_schoof(:,:)/H_float_test(:,:) elsewhere Ux_schoof(:,:) = 0. Uy_schoof(:,:) = 0. end where end subroutine Schoof_flux end module toy_retreat_mod !----------------------------------------------------------------------------------