source: trunk/SOURCES/steps_time_loop_avec_iterbeta.f90 @ 170

Last change on this file since 170 was 161, checked in by dumas, 7 years ago

Attempt to debug the classical crash in diagnoshelf : determin_tache now called by steps_time_loop | isostasy after diagnoshelf | flgzmx & flgzmy .false. if there is no ice

File size: 13.6 KB
Line 
1!> \file steps_time_loop.f90
2!! TOOOOOOOOO DOOOOOOOOOOOOO
3!! Step_grisli contient
4!<
5
6!--------------------------------------------------------------------------
7! step_time_loop est a l'interieur de la boucle temps.
8! il commence par icethick, qui détermine dt, isynchro
9
10! ensuite contient les appels aux diverses sorties
11! puis un appel a step_grisli
12
13
14subroutine step_time_loop       
15
16  use module3d_phy
17  use module_choix ! module de choix du type de run
18                   !  module_choix donne acces a tous les modules
19                   !  de declaration des packages
20  use icetempmod
21  use sorties_ncdf_grisli
22  use flottab_mod
23  use diagno_mod
24  use bilan_eau_mod
25!  use track_debug
26
27
28  implicit none
29  integer :: nt_init = 0      ! number of loops for initialisation of ice thickness
30  integer :: iter_visco       ! number of iterations for ssa viscosity
31
32  if (itracebug.eq.1)  call tracebug('Entree dans step_time_loop ')
33  if (itracebug.eq.1) write(num_tracebug,*) 'nt = ',nt,'  iter_beta = ',iter_beta
34
35
36  ! ice thickness evolution
37  !=============================
38
39
40  spinup_icethick: if ((ispinup.eq.0.or.ispinup.eq.2.or.nt.lt.nt_init) &
41       .and.(iter_beta.eq.0)) then 
42
43     !-------------------------------------------------------------------
44     ! spinup_icethick : Le calcul changement d'epaisseur (icethick)
45     !----------------
46     !   - est fait dans tous les cas pour les premieres boucles temporelles
47     !   - pour ispinup = 0 (run standard) sauf si iter_beta = 0
48     !   - pour ispinup = 2 (force avec vitesses de bilan)
49     !
50     ! else
51     !   juste un increment pour le temps
52     !    -  pour ispinup = 1 et  ispinup = 3  temperature calculation
53     
54     !-------------------------------------------------------------------
55
56     if (itracebug.eq.1)  call tracebug(' Dans spinup_icethick : appelle icethick')
57
58     call icethick3
59     call flottab
60     call determin_tache
61     call calving
62     call ablation_bord
63     call bilan_eau
64     call bilan_flux_output
65     call flottab
66     if (isynchro.eq.1.or.nt.eq.1) call determin_marais
67     if (itracebug.eq.1)  call tracebug('apres marais')
68
69     ! new calculation of ice surface elevation
70     where(.not.flot(:,:))
71        S(:,:)=Bsoc(:,:)+H(:,:)
72        S(:,:)=max(S(:,:),sealevel)
73        B(:,:)=Bsoc(:,:)
74     end where
75     
76     ! calcul de Hmx et Hmy -> shift=-1, dim=1 -> H(i-1,j)
77
78     hmx(:,:)=0.5*(H(:,:)+eoshift(H(:,:),shift=-1,boundary=0.,dim=1))
79     hmy(:,:)=0.5*(H(:,:)+eoshift(H(:,:),shift=-1,boundary=0.,dim=2))
80     hmx(1,:)=0.
81     hmy(:,1)=0.
82
83     ! mise a 0 des vitesses quand epaisseur nulle
84     where (hmx(:,:).le.1.)
85        uxbar(:,:)=0.
86     endwhere
87     where (hmy(:,:).le.1.)
88        uybar(:,:)=0.
89     endwhere
90
91  else   ! time increment when no icethick : ispinup= 1 or 3    ! loop spinup_icethick
92
93     hdot(:,:)=0.
94     time=nt*dtmax
95
96     if (itracebug.eq.1) write(num_tracebug,*) ' Dans spinup_icethick:  no icethick ',ispinup, &
97          '    time=',time
98
99
100  end if spinup_icethick
101  ! end of ice thickness evolution
102  !====================================================================
103
104
105
106  ! outputs
107  !====================================================================
108
109  !   vertical plan snapshots (profiles)
110  !------------------------------------------------------------------
111
112  !if ((NT.eq.1) &    !   .or.(NT.eq.2)
113  !     .or.(mod(abs(dble(TIME)),dble(DTPROFILE)).lt.dble(dtmin)) &
114  !     .or.(abs(TIME-TEND).lt.dtmin)) then
115  !  call sortieprofile()
116  !endif
117
118  !  sorties netcdf Hassine (aout 2010)  (2D and 3D)
119  !------------------------------------------------------------------
120
121  call testsort_time_ncdf(time)
122  if (iglob_ncdf .EQ. 1) call sortie_ncdf_cat
123  if (isynchro.eq.1) then
124    call shortoutput
125    diff_H = 0.
126    Bm_dtt(:,:) = 0.
127    bmelt_dtt(:,:) = 0.
128    calv_dtt(:,:)=0.
129    ablbord_dtt(:,:)=0.
130    diff_H_2D(:,:)=0.
131    grline_dtt(:,:)=0.
132  endif
133
134  !   sortie compteur tous les dtcpt ans
135  !------------------------------------------------------------------
136  !iout == 1 sortie cptr
137  !iout == 2 sortie nc
138  if (itracebug.eq.1)  call tracebug('dans steps_time_loop avant call out_recovery ')
139
140  call out_recovery(iout)
141
142  ! end of outputs
143
144  !====================================================================
145  !
146  call step_thermomeca()
147  !
148  ! contient tout le reste de la dynamique et du forçage
149  !====================================================================
150
151end subroutine step_time_loop
152
153
154
155
156
157
158!---------------------------------------------------------------------
159! The subroutine  step_thermomeca() calls the following processes
160!
161! - climatic forcing and its impact on surface mass balance
162! - dynamic velocities (SIA)
163! - ice temperature
164
165subroutine step_thermomeca()
166
167
168  use module3d_phy
169  use module_choix ! module de choix du type de run
170  !  module_choix donne acces a tous les modules
171  !  de declaration des packages
172  use icetempmod
173  use sorties_ncdf_grisli
174  use flottab_mod
175  use diagno_mod
176  use resolmeca_SIA_L1
177!  use track_debug
178
179
180  implicit none
181
182  integer :: m
183  logical :: shelf_vitbil = .true.   ! si vrai on prend les vitesses de bilan pour le shelf
184  integer :: nt_init_tm = 0          ! number of loops for initialisation of thermo mechanical
185  integer :: iter_visco              ! number of iterations for ssa viscosity
186
187  if (ispinup.le.1) shelf_vitbil = .false.         ! general case, ice shelves velocities are computed by diagnoshelf
188
189  timemax=time
190  if (itracebug.eq.1) write(num_tracebug,*) 'entree dans step_thermomeca', &
191       '   nt=',nt,'iter_beta=',iter_beta
192
193  !--------------------------------------------------------------------
194  !     lecture dans un petit fichier "runname.limit" au cas ou les
195  !     parametres du run seraient changes 
196  !                                                   
197  !--------------------------------------------------------------------
198
199!  call limit_file(nt,real(time),dt,tend,dtsortie,dtcpt,testdiag,dtt,runname)
200
201  !if (time.lt.tgrounded) then
202  !   marine=.false.
203  !else
204  marine=.true.
205  !endif
206
207  ! update the regions (floating, ice ...)
208  !-------------------------------------------
209  call masque()
210
211  if (iter_beta.gt.0)   call init_dragging
212
213
214  !-----------------------------------------------------------------------!
215  !   debut d'un bloc appels tous les dtt (pas de temps long- asynchrone) !
216  !-------------   -------------------------------------------------------!
217
218
219
220  isync_1: if (ISYNCHRO.eq.1) then          ! debut bloc appels tous les dtt   
221
222
223     ! climatic forcing
224     !=====================
225     call forclim                             
226     call ablation
227
228
229
230     ! SIA dynamic velocities  (deformation + loi de glissement )
231     !========================
232
233     ! remark : velocity is 3D
234     !          to update the velocity for ice temperature calculation
235     !            - surface mass balance for velocity (because ablation was called)
236     !            - impact of surface slope
237
238
239     spinup_1_veloc: if ((ispinup.ne.2).or.(nt.lt.nt_init_tm)) then        ! cas general
240
241        !---------------------------------------------------------------------------------
242        ! spinup_1 : calcul des vitesses dynamiques :
243        !------------
244        !            - est fait dans tous les cas pour les premieres boucles temporelles
245        !            - est fait pour ispinup=3 pour calibration forme de vitesse
246        !            - n'est pas fait pour ispinup=2 (test balance vel sur icethick)
247        !
248        !---------------------------------------------------------------------------------
249
250        if (itracebug.eq.1)  call tracebug(' Dans spinup_1_veloc : calculer SIA_velocities')
251
252        call Neffect()                                             
253        call SIA_velocities()
254
255     else    ! ispinup = 2
256
257        call force_balance_vel   ! fait l'ensemble du champ (flgz ou pas)
258
259     endif spinup_1_veloc
260
261
262     ! ice temperature
263     !===================================================================
264
265     ! update values in the structures Geom_g, Temperature_g, ...
266
267!afq     calc_temp:     if ((nt.gt.2).and.(geoplace(1:5).ne.'mism3')) then
268     calc_temp:     if (geoplace(1:5).ne.'mism3') then
269        if (itracebug.eq.1)  call tracebug('avant appel icetemp')
270        call icetemp
271        ! subroutines pour le calcul de la fusion basale
272        call bmeltshelf
273        call bmelt_grounded
274
275     endif calc_temp
276
277
278     ! flowlaw
279     !===================================================================
280
281     spinup_2_flowlaw: if ((ispinup.ne.2).or.(nt.lt.nt_init_tm)) then    ! cas general
282
283        !---------------------------------------------------------------------------------
284        ! spinup_2_flowlaw
285        !----------------
286        ! si ispinup =   2 , on ne fait pas le calcul de flowlaw, ni de diffusiv
287        !                    puisque les vitesses de bilan sont imposees
288        !                                             (initialisee nt < nt_init_tm)
289        !                   
290        ! Flowlaw et diffusiv : servent soit directement (ispinup = 0,1) soit pour le
291        !                       rescaling (ispinup = 3)
292        !---------------------------------------------------------------------------------
293
294        if (itracebug.eq.1)  call tracebug(' Dans spinup_2_flowlaw')
295
296           call flow_general
297           do iglen=n1poly,n2poly
298              call flowlaw(iglen)
299           end do
300
301
302     end if spinup_2_flowlaw
303
304     if (ispinup.eq.3) then
305        call calc_coef_vitbil                    ! recalibre sa_mx et s2a_mx et ddbx pour toute la zone posee
306     end if                                      ! Dans le shelf flgz = Ubar = Vbil et udef=0
307
308     call tracer      ! aurel neem (attention, position a verifier)   
309
310  endif isync_1
311
312
313  ! fin du bloc tous les dtt isync_1
314  !------------------------------------
315
316
317  spinup_4_vitdyn: if  ((ispinup.le.1).or.(nt.lt.nt_init_tm)) then 
318
319     !---------------------------------------------------------------------------------
320     ! spinup_4_vitdyn
321     !----------------
322     ! Passe dans ce bloc que quand on veut calculer les vitesses dynamiques :
323     ! - initialisations (nt.lt.nt_init_tm)
324     ! - run standard  (ispinup=0) (y compris iter_beta.ne.0)
325     ! - spinup temperature avec vitesses dynamiques  (ispinup = 1)               
326     !---------------------------------------------------------------------------------
327
328
329
330     call Neffect()
331
332     call diffusiv()
333
334
335     iterbeta:     if (iter_beta.ne.0) then          ! only for iterations on beta
336
337        if (.not.shelf_vitbil) then    !  bloc si les vitesses shelves viennent de diagnoshelf
338
339           !  Vcol declare dans le module spinup_mod
340           !  ou dans le dragging si no_spinup
341
342           where (flotmx(:,:))
343              uxbar(:,:)=uxflgz(:,:)
344           elsewhere
345              uxbar(:,:)=Vcol_x(:,:)
346           end where
347
348           where (flotmy(:,:))
349              uybar(:,:)=uyflgz(:,:)
350           elsewhere
351              uybar(:,:)=Vcol_y(:,:)
352           end where
353        else
354           uxbar(:,:)=Vcol_x(:,:)
355           uybar(:,:)=Vcol_y(:,:)
356        end if
357     end if iterbeta
358
359
360     ! strain rate calculations
361     !===========================
362
363     call strain_rate()
364
365
366     if (iter_beta.ne.0) then               ! on n'a utilise les vitesses de bilan que pour strain_rate "Dirichlet ?"
367        uxbar(:,:)=uxflgz(:,:)+uxdef(:,:)
368        uybar(:,:)=uyflgz(:,:)+uydef(:,:)
369     end if
370
371
372
373     ! ice shelves and ice streams velocities : SSA
374     !==============================================
375
376     !-----------------------------------------------------------------------
377     !              debut d'un deuxieme bloc appels tous les dtt                         
378     !-------------   -------------------------------------------------------
379
380     isync_2:    if ((shelfy).and.(marine).and.(isynchro.eq.1)) then               
381
382        !les 3 lignes ci-dessous pour avoir les champs avant diagno (en cas de pb colineaire)
383        !call init_sortie_ncdf
384        !call sortie_ncdf_cat
385
386        if (icompteur.eq.1) then    ! reprise de tout le vecteur d'etat   
387           iter_visco = 1
388
389        else if ((icompteur.ne.1).and.(nt.le.1)) then  ! initialisation avec reprise partielle
390                                                       ! ou nulle du vecteur d' etat
391          iter_visco = 10
392        else
393          iter_visco = 2
394        end if
395 
396
397           do m=1,iter_visco 
398              call diagnoshelf
399              call mix_SIA_L1 
400              call strain_rate
401           end do
402
403! iterations pour affiner beta en fonction des vitesses de bilan
404 
405           if ((time.ge.time_iter).and.(time.le.time_iter+20.)) then    !iterations seulment pendant 15 ans
406              do m=1,nb_iter_vitbil
407                 call beta_iter_vitbil 
408                 call diagnoshelf
409                 call mix_SIA_L1 
410                 call strain_rate
411              end do
412           end if
413! fin des iterations beta_iter_vitbil
414
415        endif isync_2
416
417
418     call  mix_SIA_L1 
419
420  end if spinup_4_vitdyn
421 
422  ! isostasy
423  !=================
424
425  spinup_3_bed: if ((ISYNCHRO.eq.1).and.(ispinup.eq.0.or.nt.lt.nt_init_tm)) then
426  !---------------------------------------------------------------------------------
427  ! spinup_3_bed
428  !----------------
429  ! Ne passe dans ce bloc que quand on veut l'isostasie, donc si variations epaisseur
430  ! run standard  (ispinup=0) et pour les initialisations                 
431  !---------------------------------------------------------------------------------
432
433  if (itracebug.eq.1)  call tracebug(' Dans spinup_3_bed')
434
435  if ((nbed.eq.1).and.nt.ne.1.and.isynchro.eq.1.and.(mod(abs(TIME),50.).lt.dtmin)) then
436    call bedrock                 !  bedrock adjustment
437  endif
438
439  call lake_to_slv               ! pour traiter dynamiquement les shelves sur lac
440
441  end if spinup_3_bed
442   
443
444end subroutine step_thermomeca
Note: See TracBrowser for help on using the repository browser.