source: branches/GRISLIv3/SOURCES/bmelt-seuil-profondeur_mod.f90 @ 381

Last change on this file since 381 was 381, checked in by dumas, 16 months ago

use only in subroutine bmelt_seuil_prof

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