1 | !> \file bmelt-ant-regions_mod.f90 |
---|
2 | !! Module qui calcule la fusion basale (grounded ou ice shelves). |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace BMELT_ANT_REGIONS |
---|
6 | !! Module qui calcule la fusion basale (grounded ou ice shelves) |
---|
7 | !! \author CatRitz |
---|
8 | !! \date juillet |
---|
9 | !! @note pour les ice shelves Antarctique, tient compte des différentes régions |
---|
10 | !! A choisir dans le module_choix |
---|
11 | !! @note Used module |
---|
12 | !! @note - use module3D_phy |
---|
13 | !< |
---|
14 | MODULE BMELT_ANT_REGIONS ! cat juillet 2005 |
---|
15 | |
---|
16 | USE module3D_phy |
---|
17 | |
---|
18 | |
---|
19 | implicit none |
---|
20 | |
---|
21 | REAL,dimension(nx,ny) :: bmgrz !< fusion basale a la grounding zone |
---|
22 | real,dimension(nx,ny) :: bmshelf !< fusion basale sous shelf |
---|
23 | |
---|
24 | REAL,dimension(nx,ny) :: dist_talu !< distance du point au talu continental |
---|
25 | REAL,dimension(nx,ny) :: typeshelf !< Type de shelf Ronne->1 Ross ->2 .... |
---|
26 | !< utilises pour moduler la fusion sous le shelf |
---|
27 | integer, dimension(10) :: region !< pour écrire dans le fichier param |
---|
28 | character(len=30),dimension(10) :: regname !< nom des régions |
---|
29 | real :: bsupshelf |
---|
30 | integer :: grd !< pour une sortie |
---|
31 | |
---|
32 | real :: bmelt_Ross,bmgrz_Ross !< fusion basale Ross |
---|
33 | real :: bmelt_FRis,bmgrz_FRis !< fusion basale filchner-Ronne ice shelf |
---|
34 | real :: bmelt_Amery,bmgrz_Amery !< Amery ice shelves |
---|
35 | real :: bmelt_PIG,bmgrz_PIG !< Pine Island glacier |
---|
36 | real :: bmelt_Pen,bmgrz_Pen !< Petits ice shelves Peninsule |
---|
37 | real :: bmelt_other,bmgrz_other !< Autre ice shelves |
---|
38 | real :: bmelt_talus,bmgrz_talus !< Au dela du talus continental |
---|
39 | |
---|
40 | CONTAINS |
---|
41 | !------------------------------------------------------------------------------- |
---|
42 | subroutine init_bmelt |
---|
43 | |
---|
44 | |
---|
45 | |
---|
46 | ! Cette routine fait l'initialisation pour la fusion basale. |
---|
47 | ! Elle est appelée par inputfile-vec-0.5.f90 |
---|
48 | |
---|
49 | |
---|
50 | namelist/bmelt_ant_reg/bmelt_Ross,bmgrz_Ross,bmelt_FRis,bmgrz_FRis, & |
---|
51 | bmelt_Amery,bmgrz_Amery,bmelt_PIG,bmgrz_PIG, & |
---|
52 | bmelt_Pen,bmgrz_Pen,bmelt_other,bmgrz_other, & |
---|
53 | bmelt_talus,bmgrz_talus |
---|
54 | |
---|
55 | |
---|
56 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
57 | read(num_param,bmelt_ant_reg) |
---|
58 | |
---|
59 | |
---|
60 | ! ecriture dans le fichier parametres |
---|
61 | write(num_rep_42,*)'fusion basale sous les ice shelves : bmelt-ant-regions_mod' |
---|
62 | write(num_rep_42,*)'----------------------------------------------------------' |
---|
63 | |
---|
64 | ! lecture du fichier contenant les distances au talu continental (m) |
---|
65 | open(88,file=TRIM(DIRNAMEINP)//'distance_talu-40km.xy') |
---|
66 | |
---|
67 | do j=1,ny |
---|
68 | do i=1,nx |
---|
69 | read(88,'(i3,1x,i3,1x,f10.2)') k,k,dist_talu(i,j) |
---|
70 | enddo |
---|
71 | enddo |
---|
72 | close(88) |
---|
73 | |
---|
74 | ! les bords de la grille sont mis a dist_talus=0 de force |
---|
75 | |
---|
76 | |
---|
77 | dist_talu(1,:)=0. !a gauche |
---|
78 | dist_talu(2,:)=0. |
---|
79 | |
---|
80 | dist_talu(nx-1,:)=0. ! a droite |
---|
81 | dist_talu(nx,:)=0. |
---|
82 | |
---|
83 | dist_talu(:,1)=0. ! en bas |
---|
84 | dist_talu(:,2)=0. |
---|
85 | |
---|
86 | dist_talu(:,ny-1)=0. ! en haut |
---|
87 | dist_talu(:,ny)=0. |
---|
88 | |
---|
89 | |
---|
90 | debug_3D(:,:,32)=dist_talu(:,:) |
---|
91 | |
---|
92 | ! lecture du fichier contenant les types de shelves |
---|
93 | ! Ronne-Flichner ->1, Ross -> 2 , Amery -> 3, |
---|
94 | ! PIG-> 4, les petits shelves au dessus de PIG -> 5 |
---|
95 | |
---|
96 | open(88,file=TRIM(DIRNAMEINP)//'numer-ice-shelves-juil07.dat',status='OLD') |
---|
97 | typeshelf(:,:)=100 |
---|
98 | do k=1,nx*ny |
---|
99 | read(88,*,end=36) i,j,typeshelf(i,j) |
---|
100 | end do |
---|
101 | 36 close(88) |
---|
102 | |
---|
103 | |
---|
104 | |
---|
105 | region(:)=0 |
---|
106 | regname(1)='Ronne-Fichner' |
---|
107 | regname(2)='Ross' |
---|
108 | regname(3)='Amery' |
---|
109 | regname(4)='Pig' |
---|
110 | regname(5)='Petits ice shelves peninsule' |
---|
111 | regname(6)='autres ice shelves' |
---|
112 | regname(7)='Au dela du talus continental' |
---|
113 | |
---|
114 | bms_init: do j=1,ny |
---|
115 | do i=1,nx |
---|
116 | |
---|
117 | talus: if (dist_talu(i,j).GT.5000.) then |
---|
118 | if (nint(typeshelf(i,j)).eq.1) then ! Ronne-Filchner FRis |
---|
119 | bmshelf(i,j)=bmelt_FRis |
---|
120 | bmgrz(i,j)=bmgrz_FRis |
---|
121 | |
---|
122 | |
---|
123 | if (region(1).eq.0) then ! pour l'impression dans 42 |
---|
124 | region(1)=1 |
---|
125 | write(num_rep_42,80)regname(1),bmshelf(i,j),bmgrz(i,j) |
---|
126 | endif |
---|
127 | |
---|
128 | 80 format(a32,' bmshelf=',f8.4,' bmgrz=',f8.4) |
---|
129 | |
---|
130 | else if (nint(typeshelf(i,j)).eq.2) then ! Ross |
---|
131 | bmshelf(i,j)=bmelt_Ross |
---|
132 | bmgrz(i,j)=bmgrz_Ross |
---|
133 | |
---|
134 | if (region(2).eq.0) then |
---|
135 | region(2)=1 |
---|
136 | write(num_rep_42,80)regname(2),bmshelf(i,j),bmgrz(i,j) |
---|
137 | endif |
---|
138 | |
---|
139 | |
---|
140 | else if (nint(typeshelf(i,j)).eq.3) then ! Amery |
---|
141 | bmshelf(i,j)=bmelt_Amery |
---|
142 | bmgrz(i,j)=bmgrz_Amery |
---|
143 | if (region(3).eq.0) then |
---|
144 | region(3)=1 |
---|
145 | write(num_rep_42,80)regname(3),bmshelf(i,j),bmgrz(i,j) |
---|
146 | endif |
---|
147 | |
---|
148 | |
---|
149 | else if (nint(typeshelf(i,j)).eq.4) then ! Pig |
---|
150 | bmshelf(i,j)=bmelt_PIG |
---|
151 | bmgrz(i,j)=bmgrz_PIG |
---|
152 | |
---|
153 | if (region(4).eq.0) then |
---|
154 | region(4)=1 |
---|
155 | write(num_rep_42,80)regname(4),bmshelf(i,j),bmgrz(i,j) |
---|
156 | endif |
---|
157 | |
---|
158 | else if (nint(typeshelf(i,j)).eq.5) then ! petits shelves peninsule |
---|
159 | bmshelf(i,j)=bmelt_Pen |
---|
160 | bmgrz(i,j)=bmgrz_Pen |
---|
161 | |
---|
162 | if (region(5).eq.0) then |
---|
163 | region(5)=1 |
---|
164 | write(num_rep_42,80)regname(5),bmshelf(i,j),bmgrz(i,j) |
---|
165 | endif |
---|
166 | |
---|
167 | |
---|
168 | else ! |
---|
169 | typeshelf(i,j)=6 |
---|
170 | bmshelf(i,j)=bmelt_other |
---|
171 | bmgrz(i,j)=bmgrz_other |
---|
172 | |
---|
173 | if (region(6).eq.0) then |
---|
174 | region(6)=1 |
---|
175 | write(num_rep_42,80)regname(6),bmshelf(i,j),bmgrz(i,j) |
---|
176 | endif |
---|
177 | |
---|
178 | |
---|
179 | endif |
---|
180 | |
---|
181 | |
---|
182 | else ! au dela du talus continental |
---|
183 | |
---|
184 | bmshelf(i,j)=bmelt_talus |
---|
185 | bmgrz(i,j)=bmgrz_talus |
---|
186 | |
---|
187 | if (region(7).eq.0) then |
---|
188 | region(7)=1 |
---|
189 | write(num_rep_42,80)regname(7),bmshelf(i,j),bmgrz(i,j) |
---|
190 | endif |
---|
191 | |
---|
192 | |
---|
193 | endif talus ! fin du test sur dist_talus |
---|
194 | |
---|
195 | enddo |
---|
196 | enddo bms_init |
---|
197 | |
---|
198 | debug_3D(:,:,33)=typeshelf(:,:) |
---|
199 | debug_3D(:,:,34)=bmshelf(:,:) |
---|
200 | return |
---|
201 | end subroutine init_bmelt |
---|
202 | |
---|
203 | !________________________________________________________________________________ |
---|
204 | subroutine bmeltshelf |
---|
205 | |
---|
206 | |
---|
207 | ! cette routine calcule la fusion basale proprement dite |
---|
208 | |
---|
209 | integer :: ngr ! nombre de voisins flottants |
---|
210 | real :: coef_talus ! pour ne pas changer la fusion au dessus de l'ocean profond |
---|
211 | |
---|
212 | do j=2,ny-1 |
---|
213 | do i=2,nx-1 |
---|
214 | |
---|
215 | talus_nochange : if (typeshelf(i,j).lt.10) then |
---|
216 | coef_talus = coefbmshelf |
---|
217 | else |
---|
218 | coef_talus = 1 ! pas de changement au dela du talus continental |
---|
219 | endif talus_nochange |
---|
220 | |
---|
221 | |
---|
222 | shelf: if (flot(i,j)) then ! partie flottante |
---|
223 | |
---|
224 | bmelt(i,j)=coef_talus*bmshelf(i,j) |
---|
225 | |
---|
226 | ! bloc bsupshelf à enlever |
---|
227 | ! if (time.gt.5000.) then |
---|
228 | ! bmelt(i,j)=bmelt(i,j)+bsupshelf |
---|
229 | ! endif |
---|
230 | |
---|
231 | |
---|
232 | |
---|
233 | if (fbm(i,j)) then |
---|
234 | bmelt(i,j)=coef_talus*bmgrz(i,j) |
---|
235 | endif |
---|
236 | |
---|
237 | ! if (time.gt.5000.) then |
---|
238 | ! bmelt(i,j)=bmelt(i,j)+bsupshelf |
---|
239 | ! endif |
---|
240 | |
---|
241 | ! ATTENTION le bloc suivant sert a determiner la fusion basale d'equilibre |
---|
242 | ! pour les shelves stationnaires |
---|
243 | |
---|
244 | if (igrdline.eq.1) then |
---|
245 | corrbmelt(i,j)=hdot(i,j)+bmelt(i,j) ! le bmelt d'equilibre |
---|
246 | debug_3D(i,j,28)=corrbmelt(i,j) |
---|
247 | endif |
---|
248 | |
---|
249 | |
---|
250 | else ! point posé, on compte le nombre de voisins flottants |
---|
251 | ngr=0 |
---|
252 | if (flot(i+1,j)) ngr=ngr+1 |
---|
253 | if (flot(i-1,j)) ngr=ngr+1 |
---|
254 | if (flot(i,j+1)) ngr=ngr+1 |
---|
255 | if (flot(i,j-1)) ngr=ngr+1 |
---|
256 | |
---|
257 | ! la fusion des points limites est une combinaison entre valeur posée et valeur flottante |
---|
258 | ! en fonction du nombre de points flottants |
---|
259 | |
---|
260 | bmelt(i,j)= ngr/4.*bmgrz(i,j)*coef_talus+(1.-ngr/4.)*bmelt(i,j) |
---|
261 | |
---|
262 | |
---|
263 | endif shelf |
---|
264 | |
---|
265 | end do |
---|
266 | end do |
---|
267 | |
---|
268 | ! les bords de la grille sont mis a dist_talus=0 de force |
---|
269 | ! attention, pourrait poser des problemes avec AGRIF |
---|
270 | |
---|
271 | bmelt(1,:)=bmelt_talus !a gauche |
---|
272 | bmelt(2,:)=bmelt_talus |
---|
273 | |
---|
274 | bmelt(nx-1,:)=bmelt_talus ! a droite |
---|
275 | bmelt(nx,:)=bmelt_talus |
---|
276 | |
---|
277 | bmelt(:,1)=bmelt_talus ! en bas |
---|
278 | bmelt(:,2)=bmelt_talus |
---|
279 | |
---|
280 | bmelt(:,ny-1)=bmelt_talus ! en haut |
---|
281 | bmelt(:,ny)=bmelt_talus |
---|
282 | |
---|
283 | |
---|
284 | |
---|
285 | if (igrdline.eq.1) then |
---|
286 | bmelt(:,:)=0. ! hdot donne alors le -bmelt |
---|
287 | endif |
---|
288 | |
---|
289 | ! ecriture des bmelt pour en faire des statistiques |
---|
290 | !-------------------------------------------------- |
---|
291 | !open(144,file='bmelt-test') |
---|
292 | !do j=1,ny |
---|
293 | ! do i=1,nx |
---|
294 | ! if (fbm(i,j)) then |
---|
295 | ! grd=1 |
---|
296 | ! else |
---|
297 | ! grd=0 |
---|
298 | ! endif |
---|
299 | ! if ((flot(i,j)).and.(H(i,j).gt.1.)) then |
---|
300 | ! write(144,999) i,j,grd,nint(typeshelf(i,j)), & |
---|
301 | ! bmelt(i,j),bmgrz(i,j),bmshelf(i,j),hdot(i,j),dist_talu(i,j) |
---|
302 | ! endif |
---|
303 | ! end do |
---|
304 | !end do |
---|
305 | !999 format(4(i0,2x),5(e12.3,2x)) |
---|
306 | !close(144) |
---|
307 | |
---|
308 | return |
---|
309 | end subroutine bmeltshelf |
---|
310 | |
---|
311 | |
---|
312 | |
---|
313 | END MODULE BMELT_ANT_REGIONS |
---|