source: branches/iLoveclim/SOURCES/steps_time_loop_avec_iterbeta.f90 @ 186

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

Merged branch iLOVECLIM to trunk at rev 121

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