source: branches/iLoveclim/SOURCES/main3D-0.4-40km.f90 @ 36

Last change on this file since 36 was 36, checked in by roche, 8 years ago

Tweaked the makefile and the main for library generation with gfortran.

File size: 10.7 KB
Line 
1!     **********************************************************************
2!     *       GRISLI      Grenoble Ice Shelves-Land Ice
3!     **********************************************************************
4
5
6!     Ont participe a l'ecriture de ce modele :
7!
8!                   Catherine Ritz                           (tout du long)
9!                   Adeline Fabre                      (la partie Gremlins)
10!                   Vincent Rommelaere         (ice shelves et ice streams)
11!                   Christophe Dumas (debut f90,              (Antarctique)
12!                   Vincent Peyaud      (portage HN,calving, front, hydrol)
13!                   Cyril Mazauric                                  (AGRIF)
14!                   Hassine Baya             (netcdf, doxygen, icetemp,...)
15!
16!     catritz@lgge.obs.ujf-grenoble.fr
17!
18!     **********************************************************************
19
20
21
22!> \mainpage GRISLI Modele 3D De Calotte Glaciaire
23!!
24!! \section start Pour commencer
25!! Le programme principal est dans le module main3D.
26!! Ce module est dans le fichier  main3D-0.4-40km.f90.
27!!
28!! \section tree Arbre d'appel
29!!
30!! - call grisli_init()
31!!   - step_grisli()
32!!   - sortie_ncdf_cat()
33!!   - testsort_time_ncdf()
34!!   - initial()
35!!   - sortie_hz_multi()
36!! - call step_grisli1()
37!!
38!<
39
40
41!> \file main3D-0.4-40km.f90 GRISLI Modele 3D De Calotte Glaciaire
42!! programme principal
43!! (voir l'\ref tree)
44!!
45!! @brief   modele flow line d'evolution de calotte
46!! @authors         Catherine Ritz     catritz@lgge.obs.ujf-grenoble.fr  (tout du long)
47!! @authors         Adeline Fabre                      (la partie Gremlins)
48!! @authors         Vincent Rommelaere         (ice shelves et ice streams)
49!! @authors         Christophe Dumas (debut f90,              (Antarctique)
50!! @authors         Vincent Peyaud      (portage HN,calving, front, hydrol)
51!! @authors         Cyril Mazauric                                  (AGRIF)
52!!
53!!
54!! @note use module3D_phy
55!! @note use module_choix
56!! @note use flottab_mod
57!! @note use icetempmod
58!! @note use sorties_ncdf_grisli
59!! @note use diagno_mod
60!! @note use resolmeca_SIA_L1
61!!
62!!
63!! Ce module appelle les routines suivantes :
64!! - grisli_init()
65!! - step_grisli1()
66!! - step_output()
67!!
68!<
69
70!> \namespace main3D GRISLI Modele 3D De Calotte Glaciaire
71!! programme principal
72!! (voir l'\ref tree)
73!!
74!!
75!! @brief   modele flow line d'evolution de calotte
76!! @authors         Catherine Ritz     catritz@lgge.obs.ujf-grenoble.fr  (tout du long)
77!! @authors         Adeline Fabre                      (la partie Gremlins)
78!! @authors         Vincent Rommelaere         (ice shelves et ice streams)
79!! @authors         Christophe Dumas (debut f90,              (Antarctique)
80!! @authors         Vincent Peyaud      (portage HN,calving, front, hydrol)
81!! @authors         Cyril Mazauric                                  (AGRIF)
82!!
83!!
84!! @note use module3D_phy
85!! @note use module_choix
86!! @note use flottab_mod
87!! @note use icetempmod
88!! @note use sorties_ncdf_grisli
89!! @note use diagno_mod
90!! @note use resolmeca_SIA_L1
91!!
92!!
93!! @todo itracebug : faire une routine
94!!
95!! Ce module appelle les routines suivantes :
96!! - grisli_init()
97!! - step_grisli1()
98!!
99!! Defined in file main3D-0.4-40km.f90
100!<
101
102program main3D
103
104! dmr --- void, just added for compilation
105
106end program main3D
107
108
109subroutine ISM_NORD(timCplGRISyr)
110
111  USE module3D_phy
112  USE module_choix !   module de choix du type de run
113  !  module_choix donne acces a tous les modules
114  !  de declaration des packages
115  use flottab_mod
116  use icetempmod
117  use sorties_ncdf_grisli
118  use diagno_mod 
119  use resolmeca_SIA_L1
120!  use track_debug
121!dcdmr --- GRISLI - LOVECLIM
122!  use input_timerCplGRIS
123!dcdmr --- GRISLI - LOVECLIM
124
125  implicit none
126
127  integer, intent(in) :: timCplGRISyr
128
129!dcdmr --- GRISLI - LOVECLIM
130!cdc appel d'initial au premier passage dans grisli
131! pour demarrer avec climat initialise
132  if (time.eq.tbegin) call grisli_init  ! Initializations
133! mab: timCplGRISyr corresponds to a certain amout of DAYS 
134  TEND = real(TIME) + real(timCplGRISyr)
135
136  PRINT*,'******* Appel de GRISLI Nord *******'
137  PRINT*,'TIME = ',TIME,' TEND = ',TEND
138
139  time_loop: DO WHILE (time.LT.tend)  !____________________________ debut boucle temporelle
140     call step_time_loop
141     nt= nt+1   !cdc ajoute pour incrementer nt
142     IF (nt.gt.ntmax) exit
143  end do time_loop
144
145  if (itracebug.eq.1)  call tracebug('dans main avant call out_recovery ')
146  call out_recovery(iout)
147
148!  write(6,*) "end of the run at time = ",time
149!  write(6,*) "_____________________________________________________________________"
150
151end subroutine ISM_NORD
152
153
154!---------------------------------------------------------------------------------------
155subroutine grisli_init
156
157  USE module3D_phy
158  USE module_choix ! module de choix du type de run
159  !  module_choix donne acces a tous les modules
160  !  de declaration des packages
161  use flottab_mod
162  use icetempmod
163  use sorties_ncdf_grisli
164  use util_recovery
165  use diagno_mod 
166!  use track_debug
167
168  implicit none
169
170  if (itracebug.eq.1)  call tracebug(' Entree dans routine grisli_init')
171  !      switch pour passer ou non par T lliboutry calcule => 0, ne passe pas,
172  !      1 ou 2 passe (se met a 0 tout seul si on prend un fichier .cptr)
173
174  ITEMP=0
175
176  !      switch couple physique faible =>  CP et CT independant T
177  !               0     pas de trait. vert.  A FAIRE           niveau L0
178  !               1     pas de couplage , faible physique      niveau L1
179  !               2     couplage, faible physique              niveau L2
180  !               3     couplage, physique complete sans CBT   niveau L3
181  !               4     idem 3 mais loi de def. Duval          niveau L4
182  ICOUPLE=4
183  !     switch margin IMARGIN=0 fixed, IMARGIN=1 moving
184  IMARGIN=1
185
186  TIMECG=TBEGIN
187  nt=-1   ! utilisee dans initialisation flottab
188  !     sortie profile tous les dtprofile
189  DTPROFILE=50000.
190  marine=.true.
191  !     ----------------------------------fin des modifs run les plus usuelles
192  ! DIRNAMEOUT='./'
193  DIRNAMEOUT='outputdata/ism/'
194
195  call initial  ! routine qui appel toutes les routines d'initialisation
196
197
198  !      call init_sortie_ncdf
199  !      call sortie_ncdf_cat
200  !      STOP
201
202  !     compteur tous les DTCPT
203  DTCPT=dtout
204
205  dtsortie=minval(dtsortie_hz(:))
206
207  !     ************ OPEN FILES.RITZ ****************
208
209  if ((geoplace.eq.'anteis1').or.(geoplace.eq.'ant20km')) then
210     ! fichier de reference pour le niveau des mers
211     open(num_sealevel,file=TRIM(DIRNAMEOUT)//'sealevel'//runname//'.ritz',position="append")
212     open(num_ts_ritz,file=TRIM(DIRNAMEOUT)//'ts_'//runname//'.ritz',position="append")
213     open(num_ic_vo,file=TRIM(DIRNAMEOUT)//'ic_'//runname//'vo.ritz',position="append")
214     open(num_ic_by,file=TRIM(DIRNAMEOUT)//'ic_'//runname//'by.ritz',position="append")
215     open(num_ic_dm,file=TRIM(DIRNAMEOUT)//'ic_'//runname//'dm.ritz',position="append")
216     open(num_ic_dc,file=TRIM(DIRNAMEOUT)//'ic_'//runname//'dc.ritz',position="append")
217     open(num_ic_df,file=TRIM(DIRNAMEOUT)//'ic_'//runname//'df.ritz',position="append")
218  endif
219
220  !------------------------------ INITIALISATION ----------------------------
221  !
222! ecriture netcdf apres initialisation
223
224
225
226  call testsort_time_ncdf(dble(tbegin))
227  if (iglob_ncdf .EQ. 1) call sortie_ncdf_cat
228
229
230
231  if (iter_beta.eq.0) then
232
233     if (itracebug.eq.1)  call tracebug(' Avant appel routine icethick3')
234     call icethick3
235     debug_3D(:,:,88) = S(:,:)
236     if (itracebug.eq.1)  call tracebug(' Apres appel routine icethick3')
237  end if
238
239
240  !     Tgrounded, temps pendant lequel la calotte est terrestre
241  tgrounded=tbegin-10.
242  !if (tgrounded.le.tbegin) then
243  marine=.true. ! Cas la calotte est terrestre
244  !end if
245
246  ! test vincent car certains H(i,j)=0 dans fichier de reprise
247  do j=1,ny
248     do i=1,nx
249        H(i,j)=max(1.,H(i,j))
250     enddo
251  enddo
252
253
254  ! call firstoutput()           ! ouverture fichier temporel et premieres ecritures
255
256  call forclim                   !  initialisation BM et TS         
257  call ablation
258
259
260
261  !  -----------                  CALCULATION OF INITIAL TEMPERATURES
262
263  tcpt:if (ICOMPTEUR.eq.0) then
264
265
266     if ((GEOPLACE.ne.'eismint').and.(GEOPLACE(1:6).ne.'marine')) then
267        !       ITEMP=1 => calcul de T lliboutry; ITEMP=2 => reprise d'un fichier cptr
268        !       ITEMP=0 => on ne prend pas en compte T Lliboutry
269        !       ITEMP=3 => on prend les temperatures d'un fichier cptr
270
271
272
273        if ((ITEMP.eq.0).or.(ITEMP.eq.3)) then
274           call masque()
275
276           call Neffect()
277
278           call flottab()
279
280           call Neffect()
281
282
283           !          call vitbilan_lect   ! routine de lecture des vitesses de bilan
284           !       ========================================================
285
286           if (ITEMP.eq.0) call lineartemp()
287
288           call bmelt_grounded 
289           call  bmeltshelf
290
291
292           call flow_general
293
294           do iglen=n1poly,n2poly
295              call flowlaw(iglen)
296           end do
297
298
299           call Neffect()
300           call flottab()
301
302
303           call Neffect()
304
305
306           call diffusiv()
307           call SIA_velocities()
308
309        endif
310
311
312     endif
313     !     fin du test geoplace
314
315  else  ! tcpt     on reprend un fichier compteur (ICOMPTEUR.eq.1)
316
317     time=tbegin       ! prend le temps du compteur
318
319
320     call masque()
321     call flottab()
322     call neffect()
323     call flottab()
324     call masque()
325
326     do i=1,nx
327        do j=1,ny
328           if (S(i,j).lt.0) then
329              print*,i,j,S(i,j)
330              goto 11115
331           endif
332        enddo
333     enddo
33411115 continue
335
336     !       ========================================================
337     call flow_general
338
339     do iglen=n1poly,n2poly
340        call flowlaw(iglen)
341     end do
342
343
344     call Neffect()
345     call flottab
346     call diffusiv()
347     call SIA_velocities()
348     call strain_rate
349
350  endif tcpt
351  !     fin du test sur icompteur
352
353  !      call init_sortie_ncdf
354  !      call sortie_ncdf_cat
355
356  call flottab()
357  call Neffect()
358  call flottab()
359
360  if (icompteur.eq.0) then
361     do i=1,nx
362        do j=1,ny
363           if (.not.flot(i,j)) then
364              B(i,j) = Bsoc(i,j)
365              Uxbar(i,j) = 0.
366              Uybar(i,j) = 0.
367           end if
368        end do
369     end do
370  endif
371
372  boost = .false.
373
374  do i=2,nx-1
375     do j=2,ny-1
376        hwater(i,j)=max(hwater(i,j),0.)
377     enddo
378  enddo
379  timemax=time
380
381  ! sortie des champs initiaux
382
383  call testsort_time(tbegin)
384  call sortie_hz_multi   ! pour les variables declarees dans 3D
385!  call hz_output(tbegin)
386
387!  call limit_file(1,real(time),dt,tend,dtsortie,dtcpt,testdiag,dtt,runname)
388
389  isynchro=1
390  ndebug=0
391  ndebug_max=9
392
393  call step_thermomeca()     ! un tour dans la boucle temporelle, partie avant icethick
394  call init_sortie_ncdf
395  if (itracebug.eq.1)  call tracebug(' fin routine grisli_init')
396  call testsort_time_ncdf(dble(tbegin))
397
398  if (iglob_ncdf .EQ. 1) call sortie_ncdf_cat
399
400  return
401end subroutine grisli_init
Note: See TracBrowser for help on using the repository browser.