source: branches/GRISLIv3/SOURCES/Ant40_files/bmelt-ant-regions_mod.f90 @ 446

Last change on this file since 446 was 446, checked in by aquiquet, 8 months ago

Cleaning branch: use only statement in module3D

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