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

Last change on this file since 77 was 77, checked in by dumas, 8 years ago

Merge branche iLOVECLIM sur rev 76

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
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    tostick=1.e5   ! valeurs pour les points non flgzmx
62
63    !-------------------------------------------------------------------
64    ! masque stream
65
66    mstream_mx(:,:)=1
67    mstream_my(:,:)=1
68
69
70    ! coefficient permettant de modifier le basal drag.
71    drag_mx(:,:)=1.
72    drag_my(:,:)=1.
73
74
75    return
76  end subroutine init_dragging
77  !________________________________________________________________________________
78
79  !> SUBROUTINE: dragging
80  !!Defini la localisation des streams et le frottement basal
81  !>
82  !-------------------------------------------------------------------------
83  subroutine dragging   ! defini la localisation des streams et le frottement basal
84
85
86    !         les iles n'ont pas de condition neff mais ont la condition toblim
87    !         (ce bloc etait avant dans flottab)
88
89    do j=2,ny
90       do i=2,nx
91          ilemx(i,j)=ilemx(i,j).and.(abs(tobmx(i,j)).lt.toblim)
92          ilemy(i,j)=ilemy(i,j).and.(abs(tobmy(i,j)).lt.toblim)
93       end do
94    end do
95
96
97    gzm_x: do j=1,ny  ! 2,ny
98       do i=2,nx ! -1
99
100          ! gzmx pour les hauteurs d'eau superieures a un seuil
101          ! attention, il semble qu'il y avait une erreur dans la precedente
102          ! formulation car le test etait seulement sur hwater(i,j).
103          ! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
104
105          gzmx(i,j)=gzmx(i,j).or. & 
106               (((hwater(i,j)+hwater(i-1,j))*0.5.gt.hwatstream) &
107               .and..not.flotmx(i,j))
108                  fleuvemx(i,j)=gzmx(i,j)
109
110          ! calcul du frottement basal (ce bloc etait avant dans neffect)
111          if (gzmx(i,j)) then                          ! stream
112             betamx(i,j)=cf*neffmx(i,j) 
113          else if (flotmx(i,j).or.ilemx(i,j)) then     ! flottant ou ile
114             betamx(i,j)=0.
115          else                                         ! grounded, SIA
116             betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
117          endif
118
119       end do
120    end do gzm_x
121
122    gzm_y: do j=2,ny  !-1
123       do i=1,nx  ! 2,nx
124
125          ! gzmy pour les hauteurs d'eau superieures a un seuil
126          ! attention, il semble qu'il y avait une erreur dans la precedente
127          ! formulation car le test etait seulement sur hwater(i,j).
128
129          gzmy(i,j)=gzmy(i,j).or. &
130               (((hwater(i,j)+hwater(i,j-1))*0.5.gt.hwatstream) &
131               .and..not.flotmy(i,j))
132          fleuvemy(i,j)=gzmx(i,j)
133
134          ! calcul du frottement basal (ce bloc etait avant dans neffect)
135          if (gzmy(i,j)) then                ! stream
136             betamy(i,j)=cf*neffmy(i,j) 
137          else if (flotmy(i,j).or.ilemy(i,j)) then     ! flottant ou ile
138             betamy(i,j)=0.
139          else                               ! grounded, SIA
140             betamy(i,j)=tostick             ! frottement glace posee
141          endif
142
143       end do
144    end do gzm_y
145
146    return
147  end subroutine dragging
148
149end module dragging_hwatstream
150
Note: See TracBrowser for help on using the repository browser.