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

Last change on this file since 465 was 465, checked in by aquiquet, 5 months ago

Cleaning branch: start cleaning module3D

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