source: trunk/SOURCES/bmelt-seuil-profondeur_mod.f90 @ 179

Last change on this file since 179 was 124, checked in by dumas, 7 years ago

New flag ibmelt_inv in param_list.dat : 0 standard case, 1 bmelt inversion (with igrdline=1) | bmelt-seuil-profondeur has been updated for bmelt inversion.

File size: 5.9 KB
RevLine 
[4]1!> \file bmelt-seuil-profondeur_mod.f90
2!! Prametrise la fusion basale (ice shelves)
3!<
4
5!> \namespace bmelt_seuil_prof
6!! Prametrise la fusion basale (ice shelves)
7!! \author ...
8!! \date ...
9!! @note Pour l'actuel : 2 valeurs  pour 2 domaines de profondeur
10!! + une valeur pour bmgrz
11!! A choisir dans le module_choix
12!! @note Used module
13!! @note   - use module3D_phy
14!<
15
16module  bmelt_seuil_prof         
17
[124]18  ! prametrise la fusion basale (ice shelves)
19  ! Pour l'actuel : 2 valeurs  pour 2 domaines de profondeur
20  ! + une valeur pour bmgrz
21  ! A choisir dans le module_choix
[4]22
[124]23  use module3d_phy
[4]24
25
[124]26  implicit none
[4]27
[124]28  real :: bm_grz                        !< valeur prescrite a la grounding zone
29  real,dimension(nx,ny) ::  bmgrz       !< tabelau fusion basale a la grounding zone
[4]30
[124]31  real :: bmshelf_plateau               !< valeur prescrite sur le plateau cont.
32  real :: bmshelf_abysses               !< valeur prescrite au dessus de l'ocean profond
33  real :: depth_talus                   !< profondeur de transition
34  real,dimension(nx,ny) ::  bmshelf     !< tableau fusion basale sous shelf
[4]35
36
37contains
[124]38  !-------------------------------------------------------------------------------
39  !> SUBROUTINE: init_bmelt
40  !! Cette routine fait l'initialisation pour la fusion basale.
41  !! @note Elle est appelée par inputfile-vec-0.5.f90
42  !!
43  !>
44  subroutine init_bmelt
[4]45
46
47
[124]48    ! Cette routine fait l'initialisation pour la fusion basale.
49    ! Elle est appelée par inputfile-vec-0.5.f90
[4]50
51
[124]52    namelist/bmelt_seuil/bm_grz,bmshelf_plateau,bmshelf_abysses,depth_talus
[4]53
[124]54    if (itracebug.eq.1)  call tracebug('entree dans init_bmelt de bmelt_seuil_prof')
[4]55
[124]56    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
57    read(num_param,bmelt_seuil)
58    write(num_rep_42,bmelt_seuil)
[4]59
[124]60    !    ecriture dans le fichier parametres
[4]61
62428 format(A)
63
[124]64    write(num_rep_42,428)'!___________________________________________________________' 
65    write(num_rep_42,428) '&bmelt_seuil                     ! module  bmelt_seuil_prof'
66    write(num_rep_42,*)   'bm_grz           =',bm_grz
67    write(num_rep_42,*)   'bmshelf_plateau  =',bmshelf_plateau
68    write(num_rep_42,*)   'bmshelf_abysses  =',bmshelf_abysses
69    write(num_rep_42,*)   'depth_talus      =',depth_talus
70    write(num_rep_42,*)'/'                     
71    write(num_rep_42,428) '! Pour l actuel : bm_grz a la grounding line'
72    write(num_rep_42,428) '!                 bmshelf_plateau sur le plateau continental'
73    write(num_rep_42,428) '!                 bmshelf_abysses pour les grandes profondeurs'
74    write(num_rep_42,428) '!       depth_talus, negative, separation entre les 2 domaines'
75    write(num_rep_42,*)
76    write(num_rep_42,428)'!___________________________________________________________' 
[4]77
78
79
[124]80    bmgrz(:,:) = bm_grz
[4]81
[124]82    where (Bsoc0(:,:).lt.depth_talus)
83        bmshelf(:,:)=bmshelf_abysses
84    elsewhere
85        bmshelf(:,:)=bmshelf_plateau
86    end where
[4]87
[124]88    debug_3D(:,:,34)=bmshelf(:,:)
89    return
90  end subroutine init_bmelt
[4]91
[124]92  !________________________________________________________________________________
[4]93
[124]94  !> SUBROUTINE: bmeltshelf
95  !! Cette routine calcule la fusion basale proprement dite pour les points flottants
96  !! @note coefbmshelf a ete calcule par forclim et permet de faire varier en fonction du climat
97  !>
[4]98
[124]99  subroutine bmeltshelf
[4]100
101
[124]102    ! cette routine calcule la fusion basale proprement dite pour les points flottants
103    ! coefbmshelf a ete calcule par forclim et permet de faire varier en fonction du climat
[4]104
[124]105    integer :: ngr           ! nombre de voisins flottants
106    real    :: coef_talus    ! pour ne pas changer la fusion au dessus de l'ocean profond
[4]107
[124]108    if (itracebug.eq.1)  call tracebug('entree dans bmeltshelf de bmelt_seuil_prof')
[4]109
[124]110    where (Bsoc0(:,:).lt.depth_talus)
111        bmshelf(:,:)=bmshelf_abysses
112    elsewhere
113        bmshelf(:,:)=bmshelf_plateau
114    end where
[4]115
[124]116    do j=1,ny
117      do i=1,nx
[4]118
[124]119        talus_nochange: if (Bsoc0(i,j).lt.depth_talus) then
120            coef_talus = 1         ! pas de changement au dessus de l'ocecan profond
[67]121        else
[124]122            coef_talus = coefbmshelf
[4]123        endif talus_nochange
124
125
[124]126        shelf: if (flot(i,j)) then    ! partie flottante
[4]127
128            bmelt(i,j)=coef_talus*bmshelf(i,j)
129
[124]130            ! fbm est vrai si le point est flottant mais un des voisins est pose
131            ! calcule dans flottab
[4]132
133            if (fbm(i,j)) then
[124]134                bmelt(i,j)=coef_talus*bmgrz(i,j)
135            endif
[4]136
137
138
[124]139            ! ATTENTION le bloc suivant sert a determiner la fusion basale d'equilibre
140            ! pour les shelves stationnaires
[4]141
[124]142            if ((igrdline.eq.1).and.(ibmelt_inv.eq.1)) then
143                !     corrbmelt(i,j)=hdot(i,j)+bmelt(i,j)   ! le bmelt d'equilibre
144                !     bmelt(i,j)=corrbmelt(i,j)
145                corrbmelt(i,j)=corrbmelt(i,j)+hdot(i,j)*0.85
146                bmelt(i,j)=bmelt(i,j)+corrbmelt(i,j)     
147            endif
[4]148
149
150        else                   ! point posé, on compte le nombre de voisins flottants
[124]151            ngr=0
152            if ((i.ne.1).and.(i.ne.nx).and.(j.ne.1).and.(j.ne.ny)) then
153                if (flot(i+1,j)) ngr=ngr+1
154                if (flot(i-1,j)) ngr=ngr+1
155                if (flot(i,j+1)) ngr=ngr+1
156                if (flot(i,j-1)) ngr=ngr+1
157            end if
[4]158
159
[124]160            !   la fusion des points limites est une combinaison entre valeur posée et valeur flottante
161            !   en fonction du nombre de points flottants
[4]162
[124]163            bmelt(i,j)= ngr/4.*bmgrz(i,j)*coef_talus+(1.-ngr/4.)*bmelt(i,j)
164
165            if ((igrdline.eq.1).and.(ibmelt_inv.eq.1)) then
166                corrbmelt(i,j)=corrbmelt(i,j)+hdot(i,j)  !*0.85
167                bmelt(i,j)=bmelt(i,j)+corrbmelt(i,j)
168            endif
[4]169        endif shelf
[124]170      end do
[4]171    end do
172
[124]173if ((igrdline.eq.1).and.(ibmelt_inv.eq.0)) then
174    where (i_Hp(:,:).eq.1)
175        bmelt(:,:)=0.    ! hdot donne alors le -bmelt
176    endwhere
[4]177endif
178
179return
180end subroutine  bmeltshelf
181
182
183
184end module  bmelt_seuil_prof
Note: See TracBrowser for help on using the repository browser.