source: trunk/SOURCES/furst_schoof_mod.f90 @ 94

Last change on this file since 94 was 94, checked in by aquiquet, 8 years ago

Schoof: small bug corrections + limits on U

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