source: trunk/SOURCES/steps_time_loop.f90

Last change on this file was 358, checked in by dumas, 2 years ago

For grounding line points : call bmelt_grounded before bmeltshelf in steps_time_loop.f90

File size: 13.0 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 sorties_ncdf_grisli
21  use flottab_mod
22  use icetempmod
23  use diagno_mod
24  use bilan_eau_mod
25  use bilan_flux_mod
26  !  use track_debug
27
28  implicit none
29  integer :: nt_init = 0      ! number of loops for initialisation of ice thickness
30
31  if (itracebug.eq.1)  call tracebug('Entree dans step_time_loop ')
32  if (itracebug.eq.1) write(num_tracebug,*) 'nt = ',nt
33
34
35  ! ice thickness evolution
36  !=============================
37
38
39  spinup_icethick: if (ispinup.eq.0.or.ispinup.eq.2.or.nt.lt.nt_init) then 
40
41     !-------------------------------------------------------------------
42     ! spinup_icethick : Le calcul changement d'epaisseur (icethick)
43     !----------------
44     !   - est fait dans tous les cas pour les premieres boucles temporelles
45     !   - pour ispinup = 0 (run standard)
46     !   - pour ispinup = 2 (force avec vitesses de bilan)
47     !
48     ! else
49     !   juste un increment pour le temps
50     !    -  pour ispinup = 1 et  ispinup = 3  temperature calculation
51     
52     !-------------------------------------------------------------------
53
54     if (itracebug.eq.1)  call tracebug(' Dans spinup_icethick : appelle icethick')
55
56     call icethick3
57     call flottab
58     call determin_tache
59     call calving
60     call ablation_bord
61     call bilan_eau
62     call bilan_flux_output     
63     call flottab
64     if (isynchro.eq.1.or.nt.eq.1) call determin_marais
65     if (itracebug.eq.1)  call tracebug('apres marais')
66
67     ! new calculation of ice surface elevation
68     where(.not.flot(:,:))
69        S(:,:)=Bsoc(:,:)+H(:,:)
70        S(:,:)=max(S(:,:),sealevel_2d(:,:)) !afq -- WARNING: est-ce qu'on veut vraiment mettre S a la valeur locale du niveau marin?
71        B(:,:)=Bsoc(:,:)
72     end where
73
74     ! calcul de Hmx et Hmy -> shift=-1, dim=1 -> H(i-1,j)
75
76     hmx(:,:)=0.5*(H(:,:)+eoshift(H(:,:),shift=-1,boundary=0d0,dim=1))
77     hmy(:,:)=0.5*(H(:,:)+eoshift(H(:,:),shift=-1,boundary=0d0,dim=2))
78     hmx(1,:)=0.
79     hmy(:,1)=0.
80
81     ! mise a 0 des vitesses quand epaisseur nulle
82     where (hmx(:,:).le.1.)
83        uxbar(:,:)=0.
84     endwhere
85     where (hmy(:,:).le.1.)
86        uybar(:,:)=0.
87     endwhere
88
89  else   ! time increment when no icethick : ispinup= 1 or 3    ! loop spinup_icethick
90
91     hdot(:,:)=0.
92     time=nt*dtmax
93
94     if (itracebug.eq.1) write(num_tracebug,*) ' Dans spinup_icethick:  no icethick ',ispinup, &
95          '    time=',time
96
97
98  end if spinup_icethick
99  ! end of ice thickness evolution
100  !====================================================================
101
102
103
104  ! outputs
105  !====================================================================
106
107  !   vertical plan snapshots (profiles)
108  !------------------------------------------------------------------
109
110  !if ((NT.eq.1) &    !   .or.(NT.eq.2)
111  !     .or.(mod(abs(dble(TIME)),dble(DTPROFILE)).lt.dble(dtmin)) &
112  !     .or.(abs(TIME-TEND).lt.dtmin)) then
113  !  call sortieprofile()
114  !endif
115
116  !  sorties netcdf Hassine (aout 2010)  (2D and 3D)
117  !------------------------------------------------------------------
118
119  call testsort_time_ncdf(time)
120  if (iglob_ncdf .EQ. 1) call sortie_ncdf_cat
121!  if ((isynchro.eq.1).or.(nt.eq.1).and.(mod(abs(TIME),1.).lt.dtmin)) then
122  if ((nt.eq.1).or.((isynchro.eq.1).and.(mod(abs(TIME),1.).lt.dtmin))) then
123    call shortoutput
124    diff_H = 0.
125    Bm_dtt(:,:) = 0.
126    bmelt_dtt(:,:) = 0.
127    calv_dtt(:,:)=0.
128    ablbord_dtt(:,:)=0.
129    diff_H_2D(:,:)=0.
130    grline_dtt(:,:)=0.
131  endif
132
133  !   sortie compteur tous les dtcpt ans
134  !------------------------------------------------------------------
135  !iout == 1 sortie cptr
136  !iout == 2 sortie nc
137  if (itracebug.eq.1)  call tracebug('dans steps_time_loop avant call out_recovery ')
138
139  call out_recovery(iout)
140 
141  ! end of outputs
142
143  !====================================================================
144  !
145  call step_thermomeca()
146  !
147  ! contient tout le reste de la dynamique et du forçage
148  !====================================================================
149
150end subroutine step_time_loop
151
152
153
154
155
156
157!---------------------------------------------------------------------
158! The subroutine  step_thermomeca() calls the following processes
159!
160! - climatic forcing and its impact on surface mass balance
161! - dynamic velocities (SIA)
162! - ice temperature
163
164subroutine step_thermomeca()
165
166
167  use module3d_phy
168  use module_choix ! module de choix du type de run
169  !  module_choix donne acces a tous les modules
170  !  de declaration des packages
171  use icetempmod
172  use sorties_ncdf_grisli
173  use flottab_mod
174  use diagno_mod
175  use resolmeca_SIA_L1
176  use iso_declar, only: dt_iso
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
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
212  !-----------------------------------------------------------------------!
213  !   debut d'un bloc appels tous les dtt (pas de temps long- asynchrone) !
214  !-------------   -------------------------------------------------------!
215
216
217
218  isync_1: if (ISYNCHRO.eq.1) then          ! debut bloc appels tous les dtt   
219
220     ! climatic forcing
221     !=====================
222     call forclim                             
223     call ablation
224
225     ! because eustatic sea level is updated in forclim, we need to call the RSL here
226     call rsl         ! provide the map of the local sea level (sealevel_2d)
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 bmelt_grounded
273        call bmeltshelf
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)
325     ! - spinup temperature avec vitesses dynamiques  (ispinup = 1)               
326     !---------------------------------------------------------------------------------
327
328
329
330     call Neffect()
331
332     call diffusiv()
333
334     ! strain rate calculations
335     !===========================
336
337     call strain_rate()
338
339     ! ice shelves and ice streams velocities : SSA
340     !==============================================
341
342     !-----------------------------------------------------------------------
343     !              debut d'un deuxieme bloc appels tous les dtt                         
344     !-------------   -------------------------------------------------------
345
346     isync_2:    if ((shelfy).and.(marine).and.(isynchro.eq.1)) then               
347
348        !les 3 lignes ci-dessous pour avoir les champs avant diagno (en cas de pb colineaire)
349        !call init_sortie_ncdf
350        !call sortie_ncdf_cat
351
352        if (icompteur.eq.1) then    ! reprise de tout le vecteur d'etat   
353           iter_visco = 1
354
355        else if ((icompteur.ne.1).and.(nt.le.1)) then  ! initialisation avec reprise partielle
356                                                       ! ou nulle du vecteur d' etat
357          iter_visco = 10
358        else
359          iter_visco = 2
360        end if
361 
362           test_iter_diag=0.
363           do m=1,iter_visco             
364              call diagnoshelf
365              call mix_SIA_L1 
366              call strain_rate
367              if (test_iter_diag.lt.1.e-2) exit
368           end do
369           
370           ! iterations pour inversion du beta en fonction des vitesses de bilan
371           if ((inv_beta.eq.1).and.(time.ge.time_iter).and.(time.le.time_iter_end)) then  ! iterations pendant 20 ans
372              do m=1,nb_iter_vitbil
373                 call beta_iter_vitbil(m)
374                 call diagnoshelf
375                 call mix_SIA_L1
376                 call strain_rate
377              end do
378           end if
379           ! fin des iterations beta_iter_vitbil
380
381        endif isync_2
382
383     call  mix_SIA_L1 
384
385  end if spinup_4_vitdyn
386 
387 
388  ! isostasy
389  !=================
390
391  spinup_3_bed: if ((ISYNCHRO.eq.1).and.(ispinup.eq.0.or.nt.lt.nt_init_tm)) then
392  !---------------------------------------------------------------------------------
393  ! spinup_3_bed
394  !----------------
395  ! Ne passe dans ce bloc que quand on veut l'isostasie, donc si variations epaisseur
396  ! run standard  (ispinup=0) et pour les initialisations                 
397  !---------------------------------------------------------------------------------
398
399  if (itracebug.eq.1)  call tracebug(' Dans spinup_3_bed')
400
401  if ((nbed.eq.1).and.nt.ne.1.and.isynchro.eq.1.and.(mod(abs(TIME),dt_iso).lt.dtmin)) then
402    call bedrock                 !  bedrock adjustment
403  endif
404
405
406  end if spinup_3_bed
407 
408
409end subroutine step_thermomeca
Note: See TracBrowser for help on using the repository browser.