source: trunk/SOURCES/Ant40_files/bmelt-ant-regions_mod.f90 @ 70

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

Bug fix in Qprod_icetemp (instability in ice temperature) : slowssa points heat = deformation SIA + sliding. Use fleuvemx and fleuvemy to identifie streams SSA and slow SSA points

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