source: trunk/SOURCES/Ant16_files/bmelt-ant-regions-initmip_mod.f90 @ 173

Last change on this file since 173 was 173, checked in by aquiquet, 7 years ago

Basal melt rates depend on distance from limit of the continental shelf (tof&afq)

File size: 9.0 KB
Line 
1!> \file bmelt-ant-regions-initmip_mod.f90
2!! Module qui calcule la fusion basale (grounded ou ice shelves).
3!<
4
5!> \namespace BMELT_ANT_REGIONS_INITMIP
6!! Module qui calcule la fusion basale (grounded ou ice shelves)
7!! \author CatRitz, aquiquet
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_INITMIP         ! cat juillet 2005 - afq 2017
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) :: typeshelf    !< Type de shelf Ronne->1 Ross ->2 ....
26
27  integer, parameter :: nb_reg = 18
28
29  real,dimension (nb_reg) :: bmelt_regions, bmgrz_regions
30  real :: bmelt_talus,bmgrz_talus       !< Au dela du talus continental
31  real :: bmelt_coef                    !< Coef pour corriger les valeurs std de bmelt (lu dans namelist)
32  real, parameter :: depth_talus=-2500. !< Profondeur du talus continental
33
34  ! For initMIP, we need to read an anomaly of bmelt...:
35  integer               :: bmelt_time        !< pour selectionner le type de calcul de bmelt 
36  real,dimension(nx,ny) :: bmelt_anom        !< bilan de masse des anomalies indices : nx,ny
37  real,dimension(nx,ny) ::  bmgrz_0               !< fusion basale a la grounding zone
38  real,dimension(nx,ny) ::  bmshelf_0             !< fusion basale sous shelf
39  integer :: flag_dist
40  real,dimension(nx,ny) :: coef_bmelt_dist
41 
42CONTAINS
43  !-------------------------------------------------------------------------------
44  subroutine init_bmelt
45
46
47    ! Cette routine fait l'initialisation pour la fusion basale.
48    ! Elle est appelée par inputfile-vec-0.5.f90
49
50    real*8, dimension(:,:),   pointer    :: tab               !< tableau 2d real pointer
51    real :: bmelt_dist0
52    character(len=100) :: file_number_shelves ! fichier avec les zones ice-shelves
53    character(len=100) :: file_bmelt_anom
54    real,dimension(nx,ny) :: dist_talus  ! distance au talus
55    character(len=100) :: file_dist_talus ! fichier distance au talus
56   
57    namelist/bmelt_ant_reg_initmip/bmelt_regions,bmgrz_regions,  & 
58         bmelt_talus,bmgrz_talus,bmelt_coef,file_number_shelves, &
59         flag_dist,file_dist_talus,bmelt_dist0 
60    namelist/bmelt_anom_initMIP/file_bmelt_anom, bmelt_time 
61
62
63    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
64    read(num_param,bmelt_ant_reg_initmip)
65
66    ! formats pour les ecritures dans 42
67428 format(A)
68
69    write(num_rep_42,428)'!___________________________________________________________' 
70    write(num_rep_42,428) '&bmelt-ant-regions-initmip                                 '
71    write(num_rep_42,*)
72    write(num_rep_42,'(A,18(f0.2,","))') 'bmelt_regions      = ', bmelt_regions(:)
73    write(num_rep_42,'(A,18(f0.2,","))') 'bmgrz_regions      = ', bmgrz_regions(:)
74    write(num_rep_42,*) 'bmelt_talus     = ', bmelt_talus 
75    write(num_rep_42,*) 'bmgrz_talus     = ', bmgrz_talus
76    write(num_rep_42,*) 'bmelt_coef      = ', bmelt_coef
77    write(num_rep_42,*) 'file_numer_shelves = ', file_number_shelves
78    write(num_rep_42,*) 'flag_dist = ', flag_dist
79    write(num_rep_42,*) 'file_dist_talus = ', file_dist_talus
80    write(num_rep_42,*) 'bmelt_dist0 = ', bmelt_dist0
81    write(num_rep_42,*)'/'                           
82    write(num_rep_42,428) '! bmelt_regions    :  fonte basale sous shelves pour les 18 regions initMIP'
83    write(num_rep_42,428) '! bmgrz_regions    :  fonte basale grounding zone pour les 18 regions initMIP'
84    write(num_rep_42,428) '! bmelt_talus & bmgrz_talus  :  fonte basale apres talus cont'
85    write(num_rep_42,428) '! bmelt_coef                 :  coef fonte (1 pour conserver val)'
86    write(num_rep_42,428) '! file_numer_ice-shelves     : fichier zones ice shelves'
87    write(num_rep_42,428) '! flag_dist : flag pour bmelt fnct distance talus'
88    write(num_rep_42,428) '! file_dist_talus : fichier de distance talus'
89    write(num_rep_42,428) '! bmelt_dist0 : coef de bmelt au talus'
90    write(num_rep_42,*)
91
92
93    ! lecture du fichier contenant les types de shelves
94    ! 18 regions sur la base des bassins fournis par le projet ISMIP6-initmip antarctica
95    typeshelf(:,:)=100
96    file_number_shelves=TRIM(DIRNAMEINP)//trim(file_number_shelves)
97    call Read_Ncdf_var('z',file_number_shelves,tab)
98    typeshelf(:,:)  = tab(:,:)
99
100
101    if (flag_dist.eq.1) then
102       file_dist_talus=TRIM(DIRNAMEINP)//trim(file_dist_talus)
103       call Read_Ncdf_var('z',file_dist_talus,tab)
104       dist_talus(:,:)  = min(max(tab(:,:),0.),400.)
105       coef_bmelt_dist(:,:)= 1. / ((1/bmelt_dist0)*(1+(bmelt_dist0-1)*(dist_talus(:,:)/400.)**2))
106    else
107       coef_bmelt_dist(:,:)= 1.
108    endif
109
110    bms_init:    do j=1,ny
111       do i=1,nx
112
113          talus:    if (Bsoc(i,j).GT.depth_talus) then
114
115             bmshelf(i,j) = bmelt_regions(typeshelf(i,j))
116             bmgrz(i,j)   = bmgrz_regions(typeshelf(i,j))
117
118          else   ! au dela du talus continental
119             
120             bmshelf(i,j)=bmelt_talus 
121             bmgrz(i,j)=bmgrz_talus   
122
123          endif talus   ! fin du test sur dist_talus
124
125       enddo
126    enddo bms_init
127   
128    bmshelf(:,:)=bmshelf(:,:)*bmelt_coef
129    bmgrz(:,:)=bmgrz(:,:)*bmelt_coef
130
131    bmshelf_0(:,:)=bmshelf(:,:)
132    bmgrz_0(:,:)=bmgrz(:,:)
133
134
135  ! ______ Anomalies...
136 
137  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
138  read(num_param,bmelt_anom_initMIP)
139
140  write(num_rep_42,428)'!_______________________________________________________________________' 
141  write(num_rep_42,428)'!  module bmelt_ant_regions_initmip_mod                                 '
142  write(num_rep_42,bmelt_anom_initMIP)
143  write(num_rep_42,428)'! file_bmelt_anom   = fichier anomalie bmelt                            '
144  write(num_rep_42,428)'! bmelt_time        = 0:fixe, 1:anomalies                               '
145  write(num_rep_42,428)'!_______________________________________________________________________' 
146
147  if ( bmelt_time .eq. 1 ) then 
148     file_bmelt_anom  = trim(dirnameinp)//trim(file_bmelt_anom)
149     call Read_Ncdf_var('abmb',file_bmelt_anom,tab)
150     bmelt_anom (:,:) = Tab(:,:)
151  end if
152 
153   
154    debug_3D(:,:,33)=typeshelf(:,:)
155    debug_3D(:,:,34)=bmshelf(:,:)
156
157   
158    return
159  end subroutine init_bmelt
160
161  !________________________________________________________________________________
162  subroutine bmeltshelf
163
164
165    ! cette routine calcule la fusion basale proprement dite
166
167    integer :: ngr           ! nombre de voisins flottants
168    real :: coef_talus       ! pour ne pas changer la fusion au dessus de l'ocean profond
169
170    real             :: coefanomtime
171
172    coefanomtime = min ( real(time/40.) , 1. ) 
173    if (bmelt_time.eq.0) then
174       bmshelf(:,:) = bmshelf_0(:,:)
175       bmgrz(:,:) = bmgrz_0(:,:)
176    else if (bmelt_time.eq.1) then
177       bmshelf(:,:) = bmshelf_0(:,:)+ (coefanomtime * bmelt_anom(:,:))
178       bmgrz(:,:) = bmgrz_0(:,:)    + (coefanomtime * bmelt_anom(:,:))
179    end if
180
181   
182    do j=2,ny-1
183       do i=2,nx-1
184
185          talus_nochange : if (Bsoc(i,j).GT.depth_talus) then
186             coef_talus = coefbmshelf
187          else
188             coef_talus = 1         ! pas de changement au dela du talus continental               
189          endif talus_nochange
190
191
192          shelf:    if (flot(i,j)) then    ! partie flottante
193
194             bmelt(i,j)=coef_talus*bmshelf(i,j)*coef_bmelt_dist(i,j)
195
196             if (fbm(i,j)) then
197                bmelt(i,j)=coef_talus*bmgrz(i,j)*coef_bmelt_dist(i,j)
198             endif
199
200             ! ATTENTION le bloc suivant sert a determiner la fusion basale d'equilibre
201             ! pour les shelves stationnaires
202
203             if (igrdline.eq.1) then
204                corrbmelt(i,j)=hdot(i,j)+bmelt(i,j)   ! le bmelt d'equilibre
205                debug_3D(i,j,28)=corrbmelt(i,j)
206             endif
207
208
209          else                   ! point posé, on compte le nombre de voisins flottants
210             ngr=0
211             if (flot(i+1,j)) ngr=ngr+1
212             if (flot(i-1,j)) ngr=ngr+1
213             if (flot(i,j+1)) ngr=ngr+1
214             if (flot(i,j-1)) ngr=ngr+1
215
216             !   la fusion des points limites est une combinaison entre valeur posée et valeur flottante
217             !   en fonction du nombre de points flottants
218
219             bmelt(i,j)= (ngr/4.*bmgrz(i,j)*coef_talus)*coef_bmelt_dist(i,j)+(1.-ngr/4.)*bmelt(i,j)
220
221
222          endif shelf
223
224       end do
225    end do
226
227
228    ! les bords de la grille sont mis a dist_talus=0 de force
229    ! attention, pourrait poser des problemes avec AGRIF
230
231    bmelt(1,:)=bmelt_talus   !a gauche
232    bmelt(2,:)=bmelt_talus
233
234    bmelt(nx-1,:)=bmelt_talus !  a droite
235    bmelt(nx,:)=bmelt_talus 
236
237    bmelt(:,1)=bmelt_talus    ! en bas   
238    bmelt(:,2)=bmelt_talus
239
240    bmelt(:,ny-1)=bmelt_talus    ! en haut
241    bmelt(:,ny)=bmelt_talus
242
243
244    if (igrdline.eq.1) then
245       bmelt(:,:)=0.      ! hdot donne alors le -bmelt
246    endif
247
248    return
249  end subroutine  bmeltshelf
250
251
252END MODULE  BMELT_ANT_REGIONS_INITMIP
Note: See TracBrowser for help on using the repository browser.