source: trunk/SOURCES/Ant16_files/bmelt-ant-regions-larmip_mod.f90 @ 334

Last change on this file since 334 was 183, checked in by aquiquet, 6 years ago

Additions for LARMIP experiments

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