source: branches/iLoveclim/SOURCES/Ant16_files/bmelt-ant-regions-initmip_mod.f90 @ 146

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

Series of minor bug corrections for initMIP Antarctica configuration

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