source: branches/iLoveclim/SOURCES/bmelt-seuil-profondeur_mod.f90

Last change on this file was 244, checked in by aquiquet, 5 years ago

Grisli-iloveclim branch merged to trunk at revision 243

File size: 5.9 KB
Line 
1!> \file bmelt-seuil-profondeur_mod.f90
2!! Prametrise la fusion basale (ice shelves)
3!<
4
5!> \namespace bmelt_seuil_prof
6!! Prametrise la fusion basale (ice shelves)
7!! \author ...
8!! \date ...
9!! @note Pour l'actuel : 2 valeurs  pour 2 domaines de profondeur
10!! + une valeur pour bmgrz
11!! A choisir dans le module_choix
12!! @note Used module
13!! @note   - use module3D_phy
14!<
15
16module  bmelt_seuil_prof         
17
18  ! prametrise la fusion basale (ice shelves)
19  ! Pour l'actuel : 2 valeurs  pour 2 domaines de profondeur
20  ! + une valeur pour bmgrz
21  ! A choisir dans le module_choix
22
23  use module3d_phy
24
25
26  implicit none
27
28  real :: bm_grz                        !< valeur prescrite a la grounding zone
29  real,dimension(nx,ny) ::  bmgrz       !< tabelau fusion basale a la grounding zone
30
31  real :: bmshelf_plateau               !< valeur prescrite sur le plateau cont.
32  real :: bmshelf_abysses               !< valeur prescrite au dessus de l'ocean profond
33  real :: depth_talus                   !< profondeur de transition
34  real,dimension(nx,ny) ::  bmshelf     !< tableau fusion basale sous shelf
35
36
37contains
38  !-------------------------------------------------------------------------------
39  !> SUBROUTINE: init_bmelt
40  !! Cette routine fait l'initialisation pour la fusion basale.
41  !! @note Elle est appelée par inputfile-vec-0.5.f90
42  !!
43  !>
44  subroutine init_bmelt
45
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
52    namelist/bmelt_seuil/bm_grz,bmshelf_plateau,bmshelf_abysses,depth_talus
53
54    if (itracebug.eq.1)  call tracebug('entree dans init_bmelt de bmelt_seuil_prof')
55
56    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
57    read(num_param,bmelt_seuil)
58!    write(num_rep_42,bmelt_seuil)
59
60    !    ecriture dans le fichier parametres
61
62    write(num_rep_42,'(A)')'!___________________________________________________________' 
63    write(num_rep_42,'(A)') '&bmelt_seuil                     ! module  bmelt_seuil_prof'
64    write(num_rep_42,*)   'bm_grz           =',bm_grz
65    write(num_rep_42,*)   'bmshelf_plateau  =',bmshelf_plateau
66    write(num_rep_42,*)   'bmshelf_abysses  =',bmshelf_abysses
67    write(num_rep_42,*)   'depth_talus      =',depth_talus
68    write(num_rep_42,*)'/'                     
69    write(num_rep_42,'(A)') '! Pour l actuel : bm_grz a la grounding line'
70    write(num_rep_42,'(A)') '!                 bmshelf_plateau sur le plateau continental'
71    write(num_rep_42,'(A)') '!                 bmshelf_abysses pour les grandes profondeurs'
72    write(num_rep_42,'(A)') '!       depth_talus, negative, separation entre les 2 domaines'
73    write(num_rep_42,*)
74    write(num_rep_42,'(A)')'!___________________________________________________________' 
75
76
77
78    bmgrz(:,:) = bm_grz
79
80    where (Bsoc0(:,:).lt.depth_talus)
81        bmshelf(:,:)=bmshelf_abysses
82    elsewhere
83        bmshelf(:,:)=bmshelf_plateau
84    end where
85
86    debug_3D(:,:,34)=bmshelf(:,:)
87    return
88  end subroutine init_bmelt
89
90  !________________________________________________________________________________
91
92  !> SUBROUTINE: bmeltshelf
93  !! Cette routine calcule la fusion basale proprement dite pour les points flottants
94  !! @note coefbmshelf a ete calcule par forclim et permet de faire varier en fonction du climat
95  !>
96
97  subroutine bmeltshelf
98
99
100    ! cette routine calcule la fusion basale proprement dite pour les points flottants
101    ! coefbmshelf a ete calcule par forclim et permet de faire varier en fonction du climat
102
103    integer :: ngr           ! nombre de voisins flottants
104    real    :: coef_talus    ! pour ne pas changer la fusion au dessus de l'ocean profond
105
106    if (itracebug.eq.1)  call tracebug('entree dans bmeltshelf de bmelt_seuil_prof')
107
108    where (Bsoc0(:,:).lt.depth_talus)
109        bmshelf(:,:)=bmshelf_abysses
110    elsewhere
111        bmshelf(:,:)=bmshelf_plateau
112    end where
113
114    do j=1,ny
115      do i=1,nx
116
117        talus_nochange: if (Bsoc0(i,j).lt.depth_talus) then
118            coef_talus = 1         ! pas de changement au dessus de l'ocecan profond
119        else
120            coef_talus = coefbmshelf
121        endif talus_nochange
122
123
124        shelf: if (flot(i,j)) then    ! partie flottante
125
126            bmelt(i,j)=coef_talus*bmshelf(i,j)
127
128            ! fbm est vrai si le point est flottant mais un des voisins est pose
129            ! calcule dans flottab
130
131            if (fbm(i,j)) then
132                bmelt(i,j)=coef_talus*bmgrz(i,j)
133            endif
134
135
136
137            ! ATTENTION le bloc suivant sert a determiner la fusion basale d'equilibre
138            ! pour les shelves stationnaires
139
140            if ((igrdline.eq.1).and.(ibmelt_inv.eq.1)) then
141                !     corrbmelt(i,j)=hdot(i,j)+bmelt(i,j)   ! le bmelt d'equilibre
142                !     bmelt(i,j)=corrbmelt(i,j)
143                corrbmelt(i,j)=corrbmelt(i,j)+hdot(i,j)*0.85
144                bmelt(i,j)=bmelt(i,j)+corrbmelt(i,j)     
145            endif
146
147
148        else                   ! point posé, on compte le nombre de voisins flottants
149            ngr=0
150            if ((i.ne.1).and.(i.ne.nx).and.(j.ne.1).and.(j.ne.ny)) then
151                if (flot(i+1,j)) ngr=ngr+1
152                if (flot(i-1,j)) ngr=ngr+1
153                if (flot(i,j+1)) ngr=ngr+1
154                if (flot(i,j-1)) ngr=ngr+1
155            end if
156
157
158            !   la fusion des points limites est une combinaison entre valeur posée et valeur flottante
159            !   en fonction du nombre de points flottants
160
161            bmelt(i,j)= ngr/4.*bmgrz(i,j)*coef_talus+(1.-ngr/4.)*bmelt(i,j)
162
163            if ((igrdline.eq.1).and.(ibmelt_inv.eq.1)) then
164                corrbmelt(i,j)=corrbmelt(i,j)+hdot(i,j)  !*0.85
165                bmelt(i,j)=bmelt(i,j)+corrbmelt(i,j)
166            endif
167        endif shelf
168      end do
169    end do
170
171if ((igrdline.eq.1).and.(ibmelt_inv.eq.0)) then
172    where (i_Hp(:,:).eq.1)
173        bmelt(:,:)=0.    ! hdot donne alors le -bmelt
174    endwhere
175endif
176
177return
178end subroutine  bmeltshelf
179
180
181
182end module  bmelt_seuil_prof
Note: See TracBrowser for help on using the repository browser.