Changeset 6483
- Timestamp:
- 2016-04-19T17:11:00+02:00 (8 years ago)
- Location:
- trunk/NEMOGCM/NEMO
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMOGCM/NEMO/LIM_SRC_3/ice.F90
r6478 r6483 234 234 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: u_oce, v_oce !: surface ocean velocity used in ice dynamics 235 235 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ahiu , ahiv !: hor. diffusivity coeff. at U- and V-points [m2/s] 236 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: pahu , pahv !: ice hor. eddy diffusivity coef. at U- and V-points 236 237 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ust2s, hicol !: friction velocity, ice collection thickness accreted in leads 237 238 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: strp1, strp2 !: strength at previous time steps … … 302 303 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hfx_res !: residual heat flux due to correction of ice thickness [W.m-2] 303 304 304 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: ftr_ice !: transmitted solar radiation under ice 305 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: pahu3D , pahv3D 305 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: ftr_ice !: transmitted solar radiation under ice 306 306 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rn_amax_2d !: maximum ice concentration 2d array 307 307 … … 429 429 ALLOCATE( u_oce (jpi,jpj) , v_oce (jpi,jpj) , & 430 430 & ahiu (jpi,jpj) , ahiv (jpi,jpj) , & 431 & pahu (jpi,jpj) , pahv (jpi,jpj) , & 431 432 & ust2s (jpi,jpj) , hicol (jpi,jpj) , & 432 433 & strp1 (jpi,jpj) , strp2 (jpi,jpj) , strength (jpi,jpj) , & … … 441 442 & wfx_res(jpi,jpj) , wfx_sni(jpi,jpj) , wfx_opw(jpi,jpj) , wfx_spr(jpi,jpj) , & 442 443 & afx_tot(jpi,jpj) , afx_thd(jpi,jpj), afx_dyn(jpi,jpj) , & 443 & fhtur (jpi,jpj) , ftr_ice(jpi,jpj,jpl), pahu3D(jpi,jpj,jpl+1), pahv3D(jpi,jpj,jpl+1),&444 & rn_amax_2d (jpi,jpj), qlead (jpi, jpj),&445 & sfx_res(jpi,jpj) , sfx_bri(jpi,jpj) , sfx_dyn(jpi,jpj) , 444 & fhtur (jpi,jpj) , ftr_ice(jpi,jpj,jpl), qlead (jpi,jpj) , & 445 & rn_amax_2d(jpi,jpj), & 446 & sfx_res(jpi,jpj) , sfx_bri(jpi,jpj) , sfx_dyn(jpi,jpj) , sfx_sub(jpi,jpj) , & 446 447 & sfx_bog(jpi,jpj) , sfx_bom(jpi,jpj) , sfx_sum(jpi,jpj) , sfx_sni(jpi,jpj) , sfx_opw(jpi,jpj) , & 447 448 & hfx_res(jpi,jpj) , hfx_snw(jpi,jpj) , hfx_sub(jpi,jpj) , hfx_err(jpi,jpj) , & … … 513 514 !!====================================================================== 514 515 END MODULE ice 515 -
trunk/NEMOGCM/NEMO/LIM_SRC_3/limhdf.F90
r6478 r6483 7 7 !! - ! 2001-05 (G. Madec, R. Hordoir) opa norm 8 8 !! 1.0 ! 2002-08 (C. Ethe) F90, free form 9 !! 3.0 ! 2015-08 (O. Tintó and M. Castrillo) added lim_hdf (multiple)10 9 !!---------------------------------------------------------------------- 11 10 #if defined key_lim3 … … 28 27 PRIVATE 29 28 30 PUBLIC lim_hdf ! called by lim_trp29 PUBLIC lim_hdf ! called by lim_trp 31 30 PUBLIC lim_hdf_init ! called by sbc_lim_init 32 31 … … 44 43 CONTAINS 45 44 46 SUBROUTINE lim_hdf( ptab , ihdf_vars , jpl , nlay_i)45 SUBROUTINE lim_hdf( ptab ) 47 46 !!------------------------------------------------------------------- 48 47 !! *** ROUTINE lim_hdf *** … … 55 54 !! ** Action : update ptab with the diffusive contribution 56 55 !!------------------------------------------------------------------- 57 INTEGER :: jpl, nlay_i, isize, ihdf_vars 58 REAL(wp), DIMENSION(:,:,:), INTENT( inout ),TARGET :: ptab ! Field on which the diffusion is applied 59 ! 60 INTEGER :: ji, jj, jk, jl , jm ! dummy loop indices 56 REAL(wp), DIMENSION(jpi,jpj), INTENT( inout ) :: ptab ! Field on which the diffusion is applied 57 ! 58 INTEGER :: ji, jj ! dummy loop indices 61 59 INTEGER :: iter, ierr ! local integers 62 REAL(wp) :: zrlxint ! local scalars 63 REAL(wp), POINTER , DIMENSION ( : ) :: zconv ! local scalars 64 REAL(wp), POINTER , DIMENSION(:,:,:) :: zrlx,zdiv0, ztab0 65 REAL(wp), POINTER , DIMENSION(:,:) :: zflu, zflv, zdiv 60 REAL(wp) :: zrlxint, zconv ! local scalars 61 REAL(wp), POINTER, DIMENSION(:,:) :: zrlx, zflu, zflv, zdiv0, zdiv, ztab0 66 62 CHARACTER(lc) :: charout ! local character 67 63 REAL(wp), PARAMETER :: zrelax = 0.5_wp ! relaxation constant for iterative procedure … … 69 65 INTEGER , PARAMETER :: its = 100 ! Maximum number of iteration 70 66 !!------------------------------------------------------------------- 71 TYPE(arrayptr) , ALLOCATABLE, DIMENSION(:) :: pt2d_array, zrlx_array 72 CHARACTER(len=1) , ALLOCATABLE, DIMENSION(:) :: type_array ! define the nature of ptab array grid-points 73 ! ! = T , U , V , F , W and I points 74 REAL(wp) , ALLOCATABLE, DIMENSION(:) :: psgn_array ! =-1 the sign change across the north fold boundary 75 76 !!--------------------------------------------------------------------- 67 68 CALL wrk_alloc( jpi, jpj, zrlx, zflu, zflv, zdiv0, zdiv, ztab0 ) 77 69 78 70 ! !== Initialisation ==! 79 ! +1 open water diffusion80 isize = jpl*(ihdf_vars+nlay_i)+181 ALLOCATE( zconv (isize) )82 ALLOCATE( pt2d_array(isize) , zrlx_array(isize) )83 ALLOCATE( type_array(isize) )84 ALLOCATE( psgn_array(isize) )85 86 CALL wrk_alloc( jpi, jpj, isize, zrlx, zdiv0, ztab0 )87 CALL wrk_alloc( jpi, jpj, zflu, zflv, zdiv )88 89 DO jk= 1 , isize90 pt2d_array(jk)%pt2d=>ptab(:,:,jk)91 zrlx_array(jk)%pt2d=>zrlx(:,:,jk)92 type_array(jk)='T'93 psgn_array(jk)=1.94 END DO95 96 71 ! 97 72 IF( linit ) THEN ! Metric coefficient (compute at the first call and saved in efact) … … 99 74 IF( lk_mpp ) CALL mpp_sum( ierr ) 100 75 IF( ierr /= 0 ) CALL ctl_stop( 'STOP', 'lim_hdf : unable to allocate arrays' ) 101 DO jj = 2, jpjm1 76 DO jj = 2, jpjm1 102 77 DO ji = fs_2 , fs_jpim1 ! vector opt. 103 78 efact(ji,jj) = ( e2u(ji,jj) + e2u(ji-1,jj) + e1v(ji,jj) + e1v(ji,jj-1) ) * r1_e1e2t(ji,jj) … … 108 83 ! ! Time integration parameters 109 84 ! 110 zflu (jpi,: ) = 0._wp 111 zflv (jpi,: ) = 0._wp 112 113 DO jk=1 , isize 114 ztab0(:, : , jk ) = ptab(:,:,jk) ! Arrays initialization 115 zdiv0(:, 1 , jk ) = 0._wp 116 zdiv0(:,jpj, jk ) = 0._wp 117 zdiv0(1, :, jk ) = 0._wp 118 zdiv0(jpi,:, jk ) = 0._wp 119 END DO 85 ztab0(:, : ) = ptab(:,:) ! Arrays initialization 86 zdiv0(:, 1 ) = 0._wp 87 zdiv0(:,jpj) = 0._wp 88 zflu (jpi,:) = 0._wp 89 zflv (jpi,:) = 0._wp 90 zdiv0(1, :) = 0._wp 91 zdiv0(jpi,:) = 0._wp 120 92 121 93 zconv = 1._wp !== horizontal diffusion using a Crant-Nicholson scheme ==! 122 94 iter = 0 123 95 ! 124 DO WHILE( MAXVAL(zconv(:))> ( 2._wp * 1.e-04 ) .AND. iter <= its ) ! Sub-time step loop96 DO WHILE( zconv > ( 2._wp * 1.e-04 ) .AND. iter <= its ) ! Sub-time step loop 125 97 ! 126 98 iter = iter + 1 ! incrementation of the sub-time step number 127 99 ! 128 DO jk = 1 , isize129 jl = (jk-1) /( ihdf_vars+nlay_i)+1130 IF (zconv(jk) > ( 2._wp * 1.e-04 )) THEN131 DO jj = 1, jpjm1 ! diffusive fluxes in U- and V- direction132 DO ji = 1 , fs_jpim1 ! vector opt.133 zflu(ji,jj) = pahu3D(ji,jj,jl) * e2u(ji,jj) * r1_e1u(ji,jj) * ( ptab(ji+1,jj,jk) - ptab(ji,jj,jk) )134 zflv(ji,jj) = pahv3D(ji,jj,jl) * e1v(ji,jj) * r1_e2v(ji,jj) * ( ptab(ji,jj+1,jk) - ptab(ji,jj,jk) )135 END DO136 END DO137 !138 DO jj= 2, jpjm1 ! diffusive trend : divergence of the fluxes139 DO ji = fs_2 , fs_jpim1 ! vector opt.140 zdiv(ji,jj) = ( zflu(ji,jj) - zflu(ji-1,jj) + zflv(ji,jj) - zflv(ji,jj-1) ) * r1_e1e2t(ji,jj)141 END DO142 END DO143 !144 IF( iter == 1 ) zdiv0(:,:,jk) = zdiv(:,:) ! save the 1st evaluation of the diffusive trend in zdiv0145 !146 DO jj = 2, jpjm1 ! iterative evaluation147 DO ji = fs_2 , fs_jpim1 ! vector opt.148 zrlxint = ( ztab0(ji,jj,jk) &149 & + rdt_ice * ( zalfa * ( zdiv(ji,jj) + efact(ji,jj) * ptab(ji,jj,jk) ) &150 & + ( 1.0 - zalfa ) * zdiv0(ji,jj,jk) ) &151 & ) / ( 1.0 + zalfa * rdt_ice * efact(ji,jj) )152 zrlx(ji,jj,jk) = ptab(ji,jj,jk) + zrelax * ( zrlxint - ptab(ji,jj,jk) )153 END DO154 END DO155 END IF156 157 END DO158 159 CALL lbc_lnk_multi( zrlx_array, type_array , psgn_array , isize ) ! Multiple interchange of all the variables160 !161 IF ( MOD( iter-1 , nn_convfrq ) == 0 ) THEN !Convergence test every nn_convfrq iterations (perf. optimization )162 DO jk=1,isize163 zconv(jk) = 0._wp ! convergence test164 DO jj = 2, jpjm1165 DO ji = fs_2, fs_jpim1166 zconv(jk) = MAX( zconv(jk), ABS( zrlx(ji,jj,jk) - ptab(ji,jj,jk) ) )167 END DO168 END DO169 END DO170 IF( lk_mpp ) CALL mpp_max_multiple( zconv , isize ) ! max over the global domain for all the variables171 ENDIF172 !173 DO jk=1,isize174 ptab(:,:,jk) = zrlx(:,:,jk)175 END DO176 !177 END DO ! end of sub-time step loop178 179 ! -----------------------180 !!! final step (clem) !!!181 DO jk = 1, isize182 jl = (jk-1) /( ihdf_vars+nlay_i)+1183 100 DO jj = 1, jpjm1 ! diffusive fluxes in U- and V- direction 184 101 DO ji = 1 , fs_jpim1 ! vector opt. 185 zflu(ji,jj) = pahu 3D(ji,jj,jl) * e2u(ji,jj) * r1_e1u(ji,jj) * ( ptab(ji+1,jj,jk) - ptab(ji,jj,jk) )186 zflv(ji,jj) = pahv 3D(ji,jj,jl) * e1v(ji,jj) * r1_e2v(ji,jj) * ( ptab(ji,jj+1,jk) - ptab(ji,jj,jk) )102 zflu(ji,jj) = pahu(ji,jj) * e2u(ji,jj) * r1_e1u(ji,jj) * ( ptab(ji+1,jj) - ptab(ji,jj) ) 103 zflv(ji,jj) = pahv(ji,jj) * e1v(ji,jj) * r1_e2v(ji,jj) * ( ptab(ji,jj+1) - ptab(ji,jj) ) 187 104 END DO 188 105 END DO … … 191 108 DO ji = fs_2 , fs_jpim1 ! vector opt. 192 109 zdiv(ji,jj) = ( zflu(ji,jj) - zflu(ji-1,jj) + zflv(ji,jj) - zflv(ji,jj-1) ) * r1_e1e2t(ji,jj) 193 ptab(ji,jj,jk) = ztab0(ji,jj,jk) + 0.5 * ( zdiv(ji,jj) + zdiv0(ji,jj,jk) ) 194 END DO 110 END DO 111 END DO 112 ! 113 IF( iter == 1 ) zdiv0(:,:) = zdiv(:,:) ! save the 1st evaluation of the diffusive trend in zdiv0 114 ! 115 DO jj = 2, jpjm1 ! iterative evaluation 116 DO ji = fs_2 , fs_jpim1 ! vector opt. 117 zrlxint = ( ztab0(ji,jj) & 118 & + rdt_ice * ( zalfa * ( zdiv(ji,jj) + efact(ji,jj) * ptab(ji,jj) ) & 119 & + ( 1.0 - zalfa ) * zdiv0(ji,jj) ) & 120 & ) / ( 1.0 + zalfa * rdt_ice * efact(ji,jj) ) 121 zrlx(ji,jj) = ptab(ji,jj) + zrelax * ( zrlxint - ptab(ji,jj) ) 122 END DO 123 END DO 124 CALL lbc_lnk( zrlx, 'T', 1. ) ! lateral boundary condition 125 ! 126 IF ( MOD( iter, nn_convfrq ) == 0 ) THEN ! convergence test every nn_convfrq iterations (perf. optimization) 127 zconv = 0._wp 128 DO jj = 2, jpjm1 129 DO ji = fs_2, fs_jpim1 130 zconv = MAX( zconv, ABS( zrlx(ji,jj) - ptab(ji,jj) ) ) 131 END DO 132 END DO 133 IF( lk_mpp ) CALL mpp_max( zconv ) ! max over the global domain 134 ENDIF 135 ! 136 ptab(:,:) = zrlx(:,:) 137 ! 138 END DO ! end of sub-time step loop 139 140 ! ----------------------- 141 !!! final step (clem) !!! 142 DO jj = 1, jpjm1 ! diffusive fluxes in U- and V- direction 143 DO ji = 1 , fs_jpim1 ! vector opt. 144 zflu(ji,jj) = pahu(ji,jj) * e2u(ji,jj) * r1_e1u(ji,jj) * ( ptab(ji+1,jj) - ptab(ji,jj) ) 145 zflv(ji,jj) = pahv(ji,jj) * e1v(ji,jj) * r1_e2v(ji,jj) * ( ptab(ji,jj+1) - ptab(ji,jj) ) 195 146 END DO 196 147 END DO 197 198 CALL lbc_lnk_multi( pt2d_array, type_array , psgn_array , isize ) ! Multiple interchange of all the variables 199 148 ! 149 DO jj= 2, jpjm1 ! diffusive trend : divergence of the fluxes 150 DO ji = fs_2 , fs_jpim1 ! vector opt. 151 zdiv(ji,jj) = ( zflu(ji,jj) - zflu(ji-1,jj) + zflv(ji,jj) - zflv(ji,jj-1) ) * r1_e1e2t(ji,jj) 152 ptab(ji,jj) = ztab0(ji,jj) + 0.5 * ( zdiv(ji,jj) + zdiv0(ji,jj) ) 153 END DO 154 END DO 155 CALL lbc_lnk( ptab, 'T', 1. ) ! lateral boundary condition 200 156 !!! final step (clem) !!! 201 157 ! ----------------------- 202 158 203 159 IF(ln_ctl) THEN 204 DO jk = 1 , isize 205 zrlx(:,:,jk) = ptab(:,:,jk) - ztab0(:,:,jk) 206 WRITE(charout,FMT="(' lim_hdf : zconv =',D23.16, ' iter =',I4,2X)") zconv, iter 207 CALL prt_ctl( tab2d_1=zrlx(:,:,jk), clinfo1=charout ) 208 END DO 209 ENDIF 210 ! 211 CALL wrk_dealloc( jpi, jpj, isize, zrlx, zdiv0, ztab0 ) 212 CALL wrk_dealloc( jpi, jpj, zflu, zflv, zdiv ) 213 214 DEALLOCATE( zconv ) 215 DEALLOCATE( pt2d_array , zrlx_array ) 216 DEALLOCATE( type_array ) 217 DEALLOCATE( psgn_array ) 160 zrlx(:,:) = ptab(:,:) - ztab0(:,:) 161 WRITE(charout,FMT="(' lim_hdf : zconv =',D23.16, ' iter =',I4,2X)") zconv, iter 162 CALL prt_ctl( tab2d_1=zrlx, clinfo1=charout ) 163 ENDIF 164 ! 165 CALL wrk_dealloc( jpi, jpj, zrlx, zflu, zflv, zdiv0, zdiv, ztab0 ) 218 166 ! 219 167 END SUBROUTINE lim_hdf 220 221 168 222 169 … … 232 179 !!------------------------------------------------------------------- 233 180 INTEGER :: ios ! Local integer output status for namelist read 234 NAMELIST/namicehdf/ nn_convfrq 181 NAMELIST/namicehdf/ nn_convfrq 235 182 !!------------------------------------------------------------------- 236 183 ! … … 265 212 !!====================================================================== 266 213 END MODULE limhdf 267 -
trunk/NEMOGCM/NEMO/LIM_SRC_3/limtrp.F90
r6478 r6483 63 63 INTEGER, INTENT(in) :: kt ! number of iteration 64 64 ! 65 INTEGER :: ji, jj, jk, j m , jl, jt ! dummy loop indices65 INTEGER :: ji, jj, jk, jl, jt ! dummy loop indices 66 66 INTEGER :: initad ! number of sub-timestep for the advection 67 67 REAL(wp) :: zcfl , zusnit ! - - … … 75 75 REAL(wp), POINTER, DIMENSION(:,:,:) :: zhimax ! old ice thickness 76 76 REAL(wp), POINTER, DIMENSION(:,:) :: zatold, zeiold, zesold ! old concentration, enthalpies 77 REAL(wp), POINTER, DIMENSION(:,:,:) :: zhdfptab78 77 REAL(wp) :: zdv, zvi, zvs, zsmv, zes, zei 79 78 REAL(wp) :: zvi_b, zsmv_b, zei_b, zfs_b, zfw_b, zft_b 80 !!---------------------------------------------------------------------81 INTEGER :: ihdf_vars = 6 !!Number of variables in which we apply horizontal diffusion82 !! inside limtrp for each ice category , not counting the83 !! variables corresponding to ice_layers84 79 !!--------------------------------------------------------------------- 85 80 IF( nn_timing == 1 ) CALL timing_start('limtrp') … … 90 85 CALL wrk_alloc( jpi,jpj,nlay_i,jpl, z0ei ) 91 86 CALL wrk_alloc( jpi,jpj,jpl, zhimax, zviold, zvsold, zsmvold ) 92 CALL wrk_alloc( jpi,jpj,jpl*(ihdf_vars + nlay_i)+1,zhdfptab)93 87 94 88 IF( numit == nstart .AND. lwp ) THEN … … 176 170 z0oi (:,:,jl) = oa_i (:,:, jl) * e1e2t(:,:) ! Age content 177 171 z0es (:,:,jl) = e_s (:,:,1,jl) * e1e2t(:,:) ! Snow heat content 178 DO jk = 1, nlay_i172 DO jk = 1, nlay_i 179 173 z0ei (:,:,jk,jl) = e_i (:,:,jk,jl) * e1e2t(:,:) ! Ice heat content 180 174 END DO … … 290 284 ! Diffusion of Ice fields 291 285 !------------------------------------------------------------------------------! 292 !------------------------------------ 293 ! Diffusion of other ice variables 294 !------------------------------------ 295 jm=1 296 DO jl = 1, jpl 297 ! ! Masked eddy diffusivity coefficient at ocean U- and V-points 298 ! DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row 299 ! DO ji = 1 , fs_jpim1 ! vector opt. 300 ! pahu(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji ,jj,jl) ) ) ) & 301 ! & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji+1,jj,jl) ) ) ) * ahiu(ji,jj) 302 ! pahv(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji,jj ,jl) ) ) ) & 303 ! & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- a_i(ji,jj+1,jl) ) ) ) * ahiv(ji,jj) 304 ! END DO 305 ! END DO 306 DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row 307 DO ji = 1 , fs_jpim1 ! vector opt. 308 pahu3D(ji,jj,jl) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji ,jj, jl ) ) ) ) & 309 & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji+1,jj, jl ) ) ) ) * ahiu(ji,jj) 310 pahv3D(ji,jj,jl) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji, jj, jl ) ) ) ) & 311 & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- a_i(ji, jj+1,jl ) ) ) ) * ahiv(ji,jj) 312 END DO 313 END DO 314 315 zhdfptab(:,:,jm)= a_i (:,:, jl); jm = jm + 1 316 zhdfptab(:,:,jm)= v_i (:,:, jl); jm = jm + 1 317 zhdfptab(:,:,jm)= v_s (:,:, jl); jm = jm + 1 318 zhdfptab(:,:,jm)= smv_i(:,:, jl); jm = jm + 1 319 zhdfptab(:,:,jm)= oa_i (:,:, jl); jm = jm + 1 320 zhdfptab(:,:,jm)= e_s (:,:,1,jl); jm = jm + 1 321 ! Sample of adding more variables to apply lim_hdf using lim_hdf optimization--- 322 ! zhdfptab(:,:,jm) = variable_1 (:,:,1,jl); jm = jm + 1 323 ! zhdfptab(:,:,jm) = variable_2 (:,:,1,jl); jm = jm + 1 324 ! 325 ! and in this example the parameter ihdf_vars musb be changed to 8 (necessary for allocation) 326 !---------------------------------------------------------------------------------------- 327 DO jk = 1, nlay_i 328 zhdfptab(:,:,jm)=e_i(:,:,jk,jl); jm= jm+1 329 END DO 330 END DO 286 331 287 ! 332 288 !-------------------------------- … … 334 290 !-------------------------------- 335 291 ! ! Masked eddy diffusivity coefficient at ocean U- and V-points 336 !DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row337 ! DO ji = 1 , fs_jpim1 ! vector opt.338 ! pahu(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji ,jj) ) ) ) &339 ! & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji+1,jj) ) ) ) * ahiu(ji,jj)340 ! pahv(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji,jj ) ) ) ) &341 ! & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- at_i(ji,jj+1) ) ) ) * ahiv(ji,jj)342 ! END DO343 !END DO344 345 292 DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row 346 293 DO ji = 1 , fs_jpim1 ! vector opt. 347 pahu 3D(ji,jj,jpl+1) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji ,jj) ) ) ) &348 & 349 pahv 3D(ji,jj,jpl+1) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji,jj ) ) ) ) &350 & 294 pahu(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji ,jj) ) ) ) & 295 & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji+1,jj) ) ) ) * ahiu(ji,jj) 296 pahv(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -at_i(ji,jj ) ) ) ) & 297 & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- at_i(ji,jj+1) ) ) ) * ahiv(ji,jj) 351 298 END DO 352 299 END DO 353 300 ! 354 zhdfptab(:,:,jm)= ato_i (:,:); 355 CALL lim_hdf( zhdfptab, ihdf_vars, jpl, nlay_i) 356 357 jm=1 358 DO jl = 1, jpl 359 a_i (:,:, jl) = zhdfptab(:,:,jm); jm = jm + 1 360 v_i (:,:, jl) = zhdfptab(:,:,jm); jm = jm + 1 361 v_s (:,:, jl) = zhdfptab(:,:,jm); jm = jm + 1 362 smv_i(:,:, jl) = zhdfptab(:,:,jm); jm = jm + 1 363 oa_i (:,:, jl) = zhdfptab(:,:,jm); jm = jm + 1 364 e_s (:,:,1,jl) = zhdfptab(:,:,jm); jm = jm + 1 365 ! Sample of adding more variables to apply lim_hdf--------- 366 ! variable_1 (:,:,1,jl) = zhdfptab(:,:, jm ) ; jm + 1 367 ! variable_2 (:,:,1,jl) = zhdfptab(:,:, jm ) ; jm + 1 368 !----------------------------------------------------------- 301 CALL lim_hdf( ato_i (:,:) ) 302 303 !------------------------------------ 304 ! Diffusion of other ice variables 305 !------------------------------------ 306 DO jl = 1, jpl 307 ! ! Masked eddy diffusivity coefficient at ocean U- and V-points 308 DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row 309 DO ji = 1 , fs_jpim1 ! vector opt. 310 pahu(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji ,jj,jl) ) ) ) & 311 & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji+1,jj,jl) ) ) ) * ahiu(ji,jj) 312 pahv(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -a_i(ji,jj ,jl) ) ) ) & 313 & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- a_i(ji,jj+1,jl) ) ) ) * ahiv(ji,jj) 314 END DO 315 END DO 316 317 CALL lim_hdf( v_i (:,:, jl) ) 318 CALL lim_hdf( v_s (:,:, jl) ) 319 CALL lim_hdf( smv_i(:,:, jl) ) 320 CALL lim_hdf( oa_i (:,:, jl) ) 321 CALL lim_hdf( a_i (:,:, jl) ) 322 CALL lim_hdf( e_s (:,:,1,jl) ) 369 323 DO jk = 1, nlay_i 370 e_i(:,:,jk,jl) = zhdfptab(:,:,jm);jm= jm + 1 371 END DO 372 END DO 373 374 ato_i (:,:) = zhdfptab(:,:,jm) 324 CALL lim_hdf( e_i(:,:,jk,jl) ) 325 END DO 326 END DO 375 327 376 328 !------------------------------------------------------------------------------! … … 512 464 CALL wrk_dealloc( jpi,jpj,nlay_i,jpl, z0ei ) 513 465 CALL wrk_dealloc( jpi,jpj,jpl, zviold, zvsold, zhimax, zsmvold ) 514 CALL wrk_dealloc( jpi,jpj,jpl*(ihdf_vars+nlay_i)+1,zhdfptab)515 466 ! 516 467 IF( nn_timing == 1 ) CALL timing_stop('limtrp') … … 528 479 !!====================================================================== 529 480 END MODULE limtrp 530 -
trunk/NEMOGCM/NEMO/OPA_SRC/LBC/lbclnk.F90
r6478 r6483 9 9 !! 3.5 ! 2012 (S.Mocavero, I. Epicoco) optimization of BDY comm. via lbc_bdy_lnk and lbc_obc_lnk 10 10 !! 3.4 ! 2012-12 (R. Bourdalle-Badie, G. Reffray) add a C1D case 11 !! 3.6 ! 2015-06 (O. Tintó and M. Castrillo) add lbc_lnk_multi12 11 !!---------------------------------------------------------------------- 13 12 #if defined key_mpp_mpi … … 23 22 24 23 INTERFACE lbc_lnk_multi 25 MODULE PROCEDURE mpp_lnk_2d_9 , mpp_lnk_2d_multiple24 MODULE PROCEDURE mpp_lnk_2d_9 26 25 END INTERFACE 27 26 ! … … 91 90 END INTERFACE 92 91 ! 93 INTERFACE lbc_lnk_multi94 MODULE PROCEDURE lbc_lnk_2d_9, lbc_lnk_2d_multiple95 END INTERFACE96 97 92 INTERFACE lbc_bdy_lnk 98 93 MODULE PROCEDURE lbc_bdy_lnk_2d, lbc_bdy_lnk_3d … … 102 97 MODULE PROCEDURE lbc_lnk_2d_e 103 98 END INTERFACE 104 105 TYPE arrayptr106 REAL , DIMENSION (:,:), POINTER :: pt2d107 END TYPE arrayptr108 PUBLIC arrayptr109 99 110 100 PUBLIC lbc_lnk ! ocean/ice lateral boundary conditions 111 101 PUBLIC lbc_lnk_e ! 112 PUBLIC lbc_lnk_multi ! modified ocean lateral boundary conditions113 102 PUBLIC lbc_bdy_lnk ! ocean lateral BDY boundary conditions 114 103 PUBLIC lbc_lnk_icb ! … … 192 181 ! 193 182 END SUBROUTINE lbc_lnk_2d 194 195 SUBROUTINE lbc_lnk_2d_multiple( pt2d_array , type_array , psgn_array , num_fields )196 !!197 INTEGER :: num_fields198 TYPE( arrayptr ), DIMENSION(:) :: pt2d_array199 CHARACTER(len=1), DIMENSION(:), INTENT(in ) :: type_array ! define the nature of ptab array grid-points200 ! ! = T , U , V , F , W and I points201 REAL(wp) , DIMENSION(:), INTENT(in ) :: psgn_array ! =-1 the sign change across the north fold boundary202 ! ! = 1. , the sign is kept203 !204 INTEGER :: ii !!MULTI SEND DUMMY LOOP INDICES205 !206 DO ii = 1, num_fields207 CALL lbc_lnk_2d( pt2d_array(ii)%pt2d, type_array(ii), psgn_array(ii) )208 END DO209 !210 END SUBROUTINE lbc_lnk_2d_multiple211 212 SUBROUTINE lbc_lnk_2d_9( pt2dA, cd_typeA, psgnA, pt2dB, cd_typeB, psgnB, pt2dC, cd_typeC, psgnC &213 & , pt2dD, cd_typeD, psgnD, pt2dE, cd_typeE, psgnE, pt2dF, cd_typeF, psgnF &214 & , pt2dG, cd_typeG, psgnG, pt2dH, cd_typeH, psgnH, pt2dI, cd_typeI, psgnI, cd_mpp, pval)215 !!---------------------------------------------------------------------216 ! Second 2D array on which the boundary condition is applied217 REAL(wp), DIMENSION(jpi,jpj), TARGET , INTENT(inout) :: pt2dA218 REAL(wp), DIMENSION(jpi,jpj), TARGET, OPTIONAL, INTENT(inout) :: pt2dB , pt2dC , pt2dD , pt2dE219 REAL(wp), DIMENSION(jpi,jpj), TARGET, OPTIONAL, INTENT(inout) :: pt2dF , pt2dG , pt2dH , pt2dI220 ! define the nature of ptab array grid-points221 CHARACTER(len=1) , INTENT(in ) :: cd_typeA222 CHARACTER(len=1) , OPTIONAL, INTENT(in ) :: cd_typeB , cd_typeC , cd_typeD , cd_typeE223 CHARACTER(len=1) , OPTIONAL, INTENT(in ) :: cd_typeF , cd_typeG , cd_typeH , cd_typeI224 ! =-1 the sign change across the north fold boundary225 REAL(wp) , INTENT(in ) :: psgnA226 REAL(wp) , OPTIONAL, INTENT(in ) :: psgnB , psgnC , psgnD , psgnE227 REAL(wp) , OPTIONAL, INTENT(in ) :: psgnF , psgnG , psgnH , psgnI228 CHARACTER(len=3) , OPTIONAL, INTENT(in ) :: cd_mpp ! fill the overlap area only229 REAL(wp) , OPTIONAL, INTENT(in ) :: pval ! background value (used at closed boundaries)230 !!231 !!---------------------------------------------------------------------232 233 !!The first array234 CALL lbc_lnk( pt2dA, cd_typeA, psgnA )235 236 !! Look if more arrays to process237 IF(PRESENT (psgnB) )CALL lbc_lnk( pt2dA, cd_typeA, psgnA )238 IF(PRESENT (psgnC) )CALL lbc_lnk( pt2dC, cd_typeC, psgnC )239 IF(PRESENT (psgnD) )CALL lbc_lnk( pt2dD, cd_typeD, psgnD )240 IF(PRESENT (psgnE) )CALL lbc_lnk( pt2dE, cd_typeE, psgnE )241 IF(PRESENT (psgnF) )CALL lbc_lnk( pt2dF, cd_typeF, psgnF )242 IF(PRESENT (psgnG) )CALL lbc_lnk( pt2dG, cd_typeG, psgnG )243 IF(PRESENT (psgnH) )CALL lbc_lnk( pt2dH, cd_typeH, psgnH )244 IF(PRESENT (psgnI) )CALL lbc_lnk( pt2dI, cd_typeI, psgnI )245 246 END SUBROUTINE lbc_lnk_2d_9247 248 249 250 251 183 252 184 #else … … 447 379 ! 448 380 END SUBROUTINE lbc_lnk_2d 449 450 SUBROUTINE lbc_lnk_2d_multiple( pt2d_array , type_array , psgn_array , num_fields )451 !!452 INTEGER :: num_fields453 TYPE( arrayptr ), DIMENSION(:) :: pt2d_array454 CHARACTER(len=1), DIMENSION(:), INTENT(in ) :: type_array ! define the nature of ptab array grid-points455 ! ! = T , U , V , F , W and I points456 REAL(wp) , DIMENSION(:), INTENT(in ) :: psgn_array ! =-1 the sign change across the north fold boundary457 ! ! = 1. , the sign is kept458 !459 INTEGER :: ii !!MULTI SEND DUMMY LOOP INDICES460 !461 DO ii = 1, num_fields462 CALL lbc_lnk_2d( pt2d_array(ii)%pt2d, type_array(ii), psgn_array(ii) )463 END DO464 !465 END SUBROUTINE lbc_lnk_2d_multiple466 467 SUBROUTINE lbc_lnk_2d_9( pt2dA, cd_typeA, psgnA, pt2dB, cd_typeB, psgnB, pt2dC, cd_typeC, psgnC &468 & , pt2dD, cd_typeD, psgnD, pt2dE, cd_typeE, psgnE, pt2dF, cd_typeF, psgnF &469 & , pt2dG, cd_typeG, psgnG, pt2dH, cd_typeH, psgnH, pt2dI, cd_typeI, psgnI, cd_mpp, pval)470 !!---------------------------------------------------------------------471 ! Second 2D array on which the boundary condition is applied472 REAL(wp), DIMENSION(jpi,jpj), TARGET , INTENT(inout) :: pt2dA473 REAL(wp), DIMENSION(jpi,jpj), TARGET, OPTIONAL, INTENT(inout) :: pt2dB , pt2dC , pt2dD , pt2dE474 REAL(wp), DIMENSION(jpi,jpj), TARGET, OPTIONAL, INTENT(inout) :: pt2dF , pt2dG , pt2dH , pt2dI475 ! define the nature of ptab array grid-points476 CHARACTER(len=1) , INTENT(in ) :: cd_typeA477 CHARACTER(len=1) , OPTIONAL, INTENT(in ) :: cd_typeB , cd_typeC , cd_typeD , cd_typeE478 CHARACTER(len=1) , OPTIONAL, INTENT(in ) :: cd_typeF , cd_typeG , cd_typeH , cd_typeI479 ! =-1 the sign change across the north fold boundary480 REAL(wp) , INTENT(in ) :: psgnA481 REAL(wp) , OPTIONAL, INTENT(in ) :: psgnB , psgnC , psgnD , psgnE482 REAL(wp) , OPTIONAL, INTENT(in ) :: psgnF , psgnG , psgnH , psgnI483 CHARACTER(len=3) , OPTIONAL, INTENT(in ) :: cd_mpp ! fill the overlap area only484 REAL(wp) , OPTIONAL, INTENT(in ) :: pval ! background value (used at closed boundaries)485 !!486 !!---------------------------------------------------------------------487 488 !!The first array489 CALL lbc_lnk( pt2dA, cd_typeA, psgnA )490 491 !! Look if more arrays to process492 IF(PRESENT (psgnB) )CALL lbc_lnk( pt2dA, cd_typeA, psgnA )493 IF(PRESENT (psgnC) )CALL lbc_lnk( pt2dC, cd_typeC, psgnC )494 IF(PRESENT (psgnD) )CALL lbc_lnk( pt2dD, cd_typeD, psgnD )495 IF(PRESENT (psgnE) )CALL lbc_lnk( pt2dE, cd_typeE, psgnE )496 IF(PRESENT (psgnF) )CALL lbc_lnk( pt2dF, cd_typeF, psgnF )497 IF(PRESENT (psgnG) )CALL lbc_lnk( pt2dG, cd_typeG, psgnG )498 IF(PRESENT (psgnH) )CALL lbc_lnk( pt2dH, cd_typeH, psgnH )499 IF(PRESENT (psgnI) )CALL lbc_lnk( pt2dI, cd_typeI, psgnI )500 501 END SUBROUTINE lbc_lnk_2d_9502 503 381 504 382 #endif … … 570 448 !!====================================================================== 571 449 END MODULE lbclnk 572 -
trunk/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90
r6478 r6483 24 24 !! 3.5 ! 2013 ( C. Ethe, G. Madec ) message passing arrays as local variables 25 25 !! 3.5 ! 2013 (S.Mocavero, I.Epicoco - CMCC) north fold optimizations 26 !! 3.6 ! 2015 (O. Tintó and M. Castrillo - BSC) Added 'mpp_lnk_2d_multiple', 'mpp_lbc_north_2d_multiple', 'mpp_max_multiple'27 26 !!---------------------------------------------------------------------- 28 27 … … 63 62 USE lbcnfd ! north fold treatment 64 63 USE in_out_manager ! I/O manager 65 USE wrk_nemo ! work arrays66 64 67 65 IMPLICIT NONE … … 72 70 PUBLIC mpp_ini_north, mpp_lbc_north, mpp_lbc_north_e 73 71 PUBLIC mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc 74 PUBLIC mpp_max_multiple75 72 PUBLIC mpp_lnk_3d, mpp_lnk_3d_gather, mpp_lnk_2d, mpp_lnk_2d_e 76 PUBLIC mpp_lnk_2d_9 , mpp_lnk_2d_multiple73 PUBLIC mpp_lnk_2d_9 77 74 PUBLIC mpp_lnk_sum_3d, mpp_lnk_sum_2d 78 75 PUBLIC mppscatter, mppgather … … 82 79 PUBLIC mpp_lnk_bdy_2d, mpp_lnk_bdy_3d 83 80 PUBLIC mpp_lbc_north_icb, mpp_lnk_2d_icb 84 PUBLIC mpprank85 81 86 82 TYPE arrayptr 87 83 REAL , DIMENSION (:,:), POINTER :: pt2d 88 84 END TYPE arrayptr 89 PUBLIC arrayptr90 85 91 86 !! * Interfaces … … 111 106 INTERFACE mpp_maxloc 112 107 MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d 113 END INTERFACE114 115 INTERFACE mpp_max_multiple116 MODULE PROCEDURE mppmax_real_multiple117 108 END INTERFACE 118 109 … … 735 726 ! ----------------------- 736 727 ! 728 DO ii = 1 , num_fields 737 729 !First Array 738 IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 739 ! 740 SELECT CASE ( jpni ) 741 CASE ( 1 ) ; 742 DO ii = 1 , num_fields 743 CALL lbc_nfd ( pt2d_array(ii)%pt2d( : , : ), type_array(ii) , psgn_array(ii) ) ! only 1 northern proc, no mpp 744 END DO 745 CASE DEFAULT ; CALL mpp_lbc_north_2d_multiple( pt2d_array, type_array, psgn_array, num_fields ) ! for all northern procs. 746 END SELECT 747 ! 748 ENDIF 749 ! 730 IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 731 ! 732 SELECT CASE ( jpni ) 733 CASE ( 1 ) ; CALL lbc_nfd ( pt2d_array(ii)%pt2d( : , : ), type_array(ii) , psgn_array(ii) ) ! only 1 northern proc, no mpp 734 CASE DEFAULT ; CALL mpp_lbc_north( pt2d_array(ii)%pt2d( : , : ), type_array(ii), psgn_array(ii) ) ! for all northern procs. 735 END SELECT 736 ! 737 ENDIF 738 ! 739 END DO 750 740 ! 751 741 DEALLOCATE( zt2ns, zt2sn, zt2ew, zt2we ) … … 2029 2019 END SUBROUTINE mppmax_real 2030 2020 2031 SUBROUTINE mppmax_real_multiple( ptab, NUM , kcom )2032 !!----------------------------------------------------------------------2033 !! *** routine mppmax_real ***2034 !!2035 !! ** Purpose : Maximum2036 !!2037 !!----------------------------------------------------------------------2038 REAL(wp), DIMENSION(:) , INTENT(inout) :: ptab ! ???2039 INTEGER , INTENT(in ) :: NUM2040 INTEGER , INTENT(in ), OPTIONAL :: kcom ! ???2041 !!2042 INTEGER :: ierror, localcomm2043 REAL(wp) , POINTER , DIMENSION(:) :: zwork2044 !!----------------------------------------------------------------------2045 !2046 CALL wrk_alloc(NUM , zwork)2047 localcomm = mpi_comm_opa2048 IF( PRESENT(kcom) ) localcomm = kcom2049 !2050 CALL mpi_allreduce( ptab, zwork, NUM, mpi_double_precision, mpi_max, localcomm, ierror )2051 ptab = zwork2052 CALL wrk_dealloc(NUM , zwork)2053 !2054 END SUBROUTINE mppmax_real_multiple2055 2056 2021 2057 2022 SUBROUTINE mppmin_a_real( ptab, kdim, kcom ) … … 2947 2912 END SUBROUTINE mpp_lbc_north_2d 2948 2913 2949 SUBROUTINE mpp_lbc_north_2d_multiple( pt2d_array, cd_type, psgn, num_fields)2950 !!---------------------------------------------------------------------2951 !! *** routine mpp_lbc_north_2d ***2952 !!2953 !! ** Purpose : Ensure proper north fold horizontal bondary condition2954 !! in mpp configuration in case of jpn1 > 12955 !! (for multiple 2d arrays )2956 !!2957 !! ** Method : North fold condition and mpp with more than one proc2958 !! in i-direction require a specific treatment. We gather2959 !! the 4 northern lines of the global domain on 1 processor2960 !! and apply lbc north-fold on this sub array. Then we2961 !! scatter the north fold array back to the processors.2962 !!2963 !!----------------------------------------------------------------------2964 INTEGER , INTENT (in ) :: num_fields ! number of variables contained in pt2d2965 TYPE( arrayptr ), DIMENSION(:) :: pt2d_array2966 CHARACTER(len=1), DIMENSION(:), INTENT(in ) :: cd_type ! nature of pt2d grid-points2967 ! ! = T , U , V , F or W gridpoints2968 REAL(wp), DIMENSION(:), INTENT(in ) :: psgn ! = -1. the sign change across the north fold2969 !! ! = 1. , the sign is kept2970 INTEGER :: ji, jj, jr, jk2971 INTEGER :: ierr, itaille, ildi, ilei, iilb2972 INTEGER :: ijpj, ijpjm1, ij, iproc2973 INTEGER, DIMENSION (jpmaxngh) :: ml_req_nf !for mpi_isend when avoiding mpi_allgather2974 INTEGER :: ml_err ! for mpi_isend when avoiding mpi_allgather2975 INTEGER, DIMENSION(MPI_STATUS_SIZE):: ml_stat ! for mpi_isend when avoiding mpi_allgather2976 ! ! Workspace for message transfers avoiding mpi_allgather2977 REAL(wp), DIMENSION(:,:,:) , ALLOCATABLE :: ztab2978 REAL(wp), DIMENSION(:,:,:) , ALLOCATABLE :: znorthloc, zfoldwk2979 REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: znorthgloio2980 REAL(wp), DIMENSION(:,:,:) , ALLOCATABLE :: ztabl, ztabr2981 INTEGER :: istatus(mpi_status_size)2982 INTEGER :: iflag2983 !!----------------------------------------------------------------------2984 !2985 ALLOCATE( ztab(jpiglo,4,num_fields), znorthloc(jpi,4,num_fields), zfoldwk(jpi,4,num_fields), znorthgloio(jpi,4,num_fields,jpni) ) ! expanded to 3 dimensions2986 ALLOCATE( ztabl(jpi,4,num_fields), ztabr(jpi*jpmaxngh, 4,num_fields) )2987 !2988 ijpj = 42989 ijpjm1 = 32990 !2991 2992 DO jk = 1, num_fields2993 DO jj = nlcj-ijpj+1, nlcj ! put in znorthloc the last 4 jlines of pt2d (for every variable)2994 ij = jj - nlcj + ijpj2995 znorthloc(:,ij,jk) = pt2d_array(jk)%pt2d(:,jj)2996 END DO2997 END DO2998 ! ! Build in procs of ncomm_north the znorthgloio2999 itaille = jpi * ijpj3000 3001 IF ( l_north_nogather ) THEN3002 !3003 ! Avoid the use of mpi_allgather by exchanging only with the processes already identified3004 ! (in nemo_northcomms) as being involved in this process' northern boundary exchange3005 !3006 ztabr(:,:,:) = 03007 ztabl(:,:,:) = 03008 3009 DO jk = 1, num_fields3010 DO jj = nlcj-ijpj+1, nlcj ! First put local values into the global array3011 ij = jj - nlcj + ijpj3012 DO ji = nfsloop, nfeloop3013 ztabl(ji,ij,jk) = pt2d_array(jk)%pt2d(ji,jj)3014 END DO3015 END DO3016 END DO3017 3018 DO jr = 1,nsndto3019 IF ((nfipproc(isendto(jr),jpnj) .ne. (narea-1)) .and. (nfipproc(isendto(jr),jpnj) .ne. -1)) THEN3020 CALL mppsend(5, znorthloc, itaille*num_fields, nfipproc(isendto(jr),jpnj), ml_req_nf(jr)) ! Buffer expanded "num_fields" times3021 ENDIF3022 END DO3023 DO jr = 1,nsndto3024 iproc = nfipproc(isendto(jr),jpnj)3025 IF(iproc .ne. -1) THEN3026 ilei = nleit (iproc+1)3027 ildi = nldit (iproc+1)3028 iilb = nfiimpp(isendto(jr),jpnj) - nfiimpp(isendto(1),jpnj)3029 ENDIF3030 IF((iproc .ne. (narea-1)) .and. (iproc .ne. -1)) THEN3031 CALL mpprecv(5, zfoldwk, itaille*num_fields, iproc) ! Buffer expanded "num_fields" times3032 DO jk = 1 , num_fields3033 DO jj = 1, ijpj3034 DO ji = ildi, ilei3035 ztabr(iilb+ji,jj,jk) = zfoldwk(ji,jj,jk) ! Modified to 3D3036 END DO3037 END DO3038 END DO3039 ELSE IF (iproc .eq. (narea-1)) THEN3040 DO jk = 1, num_fields3041 DO jj = 1, ijpj3042 DO ji = ildi, ilei3043 ztabr(iilb+ji,jj,jk) = pt2d_array(jk)%pt2d(ji,nlcj-ijpj+jj) ! Modified to 3D3044 END DO3045 END DO3046 END DO3047 ENDIF3048 END DO3049 IF (l_isend) THEN3050 DO jr = 1,nsndto3051 IF ((nfipproc(isendto(jr),jpnj) .ne. (narea-1)) .and. (nfipproc(isendto(jr),jpnj) .ne. -1)) THEN3052 CALL mpi_wait(ml_req_nf(jr), ml_stat, ml_err)3053 ENDIF3054 END DO3055 ENDIF3056 !3057 DO ji = 1, num_fields ! Loop to manage 3D variables3058 CALL mpp_lbc_nfd( ztabl(:,:,ji), ztabr(:,:,ji), cd_type(ji), psgn(ji) ) ! North fold boundary condition3059 END DO3060 !3061 DO jk = 1, num_fields3062 DO jj = nlcj-ijpj+1, nlcj ! Scatter back to pt2d3063 ij = jj - nlcj + ijpj3064 DO ji = 1, nlci3065 pt2d_array(jk)%pt2d(ji,jj) = ztabl(ji,ij,jk) ! Modified to 3D3066 END DO3067 END DO3068 END DO3069 3070 !3071 ELSE3072 !3073 CALL MPI_ALLGATHER( znorthloc , itaille*num_fields, MPI_DOUBLE_PRECISION, &3074 & znorthgloio, itaille*num_fields, MPI_DOUBLE_PRECISION, ncomm_north, ierr )3075 !3076 ztab(:,:,:) = 0.e03077 DO jk = 1, num_fields3078 DO jr = 1, ndim_rank_north ! recover the global north array3079 iproc = nrank_north(jr) + 13080 ildi = nldit (iproc)3081 ilei = nleit (iproc)3082 iilb = nimppt(iproc)3083 DO jj = 1, ijpj3084 DO ji = ildi, ilei3085 ztab(ji+iilb-1,jj,jk) = znorthgloio(ji,jj,jk,jr)3086 END DO3087 END DO3088 END DO3089 END DO3090 3091 DO ji = 1, num_fields3092 CALL lbc_nfd( ztab(:,:,ji), cd_type(ji), psgn(ji) ) ! North fold boundary condition3093 END DO3094 !3095 DO jk = 1, num_fields3096 DO jj = nlcj-ijpj+1, nlcj ! Scatter back to pt2d3097 ij = jj - nlcj + ijpj3098 DO ji = 1, nlci3099 pt2d_array(jk)%pt2d(ji,jj) = ztab(ji+nimpp-1,ij,jk)3100 END DO3101 END DO3102 END DO3103 !3104 !3105 ENDIF3106 DEALLOCATE( ztab, znorthloc, zfoldwk, znorthgloio )3107 DEALLOCATE( ztabl, ztabr )3108 !3109 END SUBROUTINE mpp_lbc_north_2d_multiple3110 2914 3111 2915 SUBROUTINE mpp_lbc_north_e( pt2d, cd_type, psgn)
Note: See TracChangeset
for help on using the changeset viewer.