source: trunk/SOURCES/tracer_mod.f90 @ 23

Last change on this file since 23 was 11, checked in by dumas, 9 years ago

Grice2sea compilé et validé avec le module climat_Grice2sea_years_mod. climat_GrIce2sea_years_mod.f90 inclus massb_Ice2sea_RCM et massb_Ice2sea_fixe. pdd_declar_mod.f90 supprimé, les déclarations de variables concernant le pdd sont maintenant dans le module ablation_mod.f90.

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