source: branches/iLoveclim/SOURCES/Ant40_files/bmelt-ant-regions_mod.f90 @ 52

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

Merge branche iLOVECLIM sur rev 51

File size: 11.2 KB
Line 
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!<
14MODULE  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  real :: bmelt_coef                    !< Coef pour corriger les valeurs std de bmelt (lu dans namelist)
40
41CONTAINS
42  !-------------------------------------------------------------------------------
43  subroutine init_bmelt
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,bmelt_coef     
54
55
56    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
57    read(num_param,bmelt_ant_reg)
58
59    ! formats pour les ecritures dans 42
60428 format(A)
61
62    write(num_rep_42,428)'!___________________________________________________________' 
63    write(num_rep_42,428) '&bmelt-ant-regions                                         '
64    write(num_rep_42,*)
65    write(num_rep_42,*) 'bmelt_Ross      = ',bmelt_Ross 
66    write(num_rep_42,*) 'bmgrz_Ross      = ',bmgrz_Ross
67    write(num_rep_42,*) 'bmelt_FRis      = ',bmelt_FRis 
68    write(num_rep_42,*) 'bmgrz_FRis      = ',bmgrz_FRis
69    write(num_rep_42,*) 'bmelt_Amery     = ',bmelt_Amery 
70    write(num_rep_42,*) 'bmgrz_Amery     = ',bmgrz_Amery
71    write(num_rep_42,*) 'bmelt_PIG       = ',bmelt_PIG 
72    write(num_rep_42,*) 'bmgrz_PIG       = ',bmgrz_PIG
73    write(num_rep_42,*) 'bmelt_Pen       = ',bmelt_Pen
74    write(num_rep_42,*) 'bmgrz_Pen       = ',bmgrz_Pen
75    write(num_rep_42,*) 'bmelt_other     = ',bmelt_other 
76    write(num_rep_42,*) 'bmgrz_other     = ',bmgrz_other
77    write(num_rep_42,*) 'bmelt_talus     = ',bmelt_talus 
78    write(num_rep_42,*) 'bmgrz_talus     = ',bmgrz_talus
79    write(num_rep_42,*) 'bmelt_coef      = ',bmelt_coef
80    write(num_rep_42,*)'/'                           
81    write(num_rep_42,428) '! bmelt_Ross & bmgrz_Ross    :  fonte basale secteur Ross'
82    write(num_rep_42,428) '! bmelt_FRis & bmgrz_FRis    :  fonte basale secteur Filchner Ronne'
83    write(num_rep_42,428) '! bmelt_Amery & bmgrz_Amery  :  fonte basale secteur Amery'
84    write(num_rep_42,428) '! bmelt_PIG & bmgrz_PIG      :  fonte basale secteur Pine Island'
85    write(num_rep_42,428) '! bmelt_Pen & bmgrz_Pen      :  fonte basale secteur Peninsule'
86    write(num_rep_42,428) '! bmelt_other & bmgrz_other  :  fonte basale autres secteurs'
87    write(num_rep_42,428) '! bmelt_talus & bmgrz_talus  :  fonte basale apres talus cont'
88    write(num_rep_42,428) '! bmelt_coef                 :  coef fonte (1 pour conserver val)'
89    write(num_rep_42,*)
90
91
92    !    ecriture dans le fichier parametres
93    write(num_rep_42,*)'fusion basale sous les ice shelves : bmelt-ant-regions_mod'
94    write(num_rep_42,*)'----------------------------------------------------------'
95
96
97    ! lecture du fichier contenant les distances au talu continental (m)
98    open(88,file=TRIM(DIRNAMEINP)//'distance_talu-40km.xy')
99
100    do j=1,ny
101       do i=1,nx
102          read(88,'(i3,1x,i3,1x,f10.2)') k,k,dist_talu(i,j)
103       enddo
104    enddo
105    close(88)
106
107    ! les bords de la grille sont mis a dist_talus=0 de force
108
109
110    dist_talu(1,:)=0.   !a gauche
111    dist_talu(2,:)=0.
112
113    dist_talu(nx-1,:)=0. !  a droite
114    dist_talu(nx,:)=0. 
115
116    dist_talu(:,1)=0.    ! en bas   
117    dist_talu(:,2)=0.
118
119    dist_talu(:,ny-1)=0.    ! en haut
120    dist_talu(:,ny)=0.
121
122
123    debug_3D(:,:,32)=dist_talu(:,:)
124
125    ! lecture du fichier contenant les types de shelves Ronne-Flichner
126    !  ->1, Ross -> 2 , Amery -> 3, PIG-> 4, les petits shelves au dessus de PIG -> 5
127
128    open(88,file=TRIM(DIRNAMEINP)//'numer-ice-shelves-juil07.dat',status='OLD')
129    typeshelf(:,:)=100
130    do k=1,nx*ny
131       read(88,*,end=36) i,j,typeshelf(i,j)
132    end do
13336  close(88)
134
135
136
137    region(:)=0
138    regname(1)='Ronne-Fichner'
139    regname(2)='Ross'
140    regname(3)='Amery'
141    regname(4)='Pig'
142    regname(5)='Petits ice shelves peninsule'
143    regname(6)='autres ice shelves'
144    regname(7)='Au dela du talus continental'
145
146    bms_init:    do j=1,ny
147       do i=1,nx
148
149          talus:    if (dist_talu(i,j).GT.5000.) then
150             if (nint(typeshelf(i,j)).eq.1) then ! Ronne-Filchner FRis
151                bmshelf(i,j)=bmelt_FRis         
152                bmgrz(i,j)=bmgrz_FRis         
153
154
155                if (region(1).eq.0) then          ! pour l'impression dans 42
156                   region(1)=1
157                   write(num_rep_42,80)regname(1),bmshelf(i,j),bmgrz(i,j)
158                endif
159
16080              format(a32,' bmshelf=',f8.4,'  bmgrz=',f8.4)
161
162             else  if (nint(typeshelf(i,j)).eq.2) then ! Ross
163                bmshelf(i,j)=bmelt_Ross 
164                bmgrz(i,j)=bmgrz_Ross     
165
166                if (region(2).eq.0) then
167                   region(2)=1
168                   write(num_rep_42,80)regname(2),bmshelf(i,j),bmgrz(i,j)
169                endif
170
171
172             else  if (nint(typeshelf(i,j)).eq.3) then ! Amery
173                bmshelf(i,j)=bmelt_Amery   
174                bmgrz(i,j)=bmgrz_Amery         
175                if (region(3).eq.0) then
176                   region(3)=1
177                   write(num_rep_42,80)regname(3),bmshelf(i,j),bmgrz(i,j)
178                endif
179
180
181             else  if (nint(typeshelf(i,j)).eq.4) then ! Pig
182                bmshelf(i,j)=bmelt_PIG 
183                bmgrz(i,j)=bmgrz_PIG       
184
185                if (region(4).eq.0) then
186                   region(4)=1
187                   write(num_rep_42,80)regname(4),bmshelf(i,j),bmgrz(i,j)
188                endif
189
190             else  if (nint(typeshelf(i,j)).eq.5) then ! petits shelves peninsule
191                bmshelf(i,j)=bmelt_Pen     
192                bmgrz(i,j)=bmgrz_Pen       
193
194                if (region(5).eq.0) then
195                   region(5)=1
196                   write(num_rep_42,80)regname(5),bmshelf(i,j),bmgrz(i,j)
197                endif
198
199
200             else !
201                typeshelf(i,j)=6
202                bmshelf(i,j)=bmelt_other     
203                bmgrz(i,j)=bmgrz_other       
204
205                if (region(6).eq.0) then
206                   region(6)=1
207                   write(num_rep_42,80)regname(6),bmshelf(i,j),bmgrz(i,j)
208                endif
209
210
211             endif
212
213
214          else   ! au dela du talus continental
215
216             bmshelf(i,j)=bmelt_talus 
217             bmgrz(i,j)=bmgrz_talus   
218
219             if (region(7).eq.0) then
220                region(7)=1
221                write(num_rep_42,80)regname(7),bmshelf(i,j),bmgrz(i,j)
222             endif
223
224
225          endif talus   ! fin du test sur dist_talus
226
227       enddo
228    enddo bms_init
229
230    bmshelf(:,:)=bmshelf(:,:)*bmelt_coef
231    bmgrz(:,:)=bmgrz(:,:)*bmelt_coef
232    debug_3D(:,:,33)=typeshelf(:,:)
233    debug_3D(:,:,34)=bmshelf(:,:)
234    return
235  end subroutine init_bmelt
236
237  !________________________________________________________________________________
238  subroutine bmeltshelf
239
240
241    ! cette routine calcule la fusion basale proprement dite
242
243    integer :: ngr           ! nombre de voisins flottants
244    real :: coef_talus       ! pour ne pas changer la fusion au dessus de l'ocean profond
245
246    do j=2,ny-1
247       do i=2,nx-1
248
249          talus_nochange : if (typeshelf(i,j).lt.10) then
250             coef_talus = coefbmshelf
251          else
252             coef_talus = 1         ! pas de changement au dela du talus continental
253          endif talus_nochange
254
255
256          shelf:    if (flot(i,j)) then    ! partie flottante
257
258             bmelt(i,j)=coef_talus*bmshelf(i,j)
259
260             ! bloc bsupshelf à enlever
261             !                  if (time.gt.5000.) then
262             !                     bmelt(i,j)=bmelt(i,j)+bsupshelf
263             !                  endif   
264
265
266
267             if (fbm(i,j)) then
268                bmelt(i,j)=coef_talus*bmgrz(i,j)
269             endif
270
271             !                  if (time.gt.5000.) then
272             !                     bmelt(i,j)=bmelt(i,j)+bsupshelf
273             !                  endif   
274
275             ! ATTENTION le bloc suivant sert a determiner la fusion basale d'equilibre
276             ! pour les shelves stationnaires
277
278             if (igrdline.eq.1) then
279                corrbmelt(i,j)=hdot(i,j)+bmelt(i,j)   ! le bmelt d'equilibre
280                debug_3D(i,j,28)=corrbmelt(i,j)
281             endif
282
283
284          else                   ! point posé, on compte le nombre de voisins flottants
285             ngr=0
286             if (flot(i+1,j)) ngr=ngr+1
287             if (flot(i-1,j)) ngr=ngr+1
288             if (flot(i,j+1)) ngr=ngr+1
289             if (flot(i,j-1)) ngr=ngr+1
290
291             !   la fusion des points limites est une combinaison entre valeur posée et valeur flottante
292             !   en fonction du nombre de points flottants
293
294             bmelt(i,j)= ngr/4.*bmgrz(i,j)*coef_talus+(1.-ngr/4.)*bmelt(i,j)
295
296
297          endif shelf
298
299       end do
300    end do
301
302    ! les bords de la grille sont mis a dist_talus=0 de force
303    ! attention, pourrait poser des problemes avec AGRIF
304
305    bmelt(1,:)=bmelt_talus   !a gauche
306    bmelt(2,:)=bmelt_talus
307
308    bmelt(nx-1,:)=bmelt_talus !  a droite
309    bmelt(nx,:)=bmelt_talus 
310
311    bmelt(:,1)=bmelt_talus    ! en bas   
312    bmelt(:,2)=bmelt_talus
313
314    bmelt(:,ny-1)=bmelt_talus    ! en haut
315    bmelt(:,ny)=bmelt_talus
316
317
318
319    if (igrdline.eq.1) then
320       bmelt(:,:)=0.      ! hdot donne alors le -bmelt
321    endif
322
323    ! ecriture des bmelt pour en faire des statistiques
324    !--------------------------------------------------
325    !open(144,file='bmelt-test')
326    !do j=1,ny
327    !   do i=1,nx
328    !      if (fbm(i,j)) then
329    !         grd=1
330    !      else
331    !         grd=0
332    !      endif
333    !      if ((flot(i,j)).and.(H(i,j).gt.1.)) then
334    !         write(144,999) i,j,grd,nint(typeshelf(i,j)), &
335    !              bmelt(i,j),bmgrz(i,j),bmshelf(i,j),hdot(i,j),dist_talu(i,j)
336    !      endif
337    !   end do
338    !end do
339    !999 format(4(i0,2x),5(e12.3,2x))
340    !close(144)
341
342    return
343  end subroutine  bmeltshelf
344
345
346
347END MODULE  BMELT_ANT_REGIONS
Note: See TracBrowser for help on using the repository browser.