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 9213 – NEMO

Changeset 9213


Ignore:
Timestamp:
2018-01-12T10:38:50+01:00 (6 years ago)
Author:
gm
Message:

dev_merge_2017: nemogcm.F90 : updated in SAS & OFF + data assimilation initial calls (asm_bkg_wri , tra_asm_inc ...) moved to asm_inc_init + closed sea : restructure namcfg & its control print + set ln_closea = false if domcfg file not read (ln_domcfg=F

Location:
branches/2017/dev_merge_2017/NEMOGCM
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_merge_2017/NEMOGCM/CONFIG/ORCA2_LIM3_PISCES/EXP00/namelist_cfg

    r9198 r9213  
    1616      !                    !  (=F) user defined configuration  ==>>>  see usrdef(_...) modules 
    1717      cn_domcfg = "ORCA_R2_zps_domcfg"    ! domain configuration filename 
    18       !                    !  (=F) user defined configuration  ==>>>  see usrdef(_...) modules 
    19    ln_closea    = .false.  !  T => keep closed seas (defined by closea_mask field) in the domain and apply 
    20       !                    !      special treatment of freshwater fluxes. 
    21       !                    !  F => suppress closed seas (defined by closea_mask field) from the bathymetry 
    22       !                    !      at runtime. 
    23       !                    !  If there is no closea_mask field in the domain_cfg file or we do not use 
    24       !                    !  a domain_cfg file then this logical does nothing. 
     18      ! 
     19      ln_closea    = .false.    !  T => keep closed seas (defined by closea_mask field) in the   
     20      !                         !       domain and apply special treatment of freshwater fluxes. 
     21      !                         !  F => suppress closed seas (defined by closea_mask field)  
     22      !                         !       from the bathymetry at runtime. 
     23      !                         !  If closea_mask field doesn't exist in the domain_cfg file 
     24      !                         !      then this logical does nothing. 
    2525/ 
    2626!----------------------------------------------------------------------- 
     
    2828!----------------------------------------------------------------------- 
    2929   ln_linssh   = .false.   !  =T  linear free surface  ==>>  model level are fixed in time 
     30   rn_rdt      = 5760.     !  time step for the dynamics and tracer 
     31   rn_atfp     =    0.1    !  asselin time filter parameter 
    3032/ 
    3133!----------------------------------------------------------------------- 
     
    299301/ 
    300302!----------------------------------------------------------------------- 
    301 &namctl        !   Control prints & Benchmark 
    302 !----------------------------------------------------------------------- 
     303&namctl        !   Control prints  
     304!----------------------------------------------------------------------- 
     305   ln_ctl      = .false.   !  trends control print (expensive!) 
     306   nn_print    =    0      !  level of print (0 no extra print) 
     307   ln_timing   = .false.   !  timing by routine write out in timing.output file 
     308   ln_diacfl   = .false.   !  CFL diagnostics write out in cfl_diagnostics.ascii 
    303309/ 
    304310!----------------------------------------------------------------------- 
  • branches/2017/dev_merge_2017/NEMOGCM/CONFIG/SHARED/namelist_ref

    r9198 r9213  
    7474   ln_read_cfg = .false.   !  (=T) read the domain configuration file 
    7575      !                    !  (=F) user defined configuration  ==>>>  see usrdef(_...) modules 
    76       cn_domcfg = "domain_cfg"         ! domain configuration filename 
     76      cn_domcfg = "domain_cfg"  ! domain configuration filename 
    7777      ! 
     78      ln_closea    = .false.    !  T => keep closed seas (defined by closea_mask field) in the   
     79      !                         !       domain and apply special treatment of freshwater fluxes. 
     80      !                         !  F => suppress closed seas (defined by closea_mask field)  
     81      !                         !       from the bathymetry at runtime. 
     82      !                         !  If closea_mask field doesn't exist in the domain_cfg file 
     83      !                         !       then this logical does nothing. 
    7884   ln_write_cfg= .false.   !  (=T) create the domain configuration file 
    7985      cn_domcfg_out = "domain_cfg_out" ! newly created domain configuration filename 
     
    8187   ln_use_jattr = .false.  !  use (T) the file attribute: open_ocean_jstart, if present 
    8288   !                       !  in netcdf input files, as the start j-row for reading 
    83    ln_closea    = .true.   !  T => keep closed seas (defined by closea_mask field) in the domain and apply 
    84       !                    !       special treatment of freshwater fluxes. 
    85       !                    !  F => suppress closed seas (defined by closea_mask field) from the bathymetry 
    86       !                    !       at runtime. 
    87       !                    !  If there is no closea_mask field in the domain_cfg file or we do not use 
    88       !                    !  a domain_cfg file then this logical does nothing. 
    8989/ 
    9090!----------------------------------------------------------------------- 
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/OFF_SRC/nemogcm.F90

    r9169 r9213  
    1010 
    1111   !!---------------------------------------------------------------------- 
    12    !!   nemo_gcm        : off-line: solve ocean tracer only 
    13    !!   nemo_init       : initialization of the nemo model 
    14    !!   nemo_ctl        : initialisation of algorithm flag  
    15    !!   nemo_closefile  : close remaining files 
     12   !!   nemo_gcm      : off-line: solve ocean tracer only 
     13   !!   nemo_gcm      : solve ocean dynamics, tracer, biogeochemistry and/or sea-ice 
     14   !!   nemo_init     : initialization of the NEMO system 
     15   !!   nemo_ctl      : initialisation of the contol print 
     16   !!   nemo_closefile: close remaining open files 
     17   !!   nemo_alloc    : dynamical allocation 
     18   !!   nemo_partition: calculate MPP domain decomposition 
     19   !!   factorise     : calculate the factors of the no. of MPI processes 
     20   !!   nemo_nfdcom   : Setup for north fold exchanges with explicit point-to-point messaging 
     21   !!   istate_init   : simple initialization to zero of ocean fields 
     22   !!   stp_ctl       : reduced step control (no dynamics in off-line) 
    1623   !!---------------------------------------------------------------------- 
    17    USE dom_oce         ! ocean space domain variables 
    18    USE oce             ! dynamics and tracers variables 
    19    USE trc_oce         ! Shared ocean/passive tracers variables 
    20    USE c1d             ! 1D configuration 
    21    USE domain          ! domain initialization from coordinate & bathymetry (dom_init routine) 
    22    USE usrdef_nam      ! user defined configuration 
    23    USE eosbn2          ! equation of state            (eos bn2 routine) 
     24   USE dom_oce        ! ocean space domain variables 
     25   USE oce            ! dynamics and tracers variables 
     26   USE trc_oce        ! Shared ocean/passive tracers variables 
     27   USE c1d            ! 1D configuration 
     28   USE domain         ! domain initialization from coordinate & bathymetry (dom_init routine) 
     29   USE closea         ! treatment of closed seas (for ln_closea) 
     30   USE usrdef_nam     ! user defined configuration 
     31   USE eosbn2         ! equation of state            (eos bn2 routine) 
    2432   !              ! ocean physics 
    25    USE ldftra          ! lateral diffusivity setting    (ldf_tra_init routine) 
    26    USE ldfslp          ! slopes of neutral surfaces     (ldf_slp_init routine) 
    27    USE traqsr          ! solar radiation penetration    (tra_qsr_init routine) 
    28    USE trabbl          ! bottom boundary layer          (tra_bbl_init routine) 
    29    USE traldf          ! lateral physics                (tra_ldf_init routine) 
    30    USE sbcmod          ! surface boundary condition     (sbc_init     routine) 
    31    USE phycst          ! physical constant                   (par_cst routine) 
    32    USE dtadyn          ! Lecture and Interpolation of the dynamical fields 
    33    USE trcini          ! Initilization of the passive tracers 
    34    USE daymod          ! calendar                            (day     routine) 
    35    USE trcstp          ! passive tracer time-stepping        (trc_stp routine) 
    36    USE dtadyn          ! Lecture and interpolation of the dynamical fields 
     33   USE ldftra         ! lateral diffusivity setting    (ldf_tra_init routine) 
     34   USE ldfslp         ! slopes of neutral surfaces     (ldf_slp_init routine) 
     35   USE traqsr         ! solar radiation penetration    (tra_qsr_init routine) 
     36   USE trabbl         ! bottom boundary layer          (tra_bbl_init routine) 
     37   USE traldf         ! lateral physics                (tra_ldf_init routine) 
     38   USE sbcmod         ! surface boundary condition     (sbc_init     routine) 
     39   USE phycst         ! physical constant                   (par_cst routine) 
     40   USE dtadyn         ! Lecture and Interpolation of the dynamical fields 
     41   USE trcini         ! Initilization of the passive tracers 
     42   USE daymod         ! calendar                            (day     routine) 
     43   USE trcstp         ! passive tracer time-stepping        (trc_stp routine) 
     44   USE dtadyn         ! Lecture and interpolation of the dynamical fields 
    3745   !              ! Passive tracers needs 
    38    USE trc             ! passive tracer : variables 
    39    USE trcnam          ! passive tracer : namelist 
    40    USE trcrst          ! passive tracer restart 
    41    USE diaptr          ! Need to initialise this as some variables are used in if statements later 
    42    USE sbc_oce  , ONLY : ln_rnf 
    43    USE sbcrnf          ! surface boundary condition : runoffs 
     46   USE trc            ! passive tracer : variables 
     47   USE trcnam         ! passive tracer : namelist 
     48   USE trcrst         ! passive tracer restart 
     49   USE diaptr         ! Need to initialise this as some variables are used in if statements later 
     50   USE sbc_oce , ONLY : ln_rnf 
     51   USE sbcrnf         ! surface boundary condition : runoffs 
    4452   !              ! I/O & MPP 
    45    USE iom             ! I/O library 
    46    USE in_out_manager  ! I/O manager 
    47    USE mppini          ! shared/distributed memory setting (mpp_init routine) 
    48    USE lib_mpp         ! distributed memory computing 
     53   USE iom            ! I/O library 
     54   USE in_out_manager ! I/O manager 
     55   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
     56   USE lib_mpp        ! distributed memory computing 
    4957#if defined key_iomput 
    50    USE xios 
     58   USE xios           ! xIOserver 
    5159#endif  
    52    USE prtctl          ! Print control                    (prt_ctl_init routine) 
    53    USE timing          ! Timing 
    54    USE lib_fortran     ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    55    USE lbcnfd, ONLY: isendto, nsndto, nfsloop, nfeloop ! Setup of north fold exchanges 
    56  
    57  
     60   USE prtctl         ! Print control                    (prt_ctl_init routine) 
     61   USE timing         ! Timing 
     62   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
     63   USE lbcnfd  , ONLY : isendto, nsndto, nfsloop, nfeloop   ! Setup of north fold exchanges 
    5864 
    5965   IMPLICIT NONE 
     
    6571 
    6672   !!---------------------------------------------------------------------- 
    67    !! NEMO/OFF 3.3 , NEMO Consortium (2010) 
     73   !! NEMO/OFF 4.0 , NEMO Consortium (2018) 
    6874   !! $Id$ 
    6975   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    7581      !!                     ***  ROUTINE nemo_gcm  *** 
    7682      !! 
    77       !! ** Purpose :   nemo solves the primitive equations on an orthogonal 
    78       !!      curvilinear mesh on the sphere. 
     83      !! ** Purpose :   NEMO solves the primitive equations on an orthogonal 
     84      !!              curvilinear mesh on the sphere. 
    7985      !! 
    8086      !! ** Method  : - model general initialization 
     
    99105      istp = nit000 
    100106      ! 
    101       ! Initialize arrays of runoffs structures and read data from the namelist 
    102       IF ( ln_rnf ) CALL sbc_rnf(istp) 
     107      IF( ln_rnf )   CALL sbc_rnf(istp)   ! runoffs initialization  
    103108      !  
    104       CALL iom_init( cxios_context )            ! iom_put initialization (must be done after nemo_init for AGRIF+XIOS+OASIS) 
     109      CALL iom_init( cxios_context )      ! iom_put initialization (must be done after nemo_init for AGRIF+XIOS+OASIS) 
    105110      !  
    106       DO WHILE ( istp <= nitend .AND. nstop == 0 )    ! time stepping 
     111      DO WHILE ( istp <= nitend .AND. nstop == 0 )    !==  OFF time-stepping  ==! 
    107112         ! 
    108113         IF( istp /= nit000 )   CALL day        ( istp )         ! Calendar (day was already called at nit000 in day_init) 
     
    114119                                CALL stp_ctl    ( istp, indic )  ! Time loop: control and print 
    115120         istp = istp + 1 
    116          IF( lk_mpp )   CALL mpp_max( nstop ) 
    117121      END DO 
     122      ! 
    118123#if defined key_iomput 
    119124      CALL iom_context_finalize( cxios_context ) ! needed for XIOS+AGRIF 
     
    127132      IF( nstop /= 0 .AND. lwp ) THEN                 ! error print 
    128133         WRITE(numout,cform_err) 
    129          WRITE(numout,*) nstop, ' error have been found' 
     134         WRITE(numout,*) '   ==>>>   nemo_gcm: a total of ', nstop, ' errors have been found' 
     135         WRITE(numout,*) 
    130136      ENDIF 
    131137      ! 
     
    134140      CALL nemo_closefile 
    135141      ! 
    136 # if defined key_iomput 
    137       CALL xios_finalize             ! end mpp communications 
    138 # else 
    139       IF( lk_mpp )   CALL mppstop       ! end mpp communications 
    140 # endif 
     142#if defined key_iomput 
     143                     CALL xios_finalize   ! end mpp communications with xios 
     144#else 
     145      IF( lk_mpp )   CALL mppstop         ! end mpp communications 
     146#endif 
    141147      ! 
    142148   END SUBROUTINE nemo_gcm 
     
    145151   SUBROUTINE nemo_init 
    146152      !!---------------------------------------------------------------------- 
    147       !!                     ***  ROUTINE nemo_init *** 
     153      !!                     ***  ROUTINE nemo_init  *** 
    148154      !! 
    149155      !! ** Purpose :   initialization of the nemo model in off-line mode 
    150156      !!---------------------------------------------------------------------- 
    151157      INTEGER  ::   ji                 ! dummy loop indices 
    152       INTEGER  ::   ilocal_comm        ! local integer 
    153       INTEGER  ::   ios, inum          ! local integers 
    154       INTEGER  ::   iiarea, ijarea     ! local integers 
    155       INTEGER  ::   iirest, ijrest     ! local integers 
    156       REAL(wp) ::   ziglo, zjglo, zkglo, zperio   ! local scalars 
     158      INTEGER  ::   ios, ilocal_comm   ! local integers 
     159      INTEGER  ::   iiarea, ijarea     !   -       - 
     160      INTEGER  ::   iirest, ijrest     !   -       - 
    157161      CHARACTER(len=120), DIMENSION(30) ::   cltxt, cltxt2, clnam 
    158162      !! 
    159       NAMELIST/namctl/ ln_ctl  , nn_print, nn_ictls, nn_ictle,   & 
    160          &             nn_isplt, nn_jsplt, nn_jctls, nn_jctle,   & 
     163      NAMELIST/namctl/ ln_ctl   , nn_print, nn_ictls, nn_ictle,   & 
     164         &             nn_isplt , nn_jsplt, nn_jctls, nn_jctle,   & 
    161165         &             ln_timing, ln_diacfl 
    162  
    163       NAMELIST/namcfg/ ln_read_cfg, cn_domcfg, ln_write_cfg, cn_domcfg_out, ln_use_jattr 
    164       !!---------------------------------------------------------------------- 
     166      NAMELIST/namcfg/ ln_read_cfg, cn_domcfg, ln_closea, ln_write_cfg, cn_domcfg_out, ln_use_jattr 
     167      !!---------------------------------------------------------------------- 
     168      ! 
    165169      cltxt  = '' 
    166170      cltxt2 = '' 
     
    172176      CALL ctl_opn( numnam_cfg, 'namelist_cfg', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. ) 
    173177      ! 
    174       REWIND( numnam_ref )              ! Namelist namctl in reference namelist : Control prints & Benchmark 
     178      REWIND( numnam_ref )              ! Namelist namctl in reference namelist 
    175179      READ  ( numnam_ref, namctl, IOSTAT = ios, ERR = 901 ) 
    176180901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namctl in reference namelist', .TRUE. ) 
    177       REWIND( numnam_cfg )              ! Namelist namctl in confguration namelist : Control prints & Benchmark 
     181      REWIND( numnam_cfg )              ! Namelist namctl in confguration namelist 
    178182      READ  ( numnam_cfg, namctl, IOSTAT = ios, ERR = 902 ) 
    179183902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namctl in configuration namelist', .TRUE. ) 
    180  
    181       ! 
    182       REWIND( numnam_ref )              ! Namelist namcfg in reference namelist : Control prints & Benchmark 
     184      ! 
     185      REWIND( numnam_ref )              ! Namelist namcfg in reference namelist 
    183186      READ  ( numnam_ref, namcfg, IOSTAT = ios, ERR = 903 ) 
    184187903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namcfg in reference namelist', .TRUE. ) 
    185       REWIND( numnam_cfg )              ! Namelist namcfg in confguration namelist : Control prints & Benchmark 
     188      REWIND( numnam_cfg )              ! Namelist namcfg in confguration namelist 
    186189      READ  ( numnam_cfg, namcfg, IOSTAT = ios, ERR = 904 ) 
    187190904   IF( ios >  0 )   CALL ctl_nam ( ios , 'namcfg in configuration namelist', .TRUE. )    
    188191 
    189  
    190192      !                             !--------------------------! 
    191193      !                             !  Set global domain size  !   (control print return in cltxt2) 
    192       ! 
     194      !                             !--------------------------! 
    193195      IF( ln_read_cfg ) THEN              ! Read sizes in domain configuration file 
    194196         CALL domain_cfg ( cltxt2,        cn_cfg, nn_cfg, jpiglo, jpjglo, jpkglo, jperio ) 
     
    198200      ENDIF 
    199201      ! 
    200       ! 
    201       l_offline = .true.            !  passive tracers are run offline 
     202      l_offline = .true.                  ! passive tracers are run offline 
    202203      ! 
    203204      !                             !--------------------------------------------! 
     
    219220      lwp = (narea == 1) .OR. ln_ctl          ! control of all listing output print 
    220221 
    221       IF(lwm) THEN 
    222          ! write merged namelists from earlier to output namelist now that the 
    223          ! file has been opened in call to mynode. nammpp has already been 
    224          ! written in mynode (if lk_mpp_mpi) 
     222      IF(lwm) THEN               ! write merged namelists from earlier to output namelist  
     223         !                       ! now that the file has been opened in call to mynode.  
     224         !                       ! NB: nammpp has already been written in mynode (if lk_mpp_mpi) 
    225225         WRITE( numond, namctl ) 
    226226         WRITE( numond, namcfg ) 
    227227         IF( .NOT.ln_read_cfg ) THEN 
    228228            DO ji = 1, SIZE(clnam) 
    229                IF( TRIM(clnam (ji)) /= '' )   WRITE(numond, * ) clnam(ji)    ! namusr_def print 
     229               IF( TRIM(clnam(ji)) /= '' )   WRITE(numond, * ) clnam(ji)     ! namusr_def print 
    230230            END DO 
    231231         ENDIF 
     
    234234      ! If dimensions of processor grid weren't specified in the namelist file  
    235235      ! then we calculate them here now that we have our communicator size 
    236       IF( (jpni < 1) .OR. (jpnj < 1) )THEN 
     236      IF( jpni < 1  .OR.  jpnj < 1 ) THEN 
    237237#if   defined key_mpp_mpi 
    238          CALL nemo_partition(mppsize) 
     238         CALL nemo_partition( mppsize ) 
    239239#else 
    240          jpni = 1 
    241          jpnj = 1 
     240         jpni  = 1 
     241         jpnj  = 1 
    242242         jpnij = jpni*jpnj 
    243243#endif 
    244       END IF 
     244      ENDIF 
    245245 
    246246      iiarea = 1 + MOD( narea - 1 , jpni ) 
     
    279279         WRITE(numout,*) '   CNRS - NERC - Met OFFICE - MERCATOR-ocean - INGV - CMCC' 
    280280         WRITE(numout,*) '                       NEMO team' 
    281          WRITE(numout,*) '            Ocean General Circulation Model' 
    282          WRITE(numout,*) '                  version 3.6  (2015) ' 
    283          WRITE(numout,*) 
    284          WRITE(numout,*) 
    285          DO ji = 1, SIZE(cltxt)  
    286             IF( TRIM(cltxt(ji)) /= '' )   WRITE(numout,*) cltxt(ji)      ! control print of mynode 
     281         WRITE(numout,*) '                   Off-line TOP Model' 
     282         WRITE(numout,*) '                NEMO version 4.0  (2017) ' 
     283         WRITE(numout,*) 
     284         WRITE(numout,*) 
     285         DO ji = 1, SIZE(cltxt) 
     286            IF( TRIM(cltxt (ji)) /= '' )   WRITE(numout,*) cltxt(ji)    ! control print of mynode 
    287287         END DO 
    288          WRITE(numout,cform_aaa)                                         ! Flag AAAAAAA 
    289          ! 
    290       ENDIF 
    291  
    292       ! Now we know the dimensions of the grid and numout has been set we can  
    293       ! allocate arrays 
     288         WRITE(numout,*) 
     289         WRITE(numout,*) 
     290         DO ji = 1, SIZE(cltxt2) 
     291            IF( TRIM(cltxt2(ji)) /= '' )   WRITE(numout,*) cltxt2(ji)   ! control print of domain size 
     292         END DO 
     293         ! 
     294         WRITE(numout,cform_aaa)                                        ! Flag AAAAAAA 
     295         ! 
     296      ENDIF 
     297 
     298      ! Now we know the dimensions of the grid and numout has been set: we can allocate arrays 
    294299      CALL nemo_alloc() 
    295300 
    296       !                             !--------------------------------! 
    297       !                             !  Model general initialization  ! 
    298       !                             !--------------------------------! 
    299  
    300       CALL nemo_ctl                           ! Control prints & Benchmark 
     301      !                             !-------------------------------! 
     302      !                             !  NEMO general initialization  ! 
     303      !                             !-------------------------------! 
     304 
     305      CALL nemo_ctl                          ! Control prints 
    301306 
    302307      !                                      ! Domain decomposition 
    303       CALL mpp_init 
    304       ! 
     308      CALL mpp_init                             ! MPP 
     309      IF( ln_nnogather )   CALL nemo_nfdcom     ! northfold neighbour lists 
     310      ! 
     311      !                                      ! General initialization 
    305312      IF( ln_timing    )   CALL timing_init 
    306       ! 
    307  
    308       !                                      ! General initialization 
    309313      IF( ln_timing    )   CALL timing_start( 'nemo_init') 
    310314      ! 
    311                            CALL     phy_cst    ! Physical constants 
    312                            CALL     eos_init   ! Equation of state 
    313       IF( lk_c1d       )   CALL     c1d_init   ! 1D column configuration 
    314  
    315                            CALL     dom_init   ! Domain 
    316  
    317                            CALL  istate_init   ! ocean initial state (Dynamics and tracers) 
    318  
    319       IF( ln_nnogather )   CALL nemo_northcomms   ! Initialise the northfold neighbour lists (must be done after the masks are defined) 
    320  
    321       IF( ln_ctl       )   CALL prt_ctl_init   ! Print control 
    322  
    323                            CALL     sbc_init   ! Forcings : surface module 
    324  
    325                            CALL ldf_tra_init   ! Lateral ocean tracer physics 
    326                            CALL ldf_eiv_init   ! Eddy induced velocity param 
    327                            CALL tra_ldf_init   ! lateral mixing 
    328       IF( l_ldfslp     )   CALL ldf_slp_init   ! slope of lateral mixing 
    329  
    330                            CALL tra_qsr_init   ! penetrative solar radiation qsr 
    331       IF( ln_trabbl    )   CALL tra_bbl_init   ! advective (and/or diffusive) bottom boundary layer scheme 
    332  
     315                           CALL     phy_cst     ! Physical constants 
     316                           CALL     eos_init    ! Equation of state 
     317      IF( lk_c1d       )   CALL     c1d_init    ! 1D column configuration 
     318                           CALL     dom_init    ! Domain 
     319      IF( ln_ctl       )   CALL prt_ctl_init    ! Print control 
     320 
     321                           CALL  istate_init    ! ocean initial state (Dynamics and tracers) 
     322 
     323                           CALL     sbc_init    ! Forcings : surface module 
     324 
     325      !                                      ! Tracer physics 
     326                           CALL ldf_tra_init    ! Lateral ocean tracer physics 
     327                           CALL ldf_eiv_init    ! Eddy induced velocity param 
     328                           CALL tra_ldf_init    ! lateral mixing 
     329      IF( l_ldfslp     )   CALL ldf_slp_init    ! slope of lateral mixing 
     330                           CALL tra_qsr_init    ! penetrative solar radiation qsr 
     331      IF( ln_trabbl    )   CALL tra_bbl_init    ! advective (and/or diffusive) bottom boundary layer scheme 
     332 
     333      !                                      ! Passive tracers 
    333334                           CALL trc_nam_run    ! Needed to get restart parameters for passive tracers 
    334335                           CALL trc_rst_cal( nit000, 'READ' )   ! calendar 
     
    336337 
    337338                           CALL     trc_init   ! Passive tracers initialization 
    338                            CALL dia_ptr_init   ! Initialise diaptr as some variables are used  
    339       !                                         ! in various advection and diffusion routines 
     339                           CALL dia_ptr_init   ! Poleward TRansports initialization 
     340                            
    340341      IF(lwp) WRITE(numout,cform_aaa)           ! Flag AAAAAAA 
    341342      ! 
     
    357358         WRITE(numout,*) 
    358359         WRITE(numout,*) 'nemo_ctl: Control prints' 
    359          WRITE(numout,*) '~~~~~~~ ' 
     360         WRITE(numout,*) '~~~~~~~~' 
    360361         WRITE(numout,*) '   Namelist namctl' 
    361362         WRITE(numout,*) '      run control (for debugging)     ln_ctl     = ', ln_ctl 
     
    381382      IF(lwp) THEN                  ! control print 
    382383         WRITE(numout,*) 
    383          WRITE(numout,*) 'namcfg  : configuration initialization through namelist read' 
    384          WRITE(numout,*) '~~~~~~~ ' 
    385384         WRITE(numout,*) '   Namelist namcfg' 
    386          WRITE(numout,*) '      read domain configuration files             ln_read_cfg      = ', ln_read_cfg 
     385         WRITE(numout,*) '      read domain configuration file              ln_read_cfg      = ', ln_read_cfg 
    387386         WRITE(numout,*) '         filename to be read                         cn_domcfg     = ', TRIM(cn_domcfg) 
    388          WRITE(numout,*) '      write  configuration definition files       ln_write_cfg     = ', ln_write_cfg 
     387         WRITE(numout,*) '         keep closed seas in the domain (if exist)   ln_closea     = ', TRIM(cn_domcfg) 
     388         WRITE(numout,*) '      create a configuration definition file      ln_write_cfg     = ', ln_write_cfg 
    389389         WRITE(numout,*) '         filename to be written                      cn_domcfg_out = ', TRIM(cn_domcfg_out) 
    390390         WRITE(numout,*) '      use file attribute if exists as i/p j-start ln_use_jattr     = ', ln_use_jattr 
    391391      ENDIF 
    392  
     392      IF( .NOT.ln_read_cfg )   ln_closea = .false.   ! dealing possible only with a domcfg file 
     393      ! 
    393394      !                             ! Parameter control 
    394395      ! 
     
    430431      ENDIF 
    431432      ! 
    432       IF( 1_wp /= SIGN(1._wp,-0._wp)  )   CALL ctl_stop( 'nemo_ctl: The intrinsec SIGN function follows ',  & 
    433          &                                               'f2003 standard. '                              ,  & 
    434          &                                               'Compile with key_nosignedzero enabled' ) 
     433      IF( 1._wp /= SIGN(1._wp,-0._wp)  )   CALL ctl_stop( 'nemo_ctl: The intrinsec SIGN function follows f2003 standard.',  & 
     434         &                                                'Compile with key_nosignedzero enabled' ) 
    435435      ! 
    436436   END SUBROUTINE nemo_ctl 
     
    444444      !!---------------------------------------------------------------------- 
    445445      ! 
    446       IF ( lk_mpp ) CALL mppsync 
     446      IF( lk_mpp )  CALL mppsync 
    447447      ! 
    448448      CALL iom_close                                 ! close all input/output files managed by iom_* 
     
    453453      IF( numout     /=  6 )   CLOSE( numout     )   ! standard model output file 
    454454      IF( lwm.AND.numond  /= -1 )   CLOSE( numond          )   ! oce output namelist 
    455  
     455      ! 
    456456      numout = 6                                     ! redefine numout in case it is used after this point... 
    457457      ! 
     
    467467      !! ** Method  : 
    468468      !!---------------------------------------------------------------------- 
    469       USE diawri,       ONLY: dia_wri_alloc 
    470       USE dom_oce,      ONLY: dom_oce_alloc 
    471       USE zdf_oce,      ONLY: zdf_oce_alloc 
    472       USE trc_oce,      ONLY: trc_oce_alloc 
     469      USE diawri ,   ONLY : dia_wri_alloc 
     470      USE dom_oce,   ONLY : dom_oce_alloc 
     471      USE zdf_oce,   ONLY : zdf_oce_alloc 
     472      USE trc_oce,   ONLY : trc_oce_alloc 
    473473      ! 
    474474      INTEGER :: ierr 
    475475      !!---------------------------------------------------------------------- 
    476476      ! 
    477       ierr =        oce_alloc       ()          ! ocean  
    478       ierr = ierr + dia_wri_alloc   () 
    479       ierr = ierr + dom_oce_alloc   ()          ! ocean domain 
    480       ierr = ierr + zdf_oce_alloc   ()          ! ocean vertical physics 
    481       ! 
    482       ierr = ierr + trc_oce_alloc   ()          ! shared TRC / TRA arrays 
     477      ierr =        oce_alloc    ()          ! ocean  
     478      ierr = ierr + dia_wri_alloc() 
     479      ierr = ierr + dom_oce_alloc()          ! ocean domain 
     480      ierr = ierr + zdf_oce_alloc()          ! ocean vertical physics 
     481      ierr = ierr + trc_oce_alloc()          ! shared TRC / TRA arrays 
    483482      ! 
    484483      IF( lk_mpp    )   CALL mpp_sum( ierr ) 
     
    496495      !! ** Method  : 
    497496      !!---------------------------------------------------------------------- 
    498       INTEGER, INTENT(in) :: num_pes ! The number of MPI processes we have 
     497      INTEGER, INTENT(in) ::   num_pes  ! The number of MPI processes we have 
    499498      ! 
    500499      INTEGER, PARAMETER :: nfactmax = 20 
     
    505504      INTEGER, DIMENSION(nfactmax) :: ifact ! Array of factors 
    506505      !!---------------------------------------------------------------------- 
    507  
     506      ! 
    508507      ierr = 0 
    509  
     508      ! 
    510509      CALL factorise( ifact, nfactmax, nfact, num_pes, ierr ) 
    511  
     510      ! 
    512511      IF( nfact <= 1 ) THEN 
    513512         WRITE (numout, *) 'WARNING: factorisation of number of PEs failed' 
     
    540539      !! 
    541540      !! ** Purpose :   return the prime factors of n. 
    542       !!                knfax factors are returned in array kfax which is of  
     541      !!                knfax factors are returned in array kfax which is of 
    543542      !!                maximum dimension kmaxfax. 
    544543      !! ** Method  : 
     
    550549      INTEGER :: ifac, jl, inu 
    551550      INTEGER, PARAMETER :: ntest = 14 
    552       INTEGER :: ilfax(ntest) 
     551      INTEGER, DIMENSION(ntest) ::   ilfax 
     552      !!---------------------------------------------------------------------- 
    553553      ! 
    554554      ! lfax contains the set of allowed factors. 
    555555      ilfax(:) = (/(2**jl,jl=ntest,1,-1)/) 
    556  
     556      ! 
    557557      ! Clear the error flag and initialise output vars 
    558       kerr = 0 
    559       kfax = 1 
     558      kerr  = 0 
     559      kfax  = 1 
    560560      knfax = 0 
    561  
    562       ! Find the factors of n. 
    563       IF( kn .NE. 1 ) THEN 
    564  
     561      ! 
     562      IF( kn /= 1 ) THEN      ! Find the factors of n 
     563         ! 
    565564         ! nu holds the unfactorised part of the number. 
    566565         ! knfax holds the number of factors found. 
    567566         ! l points to the allowed factor list. 
    568567         ! ifac holds the current factor. 
    569     
     568         ! 
    570569         inu   = kn 
    571570         knfax = 0 
    572     
     571         ! 
    573572         DO jl = ntest, 1, -1 
    574573            ! 
    575574            ifac = ilfax(jl) 
    576575            IF( ifac > inu )   CYCLE 
    577     
     576            ! 
    578577            ! Test whether the factor will divide. 
    579     
     578            ! 
    580579            IF( MOD(inu,ifac) == 0 ) THEN 
    581580               ! 
     
    594593            ! 
    595594         END DO 
    596     
     595         ! 
    597596      ENDIF 
    598597      ! 
     
    600599 
    601600#if defined key_mpp_mpi 
    602    SUBROUTINE nemo_northcomms 
    603       !!====================================================================== 
    604       !!                     ***  ROUTINE  nemo_northcomms  *** 
    605       !! nemo_northcomms    :  Setup for north fold exchanges with explicit  
    606       !!                       point-to-point messaging 
    607       !!===================================================================== 
    608       !!---------------------------------------------------------------------- 
    609       !! 
    610       !! ** Purpose :   Initialization of the northern neighbours lists. 
     601 
     602   SUBROUTINE nemo_nfdcom 
     603      !!---------------------------------------------------------------------- 
     604      !!                     ***  ROUTINE  nemo_nfdcom  *** 
     605      !! ** Purpose :   Setup for north fold exchanges with explicit  
     606      !!                point-to-point messaging 
     607      !! 
     608      !! ** Method :   Initialization of the northern neighbours lists. 
    611609      !!---------------------------------------------------------------------- 
    612610      !!    1.0  ! 2011-10  (A. C. Coward, NOCS & J. Donners, PRACE) 
    613       !!    2.0  ! 2013-06 Setup avoiding MPI communication (I. Epicoco, S. 
    614       !Mocavero, CMCC)  
    615       !!---------------------------------------------------------------------- 
    616  
     611      !!    2.0  ! 2013-06 Setup avoiding MPI communication (I. Epicoco, S. Mocavero, CMCC)  
     612      !!---------------------------------------------------------------------- 
    617613      INTEGER  ::   sxM, dxM, sxT, dxT, jn 
    618614      INTEGER  ::   njmppmax 
    619  
     615      !!---------------------------------------------------------------------- 
     616      ! 
    620617      njmppmax = MAXVAL( njmppt ) 
    621  
     618      ! 
    622619      !initializes the north-fold communication variables 
    623620      isendto(:) = 0 
    624       nsndto = 0 
    625  
    626       !if I am a process in the north 
    627       IF ( njmpp == njmppmax ) THEN 
    628           !sxM is the first point (in the global domain) needed to compute the 
    629           !north-fold for the current process 
    630           sxM = jpiglo - nimppt(narea) - nlcit(narea) + 1 
    631           !dxM is the last point (in the global domain) needed to compute the 
    632           !north-fold for the current process 
    633           dxM = jpiglo - nimppt(narea) + 2 
    634  
    635           !loop over the other north-fold processes to find the processes 
    636           !managing the points belonging to the sxT-dxT range 
    637  
    638           DO jn = 1, jpni 
    639                 !sxT is the first point (in the global domain) of the jn 
    640                 !process 
    641                 sxT = nfiimpp(jn, jpnj) 
    642                 !dxT is the last point (in the global domain) of the jn 
    643                 !process 
    644                 dxT = nfiimpp(jn, jpnj) + nfilcit(jn, jpnj) - 1 
    645                 IF ((sxM .gt. sxT) .AND. (sxM .lt. dxT)) THEN 
    646                    nsndto = nsndto + 1 
    647                      isendto(nsndto) = jn 
    648                 ELSEIF ((sxM .le. sxT) .AND. (dxM .ge. dxT)) THEN 
    649                    nsndto = nsndto + 1 
    650                      isendto(nsndto) = jn 
    651                 ELSEIF ((dxM .lt. dxT) .AND. (sxT .lt. dxM)) THEN 
    652                    nsndto = nsndto + 1 
    653                      isendto(nsndto) = jn 
    654                 END IF 
    655           END DO 
    656           nfsloop = 1 
    657           nfeloop = nlci 
    658           DO jn = 2,jpni-1 
    659            IF(nfipproc(jn,jpnj) .eq. (narea - 1)) THEN 
    660               IF (nfipproc(jn - 1 ,jpnj) .eq. -1) THEN 
    661                  nfsloop = nldi 
    662               ENDIF 
    663               IF (nfipproc(jn + 1,jpnj) .eq. -1) THEN 
    664                  nfeloop = nlei 
    665               ENDIF 
    666            ENDIF 
    667         END DO 
    668  
     621      nsndto     = 0 
     622      ! 
     623      IF ( njmpp == njmppmax ) THEN      ! if I am a process in the north 
     624         ! 
     625         !sxM is the first point (in the global domain) needed to compute the north-fold for the current process 
     626         sxM = jpiglo - nimppt(narea) - nlcit(narea) + 1 
     627         !dxM is the last point (in the global domain) needed to compute the north-fold for the current process 
     628         dxM = jpiglo - nimppt(narea) + 2 
     629         ! 
     630         ! loop over the other north-fold processes to find the processes 
     631         ! managing the points belonging to the sxT-dxT range 
     632         ! 
     633         DO jn = 1, jpni 
     634            ! 
     635            sxT = nfiimpp(jn, jpnj)                            ! sxT = 1st  point (in the global domain) of the jn process 
     636            dxT = nfiimpp(jn, jpnj) + nfilcit(jn, jpnj) - 1    ! dxT = last point (in the global domain) of the jn process 
     637            ! 
     638            IF    ( sxT < sxM  .AND.  sxM < dxT ) THEN 
     639               nsndto          = nsndto + 1 
     640               isendto(nsndto) = jn 
     641            ELSEIF( sxM <= sxT  .AND.  dxM >= dxT ) THEN 
     642               nsndto          = nsndto + 1 
     643               isendto(nsndto) = jn 
     644            ELSEIF( dxM <  dxT  .AND.  sxT <  dxM ) THEN 
     645               nsndto          = nsndto + 1 
     646               isendto(nsndto) = jn 
     647            ENDIF 
     648            ! 
     649         END DO 
     650         nfsloop = 1 
     651         nfeloop = nlci 
     652         DO jn = 2,jpni-1 
     653            IF( nfipproc(jn,jpnj) == (narea - 1) ) THEN 
     654               IF( nfipproc(jn-1,jpnj) == -1 )   nfsloop = nldi 
     655               IF( nfipproc(jn+1,jpnj) == -1 )   nfeloop = nlei 
     656            ENDIF 
     657         END DO 
     658         ! 
    669659      ENDIF 
    670660      l_north_nogather = .TRUE. 
    671    END SUBROUTINE nemo_northcomms 
     661      ! 
     662   END SUBROUTINE nemo_nfdcom 
     663 
    672664#else 
    673    SUBROUTINE nemo_northcomms      ! Dummy routine 
    674       WRITE(*,*) 'nemo_northcomms: You should not have seen this print! error?' 
    675    END SUBROUTINE nemo_northcomms 
     665   SUBROUTINE nemo_nfdcom      ! Dummy routine 
     666      WRITE(*,*) 'nemo_nfdcom: You should not have seen this print! error?' 
     667   END SUBROUTINE nemo_nfdcom 
    676668#endif 
    677669 
     
    696688   END SUBROUTINE istate_init 
    697689 
     690 
    698691   SUBROUTINE stp_ctl( kt, kindic ) 
    699692      !!---------------------------------------------------------------------- 
     
    722715      ! 
    723716   END SUBROUTINE stp_ctl 
     717 
    724718   !!====================================================================== 
    725719END MODULE nemogcm 
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90

    r9168 r9213  
    7474   REAL(wp), PUBLIC, DIMENSION(:)    , ALLOCATABLE ::   wgtiau               !: IAU weights for each time step 
    7575#if defined key_asminc 
    76    REAL(wp), PUBLIC, DIMENSION(:,:), ALLOCATABLE ::   ssh_iau           !: IAU-weighted sea surface height increment 
     76   REAL(wp), PUBLIC, DIMENSION(:,:)  , ALLOCATABLE ::   ssh_iau              !: IAU-weighted sea surface height increment 
    7777#endif 
    7878   !                                !!! time steps relative to the cycle interval [0,nitend-nit000-1] 
     
    254254      !-------------------------------------------------------------------- 
    255255 
    256       IF ( ln_asmiau ) THEN 
    257  
     256      IF( ln_asmiau ) THEN 
     257         ! 
    258258         ALLOCATE( wgtiau( icycper ) ) 
    259  
     259         ! 
    260260         wgtiau(:) = 0._wp 
    261  
    262          IF ( niaufn == 0 ) THEN 
    263  
    264             !--------------------------------------------------------- 
    265             ! Constant IAU forcing  
    266             !--------------------------------------------------------- 
    267  
     261         ! 
     262         !                                !--------------------------------------------------------- 
     263         IF( niaufn == 0 ) THEN           ! Constant IAU forcing  
     264            !                             !--------------------------------------------------------- 
    268265            DO jt = 1, iiauper 
    269266               wgtiau(jt+nitiaustr-1) = 1.0 / REAL( iiauper ) 
    270267            END DO 
    271  
    272          ELSEIF ( niaufn == 1 ) THEN 
    273  
    274             !--------------------------------------------------------- 
    275             ! Linear hat-like, centred in middle of IAU interval  
    276             !--------------------------------------------------------- 
    277  
     268            !                             !--------------------------------------------------------- 
     269         ELSEIF ( niaufn == 1 ) THEN      ! Linear hat-like, centred in middle of IAU interval  
     270            !                             !--------------------------------------------------------- 
    278271            ! Compute the normalization factor 
    279             znorm = 0.0 
    280             IF ( MOD( iiauper, 2 ) == 0 ) THEN  ! Even number of time steps in IAU interval 
     272            znorm = 0._wp 
     273            IF( MOD( iiauper, 2 ) == 0 ) THEN   ! Even number of time steps in IAU interval 
    281274               imid = iiauper / 2  
    282275               DO jt = 1, imid 
     
    284277               END DO 
    285278               znorm = 2.0 * znorm 
    286             ELSE                               ! Odd number of time steps in IAU interval 
     279            ELSE                                ! Odd number of time steps in IAU interval 
    287280               imid = ( iiauper + 1 ) / 2         
    288281               DO jt = 1, imid - 1 
     
    292285            ENDIF 
    293286            znorm = 1.0 / znorm 
    294  
     287            ! 
    295288            DO jt = 1, imid - 1 
    296289               wgtiau(jt+nitiaustr-1) = REAL( jt ) * znorm 
     
    299292               wgtiau(jt+nitiaustr-1) = REAL( iiauper - jt + 1 ) * znorm 
    300293            END DO 
    301  
     294            ! 
    302295         ENDIF 
    303296 
     
    325318      !-------------------------------------------------------------------- 
    326319 
    327       ALLOCATE( t_bkginc(jpi,jpj,jpk) ) 
    328       ALLOCATE( s_bkginc(jpi,jpj,jpk) ) 
    329       ALLOCATE( u_bkginc(jpi,jpj,jpk) ) 
    330       ALLOCATE( v_bkginc(jpi,jpj,jpk) ) 
    331       ALLOCATE( ssh_bkginc(jpi,jpj)   ) 
    332       ALLOCATE( seaice_bkginc(jpi,jpj)) 
    333       t_bkginc     (:,:,:) = 0._wp 
    334       s_bkginc     (:,:,:) = 0._wp 
    335       u_bkginc     (:,:,:) = 0._wp 
    336       v_bkginc     (:,:,:) = 0._wp 
    337       ssh_bkginc   (:,:)   = 0._wp 
    338       seaice_bkginc(:,:)   = 0._wp 
     320      ALLOCATE( t_bkginc     (jpi,jpj,jpk) )   ;   t_bkginc     (:,:,:) = 0._wp 
     321      ALLOCATE( s_bkginc     (jpi,jpj,jpk) )   ;   s_bkginc     (:,:,:) = 0._wp 
     322      ALLOCATE( u_bkginc     (jpi,jpj,jpk) )   ;   u_bkginc     (:,:,:) = 0._wp 
     323      ALLOCATE( v_bkginc     (jpi,jpj,jpk) )   ;   v_bkginc     (:,:,:) = 0._wp 
     324      ALLOCATE( ssh_bkginc   (jpi,jpj)     )   ;   ssh_bkginc   (:,:)   = 0._wp 
     325      ALLOCATE( seaice_bkginc(jpi,jpj)     )   ;   seaice_bkginc(:,:)   = 0._wp 
    339326#if defined key_asminc 
    340       ALLOCATE( ssh_iau(jpi,jpj)      ) 
    341       ssh_iau      (:,:)   = 0._wp 
     327      ALLOCATE( ssh_iau      (jpi,jpj)     )   ;   ssh_iau      (:,:)   = 0._wp 
    342328#endif 
    343329#if defined key_cice && defined key_asminc 
    344       ALLOCATE( ndaice_da(jpi,jpj)    ) 
    345       ndaice_da    (:,:)   = 0._wp 
    346 #endif 
    347       IF ( ( ln_trainc ).OR.( ln_dyninc ).OR.( ln_sshinc ).OR.( ln_seaiceinc ) ) THEN 
    348  
    349          !-------------------------------------------------------------------- 
    350          ! Read the increments from file 
    351          !-------------------------------------------------------------------- 
    352  
     330      ALLOCATE( ndaice_da    (jpi,jpj)     )   ;   ndaice_da    (:,:)   = 0._wp 
     331#endif 
     332      ! 
     333      IF ( ln_trainc .OR. ln_dyninc .OR.   &       !-------------------------------------- 
     334         & ln_sshinc .OR. ln_seaiceinc   ) THEN    ! Read the increments from file 
     335         !                                         !-------------------------------------- 
    353336         CALL iom_open( c_asminc, inum ) 
    354  
    355          CALL iom_get( inum, 'time', zdate_inc )  
    356  
     337         ! 
     338         CALL iom_get( inum, 'time'       , zdate_inc   )  
    357339         CALL iom_get( inum, 'z_inc_dateb', z_inc_dateb ) 
    358340         CALL iom_get( inum, 'z_inc_datef', z_inc_datef ) 
    359341         z_inc_dateb = zdate_inc 
    360342         z_inc_datef = zdate_inc 
    361  
     343         ! 
    362344         IF(lwp) THEN 
    363345            WRITE(numout,*)  
    364             WRITE(numout,*) 'asm_inc_init : Assimilation increments valid ', & 
    365                &            ' between dates ', z_inc_dateb,' and ',  & 
    366                &            z_inc_datef 
     346            WRITE(numout,*) 'asm_inc_init : Assimilation increments valid between dates ', z_inc_dateb,' and ', z_inc_datef 
    367347            WRITE(numout,*) '~~~~~~~~~~~~' 
    368348         ENDIF 
    369  
    370          IF (     ( z_inc_dateb < ndastp + nn_time0*0.0001_wp ) & 
    371             & .OR.( z_inc_datef > ditend_date ) ) & 
    372             & CALL ctl_warn( ' Validity time of assimilation increments is ', & 
    373             &                ' outside the assimilation interval' ) 
     349         ! 
     350         IF ( ( z_inc_dateb < ndastp + nn_time0*0.0001_wp ) .OR.  & 
     351            & ( z_inc_datef > ditend_date ) ) & 
     352            &    CALL ctl_warn( ' Validity time of assimilation increments is ', & 
     353            &                   ' outside the assimilation interval' ) 
    374354 
    375355         IF ( ( ln_asmdin ).AND.( zdate_inc /= ditdin_date ) ) & 
     
    418398            WHERE( ABS( seaice_bkginc(:,:) ) > 1.0e+10 ) seaice_bkginc(:,:) = 0.0 
    419399         ENDIF 
    420  
     400         ! 
    421401         CALL iom_close( inum ) 
    422   
    423       ENDIF 
    424  
    425       !----------------------------------------------------------------------- 
    426       ! Apply divergence damping filter 
    427       !----------------------------------------------------------------------- 
    428  
    429       IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN 
    430          ! 
     402         ! 
     403      ENDIF 
     404      ! 
     405      !                                            !-------------------------------------- 
     406      IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN    ! Apply divergence damping filter 
     407         !                                         !-------------------------------------- 
    431408         ALLOCATE( zhdiv(jpi,jpj) )  
    432409         ! 
     
    460437         ! 
    461438      ENDIF 
    462  
    463       !----------------------------------------------------------------------- 
    464       ! Allocate and initialize the background state arrays 
    465       !----------------------------------------------------------------------- 
    466  
    467       IF ( ln_asmdin ) THEN 
    468          ! 
    469          ALLOCATE( t_bkg(jpi,jpj,jpk) ) 
    470          ALLOCATE( s_bkg(jpi,jpj,jpk) ) 
    471          ALLOCATE( u_bkg(jpi,jpj,jpk) ) 
    472          ALLOCATE( v_bkg(jpi,jpj,jpk) ) 
    473          ALLOCATE( ssh_bkg(jpi,jpj)   ) 
    474          ! 
    475          t_bkg(:,:,:) = 0._wp 
    476          s_bkg(:,:,:) = 0._wp 
    477          u_bkg(:,:,:) = 0._wp 
    478          v_bkg(:,:,:) = 0._wp 
    479          ssh_bkg(:,:) = 0._wp 
     439      ! 
     440      !                             !----------------------------------------------------- 
     441      IF ( ln_asmdin ) THEN         ! Allocate and initialize the background state arrays 
     442         !                          !----------------------------------------------------- 
     443         ! 
     444         ALLOCATE( t_bkg  (jpi,jpj,jpk) )   ;   t_bkg  (:,:,:) = 0._wp 
     445         ALLOCATE( s_bkg  (jpi,jpj,jpk) )   ;   s_bkg  (:,:,:) = 0._wp 
     446         ALLOCATE( u_bkg  (jpi,jpj,jpk) )   ;   u_bkg  (:,:,:) = 0._wp 
     447         ALLOCATE( v_bkg  (jpi,jpj,jpk) )   ;   v_bkg  (:,:,:) = 0._wp 
     448         ALLOCATE( ssh_bkg(jpi,jpj)     )   ;   ssh_bkg(:,:)   = 0._wp 
     449         ! 
    480450         ! 
    481451         !-------------------------------------------------------------------- 
     
    489459         IF(lwp) THEN 
    490460            WRITE(numout,*)  
    491             WRITE(numout,*) 'asm_inc_init : Assimilation background state valid at : ', & 
    492                &  zdate_bkg 
    493             WRITE(numout,*) '~~~~~~~~~~~~' 
    494          ENDIF 
    495          ! 
    496          IF ( zdate_bkg /= ditdin_date ) & 
     461            WRITE(numout,*) '   ==>>>  Assimilation background state valid at : ', zdate_bkg 
     462            WRITE(numout,*) 
     463         ENDIF 
     464         ! 
     465         IF ( zdate_bkg /= ditdin_date )   & 
    497466            & CALL ctl_warn( ' Validity time of assimilation background state does', & 
    498467            &                ' not agree with Direct Initialization time' ) 
     
    521490      ENDIF 
    522491      ! 
     492      IF(lwp) WRITE(numout,*) '   ==>>>   Euler time step switch is ', neuler 
     493      ! 
     494      IF( lk_asminc ) THEN                            !==  data assimilation  ==! 
     495         IF( ln_bkgwri )   CALL asm_bkg_wri( nit000 - 1 )      ! Output background fields 
     496         IF( ln_asmdin ) THEN                                  ! Direct initialization 
     497            IF( ln_trainc )   CALL tra_asm_inc( nit000 - 1 )      ! Tracers 
     498            IF( ln_dyninc )   CALL dyn_asm_inc( nit000 - 1 )      ! Dynamics 
     499            IF( ln_sshinc )   CALL ssh_asm_inc( nit000 - 1 )      ! SSH 
     500         ENDIF 
     501      ENDIF 
     502      ! 
    523503   END SUBROUTINE asm_inc_init 
     504    
     505    
    524506   SUBROUTINE tra_asm_inc( kt ) 
    525507      !!---------------------------------------------------------------------- 
     
    786768   END SUBROUTINE ssh_asm_inc 
    787769 
     770 
    788771   SUBROUTINE ssh_asm_div( kt, phdivn ) 
    789772      !!---------------------------------------------------------------------- 
     
    824807   END SUBROUTINE ssh_asm_div 
    825808 
     809 
    826810   SUBROUTINE seaice_asm_inc( kt, kindic ) 
    827811      !!---------------------------------------------------------------------- 
     
    886870            ! seaice salinity balancing (to add) 
    887871#endif 
    888  
     872            ! 
    889873#if defined key_cice && defined key_asminc 
    890874            ! Sea-ice : CICE case. Pass ice increment tendency into CICE 
    891875            ndaice_da(:,:) = seaice_bkginc(:,:) * zincwgt / rdt 
    892876#endif 
    893  
     877            ! 
    894878            IF ( kt == nitiaufin_r ) THEN 
    895879               DEALLOCATE( seaice_bkginc ) 
    896880            ENDIF 
    897  
     881            ! 
    898882         ELSE 
    899  
     883            ! 
    900884#if defined key_cice && defined key_asminc 
    901885            ndaice_da(:,:) = 0._wp        ! Sea-ice : CICE case. Zero ice increment tendency into CICE 
    902886#endif 
    903  
     887            ! 
    904888         ENDIF 
    905889         !                          !----------------------------------------- 
     
    949933#if defined key_cice && defined key_asminc 
    950934            ndaice_da(:,:) = 0._wp     ! Sea-ice : CICE case. Zero ice increment tendency into CICE 
    951  
    952935#endif 
    953936            ! 
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/nemogcm.F90

    r9210 r9213  
    2828   !!             -   ! 2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase 
    2929   !!            3.3.1! 2011-01  (A. R. Porter, STFC Daresbury) dynamical allocation 
    30    !!            3.4  ! 2011-10  (A. C. Coward, NOCS & J. Donners, PRACE) add nemo_northcomms 
     30   !!            3.4  ! 2011-10  (A. C. Coward, NOCS & J. Donners, PRACE) add nemo_nfdcom 
    3131   !!             -   ! 2011-11  (C. Harris) decomposition changes for running with CICE 
    3232   !!            3.6  ! 2012-05  (C. Calone, J. Simeon, G. Madec, C. Ethe) Add grid coarsening  
    33    !!             -   ! 2013-06  (I. Epicoco, S. Mocavero, CMCC) nemo_northcomms: setup avoiding MPI communication  
     33   !!             -   ! 2013-06  (I. Epicoco, S. Mocavero, CMCC) nemo_nfdcom: setup avoiding MPI communication  
    3434   !!             -   ! 2014-12  (G. Madec) remove KPP scheme and cross-land advection (cla) 
    3535   !!            4.0  ! 2016-10  (G. Madec, S. Flavoni)  domain configuration / user defined interface 
     
    4444   !!   nemo_partition: calculate MPP domain decomposition 
    4545   !!   factorise     : calculate the factors of the no. of MPI processes 
     46   !!   nemo_nfdcom   : Setup for north fold exchanges with explicit point-to-point messaging 
    4647   !!---------------------------------------------------------------------- 
    4748   USE step_oce       ! module used in the ocean time stepping module (step.F90) 
     
    8889   USE lib_mpp        ! distributed memory computing 
    8990   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
    90    USE lbcnfd , ONLY : isendto, nsndto, nfsloop, nfeloop   ! Setup of north fold exchanges  
     91   USE lbcnfd  , ONLY : isendto, nsndto, nfsloop, nfeloop   ! Setup of north fold exchanges  
    9192   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    9293#if defined key_iomput 
     
    104105 
    105106   !!---------------------------------------------------------------------- 
    106    !! NEMO/OPA 3.7 , NEMO Consortium (2016) 
     107   !! NEMO/OPA 4.0 , NEMO Consortium (2018) 
    107108   !! $Id$ 
    108109   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    130131      CALL Agrif_Init_Grids()      ! AGRIF: set the meshes 
    131132#endif 
    132       ! 
    133133      !                            !-----------------------! 
    134134      CALL nemo_init               !==  Initialisations  ==! 
     
    161161      END DO 
    162162#else 
    163  
    164 !!gm     This data assimilation calls should be part of the initialisation (i.e. put in asm_inc_init) 
    165       ! 
    166       IF( lk_asminc ) THEN                            !==  data assimilation  ==!   (done prior to time stepping) 
    167          IF( ln_bkgwri )   CALL asm_bkg_wri( nit000 - 1 )      ! Output background fields 
    168          IF( ln_asmdin ) THEN                                  ! Direct initialization 
    169             IF( ln_trainc )   CALL tra_asm_inc( nit000 - 1 )      ! Tracers 
    170             IF( ln_dyninc )   CALL dyn_asm_inc( nit000 - 1 )      ! Dynamics 
    171             IF( ln_sshinc )   CALL ssh_asm_inc( nit000 - 1 )      ! SSH 
    172          ENDIF 
    173       ENDIF 
    174 !!gm end 
    175163      ! 
    176164# if defined key_agrif 
    177165      !                                               !==  AGRIF time-stepping  ==! 
    178166      CALL Agrif_Regrid() 
     167      ! 
    179168      DO WHILE( istp <= nitend .AND. nstop == 0 ) 
    180169         CALL stp 
     
    222211      IF( nstop /= 0 .AND. lwp ) THEN        ! error print 
    223212         WRITE(numout,cform_err) 
    224          WRITE(numout,*) 'nemo_gcm: a total of ', nstop, ' errors have been found' 
     213         WRITE(numout,*) '   ==>>>   nemo_gcm: a total of ', nstop, ' errors have been found' 
    225214         WRITE(numout,*) 
    226215      ENDIF 
     
    249238      !!---------------------------------------------------------------------- 
    250239      INTEGER  ::   ji                 ! dummy loop indices 
    251       INTEGER  ::   ios, ilocal_comm   ! local integer 
    252       INTEGER  ::   iiarea, ijarea     ! local integers 
    253       INTEGER  ::   iirest, ijrest     ! local integers 
     240      INTEGER  ::   ios, ilocal_comm   ! local integers 
     241      INTEGER  ::   iiarea, ijarea     !   -       - 
     242      INTEGER  ::   iirest, ijrest     !   -       - 
    254243      CHARACTER(len=120), DIMENSION(30) ::   cltxt, cltxt2, clnam 
    255       ! 
     244      !! 
    256245      NAMELIST/namctl/ ln_ctl   , nn_print, nn_ictls, nn_ictle,   & 
    257246         &             nn_isplt , nn_jsplt, nn_jctls, nn_jctle,   & 
    258247         &             ln_timing, ln_diacfl 
    259       NAMELIST/namcfg/ ln_read_cfg, cn_domcfg, ln_write_cfg, cn_domcfg_out, ln_use_jattr, ln_closea 
     248      NAMELIST/namcfg/ ln_read_cfg, cn_domcfg, ln_closea, ln_write_cfg, cn_domcfg_out, ln_use_jattr 
    260249      !!---------------------------------------------------------------------- 
    261250      ! 
     
    269258      CALL ctl_opn( numnam_cfg, 'namelist_cfg', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. ) 
    270259      ! 
    271       REWIND( numnam_ref )              ! Namelist namctl in reference namelist : Control prints 
     260      REWIND( numnam_ref )              ! Namelist namctl in reference namelist 
    272261      READ  ( numnam_ref, namctl, IOSTAT = ios, ERR = 901 ) 
    273262901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namctl in reference namelist', .TRUE. ) 
     
    276265902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namctl in configuration namelist', .TRUE. ) 
    277266      ! 
    278       REWIND( numnam_ref )              ! Namelist namcfg in reference namelist : Control prints 
     267      REWIND( numnam_ref )              ! Namelist namcfg in reference namelist 
    279268      READ  ( numnam_ref, namcfg, IOSTAT = ios, ERR = 903 ) 
    280269903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namcfg in reference namelist', .TRUE. ) 
    281       REWIND( numnam_cfg )              ! Namelist namcfg in confguration namelist : Control prints & Benchmark 
     270      REWIND( numnam_cfg )              ! Namelist namcfg in confguration namelist 
    282271      READ  ( numnam_cfg, namcfg, IOSTAT = ios, ERR = 904 ) 
    283272904   IF( ios >  0 )   CALL ctl_nam ( ios , 'namcfg in configuration namelist', .TRUE. )    
     
    435424 
    436425      !                                      ! Domain decomposition 
    437       CALL mpp_init 
    438       IF( ln_nnogather )   CALL nemo_northcomms! northfold neighbour lists 
    439       ! 
    440       IF( ln_timing    )   CALL timing_init 
     426      CALL mpp_init                             ! MPP 
     427      IF( ln_nnogather )   CALL nemo_nfdcom     ! northfold neighbour lists 
    441428      ! 
    442429      !                                      ! General initialization 
    443                            CALL     phy_cst    ! Physical constants 
    444                            CALL     eos_init   ! Equation of state 
    445       IF( lk_c1d       )   CALL     c1d_init   ! 1D column configuration 
    446                            CALL     wad_init   ! Wetting and drying options 
    447                            CALL     dom_init   ! Domain 
    448       IF( ln_crs       )   CALL     crs_init   ! coarsened grid: domain initialization  
    449       IF( ln_ctl       )   CALL prt_ctl_init   ! Print control 
     430      IF( ln_timing    )   CALL timing_init     ! timing 
     431      IF( ln_timing    )   CALL timing_start( 'nemo_init') 
     432      ! 
     433                           CALL     phy_cst     ! Physical constants 
     434                           CALL     eos_init    ! Equation of state 
     435      IF( lk_c1d       )   CALL     c1d_init    ! 1D column configuration 
     436                           CALL     wad_init    ! Wetting and drying options 
     437                           CALL     dom_init    ! Domain 
     438      IF( ln_crs       )   CALL     crs_init    ! coarsened grid: domain initialization  
     439      IF( ln_ctl       )   CALL prt_ctl_init    ! Print control 
    450440       
    451       CALL diurnal_sst_bulk_init             ! diurnal sst 
    452       IF ( ln_diurnal ) CALL diurnal_sst_coolskin_init   ! cool skin    
     441      CALL diurnal_sst_bulk_init                ! diurnal sst 
     442      IF( ln_diurnal   )   CALL diurnal_sst_coolskin_init   ! cool skin    
     443      !                             
     444      IF( ln_diurnal_only ) THEN                   ! diurnal only: a subset of the initialisation routines 
     445         CALL  istate_init                            ! ocean initial state (Dynamics and tracers) 
     446         CALL     sbc_init                            ! Forcings : surface module 
     447         CALL tra_qsr_init                            ! penetrative solar radiation qsr 
     448         IF( ln_diaobs ) THEN                         ! Observation & model comparison 
     449            CALL dia_obs_init                            ! Initialize observational data 
     450            CALL dia_obs( nit000 - 1 )                   ! Observation operator for restart 
     451         ENDIF      
     452         IF( lk_asminc )   CALL asm_inc_init          ! Assimilation increments 
     453         ! 
     454         RETURN                                       ! end of initialization 
     455      ENDIF 
    453456       
    454       ! IF ln_diurnal_only, then we only want a subset of the initialisation routines 
    455       IF( ln_diurnal_only ) THEN 
    456          CALL  istate_init   ! ocean initial state (Dynamics and tracers) 
    457          CALL     sbc_init   ! Forcings : surface module 
    458          CALL tra_qsr_init   ! penetrative solar radiation qsr 
    459          IF( ln_diaobs ) THEN                   ! Observation & model comparison 
    460             CALL dia_obs_init                      ! Initialize observational data 
    461             CALL dia_obs( nit000 - 1 )             ! Observation operator for restart 
    462          ENDIF      
    463          !                                     ! Assimilation increments 
    464          IF( lk_asminc )   CALL asm_inc_init   ! Initialize assimilation increments 
    465                   
    466          IF(lwp) WRITE(numout,*) 'Euler time step switch is ', neuler 
    467          RETURN 
    468       ENDIF 
    469        
    470                            CALL  istate_init   ! ocean initial state (Dynamics and tracers) 
     457                           CALL  istate_init    ! ocean initial state (Dynamics and tracers) 
    471458 
    472459      !                                      ! external forcing  
    473                            CALL    tide_init   ! tidal harmonics 
    474                            CALL     sbc_init   ! surface boundary conditions (including sea-ice) 
    475                            CALL     bdy_init   ! Open boundaries initialisation 
     460                           CALL    tide_init    ! tidal harmonics 
     461                           CALL     sbc_init    ! surface boundary conditions (including sea-ice) 
     462                           CALL     bdy_init    ! Open boundaries initialisation 
    476463 
    477464      !                                      ! Ocean physics 
     
    520507                           CALL     trd_init    ! Mixed-layer/Vorticity/Integral constraints trends 
    521508                           CALL dia_obs_init    ! Initialize observational data 
    522       IF( ln_diaobs    )   CALL dia_obs( nit000 - 1 )   ! Observation operator for restart 
     509                           CALL dia_tmb_init    ! TMB outputs 
     510                           CALL dia_25h_init    ! 25h mean  outputs 
     511      IF( ln_diaobs    )   CALL dia_obs( nit000-1 )   ! Observation operator for restart 
    523512 
    524513      !                                      ! Assimilation increments 
    525514      IF( lk_asminc    )   CALL asm_inc_init    ! Initialize assimilation increments 
    526       IF(lwp) WRITE(numout,*) 'Euler time step switch is ', neuler 
    527                            CALL dia_tmb_init    ! TMB outputs 
    528                            CALL dia_25h_init    ! 25h mean  outputs 
     515      ! 
     516      IF(lwp) WRITE(numout,cform_aaa)           ! Flag AAAAAAA 
     517      ! 
     518      IF( ln_timing    )   CALL timing_stop( 'nemo_init') 
    529519      ! 
    530520   END SUBROUTINE nemo_init 
     
    543533         WRITE(numout,*) 
    544534         WRITE(numout,*) 'nemo_ctl: Control prints' 
    545          WRITE(numout,*) '~~~~~~~ ' 
     535         WRITE(numout,*) '~~~~~~~~' 
    546536         WRITE(numout,*) '   Namelist namctl' 
    547537         WRITE(numout,*) '      run control (for debugging)     ln_ctl     = ', ln_ctl 
     
    563553      njctle    = nn_jctle 
    564554      isplt     = nn_isplt 
    565       jsplt     = nn_jsplt       
     555      jsplt     = nn_jsplt 
    566556 
    567557      IF(lwp) THEN                  ! control print 
    568558         WRITE(numout,*) 
    569          WRITE(numout,*) 'namcfg : configuration initialization through namelist read' 
    570          WRITE(numout,*) '~~~~~~ ' 
    571559         WRITE(numout,*) '   Namelist namcfg' 
    572560         WRITE(numout,*) '      read domain configuration file                ln_read_cfg      = ', ln_read_cfg 
    573561         WRITE(numout,*) '         filename to be read                           cn_domcfg     = ', TRIM(cn_domcfg) 
    574          WRITE(numout,*) '      write configuration definition file           ln_write_cfg     = ', ln_write_cfg 
     562         WRITE(numout,*) '         keep closed seas in the domain (if exist)     ln_closea     = ', ln_closea 
     563         WRITE(numout,*) '      create a configuration definition file        ln_write_cfg     = ', ln_write_cfg 
    575564         WRITE(numout,*) '         filename to be written                        cn_domcfg_out = ', TRIM(cn_domcfg_out) 
    576565         WRITE(numout,*) '      use file attribute if exists as i/p j-start   ln_use_jattr     = ', ln_use_jattr 
    577566      ENDIF 
     567      IF( .NOT.ln_read_cfg )   ln_closea = .false.   ! dealing possible only with a domcfg file 
     568      ! 
    578569      !                             ! Parameter control 
    579570      ! 
     
    615606      ENDIF 
    616607      ! 
    617       IF( 1_wp /= SIGN(1._wp,-0._wp)  )   CALL ctl_stop( 'nemo_ctl: The intrinsec SIGN function follows ',  & 
    618          &                                               'f2003 standard. '                              ,  & 
    619          &                                               'Compile with key_nosignedzero enabled' ) 
     608      IF( 1._wp /= SIGN(1._wp,-0._wp)  )   CALL ctl_stop( 'nemo_ctl: The intrinsec SIGN function follows f2003 standard.',  & 
     609         &                                                'Compile with key_nosignedzero enabled' ) 
    620610      ! 
    621611#if defined key_agrif 
     
    664654      !! ** Method  : 
    665655      !!---------------------------------------------------------------------- 
    666       USE diawri    , ONLY: dia_wri_alloc 
    667       USE dom_oce   , ONLY: dom_oce_alloc 
    668       USE trc_oce   , ONLY: trc_oce_alloc 
     656      USE diawri    , ONLY : dia_wri_alloc 
     657      USE dom_oce   , ONLY : dom_oce_alloc 
     658      USE trc_oce   , ONLY : trc_oce_alloc 
     659      USE bdy_oce   , ONLY : bdy_oce_alloc 
    669660#if defined key_diadct  
    670       USE diadct    , ONLY: diadct_alloc  
     661      USE diadct    , ONLY : diadct_alloc  
    671662#endif  
    672       USE bdy_oce   , ONLY: bdy_oce_alloc 
    673663      ! 
    674664      INTEGER :: ierr 
    675665      !!---------------------------------------------------------------------- 
    676666      ! 
    677       ierr =        oce_alloc       ()          ! ocean  
    678       ierr = ierr + dia_wri_alloc   () 
    679       ierr = ierr + dom_oce_alloc   ()          ! ocean domain 
    680       ierr = ierr + zdf_oce_alloc   ()          ! ocean vertical physics 
    681       ! 
    682       ierr = ierr + trc_oce_alloc   ()          ! shared TRC / TRA arrays 
     667      ierr =        oce_alloc    ()    ! ocean  
     668      ierr = ierr + dia_wri_alloc() 
     669      ierr = ierr + dom_oce_alloc()    ! ocean domain 
     670      ierr = ierr + zdf_oce_alloc()    ! ocean vertical physics 
     671      ierr = ierr + trc_oce_alloc()    ! shared TRC / TRA arrays 
     672      ierr = ierr + bdy_oce_alloc()    ! bdy masks (incl. initialization) 
    683673      ! 
    684674#if defined key_diadct  
    685       ierr = ierr + diadct_alloc    ()          !  
     675      ierr = ierr + diadct_alloc ()    !  
    686676#endif  
    687       ierr = ierr + bdy_oce_alloc   ()          ! bdy masks (incl. initialization) 
    688677      ! 
    689678      IF( lk_mpp    )   CALL mpp_sum( ierr ) 
    690       IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'nemo_alloc : unable to allocate standard ocean arrays' ) 
     679      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'nemo_alloc: unable to allocate standard ocean arrays' ) 
    691680      ! 
    692681   END SUBROUTINE nemo_alloc 
     
    766755      knfax = 0 
    767756      ! 
    768       ! Find the factors of n. 
    769       IF( kn .NE. 1 ) THEN 
    770  
     757      IF( kn /= 1 ) THEN      ! Find the factors of n 
     758         ! 
    771759         ! nu holds the unfactorised part of the number. 
    772760         ! knfax holds the number of factors found. 
     
    781769            ifac = ilfax(jl) 
    782770            IF( ifac > inu )   CYCLE 
    783     
     771            ! 
    784772            ! Test whether the factor will divide. 
    785     
     773            ! 
    786774            IF( MOD(inu,ifac) == 0 ) THEN 
    787775               ! 
     
    807795#if defined key_mpp_mpi 
    808796 
    809    SUBROUTINE nemo_northcomms 
    810       !!---------------------------------------------------------------------- 
    811       !!                     ***  ROUTINE  nemo_northcomms  *** 
     797   SUBROUTINE nemo_nfdcom 
     798      !!---------------------------------------------------------------------- 
     799      !!                     ***  ROUTINE  nemo_nfdcom  *** 
    812800      !! ** Purpose :   Setup for north fold exchanges with explicit  
    813801      !!                point-to-point messaging 
     
    828816      nsndto     = 0 
    829817      ! 
    830       !if I am a process in the north 
    831       IF ( njmpp == njmppmax ) THEN 
    832           !sxM is the first point (in the global domain) needed to compute the 
    833           !north-fold for the current process 
    834           sxM = jpiglo - nimppt(narea) - nlcit(narea) + 1 
    835           !dxM is the last point (in the global domain) needed to compute the 
    836           !north-fold for the current process 
    837           dxM = jpiglo - nimppt(narea) + 2 
    838  
    839           !loop over the other north-fold processes to find the processes 
    840           !managing the points belonging to the sxT-dxT range 
    841    
    842           DO jn = 1, jpni 
    843                 !sxT is the first point (in the global domain) of the jn 
    844                 !process 
    845                 sxT = nfiimpp(jn, jpnj) 
    846                 !dxT is the last point (in the global domain) of the jn 
    847                 !process 
    848                 dxT = nfiimpp(jn, jpnj) + nfilcit(jn, jpnj) - 1 
    849                 IF ((sxM .gt. sxT) .AND. (sxM .lt. dxT)) THEN 
    850                    nsndto = nsndto + 1 
    851                    isendto(nsndto) = jn 
    852                 ELSEIF ((sxM .le. sxT) .AND. (dxM .ge. dxT)) THEN 
    853                    nsndto = nsndto + 1 
    854                    isendto(nsndto) = jn 
    855                 ELSEIF ((dxM .lt. dxT) .AND. (sxT .lt. dxM)) THEN 
    856                    nsndto = nsndto + 1 
    857                    isendto(nsndto) = jn 
    858                 ENDIF 
    859           END DO 
    860           nfsloop = 1 
    861           nfeloop = nlci 
    862           DO jn = 2,jpni-1 
    863            IF(nfipproc(jn,jpnj) .eq. (narea - 1)) THEN 
    864               IF (nfipproc(jn - 1 ,jpnj) .eq. -1) THEN 
    865                  nfsloop = nldi 
    866               ENDIF 
    867               IF (nfipproc(jn + 1,jpnj) .eq. -1) THEN 
    868                  nfeloop = nlei 
    869               ENDIF 
    870            ENDIF 
    871         END DO 
    872  
     818      IF ( njmpp == njmppmax ) THEN      ! if I am a process in the north 
     819         ! 
     820         !sxM is the first point (in the global domain) needed to compute the north-fold for the current process 
     821         sxM = jpiglo - nimppt(narea) - nlcit(narea) + 1 
     822         !dxM is the last point (in the global domain) needed to compute the north-fold for the current process 
     823         dxM = jpiglo - nimppt(narea) + 2 
     824         ! 
     825         ! loop over the other north-fold processes to find the processes 
     826         ! managing the points belonging to the sxT-dxT range 
     827         ! 
     828         DO jn = 1, jpni 
     829            ! 
     830            sxT = nfiimpp(jn, jpnj)                            ! sxT = 1st  point (in the global domain) of the jn process 
     831            dxT = nfiimpp(jn, jpnj) + nfilcit(jn, jpnj) - 1    ! dxT = last point (in the global domain) of the jn process 
     832            ! 
     833            IF    ( sxT < sxM  .AND.  sxM < dxT ) THEN 
     834               nsndto          = nsndto + 1 
     835               isendto(nsndto) = jn 
     836            ELSEIF( sxM <= sxT  .AND.  dxM >= dxT ) THEN 
     837               nsndto          = nsndto + 1 
     838               isendto(nsndto) = jn 
     839            ELSEIF( dxM <  dxT  .AND.  sxT <  dxM ) THEN 
     840               nsndto          = nsndto + 1 
     841               isendto(nsndto) = jn 
     842            ENDIF 
     843            ! 
     844         END DO 
     845         nfsloop = 1 
     846         nfeloop = nlci 
     847         DO jn = 2,jpni-1 
     848            IF( nfipproc(jn,jpnj) == (narea - 1) ) THEN 
     849               IF( nfipproc(jn-1,jpnj) == -1 )   nfsloop = nldi 
     850               IF( nfipproc(jn+1,jpnj) == -1 )   nfeloop = nlei 
     851            ENDIF 
     852         END DO 
     853         ! 
    873854      ENDIF 
    874855      l_north_nogather = .TRUE. 
    875    END SUBROUTINE nemo_northcomms 
     856      ! 
     857   END SUBROUTINE nemo_nfdcom 
    876858 
    877859#else 
    878    SUBROUTINE nemo_northcomms      ! Dummy routine 
    879       WRITE(*,*) 'nemo_northcomms: You should not have seen this print! error?' 
    880    END SUBROUTINE nemo_northcomms 
     860   SUBROUTINE nemo_nfdcom      ! Dummy routine 
     861      WRITE(*,*) 'nemo_nfdcom: You should not have seen this print! error?' 
     862   END SUBROUTINE nemo_nfdcom 
    881863#endif 
    882864 
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/step_oce.F90

    r9023 r9213  
    8282   USE diaharm 
    8383   USE diacfl 
     84   USE diaobs          ! Observation operator 
    8485   USE flo_oce         ! floats variables 
    8586   USE floats          ! floats computation               (flo_stp routine) 
     
    9394   USE restart         ! ocean restart                    (rst_wri routine) 
    9495   USE prtctl          ! Print control                    (prt_ctl routine) 
    95  
    96    USE diaobs          ! Observation operator 
    9796 
    9897   USE in_out_manager  ! I/O manager 
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/SAS_SRC/nemogcm.F90

    r9200 r9213  
    1818   !!   nemo_partition: calculate MPP domain decomposition 
    1919   !!   factorise     : calculate the factors of the no. of MPI processes 
     20   !!   nemo_nfdcom   : Setup for north fold exchanges with explicit point-to-point messaging 
    2021   !!---------------------------------------------------------------------- 
    2122   USE step_oce       ! module used in the ocean time stepping module 
     
    3738   USE lib_mpp        ! distributed memory computing 
    3839   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
    39    USE lbcnfd   , ONLY: isendto, nsndto, nfsloop, nfeloop ! Setup of north fold exchanges 
     40   USE lbcnfd  , ONLY : isendto, nsndto, nfsloop, nfeloop  ! Setup of north fold exchanges 
    4041   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    4142#if defined key_iomput 
     
    5253 
    5354   !!---------------------------------------------------------------------- 
    54    !! NEMO/OPA 4.0 , NEMO Consortium (2016) 
     55   !! NEMO/OPA 4.0 , NEMO Consortium (2018) 
    5556   !! $Id$ 
    5657   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    7677      ! 
    7778#if defined key_agrif 
    78      CALL Agrif_Init_Grids()      ! AGRIF: set the meshes 
    79 #endif 
    80       ! 
     79      CALL Agrif_Init_Grids()      ! AGRIF: set the meshes 
     80#endif 
    8181      !                            !-----------------------! 
    8282      CALL nemo_init               !==  Initialisations  ==! 
     
    102102      !                            !-----------------------! 
    103103      istp = nit000 
     104      ! 
    104105#if defined key_agrif 
     106      !                                               !==  AGRIF time-stepping  ==! 
    105107      CALL Agrif_Regrid() 
    106 #endif 
    107  
    108       DO WHILE ( istp <= nitend .AND. nstop == 0 ) 
    109 #if defined key_agrif 
    110          CALL stp                         ! AGRIF: time stepping 
     108      ! 
     109      DO WHILE( istp <= nitend .AND. nstop == 0 ) 
     110         CALL stp 
     111         istp = istp + 1 
     112      END DO 
     113      ! 
     114      IF( .NOT. Agrif_Root() ) THEN 
     115         CALL Agrif_ParentGrid_To_ChildGrid() 
     116         IF( ln_timing )   CALL timing_finalize 
     117         CALL Agrif_ChildGrid_To_ParentGrid() 
     118      ENDIF 
     119      ! 
    111120#else 
    112          IF ( .NOT. ln_diurnal_only ) THEN 
    113             CALL stp( istp )                 ! standard time stepping 
    114          ELSE 
    115             CALL stp_diurnal( istp )        ! time step only the diurnal SST 
    116          ENDIF  
    117 #endif 
    118          istp = istp + 1 
    119          IF( lk_mpp )   CALL mpp_max( nstop ) 
     121      ! 
     122      IF( .NOT.ln_diurnal_only ) THEN                 !==  Standard time-stepping  ==! 
     123         ! 
     124         DO WHILE( istp <= nitend .AND. nstop == 0 ) 
     125            CALL stp        ( istp )  
     126            istp = istp + 1 
    120127         END DO 
     128         ! 
     129      ELSE                                            !==  diurnal SST time-steeping only  ==! 
     130         ! 
     131         DO WHILE( istp <= nitend .AND. nstop == 0 ) 
     132            CALL stp_diurnal( istp )   ! time step only the diurnal SST  
     133            istp = istp + 1 
     134         END DO 
     135         ! 
     136      ENDIF 
     137      ! 
     138#endif 
    121139      ! 
    122140      IF( ln_icebergs )   CALL icb_end( nitend ) 
     
    129147      IF( nstop /= 0 .AND. lwp ) THEN        ! error print 
    130148         WRITE(numout,cform_err) 
    131          WRITE(numout,*) nstop, ' error have been found' 
    132       ENDIF 
    133       ! 
    134 #if defined key_agrif 
    135          CALL Agrif_ParentGrid_To_ChildGrid() 
    136          IF( ln_timing )   CALL timing_finalize 
    137          CALL Agrif_ChildGrid_To_ParentGrid() 
    138 #endif 
     149         WRITE(numout,*) '   ==>>>   nemo_gcm: a total of ', nstop, ' errors have been found' 
     150         WRITE(numout,*) 
     151      ENDIF 
     152      ! 
    139153      IF( ln_timing )   CALL timing_finalize 
    140154      ! 
     
    142156      ! 
    143157#if defined key_iomput 
    144       CALL xios_finalize                     ! end mpp communications with xios 
    145       IF( lk_oasis )   CALL cpl_finalize     ! end coupling and mpp communications with OASIS 
     158                                    CALL xios_finalize  ! end mpp communications with xios 
     159      IF( lk_oasis     )            CALL cpl_finalize   ! end coupling and mpp communications with OASIS 
    146160#else 
    147       IF( lk_oasis ) THEN  
    148          CALL cpl_finalize                   ! end coupling and mpp communications with OASIS 
    149       ELSE 
    150          IF( lk_mpp )   CALL mppstop         ! end mpp communications 
     161      IF    ( lk_oasis ) THEN   ;   CALL cpl_finalize   ! end coupling and mpp communications with OASIS 
     162      ELSEIF( lk_mpp   ) THEN   ;   CALL mppstop        ! end mpp communications 
    151163      ENDIF 
    152164#endif 
     
    161173      !! ** Purpose :   initialization of the NEMO GCM 
    162174      !!---------------------------------------------------------------------- 
    163       INTEGER  ::   ji            ! dummy loop indices 
    164       INTEGER  ::   ilocal_comm   ! local integer 
    165       INTEGER  ::   ios, inum     !   -      - 
    166       INTEGER  ::   iiarea, ijarea     ! local integers 
    167       INTEGER  ::   iirest, ijrest     ! local integers 
     175      INTEGER  ::   ji                 ! dummy loop indices 
     176      INTEGER  ::   ios, ilocal_comm   ! local integers 
     177      INTEGER  ::   iiarea, ijarea     !   -       - 
     178      INTEGER  ::   iirest, ijrest     !   -       - 
    168179      CHARACTER(len=120), DIMENSION(30) ::   cltxt, cltxt2, clnam 
    169180      CHARACTER(len=80)                 ::   clname 
    170       ! 
     181      !! 
    171182      NAMELIST/namctl/ ln_ctl   , nn_print, nn_ictls, nn_ictle,   & 
    172183         &             nn_isplt , nn_jsplt, nn_jctls, nn_jctle,   & 
    173184         &             ln_timing, ln_diacfl 
    174       NAMELIST/namcfg/ ln_read_cfg, cn_domcfg, ln_write_cfg, cn_domcfg_out, ln_use_jattr, ln_closea 
     185      NAMELIST/namcfg/ ln_read_cfg, cn_domcfg, ln_closea, ln_write_cfg, cn_domcfg_out, ln_use_jattr 
    175186      !!---------------------------------------------------------------------- 
    176187      ! 
     
    193204   ENDIF 
    194205      ! 
    195       REWIND( numnam_ref )              ! Namelist namctl in reference namelist : Control prints 
     206      REWIND( numnam_ref )              ! Namelist namctl in reference namelist 
    196207      READ  ( numnam_ref, namctl, IOSTAT = ios, ERR = 901 ) 
    197208901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namctl in reference namelist', .TRUE. ) 
     
    200211902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namctl in configuration namelist', .TRUE. ) 
    201212      ! 
    202       REWIND( numnam_ref )              ! Namelist namcfg in reference namelist : Control prints 
     213      REWIND( numnam_ref )              ! Namelist namcfg in reference namelist 
    203214      READ  ( numnam_ref, namcfg, IOSTAT = ios, ERR = 903 ) 
    204215903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namcfg in reference namelist', .TRUE. ) 
    205       REWIND( numnam_cfg )              ! Namelist namcfg in confguration namelist : Control prints & Benchmark 
     216      REWIND( numnam_cfg )              ! Namelist namcfg in confguration namelist 
    206217      READ  ( numnam_cfg, namcfg, IOSTAT = ios, ERR = 904 ) 
    207218904   IF( ios >  0 )   CALL ctl_nam ( ios , 'namcfg in configuration namelist', .TRUE. )    
     
    325336      IF(lwp) THEN                            ! open listing units 
    326337         ! 
    327          IF( lk_oasis ) THEN 
    328             CALL ctl_opn( numout,   'sas.output', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE., narea ) 
    329          ELSE 
    330             CALL ctl_opn( numout, 'ocean.output', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE., narea ) 
     338         IF( lk_oasis ) THEN   ;   CALL ctl_opn( numout,   'sas.output', 'REPLACE','FORMATTED','SEQUENTIAL', -1, 6, .FALSE., narea ) 
     339         ELSE                  ;   CALL ctl_opn( numout, 'ocean.output', 'REPLACE','FORMATTED','SEQUENTIAL', -1, 6, .FALSE., narea ) 
    331340         ENDIF 
    332341         ! 
     
    335344         WRITE(numout,*) '                       NEMO team' 
    336345         WRITE(numout,*) '            Ocean General Circulation Model' 
    337          WRITE(numout,*) '                  version 3.7  (2016) ' 
     346         WRITE(numout,*) '                NEMO version 4.0  (2017) ' 
    338347         WRITE(numout,*) '             StandAlone Surface version (SAS) ' 
    339348         WRITE(numout,*) 
     
    354363      ! Now we know the dimensions of the grid and numout has been set: we can allocate arrays 
    355364      CALL nemo_alloc() 
     365 
    356366      !                             !-------------------------------! 
    357367      !                             !  NEMO general initialization  ! 
     
    361371 
    362372      !                                      ! Domain decomposition 
    363       CALL mpp_init 
    364       ! 
    365       IF( ln_timing    )   CALL timing_init 
    366       ! 
    367       !                                     ! General initialization 
    368                            CALL phy_cst    ! Physical constants 
    369                            CALL eos_init   ! Equation of state 
    370                            CALL dom_init   ! Domain 
    371  
    372      IF( ln_nnogather )    CALL nemo_northcomms   ! Initialise the northfold neighbour lists (must be done after the masks are defined) 
    373  
    374       IF( ln_ctl      )    CALL prt_ctl_init   ! Print control 
    375                            CALL day_init   ! model calendar (using both namelist and restart infos) 
    376       IF( ln_rstart )       CALL rst_read_open 
    377  
    378                            CALL sbc_init   ! Forcings : surface module  
     373      CALL mpp_init                             ! MPP 
     374      IF( ln_nnogather )   CALL nemo_nfdcom     ! northfold neighbour lists 
     375      ! 
     376      ! 
     377      !                                      ! General initialization 
     378      IF( ln_timing    )   CALL timing_init     ! timing 
     379      IF( ln_timing    )   CALL timing_start( 'nemo_init') 
     380 
     381                           CALL phy_cst         ! Physical constants 
     382                           CALL eos_init        ! Equation of seawater 
     383                           CALL dom_init        ! Domain 
     384      IF( ln_ctl      )    CALL prt_ctl_init    ! Print control 
     385       
     386                           CALL day_init        ! model calendar (using both namelist and restart infos) 
     387      IF( ln_rstart )      CALL rst_read_open 
     388 
     389      !                                      ! external forcing  
     390                           CALL sbc_init        ! Forcings : surface module  
    379391 
    380392      ! ==> clem: open boundaries init. is mandatory for sea-ice because ice BDY is not decoupled from   
     
    384396      ! ==> 
    385397                           CALL icb_init( rdt, nit000)   ! initialise icebergs instance 
    386        
    387       IF(lwp) WRITE(numout,*) 'Euler time step switch is ', neuler 
     398      ! 
     399      IF(lwp) WRITE(numout,cform_aaa)           ! Flag AAAAAAA 
     400      ! 
     401      IF( ln_timing    )   CALL timing_stop( 'nemo_init') 
    388402      ! 
    389403   END SUBROUTINE nemo_init 
     
    394408      !!                     ***  ROUTINE nemo_ctl  *** 
    395409      !! 
    396       !! ** Purpose :   control print setting  
     410      !! ** Purpose :   control print setting 
    397411      !! 
    398412      !! ** Method  : - print namctl information and check some consistencies 
     
    402416         WRITE(numout,*) 
    403417         WRITE(numout,*) 'nemo_ctl: Control prints' 
    404          WRITE(numout,*) '~~~~~~~ ' 
     418         WRITE(numout,*) '~~~~~~~~' 
    405419         WRITE(numout,*) '   Namelist namctl' 
    406420         WRITE(numout,*) '      run control (for debugging)     ln_ctl     = ', ln_ctl 
     
    426440      IF(lwp) THEN                  ! control print 
    427441         WRITE(numout,*) 
    428          WRITE(numout,*) 'namcfg  : configuration initialization through namelist read' 
    429          WRITE(numout,*) '~~~~~~~ ' 
    430442         WRITE(numout,*) '   Namelist namcfg' 
    431          WRITE(numout,*) '      read domain configuration files               ln_read_cfg      = ', ln_read_cfg 
     443         WRITE(numout,*) '      read domain configuration file                ln_read_cfg      = ', ln_read_cfg 
    432444         WRITE(numout,*) '         filename to be read                           cn_domcfg     = ', TRIM(cn_domcfg) 
    433          WRITE(numout,*) '      write  configuration definition files         ln_write_cfg     = ', ln_write_cfg 
     445         WRITE(numout,*) '         keep closed seas in the domain (if exist)     ln_closea     = ', ln_closea 
     446         WRITE(numout,*) '      create a configuration definition file        ln_write_cfg     = ', ln_write_cfg 
    434447         WRITE(numout,*) '         filename to be written                        cn_domcfg_out = ', TRIM(cn_domcfg_out) 
    435448         WRITE(numout,*) '      use file attribute if exists as i/p j-start   ln_use_jattr     = ', ln_use_jattr 
    436449      ENDIF 
     450      IF( .NOT.ln_read_cfg )   ln_closea = .false.   ! dealing possible only with a domcfg file 
     451      ! 
    437452      !                             ! Parameter control 
    438453      ! 
     
    474489      ENDIF 
    475490      ! 
    476       IF( 1_wp /= SIGN(1._wp,-0._wp)  )   CALL ctl_stop( 'nemo_ctl: The intrinsec SIGN function follows ',  & 
    477          &                                               'f2003 standard. '                              ,  & 
    478          &                                               'Compile with key_nosignedzero enabled' ) 
     491      IF( 1._wp /= SIGN(1._wp,-0._wp)  )   CALL ctl_stop( 'nemo_ctl: The intrinsec SIGN function follows f2003 standard.',  & 
     492         &                                                'Compile with key_nosignedzero enabled' ) 
     493      ! 
     494#if defined key_agrif 
     495      IF( ln_timing )   CALL ctl_stop( 'AGRIF not implemented with ln_timing = true') 
     496#endif 
    479497      ! 
    480498   END SUBROUTINE nemo_ctl 
     
    492510      CALL iom_close                                 ! close all input/output files managed by iom_* 
    493511      ! 
    494       IF( numstp          /= -1 )   CLOSE( numstp      )   ! time-step file       
     512      IF( numstp          /= -1 )   CLOSE( numstp          )   ! time-step file       
    495513      IF( numnam_ref      /= -1 )   CLOSE( numnam_ref      )   ! oce reference namelist 
    496514      IF( numnam_cfg      /= -1 )   CLOSE( numnam_cfg      )   ! oce configuration namelist 
     
    499517      IF( numnam_ice_cfg  /= -1 )   CLOSE( numnam_ice_cfg  )   ! ice configuration namelist 
    500518      IF( lwm.AND.numoni  /= -1 )   CLOSE( numoni          )   ! ice output namelist 
    501       IF( numevo_ice      /= -1 )   CLOSE( numevo_ice  )   ! ice variables (temp. evolution) 
    502       IF( numout          /=  6 )   CLOSE( numout      )   ! standard model output file 
     519      IF( numevo_ice      /= -1 )   CLOSE( numevo_ice      )   ! ice variables (temp. evolution) 
     520      IF( numout          /=  6 )   CLOSE( numout          )   ! standard model output file 
    503521      ! 
    504522      numout = 6                                     ! redefine numout in case it is used after this point... 
     
    515533      !! ** Method  : 
    516534      !!---------------------------------------------------------------------- 
    517       USE diawri    , ONLY: dia_wri_alloc 
    518       USE dom_oce   , ONLY: dom_oce_alloc 
    519       USE bdy_oce   , ONLY: ln_bdy, bdy_oce_alloc 
     535      USE diawri    , ONLY : dia_wri_alloc 
     536      USE dom_oce   , ONLY : dom_oce_alloc 
     537      USE bdy_oce   , ONLY : ln_bdy, bdy_oce_alloc 
    520538      USE oce       ! mandatory for sea-ice because needed for bdy arrays 
    521539      ! 
     
    523541      !!---------------------------------------------------------------------- 
    524542      ! 
    525       ierr =        dia_wri_alloc   () 
    526       ierr = ierr + dom_oce_alloc   ()          ! ocean domain 
    527       ierr = ierr + oce_alloc       ()          ! (tsn...) needed for agrif and/or ESIM and bdy 
    528       ierr = ierr + bdy_oce_alloc   ()          ! bdy masks (incl. initialization) 
     543      ierr =        dia_wri_alloc() 
     544      ierr = ierr + dom_oce_alloc()          ! ocean domain 
     545      ierr = ierr + oce_alloc    ()          ! (tsn...) needed for agrif and/or ESIM and bdy 
     546      ierr = ierr + bdy_oce_alloc()          ! bdy masks (incl. initialization) 
    529547      ! 
    530548      IF( lk_mpp    )   CALL mpp_sum( ierr ) 
    531       IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'nemo_alloc : unable to allocate standard ocean arrays' ) 
     549      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'nemo_alloc: unable to allocate standard ocean arrays' ) 
    532550      ! 
    533551   END SUBROUTINE nemo_alloc 
     
    538556      !!                 ***  ROUTINE nemo_partition  *** 
    539557      !! 
    540       !! ** Purpose :    
     558      !! ** Purpose : 
    541559      !! 
    542560      !! ** Method  : 
     
    607625      knfax = 0 
    608626      ! 
    609       ! Find the factors of n. 
    610       IF( kn .NE. 1 ) THEN 
    611  
     627      IF( kn /= 1 ) THEN      ! Find the factors of n 
     628         ! 
    612629         ! nu holds the unfactorised part of the number. 
    613630         ! knfax holds the number of factors found. 
     
    622639            ifac = ilfax(jl) 
    623640            IF( ifac > inu )   CYCLE 
    624     
     641            ! 
    625642            ! Test whether the factor will divide. 
    626     
     643            ! 
    627644            IF( MOD(inu,ifac) == 0 ) THEN 
    628645               ! 
     
    648665#if defined key_mpp_mpi 
    649666 
    650    SUBROUTINE nemo_northcomms 
    651       !!---------------------------------------------------------------------- 
    652       !!                     ***  ROUTINE  nemo_northcomms  *** 
     667   SUBROUTINE nemo_nfdcom 
     668      !!---------------------------------------------------------------------- 
     669      !!                     ***  ROUTINE  nemo_nfdcom  *** 
    653670      !! ** Purpose :   Setup for north fold exchanges with explicit  
    654671      !!                point-to-point messaging 
     
    669686      nsndto     = 0 
    670687      ! 
    671       !if I am a process in the north 
    672       IF ( njmpp == njmppmax ) THEN 
    673           !sxM is the first point (in the global domain) needed to compute the 
    674           !north-fold for the current process 
    675           sxM = jpiglo - nimppt(narea) - nlcit(narea) + 1 
    676           !dxM is the last point (in the global domain) needed to compute the 
    677           !north-fold for the current process 
    678           dxM = jpiglo - nimppt(narea) + 2 
    679  
    680           !loop over the other north-fold processes to find the processes 
    681           !managing the points belonging to the sxT-dxT range 
    682    
    683           DO jn = 1, jpni 
    684                 !sxT is the first point (in the global domain) of the jn 
    685                 !process 
    686                 sxT = nfiimpp(jn, jpnj) 
    687                 !dxT is the last point (in the global domain) of the jn 
    688                 !process 
    689                 dxT = nfiimpp(jn, jpnj) + nfilcit(jn, jpnj) - 1 
    690                 IF ((sxM .gt. sxT) .AND. (sxM .lt. dxT)) THEN 
    691                    nsndto = nsndto + 1 
    692                      isendto(nsndto) = jn 
    693                 ELSEIF ((sxM .le. sxT) .AND. (dxM .ge. dxT)) THEN 
    694                    nsndto = nsndto + 1 
    695                    isendto(nsndto) = jn 
    696                 ELSEIF ((dxM .lt. dxT) .AND. (sxT .lt. dxM)) THEN 
    697                    nsndto = nsndto + 1 
    698                    isendto(nsndto) = jn 
    699                 ENDIF 
    700           END DO 
    701           nfsloop = 1 
    702           nfeloop = nlci 
    703           DO jn = 2,jpni-1 
    704            IF(nfipproc(jn,jpnj) .eq. (narea - 1)) THEN 
    705               IF (nfipproc(jn - 1 ,jpnj) .eq. -1) THEN 
    706                  nfsloop = nldi 
    707               ENDIF 
    708               IF (nfipproc(jn + 1,jpnj) .eq. -1) THEN 
    709                  nfeloop = nlei 
    710               ENDIF 
    711            ENDIF 
    712         END DO 
    713  
     688      IF ( njmpp == njmppmax ) THEN      ! if I am a process in the north 
     689         ! 
     690         !sxM is the first point (in the global domain) needed to compute the north-fold for the current process 
     691         sxM = jpiglo - nimppt(narea) - nlcit(narea) + 1 
     692         !dxM is the last point (in the global domain) needed to compute the north-fold for the current process 
     693         dxM = jpiglo - nimppt(narea) + 2 
     694         ! 
     695         ! loop over the other north-fold processes to find the processes 
     696         ! managing the points belonging to the sxT-dxT range 
     697         ! 
     698         DO jn = 1, jpni 
     699            ! 
     700            sxT = nfiimpp(jn, jpnj)                            ! sxT = 1st  point (in the global domain) of the jn process 
     701            dxT = nfiimpp(jn, jpnj) + nfilcit(jn, jpnj) - 1    ! dxT = last point (in the global domain) of the jn process 
     702            ! 
     703            IF    ( sxT < sxM  .AND.  sxM < dxT ) THEN 
     704               nsndto          = nsndto + 1 
     705               isendto(nsndto) = jn 
     706            ELSEIF( sxM <= sxT  .AND.  dxM >= dxT ) THEN 
     707               nsndto          = nsndto + 1 
     708               isendto(nsndto) = jn 
     709            ELSEIF( dxM <  dxT  .AND.  sxT <  dxM ) THEN 
     710               nsndto          = nsndto + 1 
     711               isendto(nsndto) = jn 
     712            ENDIF 
     713            ! 
     714         END DO 
     715         nfsloop = 1 
     716         nfeloop = nlci 
     717         DO jn = 2,jpni-1 
     718            IF( nfipproc(jn,jpnj) == (narea - 1) ) THEN 
     719               IF( nfipproc(jn-1,jpnj) == -1 )   nfsloop = nldi 
     720               IF( nfipproc(jn+1,jpnj) == -1 )   nfeloop = nlei 
     721            ENDIF 
     722         END DO 
     723         ! 
    714724      ENDIF 
    715725      l_north_nogather = .TRUE. 
    716    END SUBROUTINE nemo_northcomms 
     726      ! 
     727   END SUBROUTINE nemo_nfdcom 
    717728 
    718729#else 
    719    SUBROUTINE nemo_northcomms      ! Dummy routine 
    720       WRITE(*,*) 'nemo_northcomms: You should not have seen this print! error?' 
    721    END SUBROUTINE nemo_northcomms 
     730   SUBROUTINE nemo_nfdcom      ! Dummy routine 
     731      WRITE(*,*) 'nemo_nfdcom: You should not have seen this print! error?' 
     732   END SUBROUTINE nemo_nfdcom 
    722733#endif 
    723734 
    724735   !!====================================================================== 
    725736END MODULE nemogcm 
     737 
Note: See TracChangeset for help on using the changeset viewer.