source: branches/iLoveclim/SOURCES/dragging_hwatermax_0.2_mod.f90 @ 187

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

Grisli-iloveclim branch merged to trunk at revision 185

File size: 4.9 KB
Line 
1!> \file dragging_hwatermax_0.2_mod.f90
2!! Defini les zones de stream avec un critere
3!<
4
5!> \namespace  dragging_hwatstream
6!! Defini les zones de stream avec 1 critere
7!! \author ...
8!! \date ...
9!! @note Le critere est:
10!! @note   * un critere sur la hauteur d'eau
11!! @note Used module
12!! @note   - use module3D_phy
13!<
14module dragging_hwatstream
15  ! Defini les zones de stream avec un critere sur la hauteur d'eau
16
17  use module3D_phy
18  use fake_beta_iter_vitbil_mod
19
20
21  implicit none
22  real :: tostick
23
24contains
25  !-------------------------------------------------------------------------------
26  !> SUBROUTINE: init_dragging
27  !!Cette routine fait l'initialisation du dragging.
28  !>
29  subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
30
31    implicit none
32
33    namelist/drag_hwat/hwatstream,cf,betamax,toblim,tostick   
34
35    if (itracebug.eq.1)  call tracebug(' Dragging avec hwatermax')
36
37
38    ! formats pour les ecritures dans 42
39428 format(A)
40
41    ! lecture des parametres du run                      block draghwat
42    !--------------------------------------------------------------------
43
44    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
45    read(num_param,drag_hwat)
46
47    write(num_rep_42,428)'!___________________________________________________________' 
48    write(num_rep_42,428) '&drag_hwat              ! nom du bloc hwatermax '
49    write(num_rep_42,*)
50    write(num_rep_42,*) 'hwatsream      = ',hwatstream 
51    write(num_rep_42,*) 'cf             = ',cf
52    write(num_rep_42,*) 'betamax         = ', betamax
53    write(num_rep_42,*) 'toblim         = ', toblim
54    write(num_rep_42,*)'/'                           
55    write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream si hwater > hwatstream '
56    write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
57    write(num_rep_42,428) '! betamax : (Pa) frottement maxi '
58    write(num_rep_42,428) '! toblim : (Pa) tes pour les iles '
59    write(num_rep_42,*)
60
61    inv_beta=0
62    tostick=1.e5   ! valeurs pour les points non flgzmx
63
64    !-------------------------------------------------------------------
65    ! masque stream
66
67    mstream_mx(:,:)=1
68    mstream_my(:,:)=1
69
70
71    ! coefficient permettant de modifier le basal drag.
72    drag_mx(:,:)=1.
73    drag_my(:,:)=1.
74
75
76    return
77  end subroutine init_dragging
78  !________________________________________________________________________________
79
80  !> SUBROUTINE: dragging
81  !!Defini la localisation des streams et le frottement basal
82  !>
83  !-------------------------------------------------------------------------
84  subroutine dragging   ! defini la localisation des streams et le frottement basal
85
86
87    !         les iles n'ont pas de condition neff mais ont la condition toblim
88    !         (ce bloc etait avant dans flottab)
89
90    do j=2,ny
91       do i=2,nx
92          ilemx(i,j)=ilemx(i,j).and.(abs(tobmx(i,j)).lt.toblim)
93          ilemy(i,j)=ilemy(i,j).and.(abs(tobmy(i,j)).lt.toblim)
94       end do
95    end do
96
97
98    gzm_x: do j=1,ny  ! 2,ny
99       do i=2,nx ! -1
100
101          ! gzmx pour les hauteurs d'eau superieures a un seuil
102          ! attention, il semble qu'il y avait une erreur dans la precedente
103          ! formulation car le test etait seulement sur hwater(i,j).
104          ! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
105
106          gzmx(i,j)=gzmx(i,j).or. & 
107               (((hwater(i,j)+hwater(i-1,j))*0.5.gt.hwatstream) &
108               .and..not.flotmx(i,j))
109                  fleuvemx(i,j)=gzmx(i,j)
110
111          ! calcul du frottement basal (ce bloc etait avant dans neffect)
112          if (gzmx(i,j)) then                          ! stream
113             betamx(i,j)=cf*neffmx(i,j) 
114          else if (flotmx(i,j).or.ilemx(i,j)) then     ! flottant ou ile
115             betamx(i,j)=0.
116          else                                         ! grounded, SIA
117             betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
118          endif
119
120       end do
121    end do gzm_x
122
123    gzm_y: do j=2,ny  !-1
124       do i=1,nx  ! 2,nx
125
126          ! gzmy pour les hauteurs d'eau superieures a un seuil
127          ! attention, il semble qu'il y avait une erreur dans la precedente
128          ! formulation car le test etait seulement sur hwater(i,j).
129
130          gzmy(i,j)=gzmy(i,j).or. &
131               (((hwater(i,j)+hwater(i,j-1))*0.5.gt.hwatstream) &
132               .and..not.flotmy(i,j))
133          fleuvemy(i,j)=gzmx(i,j)
134
135          ! calcul du frottement basal (ce bloc etait avant dans neffect)
136          if (gzmy(i,j)) then                ! stream
137             betamy(i,j)=cf*neffmy(i,j) 
138          else if (flotmy(i,j).or.ilemy(i,j)) then     ! flottant ou ile
139             betamy(i,j)=0.
140          else                               ! grounded, SIA
141             betamy(i,j)=tostick             ! frottement glace posee
142          endif
143
144       end do
145    end do gzm_y
146
147    return
148  end subroutine dragging
149
150end module dragging_hwatstream
151
Note: See TracBrowser for help on using the repository browser.