Changeset 5836 for trunk/NEMOGCM/NEMO/OPA_SRC/DOM/dom_oce.F90
- Timestamp:
- 2015-10-26T15:49:40+01:00 (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMOGCM/NEMO/OPA_SRC/DOM/dom_oce.F90
r5123 r5836 7 7 !! History : 1.0 ! 2005-10 (A. Beckmann, G. Madec) reactivate s-coordinate 8 8 !! 3.3 ! 2010-11 (G. Madec) add mbk. arrays associated to the deepest ocean level 9 !! 4.0! 2011-01 (A. R. Porter, STFC Daresbury) dynamical allocation9 !! 3.4 ! 2011-01 (A. R. Porter, STFC Daresbury) dynamical allocation 10 10 !! 3.5 ! 2012 (S. Mocavero, I. Epicoco) Add arrays associated 11 11 !! to the optimization of BDY communications 12 !! 3.7 ! 2015-11 (G. Madec) introduce surface and scale factor ratio 12 13 !!---------------------------------------------------------------------- 13 14 … … 20 21 21 22 IMPLICIT NONE 22 PUBLIC ! allows the acces to par_oce when dom_oce is used 23 ! ! exception to coding rules... to be suppressed ??? 23 PUBLIC ! allows the acces to par_oce when dom_oce is used (exception to coding rules) 24 24 25 25 PUBLIC dom_oce_alloc ! Called from nemogcm.F90 … … 107 107 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) :: r2dtra !: = 2*rdttra except at nit000 (=rdttra) if neuler=0 108 108 109 ! !!* Namelist namcla : cross land advection110 INTEGER, PUBLIC :: nn_cla !: =1 cross land advection for exchanges through some straits (ORCA2)111 112 109 !!---------------------------------------------------------------------- 113 110 !! space domain parameters … … 158 155 !! horizontal curvilinear coordinate and scale factors 159 156 !! --------------------------------------------------------------------- 160 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: glamt, glamu !: longitude of t-, u-, v- and f-points (degre) 161 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: glamv, glamf !: 162 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: gphit, gphiu !: latitude of t-, u-, v- and f-points (degre) 163 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: gphiv, gphif !: 164 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1t, e2t, r1_e1t, r1_e2t !: horizontal scale factors and inverse at t-point (m) 165 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1u, e2u, r1_e1u, r1_e2u !: horizontal scale factors and inverse at u-point (m) 166 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1v, e2v, r1_e1v, r1_e2v !: horizontal scale factors and inverse at v-point (m) 167 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1f, e2f, r1_e1f, r1_e2f !: horizontal scale factors and inverse at f-point (m) 168 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e1e2t !: surface at t-point (m2) 169 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ff !: coriolis factor (2.*omega*sin(yphi) ) (s-1) 157 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: glamt , glamu, glamv , glamf !: longitude at t, u, v, f-points [degree] 158 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: gphit , gphiu, gphiv , gphif !: latitude at t, u, v, f-points [degree] 159 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1t , e2t , r1_e1t, r1_e2t !: t-point horizontal scale factors [m] 160 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1u , e2u , r1_e1u, r1_e2u !: horizontal scale factors at u-point [m] 161 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1v , e2v , r1_e1v, r1_e2v !: horizontal scale factors at v-point [m] 162 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: e1f , e2f , r1_e1f, r1_e2f !: horizontal scale factors at f-point [m] 163 ! 164 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e1e2t , r1_e1e2t !: associated metrics at t-point 165 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e1e2u , r1_e1e2u , e2_e1u !: associated metrics at u-point 166 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e1e2v , r1_e1e2v , e1_e2v !: associated metrics at v-point 167 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e1e2f , r1_e1e2f !: associated metrics at f-point 168 ! 169 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ff !: coriolis factor [1/s] 170 170 171 171 !!---------------------------------------------------------------------- … … 216 216 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ht_0 !: reference depth at t- points (meters) 217 217 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hu_0 , hv_0 !: reference depth at u- and v-points (meters) 218 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: re2u_e1u !: scale factor coeffs at u points (e2u/e1u)219 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: re1v_e2v !: scale factor coeffs at v points (e1v/e2v)220 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e12t , r1_e12t !: horizontal cell surface and inverse at t points221 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e12u , r1_e12u !: horizontal cell surface and inverse at u points222 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e12v , r1_e12v !: horizontal cell surface and inverse at v points223 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e12f , r1_e12f !: horizontal cell surface and inverse at f points224 218 225 219 INTEGER, PUBLIC :: nla10 !: deepest W level Above ~10m (nlb10 - 1) … … 265 259 266 260 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) :: tpol, fpol !: north fold mask (jperio= 3 or 4) 267 268 #if defined key_noslip_accurate269 INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,: ) :: npcoa !: ???270 INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: nicoa, njcoa !: ???271 #endif272 261 273 262 !!---------------------------------------------------------------------- … … 333 322 INTEGER FUNCTION dom_oce_alloc() 334 323 !!---------------------------------------------------------------------- 335 INTEGER, DIMENSION(1 2) :: ierr324 INTEGER, DIMENSION(13) :: ierr 336 325 !!---------------------------------------------------------------------- 337 326 ierr(:) = 0 … … 346 335 & tpol(jpiglo) , fpol(jpiglo) , STAT=ierr(2) ) 347 336 ! 348 ALLOCATE( glamt(jpi,jpj) , gphit(jpi,jpj) , e1t(jpi,jpj) , e2t(jpi,jpj) , r1_e1t(jpi,jpj) , r1_e2t(jpi,jpj) , & 349 & glamu(jpi,jpj) , gphiu(jpi,jpj) , e1u(jpi,jpj) , e2u(jpi,jpj) , r1_e1u(jpi,jpj) , r1_e2u(jpi,jpj) , & 350 & glamv(jpi,jpj) , gphiv(jpi,jpj) , e1v(jpi,jpj) , e2v(jpi,jpj) , r1_e1v(jpi,jpj) , r1_e2v(jpi,jpj) , & 351 & glamf(jpi,jpj) , gphif(jpi,jpj) , e1f(jpi,jpj) , e2f(jpi,jpj) , r1_e1f(jpi,jpj) , r1_e2f(jpi,jpj) , & 352 & e1e2t(jpi,jpj) , ff (jpi,jpj) , STAT=ierr(3) ) 337 ALLOCATE( glamt(jpi,jpj) , glamu(jpi,jpj) , glamv(jpi,jpj) , glamf(jpi,jpj) , & 338 & gphit(jpi,jpj) , gphiu(jpi,jpj) , gphiv(jpi,jpj) , gphif(jpi,jpj) , & 339 & e1t (jpi,jpj) , e2t (jpi,jpj) , r1_e1t(jpi,jpj) , r1_e2t(jpi,jpj) , & 340 & e1u (jpi,jpj) , e2u (jpi,jpj) , r1_e1u(jpi,jpj) , r1_e2u(jpi,jpj) , & 341 & e1v (jpi,jpj) , e2v (jpi,jpj) , r1_e1v(jpi,jpj) , r1_e2v(jpi,jpj) , & 342 & e1f (jpi,jpj) , e2f (jpi,jpj) , r1_e1f(jpi,jpj) , r1_e2f(jpi,jpj) , & 343 & e1e2t(jpi,jpj) , r1_e1e2t(jpi,jpj) , & 344 & e1e2u(jpi,jpj) , r1_e1e2u(jpi,jpj) , e2_e1u(jpi,jpj) , & 345 & e1e2v(jpi,jpj) , r1_e1e2v(jpi,jpj) , e1_e2v(jpi,jpj) , & 346 & e1e2f(jpi,jpj) , r1_e1e2f(jpi,jpj) , & 347 & ff (jpi,jpj) , STAT=ierr(3) ) 353 348 ! 354 349 ALLOCATE( gdep3w_0(jpi,jpj,jpk) , e3v_0(jpi,jpj,jpk) , e3f_0 (jpi,jpj,jpk) , & … … 364 359 & gdept_b (jpi,jpj,jpk) ,gdepw_b(jpi,jpj,jpk) , e3w_b (jpi,jpj,jpk) , & 365 360 & e3t_a (jpi,jpj,jpk) , e3u_a (jpi,jpj,jpk) , e3v_a (jpi,jpj,jpk) , & 366 & ehu_a (jpi,jpj) , ehv_a(jpi,jpj), &367 & ehur_a (jpi,jpj) , ehvr_a(jpi,jpj), &368 & ehu_b (jpi,jpj) , ehv_b(jpi,jpj), &369 & ehur_b (jpi,jpj) , ehvr_b(jpi,jpj), STAT=ierr(5) )361 & ehu_a (jpi,jpj) , ehv_a (jpi,jpj), & 362 & ehur_a (jpi,jpj) , ehvr_a(jpi,jpj), & 363 & ehu_b (jpi,jpj) , ehv_b (jpi,jpj), & 364 & ehur_b (jpi,jpj) , ehvr_b(jpi,jpj), STAT=ierr(5) ) 370 365 #endif 371 366 ! 372 ALLOCATE( hu (jpi,jpj) , hur (jpi,jpj) , hu_0(jpi,jpj) , ht_0 (jpi,jpj) , & 373 & hv (jpi,jpj) , hvr (jpi,jpj) , hv_0(jpi,jpj) , ht (jpi,jpj) , & 374 & re2u_e1u(jpi,jpj) , re1v_e2v(jpi,jpj) , & 375 & e12t (jpi,jpj) , r1_e12t (jpi,jpj) , & 376 & e12u (jpi,jpj) , r1_e12u (jpi,jpj) , & 377 & e12v (jpi,jpj) , r1_e12v (jpi,jpj) , & 378 & e12f (jpi,jpj) , r1_e12f (jpi,jpj) , STAT=ierr(6) ) 367 ALLOCATE( hu(jpi,jpj) , hur(jpi,jpj) , hu_0(jpi,jpj) , ht_0(jpi,jpj) , & 368 & hv(jpi,jpj) , hvr(jpi,jpj) , hv_0(jpi,jpj) , ht (jpi,jpj) , STAT=ierr(6) ) 379 369 ! 380 370 ALLOCATE( gdept_1d(jpk) , gdepw_1d(jpk) , & … … 387 377 & scosrf(jpi,jpj) , scobot(jpi,jpj) , & 388 378 & hifv (jpi,jpj) , hiff (jpi,jpj) , & 389 & hift (jpi,jpj) , hifu (jpi,jpj) , rx1 379 & hift (jpi,jpj) , hifu (jpi,jpj) , rx1(jpi,jpj) , STAT=ierr(8) ) 390 380 391 381 ALLOCATE( mbathy(jpi,jpj) , bathy(jpi,jpj) , & 392 382 & tmask_i(jpi,jpj) , umask_i(jpi,jpj), vmask_i(jpi,jpj), fmask_i(jpi,jpj), & 393 & bmask (jpi,jpj), &383 & bmask (jpi,jpj) , & 394 384 & mbkt (jpi,jpj) , mbku (jpi,jpj) , mbkv(jpi,jpj) , STAT=ierr(9) ) 395 385 396 386 ! (ISF) Allocation of basic array 397 ALLOCATE( misfdep(jpi,jpj) , risfdep(jpi,jpj), &398 & mikt(jpi,jpj), miku(jpi,jpj), mikv(jpi,jpj) , &399 & mikf(jpi,jpj), ssmask(jpi,jpj), STAT=ierr(10) )387 ALLOCATE( misfdep(jpi,jpj) , risfdep(jpi,jpj), & 388 & mikt(jpi,jpj), miku(jpi,jpj), mikv(jpi,jpj) , & 389 & mikf(jpi,jpj), ssmask(jpi,jpj), STAT=ierr(10) ) 400 390 401 391 ALLOCATE( tmask(jpi,jpj,jpk) , umask(jpi,jpj,jpk), & … … 403 393 404 394 ALLOCATE( wmask(jpi,jpj,jpk) , wumask(jpi,jpj,jpk), wvmask(jpi,jpj,jpk) , STAT=ierr(12) ) 405 406 #if defined key_noslip_accurate407 ALLOCATE( npcoa(4,jpk), nicoa(2*(jpi+jpj),4,jpk), njcoa(2*(jpi+jpj),4,jpk), STAT=ierr(12) )408 #endif409 395 ! 410 396 dom_oce_alloc = MAXVAL(ierr)
Note: See TracChangeset
for help on using the changeset viewer.