Changeset 7557
- Timestamp:
- 2017-01-13T17:43:55+01:00 (7 years ago)
- Location:
- branches/2016/dev_INGV_UKMO_2016/NEMOGCM/NEMO/OPA_SRC
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2016/dev_INGV_UKMO_2016/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90
r7351 r7557 41 41 USE zdfddm ! vertical physics: double diffusion 42 42 USE diahth ! thermocline diagnostics 43 USE sbcwave ! wave parameters 43 44 ! 44 45 USE lbclnk ! ocean lateral boundary conditions (or mpp link) … … 682 683 CALL histdef( nid_U, "vozocrtx", "Zonal Current" , "m/s" , & ! un 683 684 & jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) 685 IF( ln_wave .AND. ln_sdw) THEN 686 CALL histdef( nid_U, "sdzocrtx", "Stokes Drift Zonal Current" , "m/s" , & ! usd 687 & jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) 688 ENDIF 684 689 ! !!! nid_U : 2D 685 690 CALL histdef( nid_U, "sozotaux", "Wind Stress along i-axis" , "N/m2" , & ! utau … … 691 696 CALL histdef( nid_V, "vomecrty", "Meridional Current" , "m/s" , & ! vn 692 697 & jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) 698 IF( ln_wave .AND. ln_sdw) THEN 699 CALL histdef( nid_V, "sdmecrty", "Stokes Drift Meridional Current" , "m/s" , & ! vsd 700 & jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) 701 ENDIF 693 702 ! !!! nid_V : 2D 694 703 CALL histdef( nid_V, "sometauy", "Wind Stress along j-axis" , "N/m2" , & ! vtau … … 707 716 IF( lk_zdfddm ) THEN 708 717 CALL histdef( nid_W,"voddmavs","Salt Vertical Eddy Diffusivity" , "m2/s" , & ! avs 718 & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 719 ENDIF 720 721 IF( ln_wave .AND. ln_sdw) THEN 722 CALL histdef( nid_W, "sdvecrtz", "Stokes Drift Vertical Current" , "m/s" , & ! wsd 709 723 & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 710 724 ENDIF … … 825 839 CALL histwrite( nid_W, "votkeavt", it, avt , ndim_T, ndex_T ) ! T vert. eddy diff. coef. 826 840 CALL histwrite( nid_W, "votkeavm", it, avmu , ndim_T, ndex_T ) ! T vert. eddy visc. coef. 841 827 842 IF( lk_zdfddm ) THEN 828 843 CALL histwrite( nid_W, "voddmavs", it, fsavs(:,:,:), ndim_T, ndex_T ) ! S vert. eddy diff. coef. 844 ENDIF 845 846 IF( ln_wave .AND. ln_sdw ) THEN 847 CALL histwrite( nid_U, "sdzocrtx", it, usd , ndim_U , ndex_U ) ! i-StokesDrift-current 848 CALL histwrite( nid_V, "sdmecrty", it, vsd , ndim_V , ndex_V ) ! j-StokesDrift-current 849 CALL histwrite( nid_W, "sdvecrtz", it, wsd , ndim_T , ndex_T ) ! StokesDrift vert. current 829 850 ENDIF 830 851 … … 939 960 & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 940 961 ENDIF 962 ! 963 IF( ln_wave .AND. ln_sdw ) THEN 964 CALL histdef( id_i, "sdzocrtx", "Stokes Drift Zonal" , "m/s" , & ! StokesDrift zonal current 965 & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 966 CALL histdef( id_i, "sdmecrty", "Stokes Drift Merid" , "m/s" , & ! StokesDrift meridonal current 967 & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 968 CALL histdef( id_i, "sdvecrtz", "Stokes Drift Vert" , "m/s" , & ! StokesDrift vertical current 969 & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 970 ENDIF 941 971 942 972 #if defined key_lim2 … … 979 1009 CALL histwrite( id_i, "vovvle3t", kt, e3t_n (:,:,:) , jpi*jpj*jpk, idex )! T-cell thickness 980 1010 END IF 1011 1012 IF( ln_wave .AND. ln_sdw ) THEN 1013 CALL histwrite( id_i, "sdzocrtx", kt, usd , jpi*jpj*jpk, idex) ! now StokesDrift i-velocity 1014 CALL histwrite( id_i, "sdmecrty", kt, vsd , jpi*jpj*jpk, idex) ! now StokesDrift j-velocity 1015 CALL histwrite( id_i, "sdvecrtz", kt, wsd , jpi*jpj*jpk, idex) ! now StokesDrift k-velocity 1016 ENDIF 1017 981 1018 ! 3. Close the file 982 1019 ! ----------------- -
branches/2016/dev_INGV_UKMO_2016/NEMOGCM/NEMO/OPA_SRC/SBC/sbcwave.F90
r7451 r7557 180 180 END DO 181 181 ! 182 CALL iom_put( "ustokes", usd ) 183 CALL iom_put( "vstokes", vsd ) 184 CALL iom_put( "wstokes", wsd ) 185 ! 182 186 CALL wrk_dealloc( jpi,jpj,jpk, ze3divh ) 183 187 CALL wrk_dealloc( jpi,jpj, zk_t, zk_u, zk_v, zu0_sd, zv0_sd )
Note: See TracChangeset
for help on using the changeset viewer.