source: trunk/SOURCES/furst_schoof_mod.f90 @ 97

Last change on this file since 97 was 97, checked in by aquiquet, 7 years ago

Small bug correction in Schoof grounding line flux

File size: 20.7 KB
Line 
1! Module Furst Schoof
2
3module furst_schoof_mod
4
5use module3d_phy
6use deformation_mod_2lois
7
8implicit none
9
10real :: frot_coef
11real :: alpha_flot
12integer :: gr_select
13real,dimension(nx,ny) :: back_force
14
15contains
16
17subroutine init_furst_schoof
18
19namelist/furst_schoof/frot_coef,gr_select
20
21rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
22read(num_param,furst_schoof)
23
24write(num_rep_42,*)'!___________________________________________________________' 
25write(num_rep_42,*)'! module furst_schoof'
26write(num_rep_42,*)
27write(num_rep_42,*)'! frot_coef  = ',frot_coef
28write(num_rep_42,*)'! gr_select  = ',gr_select
29write(num_rep_42,*)'! frot_coef : solid friction law coef 0.035 in Ritz et al 2001'
30write(num_rep_42,*)'! gr_select = 1 : Tsai , 2 : Schoof'
31write(num_rep_42,*)'!___________________________________________________________' 
32
33!  frot = 0.035 ! f in the solid friction law, 0.035 in Ritz et al. 2001       
34 
35  ! choix Tsai (loi de Coulomb) ou Schoof (loi de Weertman)
36  alpha_flot= ro/row
37  back_force(:,:)=0.1
38 
39end subroutine init_furst_schoof 
40 
41subroutine interpol_glflux
42
43  implicit none
44  real :: xpos,ypos
45  real :: Hgl
46  !  integer :: igl
47
48  real :: back_force_coef=0.1 !!!!!! TO DO !!!!!!!!
49
50  real :: phi_im2, phi_im1, phi_i, phi_ip1, phi_ip2, phi_ip3
51  real :: phi_jm2, phi_jm1, phi_j, phi_jp1, phi_jp2, phi_jp3
52  real,dimension(nx,ny) :: phi_xgl, phi_ygl ! flux gr line sur maille stag
53  integer,dimension(nx,ny) :: countx, county   ! how often do we modify ux/uy
54  real :: phi_prescr
55  real :: archim
56  real :: toutpetit = 1e-6
57  real :: denom, prodscal
58
59  !debug
60  real,dimension(nx,ny) :: xpos_tab=9999.
61  real,dimension(nx,ny) :: ypos_tab=9999.
62  real,dimension(nx,ny) :: Hglx_tab=9999.
63  real,dimension(nx,ny) :: Hgly_tab=9999.
64  real,dimension(nx,ny) :: Hmxloc, Hmyloc ! pour eviter division par 0
65  real,dimension(nx,ny) :: phi_prescr_tabx, phi_prescr_taby
66
67  real,dimension(nx,ny) :: uxbartemp, uybartemp
68 
69  real,dimension(nx,ny) :: archimtab
70
71  phi_prescr_tabx=0.
72  phi_prescr_taby=0.
73
74  Hmxloc(:,:)= max(Hmx(:,:),1.)
75  Hmyloc(:,:)=max(Hmy(:,:),1.)
76
77  uxbartemp(:,:)=0.
78  uybartemp(:,:)=0.
79  countx(:,:)=0.
80  county(:,:)=0.
81
82  archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot -sealevel
83
84  ! detection des noeuds grounding line
85  do j= 3,ny-3
86     do i=3,nx-3
87        !archim = Bsoc(i,j)+H(i,j)*ro/row -sealevel
88        if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel)) then    ! grounded and with ice
89           !           print*,'schoof gr line ij',i,j
90           ! selon x (indice i)
91           
92           if (archimtab(i-1,j).LT.0..and.Uxbar(i,j).LT.0.) then
93           !if (flot(i-1,j).and.Uxbar(i,j).LT.0.) then    ! grounding line between i-1,i    CAS WEST (ecoulement vers grille West)
94
95              call  calc_xgl4schoof(-dx,alpha_flot,H(i,j),Bsoc(i,j),H(i-1,j),Bsoc(i-1,j),Sealevel,xpos,Hgl)
96              xpos_tab(i,j)=xpos
97              Hglx_tab(i,j)=Hgl
98
99              if (gr_select.eq.1) then ! flux de Tsai
100                 call flux_Tsai4Schoof (Hgl,Abar(i,j),frot_coef,alpha_flot,glen(1),back_force(i-1,j),phi_prescr)
101              else
102                 print*,'ATTENTION FLUX AUTRE QUE TSAI NON IMPLEMENTE'
103                 stop
104              endif
105              phi_prescr=-phi_prescr ! Doit être negatif dans le cas West
106              phi_prescr_tabx(i,j)=phi_prescr
107
108              if (xpos .lt. -dx/2.) then  !  GL a west du i staggered,                          o centre, x stag
109
110                 !  grille    ! .....x......o......x......o...G..x......O......x......o......x......o
111                 !            !            i-2           i-1            i            i+1           i+2
112                 !flux        !     in            out         G out            in
113
114
115
116                 ! extrapolation pour avoir uxbar(i-1,j) qui sera ensuite prescrit dans call rempli_L2
117
118                 phi_im2 = Uxbar(i-2,j) * Hmxloc(i-2,j)
119                 phi_xgl(i-1,j)  = phi_im2 + dx * (phi_prescr-phi_im2)/(2.5 * dx + xpos)
120                 uxbartemp(i-1,j)   = uxbartemp(i-1,j)+ phi_xgl(i-1,j) / Hmxloc(i-1,j) ! attention division par 0 possible ?
121                 countx   (i-1,j)   = countx(i-1,j)+1
122                 imx_diag (i-1,j)   = 1
123
124                 ! extrapolation pour avoir uxbar(i,j)
125
126                 phi_ip1 = Uxbar(i+1,j) * Hmxloc(i+1,j)
127                 phi_xgl(i,j) = phi_ip1 + dx * (phi_prescr - phi_ip1)/(0.5 * dx - xpos)
128                 !                 print*,'uxbar(i,j)=',uxbar(i,j)
129                 uxbartemp(i,j)   = uxbartemp(i,j) + phi_xgl(i,j) / Hmxloc(i,j)
130                 countx   (i,j)   = countx(i,j)+1
131                 imx_diag (i,j)   = 1
132                 !                 print*,'schoof uxbar(i,j)=',uxbar(i,j)
133                 !                 print*,'Hgl',Hgl
134                 !                 print*,'phi_prescr',phi_prescr, phi_prescr/Hgl
135                 !                 print*
136
137              else            !  GL a Est du i staggered,                          o centre, x stag
138
139                 !  grille    ! .....x......o......x......o......x...G..O......x......o......x......o
140                 !            i-2           i-1            i            i+1           i+2
141                 !flux        !     in             in           out  G         out           in
142
143
144
145                 ! extrapolation pour avoir uxbar(i,j) qui sera ensuite prescrit dans call rempli_L2
146
147                 phi_im1 = Uxbar(i-1,j) * Hmxloc(i-1,j)
148                 phi_xgl(i,j)  = phi_im1 + dx * (phi_prescr-phi_im1)/(1.5 * dx + xpos)
149                 uxbartemp(i,j)   = uxbartemp(i,j) + phi_xgl(i,j) / Hmxloc(i,j) ! attention division par 0 possible ?
150                 countx   (i,j)   = countx(i,j)+1
151                 imx_diag (i,j) = 1
152
153                 ! extrapolation pour avoir uxbar(i+1,j)
154
155                 phi_ip2 = Uxbar(i+2,j) * Hmxloc(i+2,j)
156                 phi_xgl(i+1,j) = phi_ip2 + dx * (phi_prescr - phi_ip2)/(1.5 * dx  - xpos)
157                 uxbartemp(i+1,j)   = uxbartemp(i+1,j) +  phi_xgl(i+1,j) / Hmxloc(i+1,j)
158                 countx   (i+1,j)   = countx(i+1,j)+1
159                 imx_diag (i+1,j)   = 1
160              end if
161                                         
162           else if (archimtab(i+1,j).LT.0..and.Uxbar(i+1,j).GT.0.) then
163           !else if (flot(i+1,j).and.Uxbar(i+1,j).GT.0.) then    ! grounding line between i,i+1    CAS EST (ecoulement vers grille Est)
164
165              call  calc_xgl4schoof(dx,alpha_flot,H(i,j),Bsoc(i,j),H(i+1,j),Bsoc(i+1,j),Sealevel,xpos,Hgl)
166              xpos_tab(i,j)=xpos
167              Hglx_tab(i,j)=Hgl
168              if (gr_select.eq.1) then ! flux de Tsai
169                 call flux_Tsai4Schoof (Hgl,Abar(i,j),frot_coef,alpha_flot,glen(1),back_force(i+1,j),phi_prescr)
170              else
171                 print*,'ATTENTION FLUX AUTRE QUE TSAI NON IMPLEMENTE'
172                 stop
173              endif
174              phi_prescr_tabx(i,j)=phi_prescr
175              if (xpos .lt. -dx/2.) then  !  GL a west du i staggered,                          o centre, x stag
176
177                 !  grille    ! .....x......o......x......o......x......O..G...x......o......x......o
178                 !            !            i-2           i-1            i            i+1           i+2
179                 !flux        !                    in           out        G  out            in
180
181
182                 ! extrapolation pour avoir uxbar(i,j) qui sera ensuite prescrit dans call rempli_L2
183
184                 phi_im1 = Uxbar(i-1,j) * Hmxloc(i-1,j)
185                 phi_xgl(i,j)  = phi_im1 + dx * (phi_prescr-phi_im1)/(1.5 * dx + xpos)
186                 uxbartemp(i,j)   =  uxbartemp(i,j) + phi_xgl(i,j) / Hmxloc(i,j) ! attention division par 0 possible ?
187                 countx   (i,j)   = countx(i,j)+1
188                 imx_diag (i,j) = 1
189
190                 ! extrapolation pour avoir uxbar(i+1,j)
191
192                 phi_ip2 = Uxbar(i+2,j) * Hmxloc(i+2,j)
193                 phi_xgl(i+1,j) = phi_ip2 + dx * (phi_prescr - phi_ip2)/(1.5* dx - xpos)
194                 uxbartemp(i+1,j)   = uxbartemp(i+1,j) + phi_xgl(i+1,j) / Hmxloc(i+1,j)
195                 countx   (i+1,j)   = countx(i+1,j)+1
196                 imx_diag (i+1,j) = 1
197
198              else   !  GL a west du i staggered,                          o centre, x stag
199
200                 !  grille    ! ......x......o......x......O......x...G..o......x......o......x.....o
201                 !            !             i-1            i            i+1           i+2          i+3
202                 !flux        !                     in           out  G        out            in
203
204
205                 ! extrapolation pour avoir uxbar(i,j) qui sera ensuite prescrit dans call rempli_L2
206
207                 phi_i = Uxbar(i,j) * Hmxloc(i,j)
208                 phi_xgl(i+1,j)  = phi_i + dx * (phi_prescr-phi_i)/(0.5 * dx + xpos)
209                 uxbartemp(i+1,j)   = uxbartemp(i+1,j) + phi_xgl(i+1,j) / Hmxloc(i+1,j) ! attention division par 0 possible ?
210                 countx   (i+1,j)   = countx(i+1,j)+1
211                 imx_diag (i+1,j) = 1
212
213                 ! extrapolation pour avoir uxbar(i+1,j)
214
215                 phi_ip3 = Uxbar(i+3,j) * Hmxloc(i+3,j)
216                 phi_xgl(i+2,j) = phi_ip3 + dx * (phi_prescr - phi_ip3)/(2.5 * dx  - xpos)
217                 uxbartemp(i+2,j)   = uxbartemp(i+2,j) + phi_xgl(i+2,j) / Hmxloc(i+2,j)
218                 countx   (i+2,j)   = countx(i+2,j)+1
219                 imx_diag (i+2,j) = 1
220              end if
221
222           end if
223
224           !*******************************************************************************
225           ! selon y (indice j)
226           if (archimtab(i,j-1).LT.0..and.Uybar(i,j).LT.0.) then 
227           !if (flot(i,j-1).and.Uybar(i,j).LT.0.) then    ! grounding line between j-1,j    CAS SUD (ecoulement vers grille SUD)
228
229
230              call  calc_xgl4schoof(-dy,alpha_flot,H(i,j),Bsoc(i,j),H(i,j-1),Bsoc(i,j-1),Sealevel,ypos,Hgl)
231              ypos_tab(i,j)=ypos
232              Hgly_tab(i,j)=Hgl
233              if (gr_select.eq.1) then ! flux de Tsai
234                 call flux_Tsai4Schoof (Hgl,Abar(i,j),frot_coef,alpha_flot,glen(1),back_force(i,j-1),phi_prescr)
235              else
236                 print*,'ATTENTION FLUX AUTRE QUE TSAI NON IMPLEMENTE'
237                 stop
238              endif
239              phi_prescr=-phi_prescr ! Doit être negatif dans le cas Sud
240              phi_prescr_taby(i,j)=phi_prescr
241              if (ypos .lt. -dy/2.) then  !  GL au sud du j staggered,                          o centre, x stag
242
243                 !  grille    ! .....x......o......x......o...G..x......O......x......o......x......o
244                 !            !            j-2           j-1            j            j+1           j+2
245                 !flux        !     in            out         G out            in
246
247
248                 ! extrapolation pour avoir uybar(i,j-1) qui sera ensuite prescrit dans call rempli_L2
249
250                 phi_jm2 = Uybar(i,j-2) * Hmyloc(i,j-2)
251                 phi_ygl(i,j-1)  = phi_jm2 + dy * (phi_prescr-phi_jm2)/(2.5 * dy + ypos)
252                 uybartemp(i,j-1)   = uybartemp(i,j-1) + phi_ygl(i,j-1) / Hmyloc(i,j-1) ! attention division par 0 possible ?
253                 county   (i,j-1)   = county(i,j-1)+1
254                 imy_diag (i,j-1) = 1
255
256                 ! extrapolation pour avoir uybar(i,j)
257
258                 phi_jp1 = Uybar(i,j+1) * Hmyloc(i,j+1)
259                 phi_ygl(i,j) = phi_jp1 + dy * (phi_prescr - phi_jp1)/(0.5 * dy - ypos)
260                 uybartemp(i,j)   = uybartemp(i,j) + phi_ygl(i,j) / Hmyloc(i,j)
261                 county   (i,j)   = county(i,j)+1
262                 imy_diag (i,j) = 1
263
264              else            !  GL au nord du j staggered,                          o centre, x stag
265
266                 !  grille    ! .....x......o......x......o......x...G..O......x......o......x......o
267                 !            j-2           j-1            j            j+1           j+2
268                 !flux        !     in             in           out  G         out           in
269
270
271                 ! extrapolation pour avoir uybar(i,j) qui sera ensuite prescrit dans call rempli_L2
272
273                 phi_jm1 = Uybar(i,j-1) * Hmyloc(i,j-1)
274                 phi_ygl(i,j)  = phi_jm1 + dy * (phi_prescr-phi_jm1)/(1.5 * dy + ypos)
275                 uybartemp(i,j)   = uybartemp(i,j) + phi_ygl(i,j) / Hmyloc(i,j) ! attention division par 0 possible ?
276                 county   (i,j)   = county(i,j)+1
277                 imy_diag(i,j) = 1
278
279                 ! extrapolation pour avoir uybar(i,j+1)
280
281                 phi_jp2 = Uybar(i,j+2) * Hmyloc(i,j+2)
282                 phi_ygl(i,j+1) = phi_jp2 + dy * (phi_prescr - phi_jp2)/(1.5 * dy  - ypos)
283                 uybartemp(i,j+1)   = uybartemp(i,j+1) + phi_ygl(i,j+1) / Hmyloc(i,j+1)
284                 county   (i,j+1)   = county(i,j+1)+1
285                 imy_diag (i,j+1) = 1
286              end if
287             
288           else if (archimtab(i,j+1).LT.0..and.Uybar(i,j+1).GT.0.) then
289           !else if (flot(i,j+1).and.Uybar(i,j+1).GT.0.) then    ! grounding line between j,j+1    CAS NORD (ecoulement vers grille Nord)
290
291              call  calc_xgl4schoof(dy,alpha_flot,H(i,j),Bsoc(i,j),H(i,j+1),Bsoc(i,j+1),Sealevel,ypos,Hgl)
292              ypos_tab(i,j)=ypos
293              Hgly_tab(i,j)=Hgl
294              if (gr_select.eq.1) then ! flux de Tsai
295                 call flux_Tsai4Schoof (Hgl,Abar(i,j),frot_coef,alpha_flot,glen(1),back_force(i,j+1),phi_prescr)
296              else
297                 print*,'ATTENTION FLUX AUTRE QUE TSAI NON IMPLEMENTE'
298                 stop
299              endif
300              phi_prescr_taby(i,j)=phi_prescr
301              if (ypos .lt. -dy/2.) then  !  GL au sud du j staggered,                          o centre, x stag
302
303                 !  grille    ! .....x......o......x......o......x......O..G...x......o......x......o
304                 !            !            j-2           j-1            j            j+1           j+2
305                 !flux        !                    in           out        G  out            in
306
307
308                 ! extrapolation pour avoir uybar(i,j) qui sera ensuite prescrit dans call rempli_L2
309
310                 phi_jm1 = Uybar(i,j-1) * Hmyloc(i,j-1)
311                 phi_ygl(i,j)  = phi_jm1 + dy * (phi_prescr-phi_jm1)/(1.5 * dy + ypos)
312                 uybartemp(i,j)   = uybartemp(i,j) + phi_ygl(i,j) / Hmyloc(i,j) ! attention division par 0 possible ?
313                 county   (i,j)   = county(i,j)+1
314                 imy_diag (i,j) = 1
315
316                 ! extrapolation pour avoir uybar(i,j+1)
317
318                 phi_jp2 = Uybar(i,j+2) * Hmyloc(i,j+2)
319                 phi_ygl(i,j+1) = phi_jp2 + dy * (phi_prescr - phi_jp2)/(1.5* dy - ypos)
320                 uybartemp(i,j+1)   = uybartemp(i,j+1) + phi_ygl(i,j+1) / Hmyloc(i,j+1)
321                 county   (i,j+1)   = county(i,j+1)+1
322                 imy_diag (i,j+1) = 1
323
324              else   !  GL au nord du j staggered,                          o centre, x stag
325
326                 !  grille    ! ......x......o......x......O......x...G..o......x......o......x.....o
327                 !            !             j-1            j            j+1           j+2          j+3
328                 !flux        !                     in           out  G        out            in
329
330
331                 ! extrapolation pour avoir uybar(i,j) qui sera ensuite prescrit dans call rempli_L2
332
333                 phi_j = Uybar(i,j) * Hmyloc(i,j)
334                 phi_ygl(i,j+1)  = phi_j + dy * (phi_prescr-phi_j)/(0.5 * dy + ypos)
335                 uybartemp(i,j+1)   = uybartemp(i,j+1) +  phi_ygl(i,j+1) / Hmyloc(i,j+1) ! attention division par 0 possible ?
336                 county   (i,j+1)   = county(i,j+1)+1
337                 imy_diag (i,j+1) = 1
338
339                 ! extrapolation pour avoir uybar(i,j+1)
340
341                 phi_jp3 = Uybar(i,j+3) * Hmyloc(i,j+3)
342                 phi_ygl(i,j+2) = phi_jp3 + dy * (phi_prescr - phi_jp3)/(2.5 * dy  - ypos)
343                 uybartemp(i,j+2)   = uybartemp(i,j+2) +  phi_ygl(i,j+2) / Hmyloc(i,j+2)
344                 county   (i,j+2)   = county(i,j+2)+1
345                 imy_diag (i,j+2) = 1
346              end if
347
348           end if
349        end if
350     end do
351  end do
352  where (countx(:,:).ne.0) uxbar(:,:)=uxbartemp(:,:)/countx(:,:)
353  where (county(:,:).ne.0) uybar(:,:)=uybartemp(:,:)/county(:,:)
354 
355  uxbar(:,:)= max(min(uxbar(:,:),V_limit),-V_limit)
356  uybar(:,:)= max(min(uybar(:,:),V_limit),-V_limit)
357 
358  do j= 3,ny-3
359     do i=3,nx-3
360        denom = sqrt((uxbartemp(i,j)**2 + uybartemp(i,j)**2))*sqrt((uxbar(i,j)**2 + uybar(i,j)**2))
361        if (denom.gt.toutpetit) then
362           ! prodscal is cos theta in : U . V = u v cos theta
363           ! this number is between -1 and 1, 1 = same direction & -1 = opposite direction
364           ! we could use this as an artificial back force...
365           prodscal = (uxbartemp(i,j)*uxbar(i,j) + uybartemp(i,j)*uybar(i,j))/  &
366                (sqrt((uxbartemp(i,j)**2 + uybartemp(i,j)**2))*sqrt((uxbar(i,j)**2 + uybar(i,j)**2)))
367        endif
368     end do
369  end do
370 
371!!$  do j=1,ny
372!!$     do i=1,nx
373!!$        write(598,*) xpos_tab(i,j)
374!!$        write(599,*) ypos_tab(i,j)
375!!$        write(600,*) Hglx_tab(i,j)
376!!$        write(601,*) Hgly_tab(i,j)
377!!$        write(602,*) Abar(i,j)
378!!$        write(603,*) phi_prescr_tabx(i,j)/Hglx_tab(i,j)
379!!$        write(604,*) phi_prescr_taby(i,j)/Hgly_tab(i,j)
380!!$     enddo
381!!$  enddo
382end subroutine interpol_glflux
383 
384!-------------------------------------------------------------------------
385! calcule la position sous maille de la ligne dechouage
386 subroutine calc_xgl4schoof(dy,alpha,H_0,B_0,H_1,B_1,SL,xpos,Hgl)
387               
388
389  implicit none
390
391  real,intent(in)        ::    dy        !< longueur orientee de la maille
392  real,intent(in)        ::    alpha     !< coefficient de flottaison = rho/rhow
393  real,intent(in)        ::    H_0       !< epaisseur au point pose
394  real,intent(in)        ::    B_0       !< altitude socle au point pose
395  real,intent(in)        ::    H_1       !< epaisseur au point flottant
396  real,intent(in)        ::    B_1       !< altitude socle au point flottant
397  real,intent(in)        ::    SL        !< Sea level
398  real,intent(out)       ::    xpos      !< position de la ligne (en distance depuis le point pose)
399  real,intent(out)       ::    Hgl       !< Epaisseur de glace a la grounding line
400   
401 
402
403  real                   ::    Cgl       ! variable de travail 
404
405!~   Cgl   = (-alpha * (H_1-H_0) - (B_1-B_0)) ! -(rho/rhow (H1-H0) + (B1-B0))
406
407!~   if (abs(Cgl).gt.1.e-5) then   
408!~      xpos    = (B_0 + alpha * H_0 - SL) / Cgl
409!~   else
410!~      xpos    = 1. - 1./abs(dy)   ! verifier
411!~   end if
412 
413!~   if (xpos.LT.0..OR.xpos.GT.1.) then
414!~              print*,'calc_xgl4schoof : xpos value error i,j=',i,j
415!~              print*, 'xpos,Cgl', xpos,Cgl
416!~              print*,'B_0, alpha, H_0, SL', B_0, alpha, H_0, SL
417!~              stop
418!~      endif
419               
420  Cgl   = (B_1-B_0) + alpha * (H_1-H_0)
421
422  if (abs(Cgl).gt.1.e-5) then   
423     xpos    = (SL - (B_0 + alpha * H_0)) / Cgl 
424  else
425     xpos    = 0.5  ! verifier
426  end if
427 
428!  print*
429!  print*, 'xpos', xpos
430!       print*,'B_0, alpha, H_0, SL', B_0, alpha, H_0, SL
431!       print*,'B_1, H_1', B_1, H_1
432
433       
434  if (xpos.LT.0..OR.xpos.GT.1.) then
435                print*,'calc_xgl4schoof : xpos value error i,j=',i,j
436                print*, 'xpos,Cgl', xpos,Cgl
437                print*,'B_0, alpha, H_0, SL', B_0, alpha, H_0, SL
438                print*,'archim=',B_1+H_1*alpha - SL
439                !stop
440                xpos = min(max(0.,xpos),1.)
441  endif
442       
443  xpos = xpos * dy
444 
445  Hgl= H_0 + (H_1 - H_0) * xpos/dy
446
447end subroutine calc_xgl4schoof
448
449
450!--------------------------------------------------------------------
451! Schoof_flux : calcule le coefficient et le flux de Schoof
452!--------------------------------------------------------------------
453
454! Dans un premier temps, pas de terme horizontal (regarder dans toy_recul une tentative de confinement)
455
456subroutine flux_Schoof4Schoof (Hgl,A_mean,C_frot,alpha,n_glen,m_weert,back_force_coef,phi_schoof)
457
458 
459implicit none
460
461real,intent(in) :: Hgl
462real,intent(in) :: A_mean
463real,intent(in) :: C_frot
464real,intent(in) :: alpha
465real,intent(in) :: n_glen
466real,intent(in) :: m_weert
467real,intent(in) :: back_force_coef
468real,intent(out) :: phi_schoof
469
470
471phi_schoof = ((A_mean * rog**(n_glen+1) * (1 - alpha)**n_glen) / (4**n_glen *C_frot)) **(1./(1.+1./m_weert))
472phi_schoof = phi_schoof * Hgl**((n_glen+3.+1/m_weert)/(1.+1./m_weert))
473
474
475
476end subroutine flux_Schoof4Schoof
477
478!--------------------------------------------------------------------
479! Tsai_flux : calcule le coefficient et le flux de Schoof
480!--------------------------------------------------------------------
481
482! Dans un premier temps, pas de terme horizontal (regarder dans toy_recul une tentative de confinement)
483
484subroutine flux_Tsai4Schoof (Hgl,A_mean,f_frot,alpha,n_glen,back_force_coef,phi_Tsai)
485
486implicit none
487
488real,intent(in) :: Hgl
489real,intent(in) :: A_mean
490real,intent(in) :: f_frot
491real,intent(in) :: alpha
492real,intent(in) :: n_glen
493real,intent(in) :: back_force_coef
494real,intent(out) :: phi_Tsai
495
496real, parameter :: Q0 = 0.61
497
498
499phi_Tsai = Q0 * ((8 * A_mean * rog**(n_glen) * (1 - alpha)**(n_glen-1.)) / (4**n_glen *f_frot))*back_force_coef
500phi_Tsai = phi_Tsai * Hgl**(n_glen+2)
501
502end subroutine flux_Tsai4Schoof
503
504end module furst_schoof_mod
Note: See TracBrowser for help on using the repository browser.