source: trunk/SOURCES/prescribe-H-i2s_mod.f90 @ 10

Last change on this file since 10 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

File size: 9.3 KB
Line 
1!> \file prescribe-H-i2s_mod.f90
2!!
3!!  Determine mask and thickness for prescription of H
4!!  include ice2sea prescribed retreat  (Cat april 2012)
5!<
6
7!> \namespace prescribe_H
8!!
9!! Determine mask and thickness for prescription of H
10!! 
11!! \author CatRitz
12!! \date june 2009
13!!
14!! @note use  module3D_phy
15!! @note use param_phy_mod
16!<
17
18module prescribe_H
19
20  use module3D_phy
21  use param_phy_mod
22  use interface_input
23
24! pour les reculs de type ice2sea
25  use toy_retreat_mod
26
27  implicit none
28!  real,dimension(nx,ny)    :: Hp           !< H value if prescribed
29!  real,dimension(nx,ny)    :: Hp0          !< H value if prescribed (reference value)
30!  real,dimension(nx,ny)    :: Delta_H      !< Delta_H value if prescribed
31!  integer,dimension(nx,ny) :: i_delta_H    !< 1 if Delta_H is prescribed on this node, else 0
32!  integer,dimension(nx,ny) :: i_Hp         !< 1 if H is prescribed on this node, else 0
33!  integer,dimension(nx,ny) :: i_Hp0        !< i_hp mask reference value does not change with time
34!  integer,dimension(nx,ny) :: MK_gl0       !< mask grounding line initial
35!  integer,dimension(nx,ny) :: MK_flot0     !< mask float initial
36
37! pour grounding line retreat, ice2sea
38! nouvelle version interactive pour imposer des contraintes variees.
39! sanity check, temps de depart en fonction des bassins, ....
40
41
42! nouvelle version recul grounding line pour ice2sea
43! ------------------------------------------------------
44
45
46
47contains
48
49  !______________________________________________________________________________________
50  !> initialise prescribe H
51  !! calculate the initial grounding line position
52  !! @return MK_flot0 and Mk_gl0
53
54  subroutine init_prescribe_H
55
56    implicit none
57    integer    :: voisin   !< pour test sur des masques
58
59    if (itracebug.eq.1)  call tracebug(' Entree dans routine init_prescribe_H')
60
61    ! keep the mask of floating points
62    where (flot(:,:))
63       MK_flot0(:,:)= 1
64    elsewhere 
65       MK_flot0(:,:)= 0
66    end where
67
68
69    ! determine the initial grounding line
70
71    MK_gl0(:,:)=0
72    HP0(:,:)=H0(:,:)
73
74
75    do j=2,ny-1
76       do i=2,nx-1
77          voisin = MK_flot0(i-1,j) + MK_flot0(i+1,j) + MK_flot0(i,j-1) + MK_flot0(i,j+1)
78          if ((MK_flot0(i,j).eq.0).and.(voisin.gt.0))  then  ! grounded and at least one neighbour floating
79             MK_gl0(i,j) = 1                                ! tagged mask grounding line
80
81          end if
82       end do
83    end do
84
85
86    call prescribe_grid_boundary_0
87
88
89! Dans le cas du Groenland ice2sea, mk_init = 3 -> noeuds qui ne bougent pas
90
91    if (geoplace(1:4).eq.'GI2S') then
92
93       do j = 1,ny
94          do i = 1,nx
95             if (mk_init(i,j).eq.3) then
96                i_hp0 (i,j) = 1
97                Hp0  (i,j)  = H0(i,j)
98             endif
99          end do
100       end do
101    end if
102
103
104    if (igrdline.eq.3) then
105      where (flot(:,:))
106         mk_gr(:,:) = 0
107      elsewhere
108         mk_gr(:,:) = 1
109      end where
110      call init_proto_recul
111   end if
112
113  end subroutine init_prescribe_H
114  !______________________________________________________________________________________
115  !>  function prescribe_present_H_gl
116  !! calculate the initial grounding line position
117  !! @return i_hp(:,:) and hp(:,:)
118
119  subroutine prescribe_present_H_gl
120
121    implicit none
122
123    if (itracebug.eq.1)  call tracebug(' Entree dans routine  prescribe_present_H_gl')
124
125    where ((MK_flot0(:,:).eq. 1).or.(MK_gl0(:,:) .eq. 1)) ! floating or grounding line
126       i_hp(:,:) = 1                                      ! thickness prescribed to present value
127       hp(:,:)   = Hp0(:,:)
128    end where
129if (itracebug.eq.1)  call tracebug(' fin prescribe_present_H_gl ')
130
131
132  end subroutine prescribe_present_H_gl
133
134  !______________________________________________________________________________________
135  !>  function prescribe_fixed_points
136  !!  fixes  thickness for some points
137  !! @return i_hp(:,:) and hp(:,:)
138
139  subroutine prescribe_fixed_points
140
141if (itracebug.eq.1)  call tracebug(' Entree dans routine  prescribe_fixed_points')
142    where (i_hp0(:,:).eq.1)             ! les points i_hp0 le sont pour tout le run
143       i_hp(:,:) = i_hp0(:,:)
144       Hp(:,:)   = Hp0(:,:)
145    end where
146
147if (itracebug.eq.1)  call tracebug(' fin prescribe_fixed_points')
148
149
150  end subroutine prescribe_fixed_points
151
152  !______________________________________________________________________________________
153  !>  function prescribe_paleo_gl_shelf
154  !! calculate the  grounding line thickness
155  !! no ice shelves
156  !! @return i_hp(:,:) and hp(:,:)
157
158  subroutine prescribe_paleo_gl_shelf
159
160    implicit none
161
162    if (itracebug.eq.1)  call tracebug(' Entree dans routine prescribe_paleo_gl_shelf')
163
164
165!   noeuds qui doivent être imposés
166
167    where ((MK_flot0(:,:).eq. 1).or.(MK_gl0(:,:) .eq. 1)) ! floating or grounding line
168       i_hp(:,:) = 1                                      ! thickness prescribed
169    end where
170
171! valeur imposee
172    where (MK_flot0(:,:).eq. 1)                          ! paleo shelf epaisseur a 1
173       hp(:,:)   = 1.
174    end where
175
176
177    where (MK_gl0(:,:) .eq. 1)
178       hp(:,:) = (sealevel-Bsoc(:,:))*row/ro +1.
179    end where
180
181
182  end subroutine prescribe_paleo_gl_shelf
183
184  !______________________________________________________________________________________
185  !>  function prescribe_grid_boundary
186  !! prescribe the grid boundary thickness to 0
187  !! @return i_hp(:,:) and hp(:,:)
188
189  subroutine prescribe_grid_boundary_0
190
191    implicit none
192    if (itracebug.eq.1)  call tracebug(' Entree dans routine prescribe_grid_boundary_0')
193
194    ! prescribe thickness on the grid boundaries
195    i_hp(:,:) = 0
196
197    i_hp(1,:) = 1
198    i_hp(nx,:)= 1
199    hp(1,:)   = 0
200    hp(nx,:)  = 0
201
202    i_hp(:,1) = 1
203    i_hp(:,ny)= 1
204    hp(:,1)   = 0
205    hp(:,ny)  = 0
206
207! valeurs de reference
208    i_hp0(1,:)  = 1
209    i_hp0(nx,:) = 1
210    hp0(1,:)    = 0
211    hp0(nx,:)   = 0
212
213    i_hp0(:,1)  = 1
214    i_hp0(:,ny) = 1
215    hp0(:,1)    = 0
216    hp0(:,ny)   = 0
217
218
219  end subroutine prescribe_grid_boundary_0
220  !______________________________________________________________________________________
221
222  !______________________________________________________________________________________
223  !> break ice shelves
224  !!
225  !!
226  subroutine break_all_ice_shelves
227
228    implicit none
229    if (itracebug.eq.1)  call tracebug(' Entree dans routine break_all_ice_shelves ')
230
231
232    debug_3D(:,:,89)=0.
233
234    where (flot(:,:))                           ! floating from the beginning and eventually from grounding line retreat
235       i_hp(:,:) = 1                                      ! thickness prescribed to 1 m
236       hp (:,:)  = 1.     
237       hp0 (:,:)  = 1.   
238       H(:,:)     =1.
239       debug_3D(:,:,89)=1.
240
241    end where
242    where (MK_flot0(:,:).eq. 1)
243       i_hp(:,:) = 1                                      ! thickness prescribed to 1 m
244       hp (:,:)  = 1. 
245       hp0 (:,:)  = 1.   
246       H(:,:)     =1.
247    end where
248
249
250  end subroutine break_all_ice_shelves
251
252  !
253  !______________________________________________________________________________________
254  !> break ice shelves
255  !!
256  !!
257  subroutine melt_ice_shelves
258
259    implicit none
260    integer :: nbvoisins,nbdeglac,iv,jv
261    real    :: hmoy,hmin,hmax
262
263    if (itracebug.eq.1)  call tracebug(' Entree dans routine melt_ice_shelves')
264
265    do j = 2,ny-1
266       do i=2, nx-1
267          if  (flot(i,j)) then                                 ! floating
268             i_hp(i,j) = 1                                      ! thickness prescribed
269             hmoy=0.
270             hmin=0.
271             hmax=0.
272             nbvoisins=0
273             nbdeglac=0
274             do jv=-1,1
275                do iv=1,1
276                   if (flot(i+iv,j+jv)) then
277                      nbvoisins=nbvoisins+1
278                      hmoy=hmoy+h(i+iv,j+jv)
279                      hmin=min(hmin,h(i+iv,j+jv))
280                      hmax=max(hmax,h(i+iv,j+jv))
281                      if (h(i+jv,j+jv).lt.10.) nbdeglac=nbdeglac+1
282                   end if
283                end do
284             end do
285             hmoy=hmoy/(max(nbvoisins,1))
286
287             if(H(i,j).gt.hmoy) then  ! melt only if smooth
288                hp(i,j)=H(i,j)-1.*dt
289             else if (hmax-hmin.lt.10.) then
290                hp(i,j)=H(i,j)-1.*dt
291             else if ((H(i,j).lt.50.).and.(Hmin.lt.30.)) then
292                hp(i,j)=1.
293             else if ((H(i,j).lt.50.).and.(nbdeglac.ge.4)) then
294                hp(i,j)=1.
295             endif
296             hp(i,j)   = max(hp(i,j),1.)
297          end if
298       end do
299    end do
300
301  end subroutine melt_ice_shelves
302
303
304   !______________________________________________________________________________________
305  !>  function prescribe_present_H_gl
306  !! calculate the initial grounding line position
307  !! @return i_hp(:,:) and hp(:,:)
308
309  subroutine prescribe_present_H_gl_copy
310
311    implicit none
312    if (itracebug.eq.1)  call tracebug(' Entree dans routine  prescribe_present_H_gl-copy')
313
314    do j=1,ny
315       do i=1,nx
316          if ((MK_flot0(i,j).eq. 1).or.(MK_gl0(i,j) .eq. 1)) then
317             i_hp(i,j) = 1                                      ! thickness prescribed to present value
318             hp(i,j)   = Hp0(i,j)
319          end if
320       end do
321    end do
322
323!!$where ((MK_flot0(:,:).eq. 1).or.(MK_gl0(:,:) .eq. 1)) ! floating or grounding line
324!!$   i_hp(:,:) = 1                                      ! thickness prescribed to present value
325!!$!   hp(:,:)   = 1   !Hp0(:,:)
326!!$end where
327
328
329  end subroutine prescribe_present_H_gl_copy
330!______________________________________________________________________________________
331
332
333end module prescribe_H
Note: See TracBrowser for help on using the repository browser.