source: trunk/SOURCES/tracer_mod.f90 @ 4

Last change on this file since 4 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

File size: 15.9 KB
Line 
1!> module appele par le main (step) qui calcule les traceurs
2!!  @todo depuis tracer_mod, decider e() en local ou global
3!!  @todo attention l'appel au module climat est en dur
4!!    et peut poser des pb si different du module choix.
5
6!--------------------------------------------------------------
7! Position des forages dans la grille grisli.
8!
9! Ordre GMT :
10! mapproject  -Js-39/90/71/1:50000000 -Fk -R/-150/50/50/85 -C
11!
12! NEEM : -50.9E 77.5N  -> -281.129 / -1334.05
13!  5km:  ineem=143 /jneem=437
14!  15km: ineem=48 / jneem=146
15!  45km: ineem=16 / jneem=49
16!
17! GRIP : -37.633333E 72.583333N -> 45.4707 / -1905.94
18!  5km:  igrip=208 /jgrip=323
19!  15km: igrip=70 / jgrip=108 ! grosso modo
20!  45km: igrip=24 / jgrip=36  ! grosso modo
21!
22!--------------------------------------------------------------
23
24
25
26module tracer_mod
27
28  use module3d_phy
29  use tracer_vars             ! module de declaration de variables pour les traceurs
30  use climat_perturb_mois_mod ! on en a besoin, notamment pour NFT
31
32  implicit none
33
34! variables de travail
35  integer :: i_ij, i_jj, ip, jp, im, jm, km      ! indice dans les matrices des traceurs
36  real    :: v_xij, v_yij, v_zij                 ! vitesses 'cumulees'
37
38! pour chaque point de grille on recherche ou etait la particule au pas de temps precedent
39! Tout ce qui est mentionne traceur est lie aux proprietes de ce point
40
41  real    :: x_tra                               ! position x du traceur 
42  real    :: y_tra                               ! position y du traceur
43  real    :: e_tra                               ! position verticale du traceur
44
45! ecart entre la position du traceur et le point de grille grisli le plus proche
46  real    :: fx                                  ! ecart en x
47  real    :: fy                                  ! ecart en y
48  real    :: fz                                  ! ecart vertical ? (z ou sigma ?)
49
50  integer :: deltracer   ! pas de temps du traceur ?  -> dans la version d'origine tracer est appele tous les dtt
51  integer :: time_tra    ! temps traceur
52
53!===============================================================
54
55 
56contains
57
58subroutine init_tracer
59
60! formats pour les ecritures dans 42
61428 format(A)
62
63  real :: lambdab  ! fusion basale 'artificielle' pour age des couches
64  real :: agemax   ! age maximum de la glace (pour init_tracer)
65 
66  integer :: kk   ! indice vertical pour definition de E, mais on veut conserver la valeur de k
67!  real,dimension(nz) :: E   ! vertical coordinate in ice, scaled to H zeta
68
69
70if (itracebug.eq.1) write(num_tracebug,*)'init_tracer dans tracer_mod'
71
72! itracer pour les ecritures et lectures recovery
73 itracer = 1
74
75! profondeur reduite
76  E(1)=0.
77  E(NZ)=1.
78  do KK=1,NZ
79     if ((KK.ne.1).and.(KK.ne.NZ)) E(KK)=(KK-1.)/(NZ-1.)
80  end do
81
82
83if (itracebug.eq.1) write(num_tracebug,*)'init_tracer dans tracer_mod'
84
85  time_tra = floor(time) 
86
87  deltracer = dtt           ! sorte de pas de temps du tracer -  multiple de dtt !
88
89
90! on lit les parametres par namelist dorenavant
91  namelist/tracer/file_tr_dat,file_tr_out,file_tr_dep,lambdab,agemax !aurel : a ne pas oublier paramlist
92
93  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
94  read(num_param,tracer)
95  print*
96  print*
97  print*, 'tracer .dat file: ', file_tr_dat
98  print*, 'tracer .out file: ', file_tr_out
99  print*, 'tracer dep. file: ', file_tr_dep
100
101
102  write(num_rep_42,428)'!___________________________________________________________' 
103  write(num_rep_42,428) '&tracer                                    ! module tracer_mod'
104  write(num_rep_42,'(A,A)')   'tracer .dat file : ', file_tr_dat
105  write(num_rep_42,'(A,A)')   'tracer .out file : ', file_tr_out
106  write(num_rep_42,'(A,A)')   'tracer dep. file : ', file_tr_dep
107  write(num_rep_42,*)'/'                     
108  write(num_rep_42,428) 'declaration des repertoires *traceurs*'
109  write(num_rep_42,*)
110
111! on appelle la soubroutine qui lit le fichier .dat des traceurs
112! ce fichier contient les dates de sorties des traceurs
113  call read_tr_dat ! intervalles d'output
114
115  nij = nx * ny
116
117! initialize "previous" tracer arrays as if all ice at start
118! time was created instantaneously and in-place
119
120! attention pourquoi ne pas utiliser xcc et ycc ?
121!!$  do i=1,nx
122!!$     xgrid(i) = (i- (nx+1)/2. ) * dx
123!!$  end do
124!!$  do j=1,ny
125!!$     ygrid(j) = (j - (ny+1)/2.) * dy
126!!$  end do
127
128! notre grille a nous
129  do i=1,nx
130     xgrid(i)=xmin*1000+(i-1)*dx
131  end do
132  do j=1,ny
133     ygrid(j)=ymin*1000+(j-1)*dy
134  end do
135 
136! initialisation des abscisses et ordonnees a la valeur en surface
137! initialisation de l'age
138  do j=2,ny-1
139     i_jj = j-1
140     ydepk(:,:,i_jj) = ygrid(j) 
141     do i=2,nx-1
142        i_ij = i-1
143        xdepk(:,i_ij,i_jj) = xgrid(i)
144
145! todo0710
146! remplacer par une formule logarithmique
147! 10**-4 m/an fusion basale (ou valeur de lambda a la base)
148! + agemax 1 million d'annees
149
150!!$        if (h(i,j)>10.) then
151!!$           do k=1,4
152!!$              tdepk(k,i_ij,i_jj) = - min(3500., h(i,j)) *real((k-1)**2)
153!!$           end do
154!!$           do k=5,nz
155!!$              tdepk(k,i_ij,i_jj) = - min(3500., h(i,j)) *real((k-1)*3)
156!!$           end do
157!!$           do k=12,nz
158!!$              tdepk(k,i_ij,i_jj) = tdepk(k,i_ij,i_jj) * (real(11) + real((k-7)**3)/125. )/12.
159!!$           end do
160!!$        else
161!!$           tdepk(:,i_ij,i_jj) = 0.
162!!$        end if
163!!$     end do
164!!$  end do
165!!$
166!!$   
167!!$! a enlever eventuellement si on utilise le log
168!!$  do j=1,ny-2
169!!$     do i=1,nx-2
170!!$        if (((tdepk(20,i,j)-tdepk(21,i,j))<1000.)) then
171!!$           if ((h(i,j)>100.))  then
172!!$              tdepk(21,i,j)=tdepk(20,i,j)-10*h(i,j)
173!!$           else
174!!$              tdepk(21,i,j)=tdepk(20,i,j)- 1000.
175!!$           end if
176!!$        end if
177!!$     end do
178!!$  end do
179!!$   
180        do k=1,nz
181           if(acc(i,j).gt.lambdab) then
182              tdepk(k,i_ij,i_jj)=H(i,j)*log(    &
183                   1-E(k)*(1-lambdab/ACC(i,j)))/(ACC(i,j)-lambdab)
184           else
185              tdepk(k,i_ij,i_jj)=0.
186           end if
187        end do
188
189     end do
190  end do
191
192  where (tdepk(:,:,:).ge.agemax)
193     tdepk(:,:,:)=agemax
194  end where
195  tdepk(:,:,:) = tdepk(:,:,:) + tbegin
196
197
198 2009 format(21(f12.1))
199   
200   ! prepares accu-cumul variable for accu-interp
201
202  ! note aurel neem : ATTENTION ! dans version NL ->
203  ! chaque indice d'accucumul vaut 100 ans.
204  ! 0 pour le present, 10 pour 1000 ans
205   allocate(accucumul(0:nft-1)); accucumul = 0.0
206   ! watch out its size !!!
207   ! note: tpert(1) = tpert at -422700yr
208   ! nl 12/08/02: attention au param de l'exp(): doit etre coherent avec le
209   ! forcage climatique: ici, change de 0.062 a 0.100
210   ! attention a ce que nft <=> present, sinon mon indexage est fourré
211   ! a la rigueur je pourrais modifier le fichier d'interpolation pour
212   ! l'adapter a tdate
213   accucumul(nft-1) =  0.
214
215! todo0710 : marge d'amelioration en tenant compte des accum de depot ? 
216   do k=2,nft
217      accucumul(nft-k) = accucumul(nft-k+1) + (exp(rappact  * tpert(k) ))
218   end do
219       
220   time_max_accu = tdate(1)
221   !print*, tdate(1), nft , nftnl, rappact
222   !print*, tdate(1), nft , rappact!, nftnl, rappact
223   uxsave(:,:,:) = 0.0
224   uysave(:,:,:) = 0.0
225   uzrsave(:,:,:) = 0.0
226   xdep(:,:,:)=xdepk(:,:,:)
227   ydep(:,:,:)=ydepk(:,:,:)
228   tdep(:,:,:)=tdepk(:,:,:)
229   freezeon(:,:) = 0.   ! variable added on 19/03/03
230
231!   print*, 'tpert ',tpert(:)
232!   print*, 'accucumul ', accucumul(:)
233   print*, 'tbegin, tend : ',tbegin, tend
234
235end subroutine init_tracer
236
237
238subroutine tracer
239
240  real,dimension(nz) :: E   ! vertical coordinate in ice, scaled to H zeta
241
242  E(1)=0.
243  E(NZ)=1.
244  do K=1,NZ
245     if ((K.ne.1).and.(K.ne.NZ)) E(K)=(K-1.)/(NZ-1.)
246  end do
247
248  If (Itracebug.Eq.1) Write(num_tracebug,*)' Entree Dans Routine tracer at time=',Time
249
250   time_tra = floor(time) 
251      !! start of main task
252      ! main calculation loop omits the outer buffer cells.  update here.
253
254  !===== accucumulative velocities
255 
256     uxsave(:,:,:) = uxsave(:,:,:) + ux(:,:,:)*real(dtt)  ! todo0710 mettre les tableaux
257     uysave(:,:,:) = uysave(:,:,:) + uy(:,:,:)*real(dtt)  ! attention uxsave position en maille staggered
258                                     ! ATTENTION, il faudrait deltra pas dtt, SIGNE ?
259     do j=1,ny
260        do i=1,nx
261           ! added on 16/12/03: somehow worked before too
262           if (h(i,j)>1.5) &   
263                uzrsave(i,j,:) = uzrsave(i,j,:) + uzr(i,j,:) *real(dtt) / h(i,j)
264        end do
265     end do
266
267!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
268
269! recalcul de accucumul d'après les valeurs courantes pour
270!reconstruire accucumul de façon interactive. utile pour une calotte
271! qui varie beaucoup en géométrie.
272
273!if (mod(time_tra,100)==0) then
274!       update accu, rewriting on top of the older value
275!       i=-time_tra/100
276!       accucumul(i) = accucumul(i+1) + acc(105,49)/ 0.0295
277!end if
278
279  if (mod(time_tra,deltracer)==0) then
280     
281     !============ east and west boundaries
282     do j=2,ny-1
283        i_jj = j-1
284        xdep(:,1,i_jj) = xgrid(2)
285        xdep(:,nx-2,i_jj) = xgrid(nx-1)
286        ydep(:,1,i_jj) = ygrid(j)
287        ydep(:,nx-2,i_jj) = ygrid(j)
288!        tdep(:,1,i_jj) = time_
289!        tdep(:,nx-2,i_jj) = time_
290        tdep(:,1,i_jj) = time   
291        tdep(:,nx-2,i_jj) = time
292     end do
293      !=========== north and south boundaries
294     do i=2,nx-1
295        i_ij = i-1
296        xdep(:,i_ij,1) = xgrid(i)
297        xdep(:,i_ij,ny-2) = xgrid(i)
298        ydep(:,i_ij,1) = ygrid(2)
299        ydep(:,i_ij,ny-2) = ygrid(ny-1)
300!       tdep(:,i_ij,1) =  time_
301!       tdep(:,i_ij,ny-2) =  time_
302        tdep(:,i_ij,1) =  time   
303        tdep(:,i_ij,ny-2) =  time
304     end do
305
306     
307      !============= inside grid: now watch out: (i,j) dyn-grid <=> (i_ij, i_jj) tracer-grid
308      do j=3,ny-2
309        do i=3,nx-2
310          i_ij = i-1
311          i_jj = j-1
312
313! nl 11/03/04
314          if (h(i,j)<=1.5) then
315            xdep(:,i_ij,i_jj) = xgrid(i)
316            ydep(:,i_ij,i_jj) = ygrid(j)
317            tdep(:,i_ij,i_jj) =  time 
318            cycle
319          end if
320          freezeon(i_ij,i_jj) = max(0.0, freezeon(i_ij,i_jj)-bmelt(i,j) )
321         
322          ip = i+1
323          jp = j+1
324          do k=1,nz
325            v_xij       = (uxsave(i,j,k)+uxsave(ip,j,k)) / real(2)
326            v_yij       = (uysave(i,j,k)+uysave(i,jp,k)) / real(2)
327            v_zij       = uzrsave(i,j,k)
328
329            x_tra       = xgrid(i) - v_xij      !* real(deltracer)
330            y_tra       = ygrid(j) - v_yij      !* real(deltracer)
331            e_tra       = e(k) - v_zij          !* real(deltracer)
332
333            if ((k==nz).and.(bmelt(i,j)<-1.e-4)) then   ! bottom freezeon
334               xdep(k,i_ij,i_jj) = xgrid(i)
335               ydep(k,i_ij,i_jj) = ygrid(j)
336               if (h(i,j)>100.) then
337                 tdep(k,i_ij,i_jj) =  tdepk(k,i_ij,i_jj) - &
338                        min( abs( (tdepk(k-1,i_ij,i_jj)-tdepk(k,i_ij,i_jj)) *  &
339                          bmelt(i,j)  / ( ( e(k)-e(k-1) )*h(i,j) ) ), 2.) * deltracer
340                  ! this added 26/03/04 so that not more 2.*deltracer new years of refreezing
341               else
342                 tdep(k,i_ij,i_jj) =  tdepk(k,i_ij,i_jj)
343                 
344               end if 
345                !fix:    need to track where freeze-on happens
346                ! *1.001        proven to be dumb
347!               if((tdep(k,i_ij,i_jj)<=time_  -  1500000.) .or. (tdep(k,i_ij,i_jj)>=tend+1.0)) then
348               if((tdep(k,i_ij,i_jj)<=time  -  1500000.) .or. (tdep(k,i_ij,i_jj)>=tend+1.0)) then
349                  print *,"tdep(",k,",",i_ij,",",i_jj,") range error at time=",time_tra
350                   print *,"tdep(k,i_ij,i_jj) =",tdep(k,i_ij,i_jj), v_xij, v_yij, v_zij
351                   print *, bmelt(i,j), s(i,j), h(i,j), e_tra, tdepk(k,i_ij,i_jj), tdepk(k,i,j)
352                   print *, flot(im,jm), flot(im+1,jm), flot(im+1,jm+1), flot(im,jm+1)
353                   if (k>1) tdep(k,i_ij,i_jj) = tdepk(k-1,i_ij,i_jj) - real(h(i,j)*(k-1)**2/30.)
354                   print *, tdep(k,i_ij,i_jj)
355                end if
356
357             else if (k==1 .and. bm(i,j)>0.0) then        !! surface mass balance
358              xdep(k,i_ij,i_jj) = xgrid(i)
359              ydep(k,i_ij,i_jj) = ygrid(j)
360              tdep(k,i_ij,i_jj) =  time 
361            else if ( e_tra < e(1) ) then               !! rapid thickening by
362              xdep(k,i_ij,i_jj) = xgrid(i)                              !! surface accumulation
363              ydep(k,i_ij,i_jj) = ygrid(j)
364              tdep(k,i_ij,i_jj) =  time   
365           else
366                call celltest(time_tra,i,j,k,v_xij, v_yij, v_zij, x_tra, y_tra, e_tra, &
367                                im,jm,km,fx,fy,fz)
368                !! convert seach results from (i,j) grid to (i_ij,i_jj) tracer grid
369                !print*,'i,j,k avt interp',i,j,k
370
371                call interpolate(i_ij, i_jj, k, im-1, jm-1, km, fx,fy,fz, e_tra) 
372
373                if((tdep(k,i_ij,i_jj)<=time  -  1500000.) .or. (tdep(k,i_ij,i_jj)> time   )) then
374                   print *,"tdep(",k,",",i_ij,",",i_jj,") range error at time=",time_tra
375                   print *,"tdep(k,i_ij,i_jj) =",tdep(k,i_ij,i_jj), v_xij, v_yij, v_zij
376                   print *, bmelt(i,j), s(i,j), h(i,j), e_tra, tdepk(k,i_ij,i_jj), tdepk(k,i,j)
377                   print *, flot(im,jm), flot(im+1,jm), flot(im+1,jm+1), flot(im,jm+1) 
378                   if (k>1) tdep(k,i_ij,i_jj) = tdepk(k-1,i_ij,i_jj) - real(h(i,j)*(k-1)**2/30.)
379                   print *, tdep(k,i_ij,i_jj) 
380                end if
381
382            end if
383
384          end do
385        end do
386      end do
387
388      xdepk(:,:,:) = xdep(:,:,:)
389      ydepk(:,:,:) = ydep(:,:,:)
390      tdepk(:,:,:) = tdep(:,:,:)
391      uxsave = 0.00
392      uysave = 0.00
393      uzrsave = 0.00
394
395   end if !========== if deltracer
396
397   if (abs(time-tend)<19) then   
398      open(505,file="freezeonmap.out",status="replace")
399      print*,"save freezeon file"
400      write(505,  * ) freezeon
401      close(505)
402      print *, tdep(:,70,70)
403   end if
404
405! aurel neem -> je me cree une version formatte "out" pour les sorties
406! on traite d'abord les 4 coins
407xdep_out(1,1,:)=xdep(:,1,1)
408xdep_out(nx,1,:)=xdep(:,nx-2,1)
409xdep_out(nx,ny,:)=xdep(:,nx-2,ny-2)
410xdep_out(1,ny,:)=xdep(:,1,ny-2)
411ydep_out(1,1,:)=ydep(:,1,1)
412ydep_out(nx,1,:)=ydep(:,nx-2,1)
413ydep_out(nx,ny,:)=ydep(:,nx-2,ny-2)
414ydep_out(1,ny,:)=ydep(:,1,ny-2)
415tdep_out(1,1,:)=tdep(:,1,1)
416tdep_out(nx,1,:)=tdep(:,nx-2,1)
417tdep_out(nx,ny,:)=tdep(:,nx-2,ny-2)
418tdep_out(1,ny,:)=tdep(:,1,ny-2)
419! puis les bords Ouest et Est
420do i=2,nx-1
421   xdep_out(i,1,:)=xdep(:,i-1,1)
422   xdep_out(i,ny,:)=xdep(:,i-1,ny-2)
423   ydep_out(i,1,:)=ydep(:,i-1,1)
424   ydep_out(i,ny,:)=ydep(:,i-1,ny-2)
425   tdep_out(i,1,:)=tdep(:,i-1,1)
426   tdep_out(i,ny,:)=tdep(:,i-1,ny-2)
427end do
428! et les bords Sud et Nord
429do j=2,ny-1
430   xdep_out(1,j,:)=xdep(:,1,j-1)
431   xdep_out(nx,j,:)=xdep(:,nx-2,j-1)
432   ydep_out(1,j,:)=ydep(:,1,j-1)
433   ydep_out(nx,j,:)=ydep(:,nx-2,j-1)
434   tdep_out(1,j,:)=tdep(:,1,j-1)
435   tdep_out(nx,j,:)=tdep(:,nx-2,j-1)
436end do
437! et enfin on complete l'interieur
438do k=1,nz
439   do j=2,ny-1
440      do i=2,nx-1
441         xdep_out(i,j,k)=xdep(k,i-1,j-1)
442         ydep_out(i,j,k)=ydep(k,i-1,j-1)
443         tdep_out(i,j,k)=tdep(k,i-1,j-1)
444      end do
445   end do
446end do
447
448end subroutine tracer
449
450!=============================================================================
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470!========================================================================
471!========================================================================
472subroutine get_date_time(yr, mon, iday, hr, minu, sec)
473
474  implicit none
475  integer, intent(out) :: yr, mon, iday, hr, minu, sec
476
477  integer, dimension(8) :: values
478  character :: date*8, chtime*10, zone*5
479     
480  call date_and_time(date, chtime, zone, values)
481     
482  yr   = values(1)
483  mon  = values(2)
484  iday = values(3)
485  hr   = values(5)
486  minu  = values(6)
487  sec  = values(7)
488
489end subroutine get_date_time
490
491!========================================================================
492subroutine read_tr_dat
493!      use global_param
494  use tracer_vars
495  implicit none
496
497! lit le fichier tracer.dat du coup.
498  integer :: ioerr
499  character(len=60) :: line
500
501  open(unit=unit_tr_dat, file=file_tr_dat)
502 
503  read(unit=unit_tr_dat, fmt=*) delttrout
504  print*,delttrout
50511 format(i8)
506  nspectimes = 0
507  read(unit=unit_tr_dat, fmt="(a)",iostat=ioerr) line
508  if(ioerr/=0 .or. line(1:13)/="special times") return
509
510  do
511     nspectimes = nspectimes+1
512     read(unit=line, fmt=*,iostat=ioerr) trout(nspectimes)
513     if(ioerr/=0) then
514        nspectimes=nspectimes-1
515        exit
516     end if
517  end do
518
519  close(unit=unit_tr_dat)
520 
521end subroutine read_tr_dat
522
523   
524end module tracer_mod
Note: See TracBrowser for help on using the repository browser.