Changeset 15048
- Timestamp:
- 2021-06-23T18:02:14+02:00 (3 years ago)
- Location:
- NEMO/trunk/src
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/trunk/src/ICE/icedia.F90
r14718 r15048 65 65 INTEGER, INTENT(in) :: kt ! ocean time step 66 66 !! 67 REAL(wp) :: zbg_ivol, zbg_item, zbg_area, zbg_isal 68 REAL(wp) :: zbg_svol, zbg_stem 69 REAL(wp) :: zbg_ipvol, zbg_ilvol 70 REAL(wp) :: z_frc_voltop, z_frc_temtop, z_frc_sal 71 REAL(wp) :: z_frc_volbot, z_frc_tembot 72 REAL(wp) :: zdiff_vol, zdiff_sal, zdiff_tem 67 REAL(wp), DIMENSION(jpi,jpj,16) :: ztmp 68 REAL(wp), DIMENSION(16) :: zbg 73 69 !!--------------------------------------------------------------------------- 74 70 IF( ln_timing ) CALL timing_start('ice_dia') … … 84 80 ENDIF 85 81 86 ! ----------------------- ! 87 ! 1 - Contents ! 88 ! ----------------------- ! 89 IF( iom_use('ibgvol_tot' ) .OR. iom_use('sbgvol_tot' ) .OR. iom_use('ibgarea_tot') .OR. & 90 & iom_use('ibgsalt_tot') .OR. iom_use('ibgheat_tot') .OR. iom_use('sbgheat_tot') .OR. & 91 & iom_use('ipbgvol_tot' ) .OR. iom_use('ilbgvol_tot' ) ) THEN 92 93 zbg_ivol = glob_sum( 'icedia', vt_i(:,:) * e1e2t(:,:) ) * 1.e-9 ! ice volume (km3) 94 zbg_svol = glob_sum( 'icedia', vt_s(:,:) * e1e2t(:,:) ) * 1.e-9 ! snow volume (km3) 95 zbg_area = glob_sum( 'icedia', at_i(:,:) * e1e2t(:,:) ) * 1.e-6 ! area (km2) 96 zbg_isal = glob_sum( 'icedia', st_i(:,:) * e1e2t(:,:) ) * 1.e-9 ! salt content (pss*km3) 97 zbg_item = glob_sum( 'icedia', et_i(:,:) * e1e2t(:,:) ) * 1.e-20 ! heat content (1.e20 J) 98 zbg_stem = glob_sum( 'icedia', et_s(:,:) * e1e2t(:,:) ) * 1.e-20 ! heat content (1.e20 J) 99 ! ponds 100 zbg_ipvol = glob_sum( 'icedia', vt_ip(:,:) * e1e2t(:,:) ) * 1.e-9 ! ice pond volume (km3) 101 zbg_ilvol = glob_sum( 'icedia', vt_il(:,:) * e1e2t(:,:) ) * 1.e-9 ! ice pond lid volume (km3) 102 103 CALL iom_put( 'ibgvol_tot' , zbg_ivol ) 104 CALL iom_put( 'sbgvol_tot' , zbg_svol ) 105 CALL iom_put( 'ibgarea_tot' , zbg_area ) 106 CALL iom_put( 'ibgsalt_tot' , zbg_isal ) 107 CALL iom_put( 'ibgheat_tot' , zbg_item ) 108 CALL iom_put( 'sbgheat_tot' , zbg_stem ) 109 ! ponds 110 CALL iom_put( 'ipbgvol_tot' , zbg_ipvol ) 111 CALL iom_put( 'ilbgvol_tot' , zbg_ilvol ) 112 113 ENDIF 114 82 ztmp(:,:,:) = 0._wp ! should be better coded 83 115 84 ! ---------------------------! 116 ! 2- Trends due to forcing !85 ! 1 - Trends due to forcing ! 117 86 ! ---------------------------! 118 87 ! they must be kept outside an IF(iom_use) because of the call to dia_rst below 119 z_frc_volbot = r1_rho0 * glob_sum( 'icedia', -( wfx_ice(:,:) + wfx_snw(:,:) + wfx_err_sub(:,:) ) * e1e2t(:,:) ) * 1.e-9 ! freshwater flux ice/snow-ocean 120 z_frc_voltop = r1_rho0 * glob_sum( 'icedia', -( wfx_sub(:,:) + wfx_spr(:,:) ) * e1e2t(:,:) ) * 1.e-9 ! freshwater flux ice/snow-atm 121 z_frc_sal = r1_rho0 * glob_sum( 'icedia', - sfx(:,:) * e1e2t(:,:) ) * 1.e-9 ! salt fluxes ice/snow-ocean 122 z_frc_tembot = glob_sum( 'icedia', qt_oce_ai(:,:) * e1e2t(:,:) ) * 1.e-20 ! heat on top of ocean (and below ice) 123 z_frc_temtop = glob_sum( 'icedia', qt_atm_oi(:,:) * e1e2t(:,:) ) * 1.e-20 ! heat on top of ice-coean 124 ! 125 frc_voltop = frc_voltop + z_frc_voltop * rDt_ice ! km3 126 frc_volbot = frc_volbot + z_frc_volbot * rDt_ice ! km3 127 frc_sal = frc_sal + z_frc_sal * rDt_ice ! km3*pss 128 frc_temtop = frc_temtop + z_frc_temtop * rDt_ice ! 1.e20 J 129 frc_tembot = frc_tembot + z_frc_tembot * rDt_ice ! 1.e20 J 130 131 CALL iom_put( 'ibgfrcvoltop' , frc_voltop ) ! vol forcing ice/snw-atm (km3 equivalent ocean water) 132 CALL iom_put( 'ibgfrcvolbot' , frc_volbot ) ! vol forcing ice/snw-ocean (km3 equivalent ocean water) 133 CALL iom_put( 'ibgfrcsal' , frc_sal ) ! sal - forcing (psu*km3 equivalent ocean water) 134 CALL iom_put( 'ibgfrctemtop' , frc_temtop ) ! heat on top of ice/snw/ocean (1.e20 J) 135 CALL iom_put( 'ibgfrctembot' , frc_tembot ) ! heat on top of ocean(below ice) (1.e20 J) 136 137 IF( iom_use('ibgfrchfxtop') .OR. iom_use('ibgfrchfxbot') ) THEN 138 CALL iom_put( 'ibgfrchfxtop' , frc_temtop * z1_e1e2 * 1.e-20 * kt*rn_Dt ) ! heat on top of ice/snw/ocean (W/m2) 139 CALL iom_put( 'ibgfrchfxbot' , frc_tembot * z1_e1e2 * 1.e-20 * kt*rn_Dt ) ! heat on top of ocean(below ice) (W/m2) 140 ENDIF 88 ztmp(:,:,1) = - ( wfx_ice(:,:) + wfx_snw(:,:) + wfx_err_sub(:,:) ) * e1e2t(:,:) ! freshwater flux ice/snow-ocean 89 ztmp(:,:,2) = - ( wfx_sub(:,:) + wfx_spr(:,:) ) * e1e2t(:,:) ! freshwater flux ice/snow-atm 90 ztmp(:,:,3) = - sfx (:,:) * e1e2t(:,:) ! salt fluxes ice/snow-ocean 91 ztmp(:,:,4) = qt_atm_oi(:,:) * e1e2t(:,:) ! heat on top of ice-ocean 92 ztmp(:,:,5) = qt_oce_ai(:,:) * e1e2t(:,:) ! heat on top of ocean (and below ice) 93 94 ! ----------------------- ! 95 ! 2 - Contents ! 96 ! ----------------------- ! 97 IF( iom_use('ibgvol_tot' ) ) ztmp(:,:,6 ) = vt_i (:,:) * e1e2t(:,:) ! ice volume 98 IF( iom_use('sbgvol_tot' ) ) ztmp(:,:,7 ) = vt_s (:,:) * e1e2t(:,:) ! snow volume 99 IF( iom_use('ibgarea_tot') ) ztmp(:,:,8 ) = at_i (:,:) * e1e2t(:,:) ! area 100 IF( iom_use('ibgsalt_tot') ) ztmp(:,:,9 ) = st_i (:,:) * e1e2t(:,:) ! salt content 101 IF( iom_use('ibgheat_tot') ) ztmp(:,:,10) = et_i (:,:) * e1e2t(:,:) ! heat content 102 IF( iom_use('sbgheat_tot') ) ztmp(:,:,11) = et_s (:,:) * e1e2t(:,:) ! heat content 103 IF( iom_use('ipbgvol_tot') ) ztmp(:,:,12) = vt_ip(:,:) * e1e2t(:,:) ! ice pond volume 104 IF( iom_use('ilbgvol_tot') ) ztmp(:,:,13) = vt_il(:,:) * e1e2t(:,:) ! ice pond lid volume 141 105 142 106 ! ---------------------------------- ! 143 107 ! 3 - Content variations and drifts ! 144 108 ! ---------------------------------- ! 145 IF( iom_use('ibgvolume') .OR. iom_use('ibgsaltco') .OR. iom_use('ibgheatco') .OR. iom_use('ibgheatfx') ) THEN 146 147 zdiff_vol = r1_rho0 * glob_sum( 'icedia', ( rhoi*vt_i(:,:) + rhos*vt_s(:,:) - vol_loc_ini(:,:) ) * e1e2t(:,:) ) * 1.e-9 ! freshwater trend (km3) 148 zdiff_sal = r1_rho0 * glob_sum( 'icedia', ( rhoi*st_i(:,:) - sal_loc_ini(:,:) ) * e1e2t(:,:) ) * 1.e-9 ! salt content trend (km3*pss) 149 zdiff_tem = glob_sum( 'icedia', ( et_i(:,:) + et_s(:,:) - tem_loc_ini(:,:) ) * e1e2t(:,:) ) * 1.e-20 ! heat content trend (1.e20 J) 150 ! + SUM( qevap_ice * a_i_b, dim=3 ) !! clem: I think this term should not be there (but needs a check) 151 152 zdiff_vol = zdiff_vol - ( frc_voltop + frc_volbot ) 153 zdiff_sal = zdiff_sal - frc_sal 154 zdiff_tem = zdiff_tem - ( frc_tembot - frc_temtop ) 155 156 CALL iom_put( 'ibgvolume' , zdiff_vol ) ! ice/snow volume drift (km3 equivalent ocean water) 157 CALL iom_put( 'ibgsaltco' , zdiff_sal ) ! ice salt content drift (psu*km3 equivalent ocean water) 158 CALL iom_put( 'ibgheatco' , zdiff_tem ) ! ice/snow heat content drift (1.e20 J) 159 ! 160 ENDIF 161 109 IF( iom_use('ibgvolume') ) ztmp(:,:,14) = ( rhoi*vt_i(:,:) + rhos*vt_s(:,:) - vol_loc_ini(:,:) ) * e1e2t(:,:) ! freshwater trend 110 IF( iom_use('ibgsaltco') ) ztmp(:,:,15) = ( rhoi*st_i(:,:) - sal_loc_ini(:,:) ) * e1e2t(:,:) ! salt content trend 111 IF( iom_use('ibgheatco') .OR. iom_use('ibgheatfx') ) & 112 & ztmp(:,:,16) = ( et_i(:,:) + et_s(:,:) - tem_loc_ini(:,:) ) * e1e2t(:,:) ! heat content trend 113 114 ! global sum 115 zbg(1:16) = glob_sum_vec( 'icedia', ztmp(:,:,1:16) ) 116 117 ! change units for trends 118 zbg(1) = zbg(1) * r1_rho0 * 1.e-9 * rDt_ice ! freshwater flux ice/snow-ocean (km3) 119 zbg(2) = zbg(2) * r1_rho0 * 1.e-9 * rDt_ice ! freshwater flux ice/snow-atm (km3) 120 zbg(3) = zbg(3) * r1_rho0 * 1.e-9 * rDt_ice ! salt fluxes ice/snow-ocean (km3*pss) 121 zbg(4) = zbg(4) * 1.e-20 * rDt_ice ! heat on top of ice-ocean (1.e20 J) 122 zbg(5) = zbg(5) * 1.e-20 * rDt_ice ! heat on top of ocean (and below ice) (1.e20 J) 123 ! cumulative sum 124 frc_voltop = frc_voltop + zbg(1) 125 frc_volbot = frc_volbot + zbg(2) 126 frc_sal = frc_sal + zbg(3) 127 frc_temtop = frc_temtop + zbg(4) 128 frc_tembot = frc_tembot + zbg(5) 129 130 ! change units for contents 131 zbg(6) = zbg(6) * 1.e-9 ! ice volume (km3) 132 zbg(7) = zbg(7) * 1.e-9 ! snw volume (km3) 133 zbg(8) = zbg(8) * 1.e-6 ! ice area (km2) 134 zbg(9) = zbg(9) * 1.e-9 ! salt content (km3*pss) 135 zbg(10) = zbg(10) * 1.e-20 ! ice heat content (1.e20 J) 136 zbg(11) = zbg(11) * 1.e-20 ! snw heat content (1.e20 J) 137 zbg(12) = zbg(12) * 1.e-9 ! pnd volume (km3) 138 zbg(13) = zbg(13) * 1.e-9 ! pnd lid volume (km3) 139 140 ! change units for trends 141 zbg(14) = zbg(14) * r1_rho0 * 1.e-9 ! freshwater trend (km3) 142 zbg(15) = zbg(15) * r1_rho0 * 1.e-9 ! salt content trend (km3*pss) 143 zbg(16) = zbg(16) * 1.e-20 ! heat content trend (1.e20 J) 144 ! difference 145 zbg(14) = zbg(14) - ( frc_voltop + frc_volbot ) 146 zbg(15) = zbg(15) - frc_sal 147 zbg(16) = zbg(16) - ( frc_tembot - frc_temtop ) 148 149 ! outputs 150 CALL iom_put( 'ibgfrcvoltop' , frc_voltop ) ! vol forcing ice/snw-atm (km3 equivalent ocean water) 151 CALL iom_put( 'ibgfrcvolbot' , frc_volbot ) ! vol forcing ice/snw-ocean (km3 equivalent ocean water) 152 CALL iom_put( 'ibgfrcsal' , frc_sal ) ! sal forcing (psu*km3 equivalent ocean water) 153 CALL iom_put( 'ibgfrctemtop' , frc_temtop ) ! heat on top of ice/snw/ocean (1.e20 J) 154 CALL iom_put( 'ibgfrctembot' , frc_tembot ) ! heat on top of ocean(below ice) (1.e20 J) 155 CALL iom_put( 'ibgfrchfxtop' , frc_temtop * z1_e1e2 * 1.e-20 * kt*rn_Dt ) ! heat on top of ice/snw/ocean (W/m2) 156 CALL iom_put( 'ibgfrchfxbot' , frc_tembot * z1_e1e2 * 1.e-20 * kt*rn_Dt ) ! heat on top of ocean(below ice) (W/m2) 157 158 CALL iom_put( 'ibgvol_tot' , zbg(6) ) 159 CALL iom_put( 'sbgvol_tot' , zbg(7) ) 160 CALL iom_put( 'ibgarea_tot' , zbg(8) ) 161 CALL iom_put( 'ibgsalt_tot' , zbg(9) ) 162 CALL iom_put( 'ibgheat_tot' , zbg(10) ) 163 CALL iom_put( 'sbgheat_tot' , zbg(11) ) 164 CALL iom_put( 'ipbgvol_tot' , zbg(12) ) 165 CALL iom_put( 'ilbgvol_tot' , zbg(13) ) 166 167 CALL iom_put( 'ibgvolume' , zbg(14) ) ! ice/snow volume drift (km3 equivalent ocean water) 168 CALL iom_put( 'ibgsaltco' , zbg(15) ) ! ice salt content drift (psu*km3 equivalent ocean water) 169 CALL iom_put( 'ibgheatco' , zbg(16) ) ! ice/snow heat content drift (1.e20 J) 170 ! 171 ! restarts 162 172 IF( lrst_ice ) CALL ice_dia_rst( 'WRITE', kt_ice ) 163 173 ! -
NEMO/trunk/src/OCE/DIA/diahsb.F90
r15004 r15048 83 83 REAL(wp) :: z_wn_trd_t , z_wn_trd_s ! - - 84 84 REAL(wp) :: z_ssh_hc , z_ssh_sc ! - - 85 REAL(wp), DIMENSION(jpi,jpj) :: z2d0, z2d1 ! 2D workspace 86 REAL(wp), DIMENSION(jpi,jpj,jpkm1) :: zwrk ! 3D workspace 85 REAL(wp), DIMENSION(jpi,jpj,13) :: ztmp 86 REAL(wp), DIMENSION(jpi,jpj,jpkm1,4) :: ztmpk 87 REAL(wp), DIMENSION(17) :: zbg 87 88 !!--------------------------------------------------------------------------- 88 89 IF( ln_timing ) CALL timing_start('dia_hsb') 89 90 ! 91 ztmp (:,:,:) = 0._wp ! should be better coded 92 ztmpk(:,:,:,:) = 0._wp ! should be better coded 93 ! 90 94 ts(:,:,:,1,Kmm) = ts(:,:,:,1,Kmm) * tmask(:,:,:) ; ts(:,:,:,1,Kbb) = ts(:,:,:,1,Kbb) * tmask(:,:,:) ; 91 95 ts(:,:,:,2,Kmm) = ts(:,:,:,2,Kmm) * tmask(:,:,:) ; ts(:,:,:,2,Kbb) = ts(:,:,:,2,Kbb) * tmask(:,:,:) ; 96 ! 92 97 ! ------------------------- ! 93 98 ! 1 - Trends due to forcing ! 94 99 ! ------------------------- ! 95 z_frc_trd_v = r1_rho0 * glob_sum( 'diahsb', - ( emp(:,:) - rnf(:,:) - fwfisf_cav(:,:) - fwfisf_par(:,:) ) * surf(:,:) ) ! volume fluxes 96 z_frc_trd_t = glob_sum( 'diahsb', sbc_tsc(:,:,jp_tem) * surf(:,:) ) ! heat fluxes 97 z_frc_trd_s = glob_sum( 'diahsb', sbc_tsc(:,:,jp_sal) * surf(:,:) ) ! salt fluxes 98 ! ! Add runoff heat & salt input 99 IF( ln_rnf ) z_frc_trd_t = z_frc_trd_t + glob_sum( 'diahsb', rnf_tsc(:,:,jp_tem) * surf(:,:) ) 100 IF( ln_rnf_sal) z_frc_trd_s = z_frc_trd_s + glob_sum( 'diahsb', rnf_tsc(:,:,jp_sal) * surf(:,:) ) 101 ! ! Add ice shelf heat & salt input 102 IF( ln_isf ) z_frc_trd_t = z_frc_trd_t & 103 & + glob_sum( 'diahsb', ( risf_cav_tsc(:,:,jp_tem) + risf_par_tsc(:,:,jp_tem) ) * surf(:,:) ) 104 ! ! Add penetrative solar radiation 105 IF( ln_traqsr ) z_frc_trd_t = z_frc_trd_t + r1_rho0_rcp * glob_sum( 'diahsb', qsr (:,:) * surf(:,:) ) 106 ! ! Add geothermal heat flux 107 IF( ln_trabbc ) z_frc_trd_t = z_frc_trd_t + glob_sum( 'diahsb', qgh_trd0(:,:) * surf(:,:) ) 108 ! 109 IF( ln_linssh ) THEN 100 ! prepare trends 101 ztmp(:,:,1) = - r1_rho0 * ( emp(:,:) - rnf(:,:) - fwfisf_cav(:,:) - fwfisf_par(:,:) ) * surf(:,:) ! volume 102 ztmp(:,:,2) = sbc_tsc(:,:,jp_tem) * surf(:,:) ! heat 103 ztmp(:,:,3) = sbc_tsc(:,:,jp_sal) * surf(:,:) ! salt 104 IF( ln_rnf ) ztmp(:,:,4) = rnf_tsc(:,:,jp_tem) * surf(:,:) ! runoff temp 105 IF( ln_rnf_sal ) ztmp(:,:,5) = rnf_tsc(:,:,jp_sal) * surf(:,:) ! runoff salt 106 IF( ln_isf ) ztmp(:,:,6) = ( risf_cav_tsc(:,:,jp_tem) + risf_par_tsc(:,:,jp_tem) ) * surf(:,:) ! isf temp 107 IF( ln_traqsr ) ztmp(:,:,7) = r1_rho0_rcp * qsr(:,:) * surf(:,:) ! penetrative solar radiation 108 IF( ln_trabbc ) ztmp(:,:,8) = qgh_trd0(:,:) * surf(:,:) ! geothermal heat 109 ! 110 IF( ln_linssh ) THEN ! Advection flux through fixed surface (z=0) 110 111 IF( ln_isfcav ) THEN 111 112 DO ji=1,jpi 112 113 DO jj=1,jpj 113 z 2d0(ji,jj) =surf(ji,jj) * ww(ji,jj,mikt(ji,jj)) * ts(ji,jj,mikt(ji,jj),jp_tem,Kbb)114 z 2d1(ji,jj) =surf(ji,jj) * ww(ji,jj,mikt(ji,jj)) * ts(ji,jj,mikt(ji,jj),jp_sal,Kbb)114 ztmp(ji,jj,9 ) = - surf(ji,jj) * ww(ji,jj,mikt(ji,jj)) * ts(ji,jj,mikt(ji,jj),jp_tem,Kbb) 115 ztmp(ji,jj,10) = - surf(ji,jj) * ww(ji,jj,mikt(ji,jj)) * ts(ji,jj,mikt(ji,jj),jp_sal,Kbb) 115 116 END DO 116 117 END DO 117 118 ELSE 118 z 2d0(:,:) =surf(:,:) * ww(:,:,1) * ts(:,:,1,jp_tem,Kbb)119 z 2d1(:,:) =surf(:,:) * ww(:,:,1) * ts(:,:,1,jp_sal,Kbb)119 ztmp(:,:,9 ) = - surf(:,:) * ww(:,:,1) * ts(:,:,1,jp_tem,Kbb) 120 ztmp(:,:,10) = - surf(:,:) * ww(:,:,1) * ts(:,:,1,jp_sal,Kbb) 120 121 END IF 121 z_wn_trd_t = - glob_sum( 'diahsb', z2d0 ) 122 z_wn_trd_s = - glob_sum( 'diahsb', z2d1 ) 123 ENDIF 124 122 ENDIF 123 124 ! global sum 125 zbg(1:10) = glob_sum_vec( 'dia_hsb', ztmp(:,:,1:10) ) 126 127 ! adding up 128 z_frc_trd_v = zbg(1) ! volume fluxes 129 z_frc_trd_t = zbg(2) ! heat fluxes 130 z_frc_trd_s = zbg(3) ! salt fluxes 131 IF( ln_rnf ) z_frc_trd_t = z_frc_trd_t + zbg(4) ! runoff heat 132 IF( ln_rnf_sal) z_frc_trd_s = z_frc_trd_s + zbg(5) ! runoff salt 133 IF( ln_isf ) z_frc_trd_t = z_frc_trd_t + zbg(6) ! isf heat 134 IF( ln_traqsr ) z_frc_trd_t = z_frc_trd_t + zbg(7) ! penetrative solar flux 135 IF( ln_trabbc ) z_frc_trd_t = z_frc_trd_t + zbg(8) ! geothermal heat 136 ! 125 137 frc_v = frc_v + z_frc_trd_v * rn_Dt 126 138 frc_t = frc_t + z_frc_trd_t * rn_Dt … … 128 140 ! ! Advection flux through fixed surface (z=0) 129 141 IF( ln_linssh ) THEN 142 z_wn_trd_t = zbg(9) 143 z_wn_trd_s = zbg(10) 144 ! 130 145 frc_wn_t = frc_wn_t + z_wn_trd_t * rn_Dt 131 146 frc_wn_s = frc_wn_s + z_wn_trd_s * rn_Dt 132 147 ENDIF 133 148 134 ! ------------------------ !135 ! 2 - Content variations !136 ! ------------------------ !149 ! --------------------------------- ! 150 ! 2 - Content variations with ssh ! 151 ! --------------------------------- ! 137 152 ! glob_sum_full is needed because you keep the full interior domain to compute the sum (iscpl) 138 153 ! 139 154 ! ! volume variation (calculated with ssh) 140 z diff_v1 = glob_sum_full( 'diahsb', surf(:,:)*ssh(:,:,Kmm) - surf_ini(:,:)*ssh_ini(:,:))155 ztmp(:,:,11) = surf(:,:)*ssh(:,:,Kmm) - surf_ini(:,:)*ssh_ini(:,:) 141 156 142 157 ! ! heat & salt content variation (associated with ssh) … … 145 160 DO ji = 1, jpi 146 161 DO jj = 1, jpj 147 z 2d0(ji,jj) = surf(ji,jj) * ( ts(ji,jj,mikt(ji,jj),jp_tem,Kmm) * ssh(ji,jj,Kmm) - ssh_hc_loc_ini(ji,jj) )148 z 2d1(ji,jj) = surf(ji,jj) * ( ts(ji,jj,mikt(ji,jj),jp_sal,Kmm) * ssh(ji,jj,Kmm) - ssh_sc_loc_ini(ji,jj) )162 ztmp(ji,jj,12) = surf(ji,jj) * ( ts(ji,jj,mikt(ji,jj),jp_tem,Kmm) * ssh(ji,jj,Kmm) - ssh_hc_loc_ini(ji,jj) ) 163 ztmp(ji,jj,13) = surf(ji,jj) * ( ts(ji,jj,mikt(ji,jj),jp_sal,Kmm) * ssh(ji,jj,Kmm) - ssh_sc_loc_ini(ji,jj) ) 149 164 END DO 150 165 END DO 151 166 ELSE ! no under ice-shelf seas 152 z 2d0(:,:) = surf(:,:) * ( ts(:,:,1,jp_tem,Kmm) * ssh(:,:,Kmm) - ssh_hc_loc_ini(:,:) )153 z 2d1(:,:) = surf(:,:) * ( ts(:,:,1,jp_sal,Kmm) * ssh(:,:,Kmm) - ssh_sc_loc_ini(:,:) )167 ztmp(:,:,12) = surf(:,:) * ( ts(:,:,1,jp_tem,Kmm) * ssh(:,:,Kmm) - ssh_hc_loc_ini(:,:) ) 168 ztmp(:,:,13) = surf(:,:) * ( ts(:,:,1,jp_sal,Kmm) * ssh(:,:,Kmm) - ssh_sc_loc_ini(:,:) ) 154 169 END IF 155 z_ssh_hc = glob_sum_full( 'diahsb', z2d0 ) 156 z_ssh_sc = glob_sum_full( 'diahsb', z2d1 ) 157 ENDIF 158 ! 159 DO jk = 1, jpkm1 ! volume variation (calculated with scale factors) 160 zwrk(:,:,jk) = surf (:,:) * e3t (:,:,jk,Kmm)*tmask (:,:,jk) & 161 & - surf_ini(:,:) * e3t_ini(:,:,jk )*tmask_ini(:,:,jk) 170 ENDIF 171 172 ! global sum 173 zbg(11:13) = glob_sum_full_vec( 'dia_hsb', ztmp(:,:,11:13) ) 174 175 zdiff_v1 = zbg(11) 176 ! ! heat & salt content variation (associated with ssh) 177 IF( ln_linssh ) THEN ! linear free surface case 178 z_ssh_hc = zbg(12) 179 z_ssh_sc = zbg(13) 180 ENDIF 181 ! 182 ! --------------------------------- ! 183 ! 3 - Content variations with e3t ! 184 ! --------------------------------- ! 185 ! glob_sum_full is needed because you keep the full interior domain to compute the sum (iscpl) 186 ! 187 DO jk = 1, jpkm1 ! volume 188 ztmpk(:,:,jk,1) = surf (:,:) * e3t (:,:,jk,Kmm)*tmask (:,:,jk) & 189 & - surf_ini(:,:) * e3t_ini(:,:,jk )*tmask_ini(:,:,jk) 162 190 END DO 163 zdiff_v2 = glob_sum_full( 'diahsb', zwrk(:,:,:) ) ! glob_sum_full needed as tmask and tmask_ini could be different 164 DO jk = 1, jpkm1 ! heat content variation 165 zwrk(:,:,jk) = ( surf (:,:) * e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_tem,Kmm) & 166 & - surf_ini(:,:) * hc_loc_ini(:,:,jk) ) 191 DO jk = 1, jpkm1 ! heat 192 ztmpk(:,:,jk,2) = ( surf (:,:) * e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_tem,Kmm) & 193 & - surf_ini(:,:) * hc_loc_ini(:,:,jk) ) 167 194 END DO 168 zdiff_hc = glob_sum_full( 'diahsb', zwrk(:,:,:) ) 169 DO jk = 1, jpkm1 ! salt content variation 170 zwrk(:,:,jk) = ( surf (:,:) * e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_sal,Kmm) & 171 & - surf_ini(:,:) * sc_loc_ini(:,:,jk) ) 195 DO jk = 1, jpkm1 ! salt 196 ztmpk(:,:,jk,3) = ( surf (:,:) * e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_sal,Kmm) & 197 & - surf_ini(:,:) * sc_loc_ini(:,:,jk) ) 172 198 END DO 173 zdiff_sc = glob_sum_full( 'diahsb', zwrk(:,:,:) ) 199 DO jk = 1, jpkm1 ! total ocean volume 200 ztmpk(:,:,jk,4) = surf(:,:) * e3t(:,:,jk,Kmm) * tmask(:,:,jk) 201 END DO 202 203 ! global sum 204 zbg(14:17) = glob_sum_full_vec( 'dia_hsb', ztmpk(:,:,:,1:4) ) 205 206 zdiff_v2 = zbg(14) ! glob_sum_full needed as tmask and tmask_ini could be different 207 zdiff_hc = zbg(15) 208 zdiff_sc = zbg(16) 209 zvol_tot = zbg(17) 174 210 175 211 ! ------------------------ ! 176 ! 3- Drifts !212 ! 4 - Drifts ! 177 213 ! ------------------------ ! 178 214 zdiff_v1 = zdiff_v1 - frc_v … … 186 222 zerr_sc1 = z_ssh_sc - frc_wn_s 187 223 ENDIF 188 189 ! ----------------------- !190 ! 4 - Diagnostics writing !191 ! ----------------------- !192 DO jk = 1, jpkm1 ! total ocean volume (calculated with scale factors)193 zwrk(:,:,jk) = surf(:,:) * e3t(:,:,jk,Kmm) * tmask(:,:,jk)194 END DO195 zvol_tot = glob_sum( 'diahsb', zwrk(:,:,:) )196 224 197 225 !!gm to be added ? -
NEMO/trunk/src/OCE/LBC/lbc_lnk_call_generic.h90
r14433 r15048 18 18 #endif 19 19 20 SUBROUTINE lbc_lnk_call_/**/XD/**/_/**/PRECISION( &20 SUBROUTINE lbc_lnk_call_/**/XD/**/_/**/PRECISION( & 21 21 & cdname & 22 22 & , pt1 , cdna1 , psgn1 , pt2 , cdna2 , psgn2 , pt3 , cdna3 , psgn3 , pt4 , cdna4 , psgn4 & … … 24 24 & , pt9 , cdna9 , psgn9 , pt10, cdna10, psgn10, pt11, cdna11, psgn11, pt12, cdna12, psgn12 & 25 25 & , pt13, cdna13, psgn13, pt14, cdna14, psgn14, pt15, cdna15, psgn15, pt16, cdna16, psgn16 & 26 & , pt17, cdna17, psgn17, pt18, cdna18, psgn18, pt19, cdna19, psgn19, pt20, cdna20, psgn20 & 27 & , pt21, cdna21, psgn21, pt22, cdna22, psgn22, pt23, cdna23, psgn23, pt24, cdna24, psgn24 & 28 & , pt25, cdna25, psgn25, pt26, cdna26, psgn26, pt27, cdna27, psgn27, pt28, cdna28, psgn28 & 29 & , pt29, cdna29, psgn29, pt30, cdna30, psgn30 & 26 30 & , kfillmode, pfillval, khls, lsend, lrecv, ld4only ) 27 31 !!--------------------------------------------------------------------- 28 32 CHARACTER(len=*) , INTENT(in ) :: cdname ! name of the calling subroutine 29 33 REAL(PRECISION), DIMENSION(DIMS) , TARGET, CONTIGUOUS, INTENT(inout) :: pt1 ! arrays on which the lbc is applied 30 REAL(PRECISION), DIMENSION(DIMS), OPTIONAL, TARGET, CONTIGUOUS, INTENT(inout) :: pt2 , pt3 , pt4 , pt5 , pt6 , pt7 , pt8 , pt9, & 31 & pt10, pt11, pt12, pt13, pt14, pt15, pt16 34 REAL(PRECISION), DIMENSION(DIMS), OPTIONAL, TARGET, CONTIGUOUS, INTENT(inout) :: pt2 , pt3 , pt4 , pt5 , pt6 , pt7 , pt8 , & 35 & pt9 , pt10, pt11, pt12, pt13, pt14, pt15, & 36 & pt16, pt17, pt18, pt19, pt20, pt21, pt22, & 37 & pt23, pt24, pt25, pt26, pt27, pt28, pt29, & 38 & pt30 32 39 CHARACTER(len=1) , INTENT(in ) :: cdna1 ! nature of pt2D. array grid-points 33 CHARACTER(len=1) , OPTIONAL , INTENT(in ) :: cdna2 , cdna3 , cdna4 , cdna5 , cdna6 , cdna7 , cdna8 , cdna9, & 34 & cdna10, cdna11, cdna12, cdna13, cdna14, cdna15, cdna16 40 CHARACTER(len=1) , OPTIONAL , INTENT(in ) :: cdna2 , cdna3 , cdna4 , cdna5 , cdna6 , cdna7 , cdna8 , & 41 & cdna9 , cdna10, cdna11, cdna12, cdna13, cdna14, cdna15, & 42 & cdna16, cdna17, cdna18, cdna19, cdna20, cdna21, cdna22, & 43 & cdna23, cdna24, cdna25, cdna26, cdna27, cdna28, cdna29, & 44 & cdna30 35 45 REAL(PRECISION) , INTENT(in ) :: psgn1 ! sign used across the north fold 36 REAL(PRECISION) , OPTIONAL , INTENT(in ) :: psgn2 , psgn3 , psgn4 , psgn5 , psgn6 , psgn7 , psgn8 , psgn9, & 37 & psgn10, psgn11, psgn12, psgn13, psgn14, psgn15, psgn16 46 REAL(PRECISION) , OPTIONAL , INTENT(in ) :: psgn2 , psgn3 , psgn4 , psgn5 , psgn6 , psgn7 , psgn8 , & 47 & psgn9 , psgn10, psgn11, psgn12, psgn13, psgn14, psgn15, & 48 & psgn16, psgn17, psgn18, psgn19, psgn20, psgn21, psgn22, & 49 & psgn23, psgn24, psgn25, psgn26, psgn27, psgn28, psgn29, & 50 & psgn30 38 51 INTEGER , OPTIONAL , INTENT(in ) :: kfillmode ! filling method for halo over land (default = constant) 39 52 REAL(PRECISION) , OPTIONAL , INTENT(in ) :: pfillval ! background value (used at closed boundaries) … … 43 56 !! 44 57 INTEGER :: kfld ! number of elements that will be attributed 45 TYPE(PTR_4d_/**/PRECISION), DIMENSION( 16) :: ptab_ptr ! pointer array46 CHARACTER(len=1) , DIMENSION( 16) :: cdna_ptr ! nature of ptab_ptr grid-points47 REAL(PRECISION) , DIMENSION( 16) :: psgn_ptr ! sign used across the north fold boundary58 TYPE(PTR_4d_/**/PRECISION), DIMENSION(30) :: ptab_ptr ! pointer array 59 CHARACTER(len=1) , DIMENSION(30) :: cdna_ptr ! nature of ptab_ptr grid-points 60 REAL(PRECISION) , DIMENSION(30) :: psgn_ptr ! sign used across the north fold boundary 48 61 !!--------------------------------------------------------------------- 49 62 ! … … 69 82 IF( PRESENT(psgn15) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt15, cdna15, psgn15, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 70 83 IF( PRESENT(psgn16) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt16, cdna16, psgn16, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 84 IF( PRESENT(psgn17) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt17, cdna17, psgn17, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 85 IF( PRESENT(psgn18) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt18, cdna18, psgn18, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 86 IF( PRESENT(psgn19) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt19, cdna19, psgn19, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 87 IF( PRESENT(psgn20) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt20, cdna20, psgn20, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 88 IF( PRESENT(psgn21) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt21, cdna21, psgn21, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 89 IF( PRESENT(psgn22) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt22, cdna22, psgn22, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 90 IF( PRESENT(psgn23) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt23, cdna23, psgn16, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 91 IF( PRESENT(psgn24) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt24, cdna24, psgn24, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 92 IF( PRESENT(psgn25) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt25, cdna25, psgn25, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 93 IF( PRESENT(psgn26) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt26, cdna26, psgn26, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 94 IF( PRESENT(psgn27) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt27, cdna27, psgn27, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 95 IF( PRESENT(psgn28) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt28, cdna28, psgn28, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 96 IF( PRESENT(psgn29) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt29, cdna29, psgn29, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 97 IF( PRESENT(psgn30) ) CALL load_ptr_/**/XD/**/_/**/PRECISION( pt30, cdna30, psgn30, ptab_ptr, cdna_ptr, psgn_ptr, kfld ) 71 98 ! 72 99 IF( nn_comm == 1 ) THEN -
NEMO/trunk/src/OCE/lib_fortran.F90
r14433 r15048 32 32 PUBLIC DDPDD ! also used in closea module 33 33 PUBLIC glob_min, glob_max 34 PUBLIC glob_sum_vec 35 PUBLIC glob_sum_full_vec 34 36 #if defined key_nosignedzero 35 37 PUBLIC SIGN … … 53 55 INTERFACE glob_max 54 56 MODULE PROCEDURE glob_max_2d, glob_max_3d 57 END INTERFACE 58 INTERFACE glob_sum_vec 59 MODULE PROCEDURE glob_sum_vec_3d, glob_sum_vec_4d 60 END INTERFACE 61 INTERFACE glob_sum_full_vec 62 MODULE PROCEDURE glob_sum_full_vec_3d, glob_sum_full_vec_4d 55 63 END INTERFACE 56 64 … … 313 321 314 322 323 FUNCTION glob_sum_vec_3d( cdname, ptab ) RESULT( ptmp ) 324 !!---------------------------------------------------------------------- 325 CHARACTER(len=*), INTENT(in) :: cdname ! name of the calling subroutine 326 REAL(wp), INTENT(in) :: ptab(:,:,:) ! array on which operation is applied 327 REAL(wp), DIMENSION(SIZE(ptab,3)) :: ptmp 328 ! 329 COMPLEX(dp), DIMENSION(:), ALLOCATABLE :: ctmp 330 REAL(wp) :: ztmp 331 INTEGER :: ji , jj , jk ! dummy loop indices 332 INTEGER :: ipi, ipj, ipk ! dimensions 333 !!----------------------------------------------------------------------- 334 ! 335 ipi = SIZE(ptab,1) ! 1st dimension 336 ipj = SIZE(ptab,2) ! 2nd dimension 337 ipk = SIZE(ptab,3) ! 3rd dimension 338 ! 339 ALLOCATE( ctmp(ipk) ) 340 ! 341 DO jk = 1, ipk 342 ctmp(jk) = CMPLX( 0.e0, 0.e0, dp ) ! warning ctmp is cumulated 343 DO jj = 1, ipj 344 DO ji = 1, ipi 345 ztmp = ptab(ji,jj,jk) * tmask_i(ji,jj) 346 CALL DDPDD( CMPLX( ztmp, 0.e0, dp ), ctmp(jk) ) 347 END DO 348 END DO 349 END DO 350 CALL mpp_sum( cdname, ctmp(:) ) ! sum over the global domain 351 ! 352 ptmp = REAL( ctmp(:), wp ) 353 ! 354 DEALLOCATE( ctmp ) 355 ! 356 END FUNCTION glob_sum_vec_3d 357 358 FUNCTION glob_sum_vec_4d( cdname, ptab ) RESULT( ptmp ) 359 !!---------------------------------------------------------------------- 360 CHARACTER(len=*), INTENT(in) :: cdname ! name of the calling subroutine 361 REAL(wp), INTENT(in) :: ptab(:,:,:,:) ! array on which operation is applied 362 REAL(wp), DIMENSION(SIZE(ptab,4)) :: ptmp 363 ! 364 COMPLEX(dp), DIMENSION(:), ALLOCATABLE :: ctmp 365 REAL(wp) :: ztmp 366 INTEGER :: ji , jj , jk , jl ! dummy loop indices 367 INTEGER :: ipi, ipj, ipk, ipl ! dimensions 368 !!----------------------------------------------------------------------- 369 ! 370 ipi = SIZE(ptab,1) ! 1st dimension 371 ipj = SIZE(ptab,2) ! 2nd dimension 372 ipk = SIZE(ptab,3) ! 3rd dimension 373 ipl = SIZE(ptab,4) ! 4th dimension 374 ! 375 ALLOCATE( ctmp(ipl) ) 376 ! 377 DO jl = 1, ipl 378 ctmp(jl) = CMPLX( 0.e0, 0.e0, dp ) ! warning ctmp is cumulated 379 DO jk = 1, ipk 380 DO jj = 1, ipj 381 DO ji = 1, ipi 382 ztmp = ptab(ji,jj,jk,jl) * tmask_i(ji,jj) 383 CALL DDPDD( CMPLX( ztmp, 0.e0, dp ), ctmp(jl) ) 384 END DO 385 END DO 386 END DO 387 END DO 388 CALL mpp_sum( cdname, ctmp(:) ) ! sum over the global domain 389 ! 390 ptmp = REAL( ctmp(:), wp ) 391 ! 392 DEALLOCATE( ctmp ) 393 ! 394 END FUNCTION glob_sum_vec_4d 395 396 FUNCTION glob_sum_full_vec_3d( cdname, ptab ) RESULT( ptmp ) 397 !!---------------------------------------------------------------------- 398 CHARACTER(len=*), INTENT(in) :: cdname ! name of the calling subroutine 399 REAL(wp), INTENT(in) :: ptab(:,:,:) ! array on which operation is applied 400 REAL(wp), DIMENSION(SIZE(ptab,3)) :: ptmp 401 ! 402 COMPLEX(dp), DIMENSION(:), ALLOCATABLE :: ctmp 403 REAL(wp) :: ztmp 404 INTEGER :: ji , jj , jk ! dummy loop indices 405 INTEGER :: ipi, ipj, ipk ! dimensions 406 !!----------------------------------------------------------------------- 407 ! 408 ipi = SIZE(ptab,1) ! 1st dimension 409 ipj = SIZE(ptab,2) ! 2nd dimension 410 ipk = SIZE(ptab,3) ! 3rd dimension 411 ! 412 ALLOCATE( ctmp(ipk) ) 413 ! 414 DO jk = 1, ipk 415 ctmp(jk) = CMPLX( 0.e0, 0.e0, dp ) ! warning ctmp is cumulated 416 DO jj = 1, ipj 417 DO ji = 1, ipi 418 ztmp = ptab(ji,jj,jk) * tmask_h(ji,jj) 419 CALL DDPDD( CMPLX( ztmp, 0.e0, dp ), ctmp(jk) ) 420 END DO 421 END DO 422 END DO 423 CALL mpp_sum( cdname, ctmp(:) ) ! sum over the global domain 424 ! 425 ptmp = REAL( ctmp(:), wp ) 426 ! 427 DEALLOCATE( ctmp ) 428 ! 429 END FUNCTION glob_sum_full_vec_3d 430 431 FUNCTION glob_sum_full_vec_4d( cdname, ptab ) RESULT( ptmp ) 432 !!---------------------------------------------------------------------- 433 CHARACTER(len=*), INTENT(in) :: cdname ! name of the calling subroutine 434 REAL(wp), INTENT(in) :: ptab(:,:,:,:) ! array on which operation is applied 435 REAL(wp), DIMENSION(SIZE(ptab,4)) :: ptmp 436 ! 437 COMPLEX(dp), DIMENSION(:), ALLOCATABLE :: ctmp 438 REAL(wp) :: ztmp 439 INTEGER :: ji , jj , jk , jl ! dummy loop indices 440 INTEGER :: ipi, ipj, ipk, ipl ! dimensions 441 !!----------------------------------------------------------------------- 442 ! 443 ipi = SIZE(ptab,1) ! 1st dimension 444 ipj = SIZE(ptab,2) ! 2nd dimension 445 ipk = SIZE(ptab,3) ! 3rd dimension 446 ipl = SIZE(ptab,4) ! 4th dimension 447 ! 448 ALLOCATE( ctmp(ipl) ) 449 ! 450 DO jl = 1, ipl 451 ctmp(jl) = CMPLX( 0.e0, 0.e0, dp ) ! warning ctmp is cumulated 452 DO jk = 1, ipk 453 DO jj = 1, ipj 454 DO ji = 1, ipi 455 ztmp = ptab(ji,jj,jk,jl) * tmask_h(ji,jj) 456 CALL DDPDD( CMPLX( ztmp, 0.e0, dp ), ctmp(jl) ) 457 END DO 458 END DO 459 END DO 460 END DO 461 CALL mpp_sum( cdname, ctmp(:) ) ! sum over the global domain 462 ! 463 ptmp = REAL( ctmp(:), wp ) 464 ! 465 DEALLOCATE( ctmp ) 466 ! 467 END FUNCTION glob_sum_full_vec_4d 468 315 469 SUBROUTINE DDPDD( ydda, yddb ) 316 470 !!----------------------------------------------------------------------
Note: See TracChangeset
for help on using the changeset viewer.