source: branches/GRISLIv3/SOURCES/steps_time_loop.f90 @ 443

Last change on this file since 443 was 433, checked in by dumas, 15 months ago

Use only in step_time_loop : module_choix is used only when it is required

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, only: ispinup,isynchro,timemax,time,flot,S,B,Bsoc,H,sealevel_2d,hmx,hmy,&
17       uxbar,uybar,hdot,dtmax,dtmin,iout
18  use runparam, only: itracebug,nt,num_tracebug
19  use equat_adv_diff_2D_vect, only: icethick3
20  use calving_frange, only: calving
21  use flottab_mod, only: flottab,determin_tache,determin_marais
22  use sorties_ncdf_grisli, only: iglob_ncdf,testsort_time_ncdf,sortie_ncdf_cat
23  use resolmeca_SIA_L1, only: mix_SIA_L1
24  use bilan_eau_mod, only: diff_H,Bm_dtt,bmelt_dtt,calv_dtt,ablbord_dtt,diff_H_2D,grline_dtt,&
25       bilan_eau
26  use bilan_flux_mod, only: bilan_flux_output
27  !  module_choix donne acces aux modules interchangeables
28  use module_choix, only: shortoutput
29
30  implicit none
31  integer :: nt_init = 0      ! number of loops for initialisation of ice thickness
32
33  if (itracebug.eq.1)  call tracebug('Entree dans step_time_loop ')
34  if (itracebug.eq.1) write(num_tracebug,*) 'nt = ',nt
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) 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)
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_2d(:,:)) !afq -- WARNING: est-ce qu'on veut vraiment mettre S a la valeur locale du niveau marin?
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=0d0,dim=1))
79     hmy(:,:)=0.5*(H(:,:)+eoshift(H(:,:),shift=-1,boundary=0d0,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).or.(nt.eq.1).and.(mod(abs(TIME),1.).lt.dtmin)) then
124  if ((nt.eq.1).or.((isynchro.eq.1).and.(mod(abs(TIME),1.).lt.dtmin))) then
125     call shortoutput
126     diff_H = 0.
127     Bm_dtt(:,:) = 0.
128     bmelt_dtt(:,:) = 0.
129     calv_dtt(:,:)=0.
130     ablbord_dtt(:,:)=0.
131     diff_H_2D(:,:)=0.
132     grline_dtt(:,:)=0.
133  endif
134
135  !   sortie compteur tous les dtcpt ans
136  !------------------------------------------------------------------
137  !iout == 1 sortie cptr
138  !iout == 2 sortie nc
139  if (itracebug.eq.1)  call tracebug('dans steps_time_loop avant call out_recovery ')
140
141  call out_recovery(iout)
142
143  ! end of outputs
144
145  !====================================================================
146  !
147  call step_thermomeca()
148  !
149  ! contient tout le reste de la dynamique et du forçage
150  !====================================================================
151
152end subroutine step_time_loop
153
154
155
156
157
158
159!---------------------------------------------------------------------
160! The subroutine  step_thermomeca() calls the following processes
161!
162! - climatic forcing and its impact on surface mass balance
163! - dynamic velocities (SIA)
164! - ice temperature
165
166subroutine step_thermomeca
167 
168  use module3d_phy, only: ispinup,isynchro,timemax,time,marine,iglen,shelfy,icompteur,&
169       test_iter_diag,inv_beta,dtmin
170  use runparam, only: itracebug,nt,num_tracebug
171  use geography, only: geoplace
172  use deformation_mod_2lois, only: n1poly,n2poly
173  use iso_declar, only: nbed,dt_iso
174  use sorties_ncdf_grisli, only:
175  use resolmeca_SIA_L1, only: mix_SIA_L1
176  use diagno_mod, only: diagnoshelf
177  use beta_iter_vitbil_mod, only: beta_iter_vitbil
178  use icetempmod, only: icetemp
179  use deformation_mod_2lois, only: flow_general,flowlaw
180  !  module_choix donne acces aux modules interchangeables
181  use module_choix, only: time_iter,time_iter_end,nb_iter_vitbil,bedrock,forclim,ablation,&
182       rsl,bmeltshelf,tracer,calc_coef_vitbil,force_balance_vel
183
184  implicit none
185
186  integer :: m
187  logical :: shelf_vitbil = .true.   ! si vrai on prend les vitesses de bilan pour le shelf
188  integer :: nt_init_tm = 0          ! number of loops for initialisation of thermo mechanical
189  integer :: iter_visco              ! number of iterations for ssa viscosity
190
191  if (ispinup.le.1) shelf_vitbil = .false.         ! general case, ice shelves velocities are computed by diagnoshelf
192
193  timemax=time
194  if (itracebug.eq.1) write(num_tracebug,*) 'entree dans step_thermomeca', &
195       '   nt=',nt
196
197  !--------------------------------------------------------------------
198  !     lecture dans un petit fichier "runname.limit" au cas ou les
199  !     parametres du run seraient changes 
200  !                                                   
201  !--------------------------------------------------------------------
202
203  !  call limit_file(nt,real(time),dt,tend,dtsortie,dtcpt,testdiag,dtt,runname)
204
205  !if (time.lt.tgrounded) then
206  !   marine=.false.
207  !else
208  marine=.true.
209  !endif
210
211  ! update the regions (floating, ice ...)
212  !-------------------------------------------
213  call masque
214
215
216  !-----------------------------------------------------------------------!
217  !   debut d'un bloc appels tous les dtt (pas de temps long- asynchrone) !
218  !-------------   -------------------------------------------------------!
219
220
221
222  isync_1: if (ISYNCHRO.eq.1) then          ! debut bloc appels tous les dtt   
223
224     ! climatic forcing
225     !=====================
226     call forclim                             
227     call ablation
228
229     ! because eustatic sea level is updated in forclim, we need to call the RSL here
230     call rsl         ! provide the map of the local sea level (sealevel_2d)
231
232
233
234     ! SIA dynamic velocities  (deformation + loi de glissement )
235     !========================
236
237     ! remark : velocity is 3D
238     !          to update the velocity for ice temperature calculation
239     !            - surface mass balance for velocity (because ablation was called)
240     !            - impact of surface slope
241
242
243     spinup_1_veloc: if ((ispinup.ne.2).or.(nt.lt.nt_init_tm)) then        ! cas general
244
245        !---------------------------------------------------------------------------------
246        ! spinup_1 : calcul des vitesses dynamiques :
247        !------------
248        !            - est fait dans tous les cas pour les premieres boucles temporelles
249        !            - est fait pour ispinup=3 pour calibration forme de vitesse
250        !            - n'est pas fait pour ispinup=2 (test balance vel sur icethick)
251        !
252        !---------------------------------------------------------------------------------
253
254        if (itracebug.eq.1)  call tracebug(' Dans spinup_1_veloc : calculer SIA_velocities')
255
256        call Neffect()                                             
257        call SIA_velocities()
258
259     else    ! ispinup = 2
260
261        call force_balance_vel   ! fait l'ensemble du champ (flgz ou pas)
262
263     endif spinup_1_veloc
264
265
266     ! ice temperature
267     !===================================================================
268
269     ! update values in the structures Geom_g, Temperature_g, ...
270
271     calc_temp:     if (geoplace(1:5).ne.'mism3') then
272        if (itracebug.eq.1)  call tracebug('avant appel icetemp')
273        call icetemp
274        ! subroutines pour le calcul de la fusion basale
275        call bmelt_grounded
276        call bmeltshelf
277
278     endif calc_temp
279
280
281     ! flowlaw
282     !===================================================================
283
284     spinup_2_flowlaw: if ((ispinup.ne.2).or.(nt.lt.nt_init_tm)) then    ! cas general
285
286        !---------------------------------------------------------------------------------
287        ! spinup_2_flowlaw
288        !----------------
289        ! si ispinup =   2 , on ne fait pas le calcul de flowlaw, ni de diffusiv
290        !                    puisque les vitesses de bilan sont imposees
291        !                                             (initialisee nt < nt_init_tm)
292        !                   
293        ! Flowlaw et diffusiv : servent soit directement (ispinup = 0,1) soit pour le
294        !                       rescaling (ispinup = 3)
295        !---------------------------------------------------------------------------------
296
297        if (itracebug.eq.1)  call tracebug(' Dans spinup_2_flowlaw')
298
299        call flow_general
300        do iglen=n1poly,n2poly
301           call flowlaw(iglen)
302        end do
303
304
305     end if spinup_2_flowlaw
306
307     if (ispinup.eq.3) then
308        call calc_coef_vitbil                    ! recalibre sa_mx et s2a_mx et ddbx pour toute la zone posee
309     end if                                      ! Dans le shelf flgz = Ubar = Vbil et udef=0
310
311     call tracer      ! aurel neem (attention, position a verifier)   
312
313  endif isync_1
314
315
316  ! fin du bloc tous les dtt isync_1
317  !------------------------------------
318
319
320  spinup_4_vitdyn: if  ((ispinup.le.1).or.(nt.lt.nt_init_tm)) then 
321
322     !---------------------------------------------------------------------------------
323     ! spinup_4_vitdyn
324     !----------------
325     ! Passe dans ce bloc que quand on veut calculer les vitesses dynamiques :
326     ! - initialisations (nt.lt.nt_init_tm)
327     ! - run standard  (ispinup=0)
328     ! - spinup temperature avec vitesses dynamiques  (ispinup = 1)               
329     !---------------------------------------------------------------------------------
330
331
332
333     call Neffect()
334
335     call diffusiv()
336
337     ! strain rate calculations
338     !===========================
339
340     call strain_rate()
341
342     ! ice shelves and ice streams velocities : SSA
343     !==============================================
344
345     !-----------------------------------------------------------------------
346     !              debut d'un deuxieme bloc appels tous les dtt                         
347     !-------------   -------------------------------------------------------
348
349     isync_2:    if ((shelfy).and.(marine).and.(isynchro.eq.1)) then               
350
351        !les 3 lignes ci-dessous pour avoir les champs avant diagno (en cas de pb colineaire)
352        !call init_sortie_ncdf
353        !call sortie_ncdf_cat
354
355        if (icompteur.eq.1) then    ! reprise de tout le vecteur d'etat   
356           iter_visco = 1
357
358        else if ((icompteur.ne.1).and.(nt.le.1)) then  ! initialisation avec reprise partielle
359           ! ou nulle du vecteur d' etat
360           iter_visco = 10
361        else
362           iter_visco = 2
363        end if
364
365        test_iter_diag=0.
366        do m=1,iter_visco             
367           call diagnoshelf
368           call mix_SIA_L1 
369           call strain_rate
370           if (test_iter_diag.lt.1.e-2) exit
371        end do
372
373        ! iterations pour inversion du beta en fonction des vitesses de bilan
374        if ((inv_beta.eq.1).and.(time.ge.time_iter).and.(time.le.time_iter_end)) then  ! iterations pendant 20 ans
375           do m=1,nb_iter_vitbil
376              call beta_iter_vitbil(m)
377              call diagnoshelf
378              call mix_SIA_L1
379              call strain_rate
380           end do
381        end if
382        ! fin des iterations beta_iter_vitbil
383
384     endif isync_2
385
386     call  mix_SIA_L1 
387
388  end if spinup_4_vitdyn
389
390
391  ! isostasy
392  !=================
393
394  spinup_3_bed: if ((ISYNCHRO.eq.1).and.(ispinup.eq.0.or.nt.lt.nt_init_tm)) then
395     !---------------------------------------------------------------------------------
396     ! spinup_3_bed
397     !----------------
398     ! Ne passe dans ce bloc que quand on veut l'isostasie, donc si variations epaisseur
399     ! run standard  (ispinup=0) et pour les initialisations                 
400     !---------------------------------------------------------------------------------
401
402     if (itracebug.eq.1)  call tracebug(' Dans spinup_3_bed')
403
404     if ((nbed.eq.1).and.nt.ne.1.and.isynchro.eq.1.and.(mod(abs(TIME),dt_iso).lt.dtmin)) then
405        call bedrock                 !  bedrock adjustment
406     endif
407
408
409  end if spinup_3_bed
410
411
412end subroutine step_thermomeca
Note: See TracBrowser for help on using the repository browser.