source: branches/GRISLIv3/SOURCES/Draggings_modules/dragging_param_beta_mod.f90 @ 470

Last change on this file since 470 was 470, checked in by aquiquet, 5 months ago

Cleaning branch: some cleaning in flottab - ilemx,ilemy removed

File size: 7.9 KB
Line 
1!> \file dragging_neff_contmaj_mod.f90
2!! Defini les zones de stream avec 3 criteres
3!<
4
5!> \namespace  dragging_neff_contmaj
6!! Defini les zones de stream avec 3 criteres
7!! \author ...
8!! \date ...
9!! @note Les trois criteres sont:
10!! @note   * un critere sur la pression effective
11!! @note   * un critere de continuite depuis la cote
12!! @note   * un critere sur la courbure du socle (si negatif, vallees)
13!! @note Used module
14!! @note   - use module3D_phy
15!! @note   - use param_phy_mod
16!<
17
18module dragging_param_beta
19
20  ! Defini les zones de stream avec :
21  ! * un critere sur la pression effective
22  ! * un critere de continuite depuis la cote
23  ! * un critere sur la courbure du socle (si negatif, vallees)
24
25  use geography, only: nx,ny
26  ! pour compilation de steps_time_loop :
27  use fake_beta_iter_vitbil_mod, only:time_iter,time_iter_end,nb_iter_vitbil,beta_iter_vitbil
28
29  implicit none
30
31  logical,dimension(nx,ny) :: cote
32
33  real,dimension(nx,ny) :: beta_param ! local variable
34
35  real :: betamin   ! betamin : (Pa) frottement mini sous les streams
36
37  real :: beta_slope        ! = A in:     B = A x Neff ** K
38  real :: beta_expo         ! = K in:     B = A x Neff ** K
39
40  real,dimension(nx,ny) :: neff   ! pression effective noeuds majeurs
41  real,dimension(nx,ny) :: hwatmx ! hauteur d'eau staggered grille - afq
42  real,dimension(nx,ny) :: hwatmy ! hauteur d'eau staggered grille - afq
43
44  real :: neffmin = 1.e5 ! as in neffect, but this should be the exact same variable
45
46  real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
47  real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
48  real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
49  real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
50  logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
51
52contains
53  !-------------------------------------------------------------------------------
54
55  !> SUBROUTINE: init_dragging
56  !! Cette routine fait l'initialisation du dragging.
57  !>
58  subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
59
60    use module3d_phy, only:betamax,num_param,num_rep_42,inv_beta,mstream_mx,mstream_my,&
61         drag_mx,drag_my,betamax_2d,niter_nolin
62    use runparam, only:itracebug
63
64    implicit none
65
66    namelist/drag_param_beta/beta_slope,beta_expo,betamax,betamin
67
68    if (itracebug.eq.1)  call tracebug(' dragging avec neff et slope')
69
70
71
72    ! lecture des parametres du run                      block drag neff
73    !--------------------------------------------------------------------
74
75    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
76    read(num_param,drag_param_beta)
77
78    write(num_rep_42,'(A)')'!___________________________________________________________' 
79    write(num_rep_42,'(A)') '&drag_param_beta          ! nom du bloc dragging param beta'
80    write(num_rep_42,*)
81    write(num_rep_42,*) 'beta_slope            = ', beta_slope
82    write(num_rep_42,*) 'beta_expo             = ', beta_expo
83    write(num_rep_42,*) 'betamax               = ', betamax
84    write(num_rep_42,*) 'betamin               = ', betamin
85    write(num_rep_42,*)'/'                           
86    write(num_rep_42,'(A)') '! Slope & expo of beta = - slope x Neff ** expo'
87   
88    inv_beta=0
89    !-------------------------------------------------------------------
90    ! masque stream
91    mstream_mx(:,:)=1
92    mstream_my(:,:)=1
93
94
95    ! coefficient permettant de modifier le basal drag.
96    drag_mx(:,:)=1.
97    drag_my(:,:)=1.
98
99    betamax_2d(:,:) = betamax
100
101    niter_nolin = 1
102
103    return
104  end subroutine init_dragging
105  !________________________________________________________________________________
106
107  !> SUBROUTINE: dragging
108  !! Defini la localisation des streams et le frottement basal
109  !>
110  !-------------------------------------------------------------------------
111  subroutine dragging   ! defini le frottement basal
112
113    use module3d_phy, only:neffmx,neffmy,hwater,betamx,betamy,betamax,flot,&
114         beta_centre
115   
116    !$ USE OMP_LIB
117
118
119    implicit none
120    integer :: i,j
121
122    !$OMP PARALLEL
123    ! calcul de neff (pression effective sur noeuds majeurs)
124    if (sum(neffmx(:,:)).le.0.) neffmx(:,:) =1.e8
125    if (sum(neffmy(:,:)).le.0.) neffmy(:,:) =1.e8
126
127    !$OMP DO
128    do j=1,ny-1
129       do i=1,nx-1
130          neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4
131       enddo
132    enddo
133    !$OMP END DO
134    !aurel, for the borders:
135    !$OMP WORKSHARE
136    neff(:,ny)=neffmin
137    neff(nx,:)=neffmin
138    ! calcul de hwat (staggered)
139    hwatmx(:,:)=0.
140    hwatmy(:,:)=0.
141    !$OMP END WORKSHARE
142    !$OMP DO
143    do i=2,nx
144       hwatmx(i,:)=(hwater(i-1,:)+hwater(i,:))/2.
145    enddo
146    !$OMP END DO
147    !$OMP DO
148    do j=2,ny
149       hwatmy(:,j)=(hwater(:,j-1)+hwater(:,j))/2.
150    enddo
151    !$OMP END DO
152
153    !$OMP WORKSHARE
154
155    ! new parametrisation of beta on Neff:
156    betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo)
157    betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo)
158
159    betamx(:,:)=max(betamx(:,:),betamin)
160    betamy(:,:)=max(betamy(:,:),betamin)
161    betamx(:,:)=min(betamx(:,:),betamax)
162    betamy(:,:)=min(betamy(:,:),betamax)
163
164    where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
165    where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
166
167    !$OMP END WORKSHARE
168
169    ! points flottants
170    !$OMP DO
171    do j=2,ny
172       do i=2,nx
173          if (flot(i,j).and.(flot(i-1,j))) then
174             betamx(i,j)=0.
175          end if
176          if (flot(i,j).and.(flot(i,j-1))) then
177             betamy(i,j)=0.
178          end if
179       end do
180    end do
181    !$OMP END DO
182
183    !$OMP DO
184    do j=2,ny-1
185       do i=2,nx-1
186          beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
187               + (betamy(i,j)+betamy(i,j+1)))*0.25
188       end do
189    end do
190    !$OMP END DO
191    !$OMP END PARALLEL
192
193    return
194  end subroutine dragging
195  !________________________________________________________________________________
196
197  !> SUBROUTINE: mstream_dragging
198  !! Defini la localisation des streams
199  !>
200  !-------------------------------------------------------------------------
201  subroutine mstream_dragging   ! defini la localisation des streams
202
203    use module3d_phy, only:fleuvemx,fleuvemy,gzmx,gzmy,flgzmx,flgzmy,flotmx,flotmy,flot,&
204         betamx,betamy,betamax
205    !$ USE OMP_LIB
206   
207    implicit none
208    integer :: i,j
209
210    !$OMP PARALLEL
211
212    !$OMP WORKSHARE
213    fleuvemx(:,:)=.false.
214    fleuvemy(:,:)=.false.
215    cote(:,:)=.false.
216    gzmx(:,:)=.true.
217    gzmy(:,:)=.true.
218    flgzmx(:,:)=.false.
219    flgzmy(:,:)=.false.
220    !$OMP END WORKSHARE
221
222    ! detection des cotes
223    !$OMP DO
224    do  j=2,ny-1
225       do i=2,nx-1
226          if ((.not.flot(i,j)).and. &
227               ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
228             cote(i,j)=.true.
229          endif
230       end do
231    end do
232    !$OMP END DO
233
234
235    ! calcul de gzmx
236
237    ! points flottants : flgzmx mais pas gzmx
238    !$OMP DO
239    do j=2,ny
240       do i=2,nx
241          if (flot(i,j).and.(flot(i-1,j))) then
242             flotmx(i,j)=.true.
243             flgzmx(i,j)=.true.
244             gzmx(i,j)=.false.
245             betamx(i,j)=0.
246          end if
247          if (flot(i,j).and.(flot(i,j-1))) then
248             flotmy(i,j)=.true.
249             flgzmy(i,j)=.true.
250             gzmy(i,j)=.false.
251             betamy(i,j)=0.
252          end if
253       end do
254    end do
255    !$OMP END DO
256
257    !--------- autres criteres
258
259    ! points poses
260    !$OMP WORKSHARE
261    gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
262    gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
263
264    flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
265    flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
266
267    fleuvemx(:,:)=gzmx(:,:)
268    fleuvemy(:,:)=gzmy(:,:)
269    !$OMP END WORKSHARE
270
271    !$OMP END PARALLEL
272
273    return
274  end subroutine mstream_dragging
275
276end module dragging_param_beta
277
Note: See TracBrowser for help on using the repository browser.