Changeset 54
- Timestamp:
- 03/31/16 11:01:50 (8 years ago)
- Location:
- trunk/SOURCES
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/SOURCES/Makefile.grisli-gfortran.inc
r51 r54 726 726 727 727 beta_iter_vitbil_mod.o : Draggings_modules/beta_iter_vitbil_mod.f90 728 $(F T) Draggings_modules/beta_iter_vitbil_mod.f90728 $(F_NETCDF) Draggings_modules/beta_iter_vitbil_mod.f90 729 729 730 730 dragging_stream_impose_vitbil_mod.o : Ant40_files/dragging_stream_impose_vitbil_mod.f90 … … 1587 1587 $(diagnoshelf) \ 1588 1588 $(Liste_Netcdf) \ 1589 $(routines_communes) steps_time_loop.o $(routine_elliptiques) 1589 $(routines_communes) steps_time_loop.o $(routine_elliptiques) \ 1590 $(Liste_BLAS) 1590 1591 1591 1592 $(LK) -o ../bin/Grice2sea \ … … 1599 1600 $(Liste_Netcdf) \ 1600 1601 $(routines_communes) steps_time_loop.o \ 1601 $(routine_elliptiques) $(NCDF_LIB) $(MKL_LIB) 1602 $(routine_elliptiques) \ 1603 $(Liste_BLAS) 1604 $(NCDF_LIB) $(MKL_LIB) 1602 1605 1603 1606 Grice2sea_iterbeta : $(Dim_GrIce2sea) $(mod_dim_communs) \ … … 1628 1631 $(mod_communs) \ 1629 1632 $(mod_clim_tof) \ 1630 $(mod_ no_tracers) \1633 $(mod_tracers) \ 1631 1634 $(mod_ell) $(Liste_hemin40) \ 1632 1635 $(diagnoshelf) \ … … 1641 1644 $(mod_communs) \ 1642 1645 $(mod_clim_tof) \ 1643 $(mod_ no_tracers) \1646 $(mod_tracers) \ 1644 1647 $(mod_ell) $(Liste_hemin40) \ 1645 1648 $(diagnoshelf) \ -
trunk/SOURCES/Makefile.grisli.inc
r51 r54 1635 1635 $(mod_communs) \ 1636 1636 $(mod_clim_tof) \ 1637 $(mod_ no_tracers) \1637 $(mod_tracers) \ 1638 1638 $(mod_ell) $(Liste_hemin40) \ 1639 1639 $(diagnoshelf) \ … … 1648 1648 $(mod_communs) \ 1649 1649 $(mod_clim_tof) \ 1650 $(mod_ no_tracers) \1650 $(mod_tracers) \ 1651 1651 $(mod_ell) $(Liste_hemin40) \ 1652 1652 $(diagnoshelf) \ -
trunk/SOURCES/tracer_mod.f90
r11 r54 6 6 !-------------------------------------------------------------- 7 7 ! Position des forages dans la grille grisli. 8 ! Groenland grille GreeneemXX 8 9 ! 9 10 ! Ordre GMT : … … 29 30 use tracer_vars ! module de declaration de variables pour les traceurs 30 31 !cdc use climat_perturb_mois_mod ! on en a besoin, notamment pour NFT 31 32 ! afq, new version of tracer, does not require climat_perturb 33 32 34 implicit none 33 35 … … 51 53 integer :: time_tra ! temps traceur 52 54 53 integer :: nft ! declaration locale de NFT : a corriger plus tard pour recuperer la valeur du module climat54 real :: rappact ! pour le calcul du rapport 'accumulation : a corriger plus tard pour recuperer la valeur du module climat55 real,dimension(:),allocatable :: tpert ! temperature for climate forcing : a corriger plus tard pour recuperer la valeur du module climat56 real,dimension(:),allocatable :: TDATE ! time for climate forcing : a corriger plus tard pour recuperer la valeur du module climat57 58 55 59 56 !=============================================================== … … 69 66 real :: lambdab ! fusion basale 'artificielle' pour age des couches 70 67 real :: agemax ! age maximum de la glace (pour init_tracer) 71 68 69 integer :: pert_type ! 0 we give tpert/coeft/rappact,1 we give acc_pert only / afq 70 character(len=200) :: filforc ! nom du fichier forcage 71 character(len=200) :: filin ! variable intermediaire 72 character(len=8) :: control ! label to check clim. forc. file (filin) is usable 73 integer :: err ! pour l'allocation des tableaux 74 real :: bidon ! to skip sealevel change in fileforc 75 72 76 integer :: kk ! indice vertical pour definition de E, mais on veut conserver la valeur de k 73 77 ! real,dimension(nz) :: E ! vertical coordinate in ice, scaled to H zeta 78 79 ! on lit les parametres par namelist dorenavant 80 namelist/tracer/file_tr_dat,file_tr_out,file_tr_dep,lambdab,agemax,coefT_tra,rappact_tra,filforc,pert_type 74 81 75 82 … … 93 100 deltracer = dtt ! sorte de pas de temps du tracer - multiple de dtt ! 94 101 95 96 ! on lit les parametres par namelist dorenavant97 namelist/tracer/file_tr_dat,file_tr_out,file_tr_dep,lambdab,agemax !aurel : a ne pas oublier paramlist98 102 99 103 rewind(num_param) ! pour revenir au debut du fichier param_list.dat … … 111 115 write(num_rep_42,'(A,A)') 'tracer .out file : ', file_tr_out 112 116 write(num_rep_42,'(A,A)') 'tracer dep. file : ', file_tr_dep 117 write(num_rep_42,'(A,A)') 'lambdab : ', lambdab 118 write(num_rep_42,'(A,A)') 'agemax : ', agemax 119 write(num_rep_42,'(A,A)') 'coefT_tra : ', coefT_tra 120 write(num_rep_42,'(A,A)') 'rappact_tra : ', rappact_tra 121 write(num_rep_42,'(A,A)') 'fileforc : ', filforc 122 write(num_rep_42,'(A,A)') 'pert_type : ', pert_type 113 123 write(num_rep_42,*)'/' 114 124 write(num_rep_42,428) 'declaration des repertoires *traceurs*' … … 202 212 203 213 214 ! afq marion dufresne: 215 ! to build accucumul, we read all what is in climat_perturb here now 216 ! meaning: nft, tpert, tdate 217 ! OR if pert_type = 1 we read nft, acc_pert, tdate 218 ! acc_pert should be the ratio acc_palaeo / acc_present 219 220 filin=trim(dirforcage)//trim(filforc) 221 222 open(num_forc,file=filin,status='old') 223 224 read(num_forc,*) control,nft_tra 225 226 ! Determination of file size (line nb), allocation of perturbation array 227 if (control.ne.'nb_lines') then 228 write(6,*) filin,'indiquer le nb de ligne en debut de fichier:' 229 write(6,*) 'le nb de lignes et le label de control nb_lines' 230 stop 231 endif 232 if (.not.allocated(tdate_tra)) then 233 allocate(tdate_tra(nft_tra),stat=err) 234 if (err/=0) then 235 print *,"erreur a l'allocation du tableau Tdate",err 236 stop 4 237 end if 238 end if 239 if (.not.allocated(tpert_tra)) then 240 allocate(tpert_tra(nft_tra),stat=err) 241 if (err/=0) then 242 print *,"erreur a l'allocation du tableau Tpert",err 243 stop 4 244 end if 245 end if 246 247 do i=1,nft_tra 248 read(num_forc,*) tdate_tra(i),bidon,tpert_tra(i) 249 end do 250 close(num_forc) 251 252 204 253 2009 format(21(f12.1)) 205 254 … … 209 258 ! chaque indice d'accucumul vaut 100 ans. 210 259 ! 0 pour le present, 10 pour 1000 ans 211 allocate(accucumul(0:nft -1)); accucumul = 0.0260 allocate(accucumul(0:nft_tra-1)); accucumul = 0.0 212 261 ! watch out its size !!! 213 262 ! note: tpert(1) = tpert at -422700yr … … 217 266 ! a la rigueur je pourrais modifier le fichier d'interpolation pour 218 267 ! 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 268 accucumul(nft_tra-1) = 0. 269 270 ! todo0710 : marge d'amelioration en tenant compte des accum de depot ? 271 272 if (type_accum.eq.0) then !Tpert is given, needs to be converted in accumulation ratio 273 tpert_tra(:)=tpert_tra(:)*coefT_tra 274 do k=2,nft_tra 275 accucumul(nft_tra-k) = accucumul(nft_tra-k+1) + (exp(rappact_tra * tpert_tra(k) )) 276 end do 277 else if (type_accum.eq.1) then ! accum ration is given, nothing to do 278 do k=2,nft_tra 279 accucumul(nft_tra-k) = accucumul(nft_tra-k+1) + tpert_tra(k) 280 end do 281 else 282 write (*,*) "Problem in tracer, do you provide tpert or accpert???" 283 STOP 284 end if 285 time_max_accu = tdate_tra(1) 229 286 uxsave(:,:,:) = 0.0 230 287 uysave(:,:,:) = 0.0 … … 233 290 ydep(:,:,:)=ydepk(:,:,:) 234 291 tdep(:,:,:)=tdepk(:,:,:) 235 freezeon(:,:) = 0. ! variable added on 19/03/03 236 237 ! print*, 'tpert ',tpert(:) 238 ! print*, 'accucumul ', accucumul(:) 292 freezeon(:,:) = 0. ! variable added on 19/03/03 293 239 294 print*, 'tbegin, tend : ',tbegin, tend 240 295 … … 266 321 do i=1,nx 267 322 ! added on 16/12/03: somehow worked before too 268 if (h(i,j)>1.5) & 269 323 if (h(i,j)>1.5) & 324 uzrsave(i,j,:) = uzrsave(i,j,:) + uzr(i,j,:) *real(dtt) / h(i,j) 270 325 end do 271 326 end do … … 288 343 do j=2,ny-1 289 344 i_jj = j-1 290 345 xdep(:,1,i_jj) = xgrid(2) 291 346 xdep(:,nx-2,i_jj) = xgrid(nx-1) 292 347 ydep(:,1,i_jj) = ygrid(j) 293 348 ydep(:,nx-2,i_jj) = ygrid(j) 294 349 ! tdep(:,1,i_jj) = time_ … … 300 355 do i=2,nx-1 301 356 i_ij = i-1 302 357 xdep(:,i_ij,1) = xgrid(i) 303 358 xdep(:,i_ij,ny-2) = xgrid(i) 304 305 359 ydep(:,i_ij,1) = ygrid(2) 360 ydep(:,i_ij,ny-2) = ygrid(ny-1) 306 361 ! tdep(:,i_ij,1) = time_ 307 362 ! tdep(:,i_ij,ny-2) = time_ 308 309 363 tdep(:,i_ij,1) = time 364 tdep(:,i_ij,ny-2) = time 310 365 end do 311 366 … … 318 373 319 374 ! nl 11/03/04 320 375 if (h(i,j)<=1.5) then 321 376 xdep(:,i_ij,i_jj) = xgrid(i) 322 377 ydep(:,i_ij,i_jj) = ygrid(j) 323 378 tdep(:,i_ij,i_jj) = time 324 325 379 cycle 380 end if 326 381 freezeon(i_ij,i_jj) = max(0.0, freezeon(i_ij,i_jj)-bmelt(i,j) ) 327 382 328 383 ip = i+1 329 384 jp = j+1 … … 333 388 v_zij = uzrsave(i,j,k) 334 389 335 x_tra = xgrid(i) - v_xij 336 y_tra = ygrid(j) - v_yij 337 e_tra = e(k) - v_zij 390 x_tra = xgrid(i) - v_xij !* real(deltracer) 391 y_tra = ygrid(j) - v_yij !* real(deltracer) 392 e_tra = e(k) - v_zij !* real(deltracer) 338 393 339 394 if ((k==nz).and.(bmelt(i,j)<-1.e-4)) then ! bottom freezeon 340 395 xdep(k,i_ij,i_jj) = xgrid(i) 341 396 ydep(k,i_ij,i_jj) = ygrid(j) 342 397 if (h(i,j)>100.) then 343 398 tdep(k,i_ij,i_jj) = tdepk(k,i_ij,i_jj) - & 344 345 346 347 348 399 min( abs( (tdepk(k-1,i_ij,i_jj)-tdepk(k,i_ij,i_jj)) * & 400 bmelt(i,j) / ( ( e(k)-e(k-1) )*h(i,j) ) ), 2.) * deltracer 401 ! this added 26/03/04 so that not more 2.*deltracer new years of refreezing 402 else 403 tdep(k,i_ij,i_jj) = tdepk(k,i_ij,i_jj) 349 404 350 405 end if 351 406 !fix: need to track where freeze-on happens 352 407 ! *1.001 proven to be dumb … … 375 430 !print*,'i,j,k avt interp',i,j,k 376 431 377 call interpolate(i_ij, i_jj, k, im-1, jm-1, km, fx,fy,fz, e_tra, nft )432 call interpolate(i_ij, i_jj, k, im-1, jm-1, km, fx,fy,fz, e_tra, nft_tra) 378 433 379 434 if((tdep(k,i_ij,i_jj)<=time - 1500000.) .or. (tdep(k,i_ij,i_jj)> time )) then … … 469 524 470 525 471 472 473 474 475 476 526 !======================================================================== 477 527 !======================================================================== -
trunk/SOURCES/tracer_vars_mod.f90
r4 r54 25 25 integer :: nspectimes 26 26 real :: time_max_accu 27 28 ! afq marion dufresne: tracer with no climat perturb 29 ! i.e. we read the Forcage Tpert to compute past accumulation via 30 ! rappact 31 ! OR we read directly past accum (choice in namelist) 32 ! if past accum is given it will still be read with the name Tpert 33 ! the past accum given should be the ratio acc_palaeo / acc_present 34 integer :: nft_tra ! number of snapshot Tpert 35 integer :: type_accum ! Tpert is given (=0) or accpert (1) 36 real :: rappact_tra ! accumulation ratio to convert Tpert 37 real :: coeft_tra ! coeft to convert Tpert 38 real,dimension(:),allocatable :: tpert_tra ! temperature perturbation 39 real,dimension(:),allocatable :: tdate_tra ! time for Tpert 40 27 41 28 42 integer, dimension(maxspectimes) :: trout
Note: See TracChangeset
for help on using the changeset viewer.