New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 3680 for branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/TOP_SRC/prtctl_trc.F90 – NEMO

Ignore:
Timestamp:
2012-11-27T15:42:24+01:00 (11 years ago)
Author:
rblod
Message:

First commit of the final branch for 2012 (future nemo_3_5), see ticket #1028

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/TOP_SRC/prtctl_trc.F90

    r3294 r3680  
    1717   USE par_trc          ! TOP parameters 
    1818   USE oce_trc          ! ocean space and time domain variables 
     19   USE prtctl           ! print control for OPA 
    1920 
    2021   IMPLICIT NONE 
     
    296297   END SUBROUTINE prt_ctl_trc_init 
    297298 
    298  
    299    SUBROUTINE sub_dom 
    300       !!---------------------------------------------------------------------- 
    301       !!                  ***  ROUTINE sub_dom  *** 
    302       !!                     
    303       !! ** Purpose :   Lay out the global domain over processors.  
    304       !!                CAUTION:  
    305       !!                This part has been extracted from the mpp_init 
    306       !!                subroutine and names of variables/arrays have been  
    307       !!                slightly changed to avoid confusion but the computation 
    308       !!                is exactly the same. Any modification about indices of 
    309       !!                each sub-domain in the mppini.F90 module should be reported  
    310       !!                here. 
    311       !! 
    312       !! ** Method  :   Global domain is distributed in smaller local domains. 
    313       !!                Periodic condition is a function of the local domain position 
    314       !!                (global boundary or neighbouring domain) and of the global 
    315       !!                periodic 
    316       !!                Type :         jperio global periodic condition 
    317       !!                               nperio local  periodic condition 
    318       !! 
    319       !! ** Action  : - set domain parameters 
    320       !!                    nimpp     : longitudinal index  
    321       !!                    njmpp     : latitudinal  index 
    322       !!                    nperio    : lateral condition type  
    323       !!                    narea     : number for local area 
    324       !!                    nlcil      : first dimension 
    325       !!                    nlcjl      : second dimension 
    326       !!                    nbondil    : mark for "east-west local boundary" 
    327       !!                    nbondjl    : mark for "north-south local boundary" 
    328       !!---------------------------------------------------------------------- 
    329       INTEGER ::   ji, jj, js               ! dummy loop indices 
    330       INTEGER ::   ii, ij                   ! temporary integers 
    331       INTEGER ::   irestil, irestjl         !    "          " 
    332       INTEGER ::   ijpi  , ijpj, nlcil      ! temporary logical unit 
    333       INTEGER ::   nlcjl , nbondil, nbondjl 
    334       INTEGER ::   nrecil, nrecjl, nldil, nleil, nldjl, nlejl 
    335       REAL(wp) ::   zidom, zjdom            ! temporary scalars 
    336       INTEGER, POINTER, DIMENSION(:,:) ::   iimpptl, ijmpptl, ilcitl, ilcjtl   ! temporary workspace 
    337       !!---------------------------------------------------------------------- 
    338       ! 
    339       CALL wrk_alloc( isplt, jsplt, ilcitl, ilcjtl, iimpptl, ijmpptl ) 
    340       ! 
    341       ! Dimension arrays for subdomains 
    342       ! ------------------------------- 
    343       !  Computation of local domain sizes ilcitl() ilcjtl() 
    344       !  These dimensions depend on global sizes isplt,jsplt and jpiglo,jpjglo 
    345       !  The subdomains are squares leeser than or equal to the global 
    346       !  dimensions divided by the number of processors minus the overlap 
    347       !  array (cf. par_oce.F90). 
    348  
    349       ijpi = ( jpiglo-2*jpreci + (isplt-1) ) / isplt + 2*jpreci 
    350       ijpj = ( jpjglo-2*jprecj + (jsplt-1) ) / jsplt + 2*jprecj 
    351  
    352       nrecil  = 2 * jpreci 
    353       nrecjl  = 2 * jprecj 
    354       irestil = MOD( jpiglo - nrecil , isplt ) 
    355       irestjl = MOD( jpjglo - nrecjl , jsplt ) 
    356  
    357       IF(  irestil == 0 )   irestil = isplt 
    358       DO jj = 1, jsplt 
    359          DO ji = 1, irestil 
    360             ilcitl(ji,jj) = ijpi 
    361          END DO 
    362          DO ji = irestil+1, isplt 
    363             ilcitl(ji,jj) = ijpi -1 
    364          END DO 
    365       END DO 
    366        
    367       IF( irestjl == 0 )   irestjl = jsplt 
    368       DO ji = 1, isplt 
    369          DO jj = 1, irestjl 
    370             ilcjtl(ji,jj) = ijpj 
    371          END DO 
    372          DO jj = irestjl+1, jsplt 
    373             ilcjtl(ji,jj) = ijpj -1 
    374          END DO 
    375       END DO 
    376        
    377       zidom = nrecil 
    378       DO ji = 1, isplt 
    379          zidom = zidom + ilcitl(ji,1) - nrecil 
    380       END DO 
    381        
    382       zjdom = nrecjl 
    383       DO jj = 1, jsplt 
    384          zjdom = zjdom + ilcjtl(1,jj) - nrecjl 
    385       END DO 
    386  
    387       ! Index arrays for subdomains 
    388       ! --------------------------- 
    389  
    390       iimpptl(:,:) = 1 
    391       ijmpptl(:,:) = 1 
    392        
    393       IF( isplt > 1 ) THEN 
    394          DO jj = 1, jsplt 
    395             DO ji = 2, isplt 
    396                iimpptl(ji,jj) = iimpptl(ji-1,jj) + ilcitl(ji-1,jj) - nrecil 
    397             END DO 
    398          END DO 
    399       ENDIF 
    400  
    401       IF( jsplt > 1 ) THEN 
    402          DO jj = 2, jsplt 
    403             DO ji = 1, isplt 
    404                ijmpptl(ji,jj) = ijmpptl(ji,jj-1)+ilcjtl(ji,jj-1)-nrecjl 
    405             END DO 
    406          END DO 
    407       ENDIF 
    408        
    409       ! Subdomain description 
    410       ! --------------------- 
    411  
    412       DO js = 1, ijsplt 
    413          ii = 1 + MOD( js-1, isplt ) 
    414          ij = 1 + (js-1) / isplt 
    415          nimpptl(js) = iimpptl(ii,ij) 
    416          njmpptl(js) = ijmpptl(ii,ij) 
    417          nlcitl (js) = ilcitl (ii,ij)      
    418          nlcil       = nlcitl (js)      
    419          nlcjtl (js) = ilcjtl (ii,ij)      
    420          nlcjl       = nlcjtl (js) 
    421          nbondjl = -1                                    ! general case 
    422          IF( js   >  isplt          )   nbondjl = 0      ! first row of processor 
    423          IF( js   >  (jsplt-1)*isplt )  nbondjl = 1     ! last  row of processor 
    424          IF( jsplt == 1             )   nbondjl = 2      ! one processor only in j-direction 
    425          ibonjtl(js) = nbondjl 
    426           
    427          nbondil = 0                                     !  
    428          IF( MOD( js, isplt ) == 1 )   nbondil = -1      ! 
    429          IF( MOD( js, isplt ) == 0 )   nbondil =  1      ! 
    430          IF( isplt            == 1 )   nbondil =  2      ! one processor only in i-direction 
    431          ibonitl(js) = nbondil 
    432           
    433          nldil =  1   + jpreci 
    434          nleil = nlcil - jpreci 
    435          IF( nbondil == -1 .OR. nbondil == 2 )   nldil = 1 
    436          IF( nbondil ==  1 .OR. nbondil == 2 )   nleil = nlcil 
    437          nldjl =  1   + jprecj 
    438          nlejl = nlcjl - jprecj 
    439          IF( nbondjl == -1 .OR. nbondjl == 2 )   nldjl = 1 
    440          IF( nbondjl ==  1 .OR. nbondjl == 2 )   nlejl = nlcjl 
    441          nlditl(js) = nldil 
    442          nleitl(js) = nleil 
    443          nldjtl(js) = nldjl 
    444          nlejtl(js) = nlejl 
    445       END DO 
    446       ! 
    447       CALL wrk_dealloc( isplt, jsplt, ilcitl, ilcjtl, iimpptl, ijmpptl ) 
    448       ! 
    449    END SUBROUTINE sub_dom 
    450   
    451299#else 
    452300   !!---------------------------------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.