Changeset 1334 for trunk/NEMO/OPA_SRC
- Timestamp:
- 2009-03-03T15:07:48+01:00 (15 years ago)
- Location:
- trunk/NEMO/OPA_SRC
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMO/OPA_SRC/DIA/diagap.F90
r1317 r1334 82 82 !! * local declarations 83 83 INTEGER :: ji, jj, jk ! dummy loop indices 84 INTEGER :: it 84 INTEGER :: it, itmod ! time step indices 85 85 CHARACTER (len=40) :: clhstnam 86 86 CHARACTER (len=40) :: clop … … 170 170 ! Horizontal grid : zphi() 171 171 CALL histbeg( clhstnam, 1, zfoo, 1, zloo, & 172 1, 1, 1, 1, 0, zjulian, zdt, nhoridg, numgap, domain_id=nidom )172 1, 1, 1, 1, nit000-1, zjulian, zdt, nhoridg, numgap, domain_id=nidom ) 173 173 ! Vertical grids : gdept, gdepw 174 174 CALL histvert( numgap, "deptht", "Vertical T levels", & 175 "m", jpk, gdept_0, ndepidg )175 "m", jpk, gdept_0, ndepidg, "down" ) 176 176 177 177 ! define fields to be stored … … 198 198 ! ---------------------- 199 199 200 it = kt - nit000 + 1 ! define time axis 201 IF( MOD( it, ngap ) == 0 ) THEN 200 itmod = kt - nit000 + 1 ! define time axis 201 it = kt 202 IF( MOD( itmod, ngap ) == 0 ) THEN 202 203 203 204 ! initialization … … 240 241 ! ----------------------------====== 241 242 242 IF( MOD( it , nprg ) == 0 ) THEN243 IF( MOD( itmod, nprg ) == 0 ) THEN 243 244 IF(lwp) THEN 244 245 WRITE(numout,*) 'dia_gap: time step = ', kt, 'model - data' -
trunk/NEMO/OPA_SRC/DIA/diaptr.F90
r1317 r1334 418 418 419 419 CHARACTER (len=40) :: clhstnam, clop ! temporary names 420 INTEGER :: iline, it, ji 420 INTEGER :: iline, it, ji, itmod ! 421 421 REAL(wp) :: zsto, zout, zdt, zjulian ! temporary scalars 422 422 REAL(wp), DIMENSION(jpj) :: zphi, zfoo … … 424 424 425 425 ! define time axis 426 it = kt - nit000 + 1 426 it = kt 427 itmod = kt - nit000 + 1 427 428 428 429 ! Initialization … … 481 482 ! Horizontal grid : zphi() 482 483 CALL histbeg(clhstnam, 1, zfoo, jpj, zphi, & 483 1, 1, 1, jpj, 0, zjulian, zdt, nhoridz, numptr, domain_id=nidom )484 1, 1, 1, jpj, nit000-1, zjulian, zdt, nhoridz, numptr, domain_id=nidom ) 484 485 ! Vertical grids : gdept_0, gdepw_0 485 486 CALL histvert( numptr, "deptht", "Vertical T levels", & 486 "m", jpk, gdept_0, ndepidzt )487 "m", jpk, gdept_0, ndepidzt, "down" ) 487 488 CALL histvert( numptr, "depthw", "Vertical W levels", & 488 "m", jpk, gdepw_0, ndepidzw )489 "m", jpk, gdepw_0, ndepidzw, "down" ) 489 490 490 491 ! Zonal mean T and S … … 555 556 ENDIF 556 557 557 IF( MOD( it , nf_ptr ) == 0 ) THEN558 IF( MOD( itmod, nf_ptr ) == 0 ) THEN 558 559 559 560 IF(lwp) THEN -
trunk/NEMO/OPA_SRC/DIA/diawri.F90
r1318 r1334 109 109 INTEGER :: inum = 11 ! temporary logical unit 110 110 INTEGER :: & 111 iimi, iima, ipk, it, 111 iimi, iima, ipk, it, itmod, & ! temporary integers 112 112 ijmi, ijma ! " " 113 113 REAL(wp) :: & … … 148 148 149 149 ! define time axis 150 it = kt - nit000 + 1 150 it = kt 151 itmod = kt - nit000 + 1 151 152 152 153 … … 184 185 CALL histbeg( clhstnam, jpi, glamt, jpj, gphit, & ! Horizontal grid: glamt and gphit 185 186 & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & 186 & 0, zjulian, zdt, nh_T, nid_T, domain_id=nidom )187 & nit000-1, zjulian, zdt, nh_T, nid_T, domain_id=nidom ) 187 188 CALL histvert( nid_T, "deptht", "Vertical T levels", & ! Vertical grid: gdept 188 & "m", ipk, gdept_0, nz_T )189 & "m", ipk, gdept_0, nz_T, "down" ) 189 190 ! ! Index of ocean points 190 191 CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T ) ! volume … … 197 198 CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu, & ! Horizontal grid: glamu and gphiu 198 199 & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & 199 & 0, zjulian, zdt, nh_U, nid_U, domain_id=nidom )200 & nit000-1, zjulian, zdt, nh_U, nid_U, domain_id=nidom ) 200 201 CALL histvert( nid_U, "depthu", "Vertical U levels", & ! Vertical grid: gdept 201 & "m", ipk, gdept_0, nz_U )202 & "m", ipk, gdept_0, nz_U, "down" ) 202 203 ! ! Index of ocean points 203 204 CALL wheneq( jpi*jpj*ipk, umask, 1, 1., ndex_U , ndim_U ) ! volume … … 210 211 CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv, & ! Horizontal grid: glamv and gphiv 211 212 & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & 212 & 0, zjulian, zdt, nh_V, nid_V, domain_id=nidom )213 & nit000-1, zjulian, zdt, nh_V, nid_V, domain_id=nidom ) 213 214 CALL histvert( nid_V, "depthv", "Vertical V levels", & ! Vertical grid : gdept 214 & "m", ipk, gdept_0, nz_V )215 & "m", ipk, gdept_0, nz_V, "down" ) 215 216 ! ! Index of ocean points 216 217 CALL wheneq( jpi*jpj*ipk, vmask, 1, 1., ndex_V , ndim_V ) ! volume … … 223 224 CALL histbeg( clhstnam, jpi, glamt, jpj, gphit, & ! Horizontal grid: glamt and gphit 224 225 & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & 225 & 0, zjulian, zdt, nh_W, nid_W, domain_id=nidom )226 & nit000-1, zjulian, zdt, nh_W, nid_W, domain_id=nidom ) 226 227 CALL histvert( nid_W, "depthw", "Vertical W levels", & ! Vertical grid: gdepw 227 & "m", ipk, gdepw_0, nz_W )228 & "m", ipk, gdepw_0, nz_W, "down" ) 228 229 229 230 … … 406 407 ! donne le nombre d'elements, et ndex la liste des indices a sortir 407 408 408 IF( lwp .AND. MOD( it , nwrite ) == 0 ) THEN409 IF( lwp .AND. MOD( itmod, nwrite ) == 0 ) THEN 409 410 WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step' 410 411 WRITE(numout,*) '~~~~~~ ' … … 512 513 513 514 ! Create an output files (output.abort.nc) if S < 0 or u > 20 m/s 514 IF( kindic < 0 ) CALL dia_wri_state( 'output.abort' )515 IF( kindic < 0 ) CALL dia_wri_state( 'output.abort', kt ) 515 516 516 517 ! 3. Close all files … … 526 527 527 528 528 SUBROUTINE dia_wri_state( cdfile_name )529 SUBROUTINE dia_wri_state( cdfile_name, kt ) 529 530 !!--------------------------------------------------------------------- 530 531 !! *** ROUTINE dia_wri_state *** … … 546 547 !!---------------------------------------------------------------------- 547 548 !! * Arguments 548 CHARACTER (len=* ), INTENT( in ) :: &549 cdfile_name ! name of the file created549 CHARACTER (len=* ), INTENT( in ) :: cdfile_name ! name of the file created 550 INTEGER , INTENT( in ) :: kt ! ocean time-step index 550 551 551 552 !! * Local declarations … … 588 589 zjulian = zjulian - adatrj ! set calendar origin to the beginning of the experiment 589 590 CALL histbeg( clname, jpi, glamt, jpj, gphit, & 590 1, jpi, 1, jpj, 0, zjulian, zdt, nh_i, id_i, domain_id=nidom ) ! Horizontal grid : glamt and gphit591 1, jpi, 1, jpj, nit000-1, zjulian, zdt, nh_i, id_i, domain_id=nidom ) ! Horizontal grid : glamt and gphit 591 592 CALL histvert( id_i, "deptht", "Vertical T levels", & ! Vertical grid : gdept 592 "m", jpk, gdept_0, nz_i )593 "m", jpk, gdept_0, nz_i, "down") 593 594 594 595 ! Declare all the output fields as NetCDF variables … … 634 635 635 636 ! Write all fields on T grid 636 CALL histwrite( id_i, "votemper", 1, tn , jpi*jpj*jpk, idex ) ! now temperature637 CALL histwrite( id_i, "vosaline", 1, sn , jpi*jpj*jpk, idex ) ! now salinity638 #if defined key_dynspg_rl 639 CALL histwrite( id_i, "sobarstf", 1, bsfn , jpi*jpj , idex ) ! barotropic streamfunction637 CALL histwrite( id_i, "votemper", kt, tn , jpi*jpj*jpk, idex ) ! now temperature 638 CALL histwrite( id_i, "vosaline", kt, sn , jpi*jpj*jpk, idex ) ! now salinity 639 #if defined key_dynspg_rl 640 CALL histwrite( id_i, "sobarstf", kt, bsfn , jpi*jpj , idex ) ! barotropic streamfunction 640 641 #else 641 CALL histwrite( id_i, "sossheig", 1, sshn , jpi*jpj , idex ) ! sea surface height642 #endif 643 CALL histwrite( id_i, "vozocrtx", 1, un , jpi*jpj*jpk, idex ) ! now i-velocity644 CALL histwrite( id_i, "vomecrty", 1, vn , jpi*jpj*jpk, idex ) ! now j-velocity645 CALL histwrite( id_i, "vovecrtz", 1, wn , jpi*jpj*jpk, idex ) ! now k-velocity646 CALL histwrite( id_i, "sowaflup", 1, emp , jpi*jpj , idex ) ! freshwater budget647 CALL histwrite( id_i, "sohefldo", 1, qsr + qns, jpi*jpj , idex ) ! total heat flux648 CALL histwrite( id_i, "soshfldo", 1, qsr , jpi*jpj , idex ) ! solar heat flux649 CALL histwrite( id_i, "soicecov", 1, fr_i , jpi*jpj , idex ) ! ice fraction650 CALL histwrite( id_i, "sozotaux", 1, utau , jpi*jpj , idex ) ! i-wind stress651 CALL histwrite( id_i, "sometauy", 1, vtau , jpi*jpj , idex ) ! j-wind stress642 CALL histwrite( id_i, "sossheig", kt, sshn , jpi*jpj , idex ) ! sea surface height 643 #endif 644 CALL histwrite( id_i, "vozocrtx", kt, un , jpi*jpj*jpk, idex ) ! now i-velocity 645 CALL histwrite( id_i, "vomecrty", kt, vn , jpi*jpj*jpk, idex ) ! now j-velocity 646 CALL histwrite( id_i, "vovecrtz", kt, wn , jpi*jpj*jpk, idex ) ! now k-velocity 647 CALL histwrite( id_i, "sowaflup", kt, emp , jpi*jpj , idex ) ! freshwater budget 648 CALL histwrite( id_i, "sohefldo", kt, qsr + qns, jpi*jpj , idex ) ! total heat flux 649 CALL histwrite( id_i, "soshfldo", kt, qsr , jpi*jpj , idex ) ! solar heat flux 650 CALL histwrite( id_i, "soicecov", kt, fr_i , jpi*jpj , idex ) ! ice fraction 651 CALL histwrite( id_i, "sozotaux", kt, utau , jpi*jpj , idex ) ! i-wind stress 652 CALL histwrite( id_i, "sometauy", kt, vtau , jpi*jpj , idex ) ! j-wind stress 652 653 653 654 ! 3. Close the file -
trunk/NEMO/OPA_SRC/TRD/trdmld.F90
r1317 r1334 230 230 INTEGER, INTENT( in ) :: kt ! ocean time-step index 231 231 !! 232 INTEGER :: ji, jj, jk, jl, ik, it 232 INTEGER :: ji, jj, jk, jl, ik, it, itmod 233 233 LOGICAL :: lldebug = .TRUE. 234 234 REAL(wp) :: zavt, zfn, zfn2 … … 389 389 smltrd(:,:,:) = smltrd(:,:,:) * ucf ! is no longer used, and is reset to 0. at next time step) 390 390 391 it = kt - nit000 + 1 392 393 MODULO_NTRD : IF( MOD( it, ntrd ) == 0 ) THEN ! nitend MUST be multiple of ntrd 391 ! define time axis 392 it = kt 393 itmod = kt - nit000 + 1 394 395 MODULO_NTRD : IF( MOD( itmod, ntrd ) == 0 ) THEN ! nitend MUST be multiple of ntrd 394 396 ! 395 397 ztmltot (:,:) = 0.e0 ; zsmltot (:,:) = 0.e0 ! reset arrays to zero … … 576 578 #if defined key_dimgout 577 579 578 IF( MOD( it , ntrd ) == 0 ) THEN580 IF( MOD( itmod, ntrd ) == 0 ) THEN 579 581 iyear = ndastp/10000 580 582 imon = (ndastp-iyear*10000)/100 … … 593 595 ! ---------------------------------- 594 596 595 IF( lwp .AND. MOD( it , ntrd ) == 0 ) THEN597 IF( lwp .AND. MOD( itmod , ntrd ) == 0 ) THEN 596 598 WRITE(numout,*) ' ' 597 599 WRITE(numout,*) 'trd_mld : write trends in the NetCDF file :' … … 683 685 #endif 684 686 685 IF( MOD( it , ntrd ) == 0 ) THEN687 IF( MOD( itmod, ntrd ) == 0 ) THEN 686 688 ! 687 689 ! III.5 Reset cumulative arrays to zero … … 876 878 IF(lwp) WRITE(numout,*) ' Name of NETCDF file ', clhstnam 877 879 CALL histbeg( clhstnam, jpi, glamt, jpj, gphit, & 878 & 1, jpi, 1, jpj, 0, zjulian, rdt, nh_t, nidtrd, domain_id=nidom )880 & 1, jpi, 1, jpj, nit000-1, zjulian, rdt, nh_t, nidtrd, domain_id=nidom ) 879 881 880 882 !-- Define the ML depth variable -
trunk/NEMO/OPA_SRC/TRD/trdvor.F90
r1317 r1334 312 312 INTEGER, INTENT( in ) :: kt ! ocean time-step index 313 313 !! 314 INTEGER :: ji, jj, jk, jl, it 314 INTEGER :: ji, jj, jk, jl, it, itmod 315 315 REAL(wp) :: zmean 316 316 REAL(wp), DIMENSION(jpi,jpj) :: zun, zvn … … 406 406 407 407 ! define time axis 408 it= kt - nit000 + 1 408 it = kt 409 itmod = kt - nit000 + 1 409 410 410 411 IF( MOD( it, ntrd ) == 0 ) THEN … … 455 456 IF( kt >= nit000+1 ) THEN 456 457 457 IF( lwp .AND. MOD( it , ntrd ) == 0 ) THEN458 IF( lwp .AND. MOD( itmod, ntrd ) == 0 ) THEN 458 459 WRITE(numout,*) '' 459 460 WRITE(numout,*) 'trd_vor : write trends in the NetCDF file at kt = ', kt … … 568 569 IF(lwp) WRITE(numout,*) ' Name of NETCDF file ', clhstnam 569 570 CALL histbeg( clhstnam, jpi, glamf, jpj, gphif,1, jpi, & ! Horizontal grid : glamt and gphit 570 & 1, jpj, 0, zjulian, rdt, nh_t, nidvor, domain_id=nidom )571 & 1, jpj, nit000-1, zjulian, rdt, nh_t, nidvor, domain_id=nidom ) 571 572 CALL wheneq( jpi*jpj, fmask, 1, 1., ndexvor1, ndimvor1 ) ! surface 572 573 -
trunk/NEMO/OPA_SRC/step.F90
r1246 r1334 195 195 196 196 IF( ninist == 1 ) THEN ! Output the initial state and forcings 197 CALL dia_wri_state( 'output.init' )197 CALL dia_wri_state( 'output.init', kstp ) 198 198 ninist = 0 199 199 ENDIF
Note: See TracChangeset
for help on using the changeset viewer.