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