source: branches/GRISLIv3/SOURCES/tracer_mod.f90 @ 425

Last change on this file since 425 was 387, checked in by dumas, 16 months ago

use only for tracers : compile but not tested in a simulation

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