source: branches/iLoveclim/SOURCES/bmelt-seuil-profondeur_mod.f90 @ 77

Last change on this file since 77 was 77, checked in by dumas, 8 years ago

Merge branche iLOVECLIM sur rev 76

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