[3] | 1 | MODULE lib_mpp |
---|
[13] | 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE lib_mpp *** |
---|
[1344] | 4 | !! Ocean numerics: massively parallel processing library |
---|
[13] | 5 | !!===================================================================== |
---|
[1344] | 6 | !! History : OPA ! 1994 (M. Guyon, J. Escobar, M. Imbard) Original code |
---|
| 7 | !! 7.0 ! 1997 (A.M. Treguier) SHMEM additions |
---|
| 8 | !! 8.0 ! 1998 (M. Imbard, J. Escobar, L. Colombet ) SHMEM and MPI |
---|
| 9 | !! ! 1998 (J.M. Molines) Open boundary conditions |
---|
[9019] | 10 | !! NEMO 1.0 ! 2003 (J.M. Molines, G. Madec) F90, free form |
---|
[1344] | 11 | !! ! 2003 (J.M. Molines) add mpp_ini_north(_3d,_2d) |
---|
| 12 | !! - ! 2004 (R. Bourdalle Badie) isend option in mpi |
---|
| 13 | !! ! 2004 (J.M. Molines) minloc, maxloc |
---|
| 14 | !! - ! 2005 (G. Madec, S. Masson) npolj=5,6 F-point & ice cases |
---|
| 15 | !! - ! 2005 (R. Redler) Replacement of MPI_COMM_WORLD except for MPI_Abort |
---|
| 16 | !! - ! 2005 (R. Benshila, G. Madec) add extra halo case |
---|
| 17 | !! - ! 2008 (R. Benshila) add mpp_ini_ice |
---|
| 18 | !! 3.2 ! 2009 (R. Benshila) SHMEM suppression, north fold in lbc_nfd |
---|
[3764] | 19 | !! 3.2 ! 2009 (O. Marti) add mpp_ini_znl |
---|
[2715] | 20 | !! 4.0 ! 2011 (G. Madec) move ctl_ routines from in_out_manager |
---|
[9019] | 21 | !! 3.5 ! 2012 (S.Mocavero, I. Epicoco) Add mpp_lnk_bdy_3d/2d routines to optimize the BDY comm. |
---|
| 22 | !! 3.5 ! 2013 (C. Ethe, G. Madec) message passing arrays as local variables |
---|
[6140] | 23 | !! 3.5 ! 2013 (S.Mocavero, I.Epicoco - CMCC) north fold optimizations |
---|
[9019] | 24 | !! 3.6 ! 2015 (O. Tintó and M. Castrillo - BSC) Added '_multiple' case for 2D lbc and max |
---|
| 25 | !! 4.0 ! 2017 (G. Madec) automatique allocation of array argument (use any 3rd dimension) |
---|
| 26 | !! - ! 2017 (G. Madec) create generic.h90 files to generate all lbc and north fold routines |
---|
[13] | 27 | !!---------------------------------------------------------------------- |
---|
[2715] | 28 | |
---|
| 29 | !!---------------------------------------------------------------------- |
---|
[6140] | 30 | !! ctl_stop : update momentum and tracer Kz from a tke scheme |
---|
| 31 | !! ctl_warn : initialization, namelist read, and parameters control |
---|
| 32 | !! ctl_opn : Open file and check if required file is available. |
---|
| 33 | !! ctl_nam : Prints informations when an error occurs while reading a namelist |
---|
| 34 | !! get_unit : give the index of an unused logical unit |
---|
[2715] | 35 | !!---------------------------------------------------------------------- |
---|
[3764] | 36 | #if defined key_mpp_mpi |
---|
[13] | 37 | !!---------------------------------------------------------------------- |
---|
[1344] | 38 | !! 'key_mpp_mpi' MPI massively parallel processing library |
---|
| 39 | !!---------------------------------------------------------------------- |
---|
[2715] | 40 | !! lib_mpp_alloc : allocate mpp arrays |
---|
| 41 | !! mynode : indentify the processor unit |
---|
| 42 | !! mpp_lnk : interface (defined in lbclnk) for message passing of 2d or 3d arrays (mpp_lnk_2d, mpp_lnk_3d) |
---|
[4990] | 43 | !! mpp_lnk_icb : interface for message passing of 2d arrays with extra halo for icebergs (mpp_lnk_2d_icb) |
---|
[6140] | 44 | !! mpprecv : |
---|
[9019] | 45 | !! mppsend : |
---|
[2715] | 46 | !! mppscatter : |
---|
| 47 | !! mppgather : |
---|
| 48 | !! mpp_min : generic interface for mppmin_int , mppmin_a_int , mppmin_real, mppmin_a_real |
---|
| 49 | !! mpp_max : generic interface for mppmax_int , mppmax_a_int , mppmax_real, mppmax_a_real |
---|
| 50 | !! mpp_sum : generic interface for mppsum_int , mppsum_a_int , mppsum_real, mppsum_a_real |
---|
| 51 | !! mpp_minloc : |
---|
| 52 | !! mpp_maxloc : |
---|
| 53 | !! mppsync : |
---|
| 54 | !! mppstop : |
---|
[1344] | 55 | !! mpp_ini_north : initialisation of north fold |
---|
[9019] | 56 | !! mpp_lbc_north_icb : alternative to mpp_nfd for extra outer halo with icebergs |
---|
[13] | 57 | !!---------------------------------------------------------------------- |
---|
[3764] | 58 | USE dom_oce ! ocean space and time domain |
---|
[2715] | 59 | USE lbcnfd ! north fold treatment |
---|
| 60 | USE in_out_manager ! I/O manager |
---|
[3] | 61 | |
---|
[13] | 62 | IMPLICIT NONE |
---|
[415] | 63 | PRIVATE |
---|
[9019] | 64 | |
---|
| 65 | INTERFACE mpp_nfd |
---|
[10425] | 66 | MODULE PROCEDURE mpp_nfd_2d , mpp_nfd_3d , mpp_nfd_4d |
---|
[9019] | 67 | MODULE PROCEDURE mpp_nfd_2d_ptr, mpp_nfd_3d_ptr, mpp_nfd_4d_ptr |
---|
| 68 | END INTERFACE |
---|
| 69 | |
---|
| 70 | ! Interface associated to the mpp_lnk_... routines is defined in lbclnk |
---|
[10425] | 71 | PUBLIC mpp_lnk_2d , mpp_lnk_3d , mpp_lnk_4d |
---|
[9019] | 72 | PUBLIC mpp_lnk_2d_ptr, mpp_lnk_3d_ptr, mpp_lnk_4d_ptr |
---|
| 73 | ! |
---|
| 74 | !!gm this should be useless |
---|
| 75 | PUBLIC mpp_nfd_2d , mpp_nfd_3d , mpp_nfd_4d |
---|
| 76 | PUBLIC mpp_nfd_2d_ptr, mpp_nfd_3d_ptr, mpp_nfd_4d_ptr |
---|
| 77 | !!gm end |
---|
| 78 | ! |
---|
[4147] | 79 | PUBLIC ctl_stop, ctl_warn, get_unit, ctl_opn, ctl_nam |
---|
[1344] | 80 | PUBLIC mynode, mppstop, mppsync, mpp_comm_free |
---|
[9019] | 81 | PUBLIC mpp_ini_north |
---|
| 82 | PUBLIC mpp_lnk_2d_icb |
---|
| 83 | PUBLIC mpp_lbc_north_icb |
---|
[1344] | 84 | PUBLIC mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc |
---|
[10425] | 85 | PUBLIC mpp_delay_max, mpp_delay_sum, mpp_delay_rcv |
---|
[3294] | 86 | PUBLIC mppscatter, mppgather |
---|
[10425] | 87 | PUBLIC mpp_ini_znl |
---|
[3764] | 88 | PUBLIC mppsend, mpprecv ! needed by TAM and ICB routines |
---|
[9890] | 89 | PUBLIC mpp_lnk_bdy_2d, mpp_lnk_bdy_3d, mpp_lnk_bdy_4d |
---|
[5429] | 90 | |
---|
[13] | 91 | !! * Interfaces |
---|
| 92 | !! define generic interface for these routine as they are called sometimes |
---|
[1344] | 93 | !! with scalar arguments instead of array arguments, which causes problems |
---|
| 94 | !! for the compilation on AIX system as well as NEC and SGI. Ok on COMPACQ |
---|
[13] | 95 | INTERFACE mpp_min |
---|
| 96 | MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real |
---|
| 97 | END INTERFACE |
---|
| 98 | INTERFACE mpp_max |
---|
[681] | 99 | MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real |
---|
[13] | 100 | END INTERFACE |
---|
| 101 | INTERFACE mpp_sum |
---|
[6140] | 102 | MODULE PROCEDURE mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real, & |
---|
[9019] | 103 | & mppsum_realdd, mppsum_a_realdd |
---|
[13] | 104 | END INTERFACE |
---|
[1344] | 105 | INTERFACE mpp_minloc |
---|
| 106 | MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d |
---|
| 107 | END INTERFACE |
---|
| 108 | INTERFACE mpp_maxloc |
---|
| 109 | MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d |
---|
| 110 | END INTERFACE |
---|
[6490] | 111 | |
---|
[51] | 112 | !! ========================= !! |
---|
| 113 | !! MPI variable definition !! |
---|
| 114 | !! ========================= !! |
---|
[1629] | 115 | !$AGRIF_DO_NOT_TREAT |
---|
[2004] | 116 | INCLUDE 'mpif.h' |
---|
[1629] | 117 | !$AGRIF_END_DO_NOT_TREAT |
---|
[3764] | 118 | |
---|
[1344] | 119 | LOGICAL, PUBLIC, PARAMETER :: lk_mpp = .TRUE. !: mpp flag |
---|
[3] | 120 | |
---|
[1344] | 121 | INTEGER, PARAMETER :: nprocmax = 2**10 ! maximun dimension (required to be a power of 2) |
---|
[3764] | 122 | |
---|
[10425] | 123 | INTEGER, PUBLIC :: mppsize ! number of process |
---|
| 124 | INTEGER, PUBLIC :: mpprank ! process number [ 0 - size-1 ] |
---|
[2363] | 125 | !$AGRIF_DO_NOT_TREAT |
---|
[9570] | 126 | INTEGER, PUBLIC :: mpi_comm_oce ! opa local communicator |
---|
[2363] | 127 | !$AGRIF_END_DO_NOT_TREAT |
---|
[3] | 128 | |
---|
[2480] | 129 | INTEGER :: MPI_SUMDD |
---|
[1976] | 130 | |
---|
[1345] | 131 | ! variables used for zonal integration |
---|
| 132 | INTEGER, PUBLIC :: ncomm_znl !: communicator made by the processors on the same zonal average |
---|
[9019] | 133 | LOGICAL, PUBLIC :: l_znl_root !: True on the 'left'most processor on the same row |
---|
| 134 | INTEGER :: ngrp_znl ! group ID for the znl processors |
---|
| 135 | INTEGER :: ndim_rank_znl ! number of processors on the same zonal average |
---|
[2715] | 136 | INTEGER, DIMENSION(:), ALLOCATABLE, SAVE :: nrank_znl ! dimension ndim_rank_znl, number of the procs into the same znl domain |
---|
[3] | 137 | |
---|
[3764] | 138 | ! North fold condition in mpp_mpi with jpni > 1 (PUBLIC for TAM) |
---|
[9019] | 139 | INTEGER, PUBLIC :: ngrp_world !: group ID for the world processors |
---|
| 140 | INTEGER, PUBLIC :: ngrp_opa !: group ID for the opa processors |
---|
| 141 | INTEGER, PUBLIC :: ngrp_north !: group ID for the northern processors (to be fold) |
---|
| 142 | INTEGER, PUBLIC :: ncomm_north !: communicator made by the processors belonging to ngrp_north |
---|
| 143 | INTEGER, PUBLIC :: ndim_rank_north !: number of 'sea' processor in the northern line (can be /= jpni !) |
---|
| 144 | INTEGER, PUBLIC :: njmppmax !: value of njmpp for the processors of the northern line |
---|
| 145 | INTEGER, PUBLIC :: north_root !: number (in the comm_opa) of proc 0 in the northern comm |
---|
| 146 | INTEGER, PUBLIC, DIMENSION(:), ALLOCATABLE, SAVE :: nrank_north !: dimension ndim_rank_north |
---|
[3764] | 147 | |
---|
[1344] | 148 | ! Type of send : standard, buffered, immediate |
---|
[9019] | 149 | CHARACTER(len=1), PUBLIC :: cn_mpi_send !: type od mpi send/recieve (S=standard, B=bsend, I=isend) |
---|
| 150 | LOGICAL , PUBLIC :: l_isend = .FALSE. !: isend use indicator (T if cn_mpi_send='I') |
---|
| 151 | INTEGER , PUBLIC :: nn_buffer !: size of the buffer in case of mpi_bsend |
---|
[3764] | 152 | |
---|
[10425] | 153 | ! Communications summary report |
---|
| 154 | CHARACTER(len=128), DIMENSION(:), ALLOCATABLE :: crname_lbc !: names of lbc_lnk calling routines |
---|
[10437] | 155 | CHARACTER(len=128), DIMENSION(:), ALLOCATABLE :: crname_glb !: names of global comm calling routines |
---|
| 156 | CHARACTER(len=128), DIMENSION(:), ALLOCATABLE :: crname_dlg !: names of delayed global comm calling routines |
---|
[10425] | 157 | INTEGER, PUBLIC :: ncom_stp = 0 !: copy of time step # istp |
---|
| 158 | INTEGER, PUBLIC :: ncom_fsbc = 1 !: copy of sbc time step # nn_fsbc |
---|
| 159 | INTEGER, PUBLIC :: ncom_dttrc = 1 !: copy of top time step # nn_dttrc |
---|
| 160 | INTEGER, PUBLIC :: ncom_freq !: frequency of comm diagnostic |
---|
| 161 | INTEGER, PUBLIC , DIMENSION(:,:), ALLOCATABLE :: ncomm_sequence !: size of communicated arrays (halos) |
---|
[10782] | 162 | INTEGER, PARAMETER, PUBLIC :: ncom_rec_max = 5000 !: max number of communication record |
---|
[10425] | 163 | INTEGER, PUBLIC :: n_sequence_lbc = 0 !: # of communicated arraysvia lbc |
---|
| 164 | INTEGER, PUBLIC :: n_sequence_glb = 0 !: # of global communications |
---|
[10437] | 165 | INTEGER, PUBLIC :: n_sequence_dlg = 0 !: # of delayed global communications |
---|
[10425] | 166 | INTEGER, PUBLIC :: numcom = -1 !: logical unit for communicaton report |
---|
| 167 | LOGICAL, PUBLIC :: l_full_nf_update = .TRUE. !: logical for a full (2lines) update of bc at North fold report |
---|
| 168 | INTEGER, PARAMETER, PUBLIC :: nbdelay = 2 !: number of delayed operations |
---|
| 169 | !: name (used as id) of allreduce-delayed operations |
---|
[10521] | 170 | ! Warning: we must use the same character length in an array constructor (at least for gcc compiler) |
---|
| 171 | CHARACTER(len=32), DIMENSION(nbdelay), PUBLIC :: c_delaylist = (/ 'cflice', 'fwb ' /) |
---|
[10425] | 172 | !: component name where the allreduce-delayed operation is performed |
---|
| 173 | CHARACTER(len=3), DIMENSION(nbdelay), PUBLIC :: c_delaycpnt = (/ 'ICE' , 'OCE' /) |
---|
| 174 | TYPE, PUBLIC :: DELAYARR |
---|
| 175 | REAL( wp), POINTER, DIMENSION(:) :: z1d => NULL() |
---|
| 176 | COMPLEX(wp), POINTER, DIMENSION(:) :: y1d => NULL() |
---|
| 177 | END TYPE DELAYARR |
---|
[10815] | 178 | TYPE( DELAYARR ), DIMENSION(nbdelay), PUBLIC, SAVE :: todelay !: must have SAVE for default initialization of DELAYARR |
---|
| 179 | INTEGER, DIMENSION(nbdelay), PUBLIC :: ndelayid = -1 !: mpi request id of the delayed operations |
---|
[10425] | 180 | |
---|
| 181 | ! timing summary report |
---|
| 182 | REAL(wp), DIMENSION(2), PUBLIC :: waiting_time = 0._wp |
---|
| 183 | REAL(wp) , PUBLIC :: compute_time = 0._wp, elapsed_time = 0._wp |
---|
| 184 | |
---|
[9019] | 185 | REAL(wp), DIMENSION(:), ALLOCATABLE, SAVE :: tampon ! buffer in case of bsend |
---|
[3] | 186 | |
---|
[9019] | 187 | LOGICAL, PUBLIC :: ln_nnogather !: namelist control of northfold comms |
---|
| 188 | LOGICAL, PUBLIC :: l_north_nogather = .FALSE. !: internal control of northfold comms |
---|
| 189 | |
---|
[51] | 190 | !!---------------------------------------------------------------------- |
---|
[9598] | 191 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
[888] | 192 | !! $Id$ |
---|
[10068] | 193 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
[1344] | 194 | !!---------------------------------------------------------------------- |
---|
[3] | 195 | CONTAINS |
---|
| 196 | |
---|
[9019] | 197 | FUNCTION mynode( ldtxt, ldname, kumnam_ref, kumnam_cfg, kumond, kstop, localComm ) |
---|
[2715] | 198 | !!---------------------------------------------------------------------- |
---|
[51] | 199 | !! *** routine mynode *** |
---|
[3764] | 200 | !! |
---|
[51] | 201 | !! ** Purpose : Find processor unit |
---|
| 202 | !!---------------------------------------------------------------------- |
---|
[6140] | 203 | CHARACTER(len=*),DIMENSION(:), INTENT( out) :: ldtxt ! |
---|
| 204 | CHARACTER(len=*) , INTENT(in ) :: ldname ! |
---|
| 205 | INTEGER , INTENT(in ) :: kumnam_ref ! logical unit for reference namelist |
---|
| 206 | INTEGER , INTENT(in ) :: kumnam_cfg ! logical unit for configuration namelist |
---|
| 207 | INTEGER , INTENT(inout) :: kumond ! logical unit for namelist output |
---|
| 208 | INTEGER , INTENT(inout) :: kstop ! stop indicator |
---|
| 209 | INTEGER , OPTIONAL , INTENT(in ) :: localComm ! |
---|
[2715] | 210 | ! |
---|
[4147] | 211 | INTEGER :: mynode, ierr, code, ji, ii, ios |
---|
[532] | 212 | LOGICAL :: mpi_was_called |
---|
[2715] | 213 | ! |
---|
[10428] | 214 | NAMELIST/nammpp/ cn_mpi_send, nn_buffer, jpni, jpnj, ln_nnogather |
---|
[51] | 215 | !!---------------------------------------------------------------------- |
---|
[1344] | 216 | ! |
---|
[2481] | 217 | ii = 1 |
---|
[6140] | 218 | WRITE(ldtxt(ii),*) ; ii = ii + 1 |
---|
| 219 | WRITE(ldtxt(ii),*) 'mynode : mpi initialisation' ; ii = ii + 1 |
---|
| 220 | WRITE(ldtxt(ii),*) '~~~~~~ ' ; ii = ii + 1 |
---|
[1344] | 221 | ! |
---|
[4147] | 222 | REWIND( kumnam_ref ) ! Namelist nammpp in reference namelist: mpi variables |
---|
| 223 | READ ( kumnam_ref, nammpp, IOSTAT = ios, ERR = 901) |
---|
[9168] | 224 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nammpp in reference namelist', lwp ) |
---|
[9019] | 225 | ! |
---|
[4147] | 226 | REWIND( kumnam_cfg ) ! Namelist nammpp in configuration namelist: mpi variables |
---|
| 227 | READ ( kumnam_cfg, nammpp, IOSTAT = ios, ERR = 902 ) |
---|
[9168] | 228 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'nammpp in configuration namelist', lwp ) |
---|
[9019] | 229 | ! |
---|
[1344] | 230 | ! ! control print |
---|
[6140] | 231 | WRITE(ldtxt(ii),*) ' Namelist nammpp' ; ii = ii + 1 |
---|
| 232 | WRITE(ldtxt(ii),*) ' mpi send type cn_mpi_send = ', cn_mpi_send ; ii = ii + 1 |
---|
| 233 | WRITE(ldtxt(ii),*) ' size exported buffer nn_buffer = ', nn_buffer,' bytes'; ii = ii + 1 |
---|
[9019] | 234 | ! |
---|
| 235 | IF( jpni < 1 .OR. jpnj < 1 ) THEN |
---|
[10428] | 236 | WRITE(ldtxt(ii),*) ' jpni and jpnj will be calculated automatically' ; ii = ii + 1 |
---|
[2715] | 237 | ELSE |
---|
[6140] | 238 | WRITE(ldtxt(ii),*) ' processor grid extent in i jpni = ',jpni ; ii = ii + 1 |
---|
| 239 | WRITE(ldtxt(ii),*) ' processor grid extent in j jpnj = ',jpnj ; ii = ii + 1 |
---|
[9019] | 240 | ENDIF |
---|
[2715] | 241 | |
---|
[3294] | 242 | WRITE(ldtxt(ii),*) ' avoid use of mpi_allgather at the north fold ln_nnogather = ', ln_nnogather ; ii = ii + 1 |
---|
| 243 | |
---|
[2480] | 244 | CALL mpi_initialized ( mpi_was_called, code ) |
---|
| 245 | IF( code /= MPI_SUCCESS ) THEN |
---|
[3764] | 246 | DO ji = 1, SIZE(ldtxt) |
---|
[2481] | 247 | IF( TRIM(ldtxt(ji)) /= '' ) WRITE(*,*) ldtxt(ji) ! control print of mynode |
---|
[3764] | 248 | END DO |
---|
[2480] | 249 | WRITE(*, cform_err) |
---|
| 250 | WRITE(*, *) 'lib_mpp: Error in routine mpi_initialized' |
---|
| 251 | CALL mpi_abort( mpi_comm_world, code, ierr ) |
---|
| 252 | ENDIF |
---|
[415] | 253 | |
---|
[2480] | 254 | IF( mpi_was_called ) THEN |
---|
| 255 | ! |
---|
| 256 | SELECT CASE ( cn_mpi_send ) |
---|
| 257 | CASE ( 'S' ) ! Standard mpi send (blocking) |
---|
[6140] | 258 | WRITE(ldtxt(ii),*) ' Standard blocking mpi send (send)' ; ii = ii + 1 |
---|
[2480] | 259 | CASE ( 'B' ) ! Buffer mpi send (blocking) |
---|
[6140] | 260 | WRITE(ldtxt(ii),*) ' Buffer blocking mpi send (bsend)' ; ii = ii + 1 |
---|
[9570] | 261 | IF( Agrif_Root() ) CALL mpi_init_oce( ldtxt, ii, ierr ) |
---|
[2480] | 262 | CASE ( 'I' ) ! Immediate mpi send (non-blocking send) |
---|
[6140] | 263 | WRITE(ldtxt(ii),*) ' Immediate non-blocking send (isend)' ; ii = ii + 1 |
---|
[2480] | 264 | l_isend = .TRUE. |
---|
| 265 | CASE DEFAULT |
---|
[6140] | 266 | WRITE(ldtxt(ii),cform_err) ; ii = ii + 1 |
---|
| 267 | WRITE(ldtxt(ii),*) ' bad value for cn_mpi_send = ', cn_mpi_send ; ii = ii + 1 |
---|
[2715] | 268 | kstop = kstop + 1 |
---|
[2480] | 269 | END SELECT |
---|
[9019] | 270 | ! |
---|
| 271 | ELSEIF ( PRESENT(localComm) .AND. .NOT. mpi_was_called ) THEN |
---|
[10425] | 272 | WRITE(ldtxt(ii),cform_err) ; ii = ii + 1 |
---|
[6140] | 273 | WRITE(ldtxt(ii),*) ' lib_mpp: You cannot provide a local communicator ' ; ii = ii + 1 |
---|
| 274 | WRITE(ldtxt(ii),*) ' without calling MPI_Init before ! ' ; ii = ii + 1 |
---|
[2715] | 275 | kstop = kstop + 1 |
---|
[532] | 276 | ELSE |
---|
[1601] | 277 | SELECT CASE ( cn_mpi_send ) |
---|
[524] | 278 | CASE ( 'S' ) ! Standard mpi send (blocking) |
---|
[6140] | 279 | WRITE(ldtxt(ii),*) ' Standard blocking mpi send (send)' ; ii = ii + 1 |
---|
[2480] | 280 | CALL mpi_init( ierr ) |
---|
[524] | 281 | CASE ( 'B' ) ! Buffer mpi send (blocking) |
---|
[6140] | 282 | WRITE(ldtxt(ii),*) ' Buffer blocking mpi send (bsend)' ; ii = ii + 1 |
---|
[9570] | 283 | IF( Agrif_Root() ) CALL mpi_init_oce( ldtxt, ii, ierr ) |
---|
[524] | 284 | CASE ( 'I' ) ! Immediate mpi send (non-blocking send) |
---|
[6140] | 285 | WRITE(ldtxt(ii),*) ' Immediate non-blocking send (isend)' ; ii = ii + 1 |
---|
[524] | 286 | l_isend = .TRUE. |
---|
[2480] | 287 | CALL mpi_init( ierr ) |
---|
[524] | 288 | CASE DEFAULT |
---|
[6140] | 289 | WRITE(ldtxt(ii),cform_err) ; ii = ii + 1 |
---|
| 290 | WRITE(ldtxt(ii),*) ' bad value for cn_mpi_send = ', cn_mpi_send ; ii = ii + 1 |
---|
[2715] | 291 | kstop = kstop + 1 |
---|
[524] | 292 | END SELECT |
---|
[2480] | 293 | ! |
---|
[415] | 294 | ENDIF |
---|
[570] | 295 | |
---|
[3764] | 296 | IF( PRESENT(localComm) ) THEN |
---|
[2480] | 297 | IF( Agrif_Root() ) THEN |
---|
[9570] | 298 | mpi_comm_oce = localComm |
---|
[2480] | 299 | ENDIF |
---|
| 300 | ELSE |
---|
[9570] | 301 | CALL mpi_comm_dup( mpi_comm_world, mpi_comm_oce, code) |
---|
[2480] | 302 | IF( code /= MPI_SUCCESS ) THEN |
---|
[3764] | 303 | DO ji = 1, SIZE(ldtxt) |
---|
[2481] | 304 | IF( TRIM(ldtxt(ji)) /= '' ) WRITE(*,*) ldtxt(ji) ! control print of mynode |
---|
| 305 | END DO |
---|
[2480] | 306 | WRITE(*, cform_err) |
---|
| 307 | WRITE(*, *) ' lib_mpp: Error in routine mpi_comm_dup' |
---|
| 308 | CALL mpi_abort( mpi_comm_world, code, ierr ) |
---|
| 309 | ENDIF |
---|
[3764] | 310 | ENDIF |
---|
[2480] | 311 | |
---|
[5656] | 312 | #if defined key_agrif |
---|
[9019] | 313 | IF( Agrif_Root() ) THEN |
---|
[9570] | 314 | CALL Agrif_MPI_Init(mpi_comm_oce) |
---|
[5656] | 315 | ELSE |
---|
[9570] | 316 | CALL Agrif_MPI_set_grid_comm(mpi_comm_oce) |
---|
[5656] | 317 | ENDIF |
---|
| 318 | #endif |
---|
| 319 | |
---|
[9570] | 320 | CALL mpi_comm_rank( mpi_comm_oce, mpprank, ierr ) |
---|
| 321 | CALL mpi_comm_size( mpi_comm_oce, mppsize, ierr ) |
---|
[629] | 322 | mynode = mpprank |
---|
[4624] | 323 | |
---|
| 324 | IF( mynode == 0 ) THEN |
---|
[5407] | 325 | CALL ctl_opn( kumond, TRIM(ldname), 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. , 1 ) |
---|
| 326 | WRITE(kumond, nammpp) |
---|
[4624] | 327 | ENDIF |
---|
[3764] | 328 | ! |
---|
[1976] | 329 | CALL MPI_OP_CREATE(DDPDD_MPI, .TRUE., MPI_SUMDD, ierr) |
---|
| 330 | ! |
---|
[51] | 331 | END FUNCTION mynode |
---|
[3] | 332 | |
---|
[9019] | 333 | !!---------------------------------------------------------------------- |
---|
| 334 | !! *** routine mpp_lnk_(2,3,4)d *** |
---|
| 335 | !! |
---|
| 336 | !! * Argument : dummy argument use in mpp_lnk_... routines |
---|
| 337 | !! ptab : array or pointer of arrays on which the boundary condition is applied |
---|
| 338 | !! cd_nat : nature of array grid-points |
---|
| 339 | !! psgn : sign used across the north fold boundary |
---|
| 340 | !! kfld : optional, number of pt3d arrays |
---|
| 341 | !! cd_mpp : optional, fill the overlap area only |
---|
| 342 | !! pval : optional, background value (used at closed boundaries) |
---|
| 343 | !!---------------------------------------------------------------------- |
---|
| 344 | ! |
---|
| 345 | ! !== 2D array and array of 2D pointer ==! |
---|
| 346 | ! |
---|
| 347 | # define DIM_2d |
---|
| 348 | # define ROUTINE_LNK mpp_lnk_2d |
---|
| 349 | # include "mpp_lnk_generic.h90" |
---|
| 350 | # undef ROUTINE_LNK |
---|
| 351 | # define MULTI |
---|
| 352 | # define ROUTINE_LNK mpp_lnk_2d_ptr |
---|
| 353 | # include "mpp_lnk_generic.h90" |
---|
| 354 | # undef ROUTINE_LNK |
---|
| 355 | # undef MULTI |
---|
| 356 | # undef DIM_2d |
---|
| 357 | ! |
---|
| 358 | ! !== 3D array and array of 3D pointer ==! |
---|
| 359 | ! |
---|
| 360 | # define DIM_3d |
---|
| 361 | # define ROUTINE_LNK mpp_lnk_3d |
---|
| 362 | # include "mpp_lnk_generic.h90" |
---|
| 363 | # undef ROUTINE_LNK |
---|
| 364 | # define MULTI |
---|
| 365 | # define ROUTINE_LNK mpp_lnk_3d_ptr |
---|
| 366 | # include "mpp_lnk_generic.h90" |
---|
| 367 | # undef ROUTINE_LNK |
---|
| 368 | # undef MULTI |
---|
| 369 | # undef DIM_3d |
---|
| 370 | ! |
---|
| 371 | ! !== 4D array and array of 4D pointer ==! |
---|
| 372 | ! |
---|
| 373 | # define DIM_4d |
---|
| 374 | # define ROUTINE_LNK mpp_lnk_4d |
---|
| 375 | # include "mpp_lnk_generic.h90" |
---|
| 376 | # undef ROUTINE_LNK |
---|
| 377 | # define MULTI |
---|
| 378 | # define ROUTINE_LNK mpp_lnk_4d_ptr |
---|
| 379 | # include "mpp_lnk_generic.h90" |
---|
| 380 | # undef ROUTINE_LNK |
---|
| 381 | # undef MULTI |
---|
| 382 | # undef DIM_4d |
---|
[6140] | 383 | |
---|
[9019] | 384 | !!---------------------------------------------------------------------- |
---|
| 385 | !! *** routine mpp_nfd_(2,3,4)d *** |
---|
| 386 | !! |
---|
| 387 | !! * Argument : dummy argument use in mpp_nfd_... routines |
---|
| 388 | !! ptab : array or pointer of arrays on which the boundary condition is applied |
---|
| 389 | !! cd_nat : nature of array grid-points |
---|
| 390 | !! psgn : sign used across the north fold boundary |
---|
| 391 | !! kfld : optional, number of pt3d arrays |
---|
| 392 | !! cd_mpp : optional, fill the overlap area only |
---|
| 393 | !! pval : optional, background value (used at closed boundaries) |
---|
| 394 | !!---------------------------------------------------------------------- |
---|
| 395 | ! |
---|
| 396 | ! !== 2D array and array of 2D pointer ==! |
---|
| 397 | ! |
---|
| 398 | # define DIM_2d |
---|
| 399 | # define ROUTINE_NFD mpp_nfd_2d |
---|
| 400 | # include "mpp_nfd_generic.h90" |
---|
| 401 | # undef ROUTINE_NFD |
---|
| 402 | # define MULTI |
---|
| 403 | # define ROUTINE_NFD mpp_nfd_2d_ptr |
---|
| 404 | # include "mpp_nfd_generic.h90" |
---|
| 405 | # undef ROUTINE_NFD |
---|
| 406 | # undef MULTI |
---|
| 407 | # undef DIM_2d |
---|
| 408 | ! |
---|
| 409 | ! !== 3D array and array of 3D pointer ==! |
---|
| 410 | ! |
---|
| 411 | # define DIM_3d |
---|
| 412 | # define ROUTINE_NFD mpp_nfd_3d |
---|
| 413 | # include "mpp_nfd_generic.h90" |
---|
| 414 | # undef ROUTINE_NFD |
---|
| 415 | # define MULTI |
---|
| 416 | # define ROUTINE_NFD mpp_nfd_3d_ptr |
---|
| 417 | # include "mpp_nfd_generic.h90" |
---|
| 418 | # undef ROUTINE_NFD |
---|
| 419 | # undef MULTI |
---|
| 420 | # undef DIM_3d |
---|
| 421 | ! |
---|
| 422 | ! !== 4D array and array of 4D pointer ==! |
---|
| 423 | ! |
---|
| 424 | # define DIM_4d |
---|
| 425 | # define ROUTINE_NFD mpp_nfd_4d |
---|
| 426 | # include "mpp_nfd_generic.h90" |
---|
| 427 | # undef ROUTINE_NFD |
---|
| 428 | # define MULTI |
---|
| 429 | # define ROUTINE_NFD mpp_nfd_4d_ptr |
---|
| 430 | # include "mpp_nfd_generic.h90" |
---|
| 431 | # undef ROUTINE_NFD |
---|
| 432 | # undef MULTI |
---|
| 433 | # undef DIM_4d |
---|
[3] | 434 | |
---|
[1344] | 435 | |
---|
[9019] | 436 | !!---------------------------------------------------------------------- |
---|
| 437 | !! *** routine mpp_lnk_bdy_(2,3,4)d *** |
---|
| 438 | !! |
---|
| 439 | !! * Argument : dummy argument use in mpp_lnk_... routines |
---|
| 440 | !! ptab : array or pointer of arrays on which the boundary condition is applied |
---|
| 441 | !! cd_nat : nature of array grid-points |
---|
| 442 | !! psgn : sign used across the north fold boundary |
---|
| 443 | !! kb_bdy : BDY boundary set |
---|
| 444 | !! kfld : optional, number of pt3d arrays |
---|
| 445 | !!---------------------------------------------------------------------- |
---|
| 446 | ! |
---|
| 447 | ! !== 2D array and array of 2D pointer ==! |
---|
| 448 | ! |
---|
| 449 | # define DIM_2d |
---|
| 450 | # define ROUTINE_BDY mpp_lnk_bdy_2d |
---|
| 451 | # include "mpp_bdy_generic.h90" |
---|
| 452 | # undef ROUTINE_BDY |
---|
| 453 | # undef DIM_2d |
---|
| 454 | ! |
---|
| 455 | ! !== 3D array and array of 3D pointer ==! |
---|
| 456 | ! |
---|
| 457 | # define DIM_3d |
---|
| 458 | # define ROUTINE_BDY mpp_lnk_bdy_3d |
---|
| 459 | # include "mpp_bdy_generic.h90" |
---|
| 460 | # undef ROUTINE_BDY |
---|
| 461 | # undef DIM_3d |
---|
| 462 | ! |
---|
| 463 | ! !== 4D array and array of 4D pointer ==! |
---|
| 464 | ! |
---|
[9890] | 465 | # define DIM_4d |
---|
| 466 | # define ROUTINE_BDY mpp_lnk_bdy_4d |
---|
| 467 | # include "mpp_bdy_generic.h90" |
---|
| 468 | # undef ROUTINE_BDY |
---|
| 469 | # undef DIM_4d |
---|
[3] | 470 | |
---|
[9019] | 471 | !!---------------------------------------------------------------------- |
---|
| 472 | !! |
---|
| 473 | !! load_array & mpp_lnk_2d_9 à generaliser a 3D et 4D |
---|
[5429] | 474 | |
---|
| 475 | |
---|
[9019] | 476 | !! mpp_lnk_sum_2d et 3D ====>>>>>> à virer du code !!!! |
---|
[5429] | 477 | |
---|
[9019] | 478 | |
---|
| 479 | !!---------------------------------------------------------------------- |
---|
[5429] | 480 | |
---|
| 481 | |
---|
[888] | 482 | |
---|
[1344] | 483 | SUBROUTINE mppsend( ktyp, pmess, kbytes, kdest, md_req ) |
---|
[51] | 484 | !!---------------------------------------------------------------------- |
---|
| 485 | !! *** routine mppsend *** |
---|
[3764] | 486 | !! |
---|
[51] | 487 | !! ** Purpose : Send messag passing array |
---|
| 488 | !! |
---|
| 489 | !!---------------------------------------------------------------------- |
---|
[1344] | 490 | REAL(wp), INTENT(inout) :: pmess(*) ! array of real |
---|
| 491 | INTEGER , INTENT(in ) :: kbytes ! size of the array pmess |
---|
| 492 | INTEGER , INTENT(in ) :: kdest ! receive process number |
---|
| 493 | INTEGER , INTENT(in ) :: ktyp ! tag of the message |
---|
| 494 | INTEGER , INTENT(in ) :: md_req ! argument for isend |
---|
| 495 | !! |
---|
| 496 | INTEGER :: iflag |
---|
[51] | 497 | !!---------------------------------------------------------------------- |
---|
[1344] | 498 | ! |
---|
[1601] | 499 | SELECT CASE ( cn_mpi_send ) |
---|
[300] | 500 | CASE ( 'S' ) ! Standard mpi send (blocking) |
---|
[9570] | 501 | CALL mpi_send ( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce , iflag ) |
---|
[300] | 502 | CASE ( 'B' ) ! Buffer mpi send (blocking) |
---|
[9570] | 503 | CALL mpi_bsend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce , iflag ) |
---|
[300] | 504 | CASE ( 'I' ) ! Immediate mpi send (non-blocking send) |
---|
[1344] | 505 | ! be carefull, one more argument here : the mpi request identifier.. |
---|
[9570] | 506 | CALL mpi_isend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce, md_req, iflag ) |
---|
[300] | 507 | END SELECT |
---|
[1344] | 508 | ! |
---|
[51] | 509 | END SUBROUTINE mppsend |
---|
[3] | 510 | |
---|
| 511 | |
---|
[3294] | 512 | SUBROUTINE mpprecv( ktyp, pmess, kbytes, ksource ) |
---|
[51] | 513 | !!---------------------------------------------------------------------- |
---|
| 514 | !! *** routine mpprecv *** |
---|
| 515 | !! |
---|
| 516 | !! ** Purpose : Receive messag passing array |
---|
| 517 | !! |
---|
| 518 | !!---------------------------------------------------------------------- |
---|
[1344] | 519 | REAL(wp), INTENT(inout) :: pmess(*) ! array of real |
---|
| 520 | INTEGER , INTENT(in ) :: kbytes ! suze of the array pmess |
---|
| 521 | INTEGER , INTENT(in ) :: ktyp ! Tag of the recevied message |
---|
[3764] | 522 | INTEGER, OPTIONAL, INTENT(in) :: ksource ! source process number |
---|
[1344] | 523 | !! |
---|
[51] | 524 | INTEGER :: istatus(mpi_status_size) |
---|
| 525 | INTEGER :: iflag |
---|
[3294] | 526 | INTEGER :: use_source |
---|
[1344] | 527 | !!---------------------------------------------------------------------- |
---|
| 528 | ! |
---|
[3764] | 529 | ! If a specific process number has been passed to the receive call, |
---|
[3294] | 530 | ! use that one. Default is to use mpi_any_source |
---|
[6140] | 531 | use_source = mpi_any_source |
---|
| 532 | IF( PRESENT(ksource) ) use_source = ksource |
---|
| 533 | ! |
---|
[9570] | 534 | CALL mpi_recv( pmess, kbytes, mpi_double_precision, use_source, ktyp, mpi_comm_oce, istatus, iflag ) |
---|
[1344] | 535 | ! |
---|
[51] | 536 | END SUBROUTINE mpprecv |
---|
[3] | 537 | |
---|
| 538 | |
---|
[51] | 539 | SUBROUTINE mppgather( ptab, kp, pio ) |
---|
| 540 | !!---------------------------------------------------------------------- |
---|
| 541 | !! *** routine mppgather *** |
---|
[3764] | 542 | !! |
---|
| 543 | !! ** Purpose : Transfert between a local subdomain array and a work |
---|
[51] | 544 | !! array which is distributed following the vertical level. |
---|
| 545 | !! |
---|
[1344] | 546 | !!---------------------------------------------------------------------- |
---|
[6140] | 547 | REAL(wp), DIMENSION(jpi,jpj) , INTENT(in ) :: ptab ! subdomain input array |
---|
| 548 | INTEGER , INTENT(in ) :: kp ! record length |
---|
[1344] | 549 | REAL(wp), DIMENSION(jpi,jpj,jpnij), INTENT( out) :: pio ! subdomain input array |
---|
[51] | 550 | !! |
---|
[1344] | 551 | INTEGER :: itaille, ierror ! temporary integer |
---|
[51] | 552 | !!--------------------------------------------------------------------- |
---|
[1344] | 553 | ! |
---|
| 554 | itaille = jpi * jpj |
---|
| 555 | CALL mpi_gather( ptab, itaille, mpi_double_precision, pio, itaille , & |
---|
[9570] | 556 | & mpi_double_precision, kp , mpi_comm_oce, ierror ) |
---|
[1344] | 557 | ! |
---|
[51] | 558 | END SUBROUTINE mppgather |
---|
[3] | 559 | |
---|
| 560 | |
---|
[51] | 561 | SUBROUTINE mppscatter( pio, kp, ptab ) |
---|
| 562 | !!---------------------------------------------------------------------- |
---|
| 563 | !! *** routine mppscatter *** |
---|
| 564 | !! |
---|
[3764] | 565 | !! ** Purpose : Transfert between awork array which is distributed |
---|
[51] | 566 | !! following the vertical level and the local subdomain array. |
---|
| 567 | !! |
---|
| 568 | !!---------------------------------------------------------------------- |
---|
[6140] | 569 | REAL(wp), DIMENSION(jpi,jpj,jpnij) :: pio ! output array |
---|
| 570 | INTEGER :: kp ! Tag (not used with MPI |
---|
| 571 | REAL(wp), DIMENSION(jpi,jpj) :: ptab ! subdomain array input |
---|
[1344] | 572 | !! |
---|
| 573 | INTEGER :: itaille, ierror ! temporary integer |
---|
[51] | 574 | !!--------------------------------------------------------------------- |
---|
[1344] | 575 | ! |
---|
[6140] | 576 | itaille = jpi * jpj |
---|
[1344] | 577 | ! |
---|
| 578 | CALL mpi_scatter( pio, itaille, mpi_double_precision, ptab, itaille , & |
---|
[9570] | 579 | & mpi_double_precision, kp , mpi_comm_oce, ierror ) |
---|
[1344] | 580 | ! |
---|
[51] | 581 | END SUBROUTINE mppscatter |
---|
[3] | 582 | |
---|
[10425] | 583 | |
---|
| 584 | SUBROUTINE mpp_delay_sum( cdname, cdelay, y_in, pout, ldlast, kcom ) |
---|
| 585 | !!---------------------------------------------------------------------- |
---|
| 586 | !! *** routine mpp_delay_sum *** |
---|
[1344] | 587 | !! |
---|
[10425] | 588 | !! ** Purpose : performed delayed mpp_sum, the result is received on next call |
---|
| 589 | !! |
---|
[1344] | 590 | !!---------------------------------------------------------------------- |
---|
[10425] | 591 | CHARACTER(len=*), INTENT(in ) :: cdname ! name of the calling subroutine |
---|
| 592 | CHARACTER(len=*), INTENT(in ) :: cdelay ! name (used as id) of the delayed operation |
---|
| 593 | COMPLEX(wp), INTENT(in ), DIMENSION(:) :: y_in |
---|
| 594 | REAL(wp), INTENT( out), DIMENSION(:) :: pout |
---|
| 595 | LOGICAL, INTENT(in ) :: ldlast ! true if this is the last time we call this routine |
---|
| 596 | INTEGER, INTENT(in ), OPTIONAL :: kcom |
---|
[1344] | 597 | !! |
---|
[10425] | 598 | INTEGER :: ji, isz |
---|
| 599 | INTEGER :: idvar |
---|
| 600 | INTEGER :: ierr, ilocalcomm |
---|
| 601 | COMPLEX(wp), ALLOCATABLE, DIMENSION(:) :: ytmp |
---|
[1344] | 602 | !!---------------------------------------------------------------------- |
---|
[9570] | 603 | ilocalcomm = mpi_comm_oce |
---|
[9019] | 604 | IF( PRESENT(kcom) ) ilocalcomm = kcom |
---|
[3] | 605 | |
---|
[10425] | 606 | isz = SIZE(y_in) |
---|
| 607 | |
---|
[10437] | 608 | IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ld_dlg = .TRUE. ) |
---|
[13] | 609 | |
---|
[10425] | 610 | idvar = -1 |
---|
| 611 | DO ji = 1, nbdelay |
---|
| 612 | IF( TRIM(cdelay) == TRIM(c_delaylist(ji)) ) idvar = ji |
---|
| 613 | END DO |
---|
| 614 | IF ( idvar == -1 ) CALL ctl_stop( 'STOP',' mpp_delay_sum : please add a new delayed exchange for '//TRIM(cdname) ) |
---|
| 615 | |
---|
| 616 | IF ( ndelayid(idvar) == 0 ) THEN ! first call with restart: %z1d defined in iom_delay_rst |
---|
| 617 | ! -------------------------- |
---|
| 618 | IF ( SIZE(todelay(idvar)%z1d) /= isz ) THEN ! Check dimension coherence |
---|
| 619 | IF(lwp) WRITE(numout,*) ' WARNING: the nb of delayed variables in restart file is not the model one' |
---|
| 620 | DEALLOCATE(todelay(idvar)%z1d) |
---|
| 621 | ndelayid(idvar) = -1 ! do as if we had no restart |
---|
| 622 | ELSE |
---|
| 623 | ALLOCATE(todelay(idvar)%y1d(isz)) |
---|
| 624 | todelay(idvar)%y1d(:) = CMPLX(todelay(idvar)%z1d(:), 0., wp) ! create %y1d, complex variable needed by mpi_sumdd |
---|
| 625 | END IF |
---|
| 626 | ENDIF |
---|
| 627 | |
---|
| 628 | IF( ndelayid(idvar) == -1 ) THEN ! first call without restart: define %y1d and %z1d from y_in with blocking allreduce |
---|
| 629 | ! -------------------------- |
---|
| 630 | ALLOCATE(todelay(idvar)%z1d(isz), todelay(idvar)%y1d(isz)) ! allocate also %z1d as used for the restart |
---|
| 631 | CALL mpi_allreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ierr ) ! get %y1d |
---|
| 632 | todelay(idvar)%z1d(:) = REAL(todelay(idvar)%y1d(:), wp) ! define %z1d from %y1d |
---|
| 633 | ENDIF |
---|
| 634 | |
---|
| 635 | IF( ndelayid(idvar) > 0 ) CALL mpp_delay_rcv( idvar ) ! make sure %z1d is received |
---|
| 636 | |
---|
| 637 | ! send back pout from todelay(idvar)%z1d defined at previous call |
---|
| 638 | pout(:) = todelay(idvar)%z1d(:) |
---|
| 639 | |
---|
| 640 | ! send y_in into todelay(idvar)%y1d with a non-blocking communication |
---|
[10526] | 641 | #if defined key_mpi2 |
---|
| 642 | IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.) |
---|
| 643 | CALL mpi_allreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ndelayid(idvar), ierr ) |
---|
| 644 | IF( ln_timing ) CALL tic_tac(.FALSE., ld_global = .TRUE.) |
---|
| 645 | #else |
---|
[10425] | 646 | CALL mpi_iallreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ndelayid(idvar), ierr ) |
---|
[10526] | 647 | #endif |
---|
[10425] | 648 | |
---|
| 649 | END SUBROUTINE mpp_delay_sum |
---|
| 650 | |
---|
[9019] | 651 | |
---|
[10425] | 652 | SUBROUTINE mpp_delay_max( cdname, cdelay, p_in, pout, ldlast, kcom ) |
---|
[9019] | 653 | !!---------------------------------------------------------------------- |
---|
[10425] | 654 | !! *** routine mpp_delay_max *** |
---|
[9019] | 655 | !! |
---|
[10425] | 656 | !! ** Purpose : performed delayed mpp_max, the result is received on next call |
---|
[9019] | 657 | !! |
---|
| 658 | !!---------------------------------------------------------------------- |
---|
[10425] | 659 | CHARACTER(len=*), INTENT(in ) :: cdname ! name of the calling subroutine |
---|
| 660 | CHARACTER(len=*), INTENT(in ) :: cdelay ! name (used as id) of the delayed operation |
---|
| 661 | REAL(wp), INTENT(in ), DIMENSION(:) :: p_in ! |
---|
| 662 | REAL(wp), INTENT( out), DIMENSION(:) :: pout ! |
---|
| 663 | LOGICAL, INTENT(in ) :: ldlast ! true if this is the last time we call this routine |
---|
| 664 | INTEGER, INTENT(in ), OPTIONAL :: kcom |
---|
[9019] | 665 | !! |
---|
[10425] | 666 | INTEGER :: ji, isz |
---|
| 667 | INTEGER :: idvar |
---|
| 668 | INTEGER :: ierr, ilocalcomm |
---|
[9019] | 669 | !!---------------------------------------------------------------------- |
---|
[9570] | 670 | ilocalcomm = mpi_comm_oce |
---|
[9019] | 671 | IF( PRESENT(kcom) ) ilocalcomm = kcom |
---|
[6140] | 672 | |
---|
[10425] | 673 | isz = SIZE(p_in) |
---|
[9019] | 674 | |
---|
[10437] | 675 | IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ld_dlg = .TRUE. ) |
---|
[13] | 676 | |
---|
[10425] | 677 | idvar = -1 |
---|
| 678 | DO ji = 1, nbdelay |
---|
| 679 | IF( TRIM(cdelay) == TRIM(c_delaylist(ji)) ) idvar = ji |
---|
| 680 | END DO |
---|
| 681 | IF ( idvar == -1 ) CALL ctl_stop( 'STOP',' mpp_delay_max : please add a new delayed exchange for '//TRIM(cdname) ) |
---|
[3] | 682 | |
---|
[10425] | 683 | IF ( ndelayid(idvar) == 0 ) THEN ! first call with restart: %z1d defined in iom_delay_rst |
---|
| 684 | ! -------------------------- |
---|
| 685 | IF ( SIZE(todelay(idvar)%z1d) /= isz ) THEN ! Check dimension coherence |
---|
| 686 | IF(lwp) WRITE(numout,*) ' WARNING: the nb of delayed variables in restart file is not the model one' |
---|
| 687 | DEALLOCATE(todelay(idvar)%z1d) |
---|
| 688 | ndelayid(idvar) = -1 ! do as if we had no restart |
---|
| 689 | END IF |
---|
| 690 | ENDIF |
---|
[13] | 691 | |
---|
[10425] | 692 | IF( ndelayid(idvar) == -1 ) THEN ! first call without restart: define %z1d from p_in with a blocking allreduce |
---|
| 693 | ! -------------------------- |
---|
| 694 | ALLOCATE(todelay(idvar)%z1d(isz)) |
---|
| 695 | CALL mpi_allreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ierr ) ! get %z1d |
---|
| 696 | ENDIF |
---|
[3] | 697 | |
---|
[10425] | 698 | IF( ndelayid(idvar) > 0 ) CALL mpp_delay_rcv( idvar ) ! make sure %z1d is received |
---|
[3] | 699 | |
---|
[10425] | 700 | ! send back pout from todelay(idvar)%z1d defined at previous call |
---|
| 701 | pout(:) = todelay(idvar)%z1d(:) |
---|
[13] | 702 | |
---|
[10425] | 703 | ! send p_in into todelay(idvar)%z1d with a non-blocking communication |
---|
[10526] | 704 | #if defined key_mpi2 |
---|
| 705 | IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.) |
---|
| 706 | CALL mpi_allreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ndelayid(idvar), ierr ) |
---|
| 707 | IF( ln_timing ) CALL tic_tac(.FALSE., ld_global = .TRUE.) |
---|
| 708 | #else |
---|
[10425] | 709 | CALL mpi_iallreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ndelayid(idvar), ierr ) |
---|
[10526] | 710 | #endif |
---|
[10425] | 711 | |
---|
| 712 | END SUBROUTINE mpp_delay_max |
---|
| 713 | |
---|
| 714 | |
---|
| 715 | SUBROUTINE mpp_delay_rcv( kid ) |
---|
| 716 | !!---------------------------------------------------------------------- |
---|
| 717 | !! *** routine mpp_delay_rcv *** |
---|
[1344] | 718 | !! |
---|
[10425] | 719 | !! ** Purpose : force barrier for delayed mpp (needed for restart) |
---|
[1344] | 720 | !! |
---|
[10425] | 721 | !!---------------------------------------------------------------------- |
---|
| 722 | INTEGER,INTENT(in ) :: kid |
---|
| 723 | INTEGER :: ierr |
---|
| 724 | !!---------------------------------------------------------------------- |
---|
| 725 | IF( ndelayid(kid) /= -2 ) THEN |
---|
[10526] | 726 | #if ! defined key_mpi2 |
---|
[10425] | 727 | IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.) |
---|
| 728 | CALL mpi_wait( ndelayid(kid), MPI_STATUS_IGNORE, ierr ) ! make sure todelay(kid) is received |
---|
| 729 | IF( ln_timing ) CALL tic_tac(.FALSE., ld_global = .TRUE.) |
---|
[10526] | 730 | #endif |
---|
[10425] | 731 | IF( ASSOCIATED(todelay(kid)%y1d) ) todelay(kid)%z1d(:) = REAL(todelay(kid)%y1d(:), wp) ! define %z1d from %y1d |
---|
| 732 | ndelayid(kid) = -2 ! add flag to know that mpi_wait was already called on kid |
---|
| 733 | ENDIF |
---|
| 734 | END SUBROUTINE mpp_delay_rcv |
---|
[3] | 735 | |
---|
[10425] | 736 | |
---|
| 737 | !!---------------------------------------------------------------------- |
---|
| 738 | !! *** mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real *** |
---|
| 739 | !! |
---|
| 740 | !!---------------------------------------------------------------------- |
---|
| 741 | !! |
---|
| 742 | # define OPERATION_MAX |
---|
| 743 | # define INTEGER_TYPE |
---|
| 744 | # define DIM_0d |
---|
| 745 | # define ROUTINE_ALLREDUCE mppmax_int |
---|
| 746 | # include "mpp_allreduce_generic.h90" |
---|
| 747 | # undef ROUTINE_ALLREDUCE |
---|
| 748 | # undef DIM_0d |
---|
| 749 | # define DIM_1d |
---|
| 750 | # define ROUTINE_ALLREDUCE mppmax_a_int |
---|
| 751 | # include "mpp_allreduce_generic.h90" |
---|
| 752 | # undef ROUTINE_ALLREDUCE |
---|
| 753 | # undef DIM_1d |
---|
| 754 | # undef INTEGER_TYPE |
---|
| 755 | ! |
---|
| 756 | # define REAL_TYPE |
---|
| 757 | # define DIM_0d |
---|
| 758 | # define ROUTINE_ALLREDUCE mppmax_real |
---|
| 759 | # include "mpp_allreduce_generic.h90" |
---|
| 760 | # undef ROUTINE_ALLREDUCE |
---|
| 761 | # undef DIM_0d |
---|
| 762 | # define DIM_1d |
---|
| 763 | # define ROUTINE_ALLREDUCE mppmax_a_real |
---|
| 764 | # include "mpp_allreduce_generic.h90" |
---|
| 765 | # undef ROUTINE_ALLREDUCE |
---|
| 766 | # undef DIM_1d |
---|
| 767 | # undef REAL_TYPE |
---|
| 768 | # undef OPERATION_MAX |
---|
| 769 | !!---------------------------------------------------------------------- |
---|
| 770 | !! *** mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real *** |
---|
| 771 | !! |
---|
| 772 | !!---------------------------------------------------------------------- |
---|
| 773 | !! |
---|
| 774 | # define OPERATION_MIN |
---|
| 775 | # define INTEGER_TYPE |
---|
| 776 | # define DIM_0d |
---|
| 777 | # define ROUTINE_ALLREDUCE mppmin_int |
---|
| 778 | # include "mpp_allreduce_generic.h90" |
---|
| 779 | # undef ROUTINE_ALLREDUCE |
---|
| 780 | # undef DIM_0d |
---|
| 781 | # define DIM_1d |
---|
| 782 | # define ROUTINE_ALLREDUCE mppmin_a_int |
---|
| 783 | # include "mpp_allreduce_generic.h90" |
---|
| 784 | # undef ROUTINE_ALLREDUCE |
---|
| 785 | # undef DIM_1d |
---|
| 786 | # undef INTEGER_TYPE |
---|
| 787 | ! |
---|
| 788 | # define REAL_TYPE |
---|
| 789 | # define DIM_0d |
---|
| 790 | # define ROUTINE_ALLREDUCE mppmin_real |
---|
| 791 | # include "mpp_allreduce_generic.h90" |
---|
| 792 | # undef ROUTINE_ALLREDUCE |
---|
| 793 | # undef DIM_0d |
---|
| 794 | # define DIM_1d |
---|
| 795 | # define ROUTINE_ALLREDUCE mppmin_a_real |
---|
| 796 | # include "mpp_allreduce_generic.h90" |
---|
| 797 | # undef ROUTINE_ALLREDUCE |
---|
| 798 | # undef DIM_1d |
---|
| 799 | # undef REAL_TYPE |
---|
| 800 | # undef OPERATION_MIN |
---|
[869] | 801 | |
---|
[10425] | 802 | !!---------------------------------------------------------------------- |
---|
| 803 | !! *** mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real *** |
---|
| 804 | !! |
---|
| 805 | !! Global sum of 1D array or a variable (integer, real or complex) |
---|
| 806 | !!---------------------------------------------------------------------- |
---|
| 807 | !! |
---|
| 808 | # define OPERATION_SUM |
---|
| 809 | # define INTEGER_TYPE |
---|
| 810 | # define DIM_0d |
---|
| 811 | # define ROUTINE_ALLREDUCE mppsum_int |
---|
| 812 | # include "mpp_allreduce_generic.h90" |
---|
| 813 | # undef ROUTINE_ALLREDUCE |
---|
| 814 | # undef DIM_0d |
---|
| 815 | # define DIM_1d |
---|
| 816 | # define ROUTINE_ALLREDUCE mppsum_a_int |
---|
| 817 | # include "mpp_allreduce_generic.h90" |
---|
| 818 | # undef ROUTINE_ALLREDUCE |
---|
| 819 | # undef DIM_1d |
---|
| 820 | # undef INTEGER_TYPE |
---|
| 821 | ! |
---|
| 822 | # define REAL_TYPE |
---|
| 823 | # define DIM_0d |
---|
| 824 | # define ROUTINE_ALLREDUCE mppsum_real |
---|
| 825 | # include "mpp_allreduce_generic.h90" |
---|
| 826 | # undef ROUTINE_ALLREDUCE |
---|
| 827 | # undef DIM_0d |
---|
| 828 | # define DIM_1d |
---|
| 829 | # define ROUTINE_ALLREDUCE mppsum_a_real |
---|
| 830 | # include "mpp_allreduce_generic.h90" |
---|
| 831 | # undef ROUTINE_ALLREDUCE |
---|
| 832 | # undef DIM_1d |
---|
| 833 | # undef REAL_TYPE |
---|
| 834 | # undef OPERATION_SUM |
---|
| 835 | |
---|
| 836 | # define OPERATION_SUM_DD |
---|
| 837 | # define COMPLEX_TYPE |
---|
| 838 | # define DIM_0d |
---|
| 839 | # define ROUTINE_ALLREDUCE mppsum_realdd |
---|
| 840 | # include "mpp_allreduce_generic.h90" |
---|
| 841 | # undef ROUTINE_ALLREDUCE |
---|
| 842 | # undef DIM_0d |
---|
| 843 | # define DIM_1d |
---|
| 844 | # define ROUTINE_ALLREDUCE mppsum_a_realdd |
---|
| 845 | # include "mpp_allreduce_generic.h90" |
---|
| 846 | # undef ROUTINE_ALLREDUCE |
---|
| 847 | # undef DIM_1d |
---|
| 848 | # undef COMPLEX_TYPE |
---|
| 849 | # undef OPERATION_SUM_DD |
---|
| 850 | |
---|
| 851 | !!---------------------------------------------------------------------- |
---|
| 852 | !! *** mpp_minloc2d, mpp_minloc3d, mpp_maxloc2d, mpp_maxloc3d |
---|
| 853 | !! |
---|
| 854 | !!---------------------------------------------------------------------- |
---|
| 855 | !! |
---|
| 856 | # define OPERATION_MINLOC |
---|
| 857 | # define DIM_2d |
---|
| 858 | # define ROUTINE_LOC mpp_minloc2d |
---|
| 859 | # include "mpp_loc_generic.h90" |
---|
| 860 | # undef ROUTINE_LOC |
---|
| 861 | # undef DIM_2d |
---|
| 862 | # define DIM_3d |
---|
| 863 | # define ROUTINE_LOC mpp_minloc3d |
---|
| 864 | # include "mpp_loc_generic.h90" |
---|
| 865 | # undef ROUTINE_LOC |
---|
| 866 | # undef DIM_3d |
---|
| 867 | # undef OPERATION_MINLOC |
---|
| 868 | |
---|
| 869 | # define OPERATION_MAXLOC |
---|
| 870 | # define DIM_2d |
---|
| 871 | # define ROUTINE_LOC mpp_maxloc2d |
---|
| 872 | # include "mpp_loc_generic.h90" |
---|
| 873 | # undef ROUTINE_LOC |
---|
| 874 | # undef DIM_2d |
---|
| 875 | # define DIM_3d |
---|
| 876 | # define ROUTINE_LOC mpp_maxloc3d |
---|
| 877 | # include "mpp_loc_generic.h90" |
---|
| 878 | # undef ROUTINE_LOC |
---|
| 879 | # undef DIM_3d |
---|
| 880 | # undef OPERATION_MAXLOC |
---|
| 881 | |
---|
[1344] | 882 | SUBROUTINE mppsync() |
---|
| 883 | !!---------------------------------------------------------------------- |
---|
| 884 | !! *** routine mppsync *** |
---|
[3764] | 885 | !! |
---|
[1344] | 886 | !! ** Purpose : Massively parallel processors, synchroneous |
---|
| 887 | !! |
---|
| 888 | !!----------------------------------------------------------------------- |
---|
| 889 | INTEGER :: ierror |
---|
| 890 | !!----------------------------------------------------------------------- |
---|
| 891 | ! |
---|
[9570] | 892 | CALL mpi_barrier( mpi_comm_oce, ierror ) |
---|
[1344] | 893 | ! |
---|
| 894 | END SUBROUTINE mppsync |
---|
[3] | 895 | |
---|
| 896 | |
---|
[10425] | 897 | SUBROUTINE mppstop( ldfinal, ld_force_abort ) |
---|
[1344] | 898 | !!---------------------------------------------------------------------- |
---|
| 899 | !! *** routine mppstop *** |
---|
[3764] | 900 | !! |
---|
[3294] | 901 | !! ** purpose : Stop massively parallel processors method |
---|
[1344] | 902 | !! |
---|
| 903 | !!---------------------------------------------------------------------- |
---|
[10425] | 904 | LOGICAL, OPTIONAL, INTENT(in) :: ldfinal ! source process number |
---|
| 905 | LOGICAL, OPTIONAL, INTENT(in) :: ld_force_abort ! source process number |
---|
| 906 | LOGICAL :: llfinal, ll_force_abort |
---|
[1344] | 907 | INTEGER :: info |
---|
| 908 | !!---------------------------------------------------------------------- |
---|
[10425] | 909 | llfinal = .FALSE. |
---|
| 910 | IF( PRESENT(ldfinal) ) llfinal = ldfinal |
---|
| 911 | ll_force_abort = .FALSE. |
---|
| 912 | IF( PRESENT(ld_force_abort) ) ll_force_abort = ld_force_abort |
---|
[1344] | 913 | ! |
---|
[10425] | 914 | IF(ll_force_abort) THEN |
---|
| 915 | CALL mpi_abort( MPI_COMM_WORLD ) |
---|
| 916 | ELSE |
---|
| 917 | CALL mppsync |
---|
| 918 | CALL mpi_finalize( info ) |
---|
| 919 | ENDIF |
---|
[10815] | 920 | IF( .NOT. llfinal ) STOP 123 |
---|
[1344] | 921 | ! |
---|
| 922 | END SUBROUTINE mppstop |
---|
[3] | 923 | |
---|
| 924 | |
---|
[1344] | 925 | SUBROUTINE mpp_comm_free( kcom ) |
---|
| 926 | !!---------------------------------------------------------------------- |
---|
| 927 | INTEGER, INTENT(in) :: kcom |
---|
| 928 | !! |
---|
| 929 | INTEGER :: ierr |
---|
| 930 | !!---------------------------------------------------------------------- |
---|
| 931 | ! |
---|
| 932 | CALL MPI_COMM_FREE(kcom, ierr) |
---|
| 933 | ! |
---|
| 934 | END SUBROUTINE mpp_comm_free |
---|
[3] | 935 | |
---|
[869] | 936 | |
---|
[2715] | 937 | SUBROUTINE mpp_ini_znl( kumout ) |
---|
[1345] | 938 | !!---------------------------------------------------------------------- |
---|
| 939 | !! *** routine mpp_ini_znl *** |
---|
| 940 | !! |
---|
| 941 | !! ** Purpose : Initialize special communicator for computing zonal sum |
---|
| 942 | !! |
---|
| 943 | !! ** Method : - Look for processors in the same row |
---|
| 944 | !! - Put their number in nrank_znl |
---|
| 945 | !! - Create group for the znl processors |
---|
| 946 | !! - Create a communicator for znl processors |
---|
| 947 | !! - Determine if processor should write znl files |
---|
| 948 | !! |
---|
| 949 | !! ** output |
---|
| 950 | !! ndim_rank_znl = number of processors on the same row |
---|
| 951 | !! ngrp_znl = group ID for the znl processors |
---|
| 952 | !! ncomm_znl = communicator for the ice procs. |
---|
| 953 | !! n_znl_root = number (in the world) of proc 0 in the ice comm. |
---|
| 954 | !! |
---|
| 955 | !!---------------------------------------------------------------------- |
---|
[2715] | 956 | INTEGER, INTENT(in) :: kumout ! ocean.output logical units |
---|
[1345] | 957 | ! |
---|
[2715] | 958 | INTEGER :: jproc ! dummy loop integer |
---|
| 959 | INTEGER :: ierr, ii ! local integer |
---|
| 960 | INTEGER, ALLOCATABLE, DIMENSION(:) :: kwork |
---|
| 961 | !!---------------------------------------------------------------------- |
---|
[1345] | 962 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_world : ', ngrp_world |
---|
| 963 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_world : ', mpi_comm_world |
---|
[9570] | 964 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_oce : ', mpi_comm_oce |
---|
[1345] | 965 | ! |
---|
[2715] | 966 | ALLOCATE( kwork(jpnij), STAT=ierr ) |
---|
| 967 | IF( ierr /= 0 ) THEN |
---|
| 968 | WRITE(kumout, cform_err) |
---|
| 969 | WRITE(kumout,*) 'mpp_ini_znl : failed to allocate 1D array of length jpnij' |
---|
| 970 | CALL mppstop |
---|
| 971 | ENDIF |
---|
| 972 | |
---|
| 973 | IF( jpnj == 1 ) THEN |
---|
[1345] | 974 | ngrp_znl = ngrp_world |
---|
[9570] | 975 | ncomm_znl = mpi_comm_oce |
---|
[1345] | 976 | ELSE |
---|
| 977 | ! |
---|
[9570] | 978 | CALL MPI_ALLGATHER ( njmpp, 1, mpi_integer, kwork, 1, mpi_integer, mpi_comm_oce, ierr ) |
---|
[1345] | 979 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - kwork pour njmpp : ', kwork |
---|
| 980 | !-$$ CALL flush(numout) |
---|
| 981 | ! |
---|
| 982 | ! Count number of processors on the same row |
---|
| 983 | ndim_rank_znl = 0 |
---|
| 984 | DO jproc=1,jpnij |
---|
| 985 | IF ( kwork(jproc) == njmpp ) THEN |
---|
| 986 | ndim_rank_znl = ndim_rank_znl + 1 |
---|
| 987 | ENDIF |
---|
| 988 | END DO |
---|
| 989 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ndim_rank_znl : ', ndim_rank_znl |
---|
| 990 | !-$$ CALL flush(numout) |
---|
| 991 | ! Allocate the right size to nrank_znl |
---|
[1441] | 992 | IF (ALLOCATED (nrank_znl)) DEALLOCATE(nrank_znl) |
---|
[1345] | 993 | ALLOCATE(nrank_znl(ndim_rank_znl)) |
---|
[3764] | 994 | ii = 0 |
---|
[1345] | 995 | nrank_znl (:) = 0 |
---|
| 996 | DO jproc=1,jpnij |
---|
| 997 | IF ( kwork(jproc) == njmpp) THEN |
---|
| 998 | ii = ii + 1 |
---|
[3764] | 999 | nrank_znl(ii) = jproc -1 |
---|
[1345] | 1000 | ENDIF |
---|
| 1001 | END DO |
---|
| 1002 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - nrank_znl : ', nrank_znl |
---|
| 1003 | !-$$ CALL flush(numout) |
---|
| 1004 | |
---|
| 1005 | ! Create the opa group |
---|
[9570] | 1006 | CALL MPI_COMM_GROUP(mpi_comm_oce,ngrp_opa,ierr) |
---|
[1345] | 1007 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_opa : ', ngrp_opa |
---|
| 1008 | !-$$ CALL flush(numout) |
---|
| 1009 | |
---|
| 1010 | ! Create the znl group from the opa group |
---|
| 1011 | CALL MPI_GROUP_INCL ( ngrp_opa, ndim_rank_znl, nrank_znl, ngrp_znl, ierr ) |
---|
| 1012 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_znl ', ngrp_znl |
---|
| 1013 | !-$$ CALL flush(numout) |
---|
| 1014 | |
---|
| 1015 | ! Create the znl communicator from the opa communicator, ie the pool of procs in the same row |
---|
[9570] | 1016 | CALL MPI_COMM_CREATE ( mpi_comm_oce, ngrp_znl, ncomm_znl, ierr ) |
---|
[1345] | 1017 | !-$$ WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ncomm_znl ', ncomm_znl |
---|
| 1018 | !-$$ CALL flush(numout) |
---|
| 1019 | ! |
---|
| 1020 | END IF |
---|
| 1021 | |
---|
| 1022 | ! Determines if processor if the first (starting from i=1) on the row |
---|
[3764] | 1023 | IF ( jpni == 1 ) THEN |
---|
[1345] | 1024 | l_znl_root = .TRUE. |
---|
| 1025 | ELSE |
---|
| 1026 | l_znl_root = .FALSE. |
---|
| 1027 | kwork (1) = nimpp |
---|
[10425] | 1028 | CALL mpp_min ( 'lib_mpp', kwork(1), kcom = ncomm_znl) |
---|
[1345] | 1029 | IF ( nimpp == kwork(1)) l_znl_root = .TRUE. |
---|
| 1030 | END IF |
---|
| 1031 | |
---|
[2715] | 1032 | DEALLOCATE(kwork) |
---|
| 1033 | |
---|
[1345] | 1034 | END SUBROUTINE mpp_ini_znl |
---|
| 1035 | |
---|
| 1036 | |
---|
[1344] | 1037 | SUBROUTINE mpp_ini_north |
---|
| 1038 | !!---------------------------------------------------------------------- |
---|
| 1039 | !! *** routine mpp_ini_north *** |
---|
| 1040 | !! |
---|
[3764] | 1041 | !! ** Purpose : Initialize special communicator for north folding |
---|
[1344] | 1042 | !! condition together with global variables needed in the mpp folding |
---|
| 1043 | !! |
---|
| 1044 | !! ** Method : - Look for northern processors |
---|
| 1045 | !! - Put their number in nrank_north |
---|
| 1046 | !! - Create groups for the world processors and the north processors |
---|
| 1047 | !! - Create a communicator for northern processors |
---|
| 1048 | !! |
---|
| 1049 | !! ** output |
---|
| 1050 | !! njmppmax = njmpp for northern procs |
---|
| 1051 | !! ndim_rank_north = number of processors in the northern line |
---|
| 1052 | !! nrank_north (ndim_rank_north) = number of the northern procs. |
---|
| 1053 | !! ngrp_world = group ID for the world processors |
---|
| 1054 | !! ngrp_north = group ID for the northern processors |
---|
| 1055 | !! ncomm_north = communicator for the northern procs. |
---|
| 1056 | !! north_root = number (in the world) of proc 0 in the northern comm. |
---|
| 1057 | !! |
---|
| 1058 | !!---------------------------------------------------------------------- |
---|
| 1059 | INTEGER :: ierr |
---|
| 1060 | INTEGER :: jjproc |
---|
| 1061 | INTEGER :: ii, ji |
---|
| 1062 | !!---------------------------------------------------------------------- |
---|
| 1063 | ! |
---|
| 1064 | njmppmax = MAXVAL( njmppt ) |
---|
| 1065 | ! |
---|
| 1066 | ! Look for how many procs on the northern boundary |
---|
| 1067 | ndim_rank_north = 0 |
---|
| 1068 | DO jjproc = 1, jpnij |
---|
| 1069 | IF( njmppt(jjproc) == njmppmax ) ndim_rank_north = ndim_rank_north + 1 |
---|
| 1070 | END DO |
---|
| 1071 | ! |
---|
| 1072 | ! Allocate the right size to nrank_north |
---|
[1441] | 1073 | IF (ALLOCATED (nrank_north)) DEALLOCATE(nrank_north) |
---|
[1344] | 1074 | ALLOCATE( nrank_north(ndim_rank_north) ) |
---|
[869] | 1075 | |
---|
[1344] | 1076 | ! Fill the nrank_north array with proc. number of northern procs. |
---|
| 1077 | ! Note : the rank start at 0 in MPI |
---|
| 1078 | ii = 0 |
---|
| 1079 | DO ji = 1, jpnij |
---|
| 1080 | IF ( njmppt(ji) == njmppmax ) THEN |
---|
| 1081 | ii=ii+1 |
---|
| 1082 | nrank_north(ii)=ji-1 |
---|
| 1083 | END IF |
---|
| 1084 | END DO |
---|
| 1085 | ! |
---|
| 1086 | ! create the world group |
---|
[9570] | 1087 | CALL MPI_COMM_GROUP( mpi_comm_oce, ngrp_world, ierr ) |
---|
[1344] | 1088 | ! |
---|
| 1089 | ! Create the North group from the world group |
---|
| 1090 | CALL MPI_GROUP_INCL( ngrp_world, ndim_rank_north, nrank_north, ngrp_north, ierr ) |
---|
| 1091 | ! |
---|
| 1092 | ! Create the North communicator , ie the pool of procs in the north group |
---|
[9570] | 1093 | CALL MPI_COMM_CREATE( mpi_comm_oce, ngrp_north, ncomm_north, ierr ) |
---|
[1344] | 1094 | ! |
---|
| 1095 | END SUBROUTINE mpp_ini_north |
---|
[869] | 1096 | |
---|
| 1097 | |
---|
[9570] | 1098 | SUBROUTINE mpi_init_oce( ldtxt, ksft, code ) |
---|
[1344] | 1099 | !!--------------------------------------------------------------------- |
---|
| 1100 | !! *** routine mpp_init.opa *** |
---|
| 1101 | !! |
---|
| 1102 | !! ** Purpose :: export and attach a MPI buffer for bsend |
---|
| 1103 | !! |
---|
| 1104 | !! ** Method :: define buffer size in namelist, if 0 no buffer attachment |
---|
| 1105 | !! but classical mpi_init |
---|
[3764] | 1106 | !! |
---|
| 1107 | !! History :: 01/11 :: IDRIS initial version for IBM only |
---|
[1344] | 1108 | !! 08/04 :: R. Benshila, generalisation |
---|
| 1109 | !!--------------------------------------------------------------------- |
---|
[3764] | 1110 | CHARACTER(len=*),DIMENSION(:), INTENT( out) :: ldtxt |
---|
[2481] | 1111 | INTEGER , INTENT(inout) :: ksft |
---|
| 1112 | INTEGER , INTENT( out) :: code |
---|
| 1113 | INTEGER :: ierr, ji |
---|
| 1114 | LOGICAL :: mpi_was_called |
---|
[1344] | 1115 | !!--------------------------------------------------------------------- |
---|
| 1116 | ! |
---|
| 1117 | CALL mpi_initialized( mpi_was_called, code ) ! MPI initialization |
---|
[532] | 1118 | IF ( code /= MPI_SUCCESS ) THEN |
---|
[3764] | 1119 | DO ji = 1, SIZE(ldtxt) |
---|
[2481] | 1120 | IF( TRIM(ldtxt(ji)) /= '' ) WRITE(*,*) ldtxt(ji) ! control print of mynode |
---|
[3764] | 1121 | END DO |
---|
[2481] | 1122 | WRITE(*, cform_err) |
---|
| 1123 | WRITE(*, *) ' lib_mpp: Error in routine mpi_initialized' |
---|
[1344] | 1124 | CALL mpi_abort( mpi_comm_world, code, ierr ) |
---|
[532] | 1125 | ENDIF |
---|
[1344] | 1126 | ! |
---|
| 1127 | IF( .NOT. mpi_was_called ) THEN |
---|
| 1128 | CALL mpi_init( code ) |
---|
[9570] | 1129 | CALL mpi_comm_dup( mpi_comm_world, mpi_comm_oce, code ) |
---|
[532] | 1130 | IF ( code /= MPI_SUCCESS ) THEN |
---|
[3764] | 1131 | DO ji = 1, SIZE(ldtxt) |
---|
[2481] | 1132 | IF( TRIM(ldtxt(ji)) /= '' ) WRITE(*,*) ldtxt(ji) ! control print of mynode |
---|
| 1133 | END DO |
---|
| 1134 | WRITE(*, cform_err) |
---|
| 1135 | WRITE(*, *) ' lib_mpp: Error in routine mpi_comm_dup' |
---|
[532] | 1136 | CALL mpi_abort( mpi_comm_world, code, ierr ) |
---|
| 1137 | ENDIF |
---|
| 1138 | ENDIF |
---|
[1344] | 1139 | ! |
---|
[897] | 1140 | IF( nn_buffer > 0 ) THEN |
---|
[2481] | 1141 | WRITE(ldtxt(ksft),*) 'mpi_bsend, buffer allocation of : ', nn_buffer ; ksft = ksft + 1 |
---|
[897] | 1142 | ! Buffer allocation and attachment |
---|
[2481] | 1143 | ALLOCATE( tampon(nn_buffer), stat = ierr ) |
---|
[3764] | 1144 | IF( ierr /= 0 ) THEN |
---|
| 1145 | DO ji = 1, SIZE(ldtxt) |
---|
[2481] | 1146 | IF( TRIM(ldtxt(ji)) /= '' ) WRITE(*,*) ldtxt(ji) ! control print of mynode |
---|
| 1147 | END DO |
---|
| 1148 | WRITE(*, cform_err) |
---|
| 1149 | WRITE(*, *) ' lib_mpp: Error in ALLOCATE', ierr |
---|
| 1150 | CALL mpi_abort( mpi_comm_world, code, ierr ) |
---|
| 1151 | END IF |
---|
| 1152 | CALL mpi_buffer_attach( tampon, nn_buffer, code ) |
---|
[897] | 1153 | ENDIF |
---|
[1344] | 1154 | ! |
---|
[9570] | 1155 | END SUBROUTINE mpi_init_oce |
---|
[3] | 1156 | |
---|
[9019] | 1157 | |
---|
| 1158 | SUBROUTINE DDPDD_MPI( ydda, yddb, ilen, itype ) |
---|
[1976] | 1159 | !!--------------------------------------------------------------------- |
---|
| 1160 | !! Routine DDPDD_MPI: used by reduction operator MPI_SUMDD |
---|
| 1161 | !! |
---|
| 1162 | !! Modification of original codes written by David H. Bailey |
---|
| 1163 | !! This subroutine computes yddb(i) = ydda(i)+yddb(i) |
---|
| 1164 | !!--------------------------------------------------------------------- |
---|
[9019] | 1165 | INTEGER , INTENT(in) :: ilen, itype |
---|
| 1166 | COMPLEX(wp), DIMENSION(ilen), INTENT(in) :: ydda |
---|
| 1167 | COMPLEX(wp), DIMENSION(ilen), INTENT(inout) :: yddb |
---|
[1976] | 1168 | ! |
---|
| 1169 | REAL(wp) :: zerr, zt1, zt2 ! local work variables |
---|
[9019] | 1170 | INTEGER :: ji, ztmp ! local scalar |
---|
| 1171 | !!--------------------------------------------------------------------- |
---|
| 1172 | ! |
---|
[1976] | 1173 | ztmp = itype ! avoid compilation warning |
---|
[9019] | 1174 | ! |
---|
[1976] | 1175 | DO ji=1,ilen |
---|
| 1176 | ! Compute ydda + yddb using Knuth's trick. |
---|
| 1177 | zt1 = real(ydda(ji)) + real(yddb(ji)) |
---|
| 1178 | zerr = zt1 - real(ydda(ji)) |
---|
| 1179 | zt2 = ((real(yddb(ji)) - zerr) + (real(ydda(ji)) - (zt1 - zerr))) & |
---|
| 1180 | + aimag(ydda(ji)) + aimag(yddb(ji)) |
---|
| 1181 | |
---|
| 1182 | ! The result is zt1 + zt2, after normalization. |
---|
| 1183 | yddb(ji) = cmplx ( zt1 + zt2, zt2 - ((zt1 + zt2) - zt1),wp ) |
---|
| 1184 | END DO |
---|
[9019] | 1185 | ! |
---|
[1976] | 1186 | END SUBROUTINE DDPDD_MPI |
---|
| 1187 | |
---|
[6140] | 1188 | |
---|
[9019] | 1189 | SUBROUTINE mpp_lbc_north_icb( pt2d, cd_type, psgn, kextj) |
---|
[4990] | 1190 | !!--------------------------------------------------------------------- |
---|
| 1191 | !! *** routine mpp_lbc_north_icb *** |
---|
| 1192 | !! |
---|
| 1193 | !! ** Purpose : Ensure proper north fold horizontal bondary condition |
---|
| 1194 | !! in mpp configuration in case of jpn1 > 1 and for 2d |
---|
| 1195 | !! array with outer extra halo |
---|
| 1196 | !! |
---|
| 1197 | !! ** Method : North fold condition and mpp with more than one proc |
---|
| 1198 | !! in i-direction require a specific treatment. We gather |
---|
[9019] | 1199 | !! the 4+kextj northern lines of the global domain on 1 |
---|
[4990] | 1200 | !! processor and apply lbc north-fold on this sub array. |
---|
| 1201 | !! Then we scatter the north fold array back to the processors. |
---|
[9019] | 1202 | !! This routine accounts for an extra halo with icebergs |
---|
| 1203 | !! and assumes ghost rows and columns have been suppressed. |
---|
[4990] | 1204 | !! |
---|
| 1205 | !!---------------------------------------------------------------------- |
---|
| 1206 | REAL(wp), DIMENSION(:,:), INTENT(inout) :: pt2d ! 2D array with extra halo |
---|
| 1207 | CHARACTER(len=1) , INTENT(in ) :: cd_type ! nature of pt3d grid-points |
---|
| 1208 | ! ! = T , U , V , F or W -points |
---|
| 1209 | REAL(wp) , INTENT(in ) :: psgn ! = -1. the sign change across the |
---|
| 1210 | !! ! north fold, = 1. otherwise |
---|
[9019] | 1211 | INTEGER , INTENT(in ) :: kextj ! Extra halo width at north fold |
---|
[6140] | 1212 | ! |
---|
[4990] | 1213 | INTEGER :: ji, jj, jr |
---|
| 1214 | INTEGER :: ierr, itaille, ildi, ilei, iilb |
---|
[9019] | 1215 | INTEGER :: ipj, ij, iproc |
---|
[4990] | 1216 | ! |
---|
| 1217 | REAL(wp), DIMENSION(:,:) , ALLOCATABLE :: ztab_e, znorthloc_e |
---|
| 1218 | REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: znorthgloio_e |
---|
| 1219 | !!---------------------------------------------------------------------- |
---|
| 1220 | ! |
---|
[9019] | 1221 | ipj=4 |
---|
[9467] | 1222 | ALLOCATE( ztab_e(jpiglo, 1-kextj:ipj+kextj) , & |
---|
| 1223 | & znorthloc_e(jpimax, 1-kextj:ipj+kextj) , & |
---|
| 1224 | & znorthgloio_e(jpimax, 1-kextj:ipj+kextj,jpni) ) |
---|
[4990] | 1225 | ! |
---|
[9019] | 1226 | ztab_e(:,:) = 0._wp |
---|
| 1227 | znorthloc_e(:,:) = 0._wp |
---|
[6140] | 1228 | ! |
---|
[9467] | 1229 | ij = 1 - kextj |
---|
| 1230 | ! put the last ipj+2*kextj lines of pt2d into znorthloc_e |
---|
| 1231 | DO jj = jpj - ipj + 1 - kextj , jpj + kextj |
---|
| 1232 | znorthloc_e(1:jpi,ij)=pt2d(1:jpi,jj) |
---|
[4990] | 1233 | ij = ij + 1 |
---|
| 1234 | END DO |
---|
| 1235 | ! |
---|
[9467] | 1236 | itaille = jpimax * ( ipj + 2*kextj ) |
---|
[10425] | 1237 | ! |
---|
| 1238 | IF( ln_timing ) CALL tic_tac(.TRUE.) |
---|
[9467] | 1239 | CALL MPI_ALLGATHER( znorthloc_e(1,1-kextj) , itaille, MPI_DOUBLE_PRECISION, & |
---|
| 1240 | & znorthgloio_e(1,1-kextj,1), itaille, MPI_DOUBLE_PRECISION, & |
---|
| 1241 | & ncomm_north, ierr ) |
---|
[4990] | 1242 | ! |
---|
[10425] | 1243 | IF( ln_timing ) CALL tic_tac(.FALSE.) |
---|
| 1244 | ! |
---|
[4990] | 1245 | DO jr = 1, ndim_rank_north ! recover the global north array |
---|
| 1246 | iproc = nrank_north(jr) + 1 |
---|
| 1247 | ildi = nldit (iproc) |
---|
| 1248 | ilei = nleit (iproc) |
---|
| 1249 | iilb = nimppt(iproc) |
---|
[9467] | 1250 | DO jj = 1-kextj, ipj+kextj |
---|
[4990] | 1251 | DO ji = ildi, ilei |
---|
| 1252 | ztab_e(ji+iilb-1,jj) = znorthgloio_e(ji,jj,jr) |
---|
| 1253 | END DO |
---|
| 1254 | END DO |
---|
| 1255 | END DO |
---|
| 1256 | |
---|
| 1257 | ! 2. North-Fold boundary conditions |
---|
| 1258 | ! ---------------------------------- |
---|
[9467] | 1259 | CALL lbc_nfd( ztab_e(:,1-kextj:ipj+kextj), cd_type, psgn, kextj ) |
---|
[4990] | 1260 | |
---|
[9467] | 1261 | ij = 1 - kextj |
---|
[4990] | 1262 | !! Scatter back to pt2d |
---|
[9467] | 1263 | DO jj = jpj - ipj + 1 - kextj , jpj + kextj |
---|
[9019] | 1264 | DO ji= 1, jpi |
---|
[4990] | 1265 | pt2d(ji,jj) = ztab_e(ji+nimpp-1,ij) |
---|
| 1266 | END DO |
---|
[9467] | 1267 | ij = ij +1 |
---|
[4990] | 1268 | END DO |
---|
| 1269 | ! |
---|
| 1270 | DEALLOCATE( ztab_e, znorthloc_e, znorthgloio_e ) |
---|
| 1271 | ! |
---|
| 1272 | END SUBROUTINE mpp_lbc_north_icb |
---|
| 1273 | |
---|
[6140] | 1274 | |
---|
[10425] | 1275 | SUBROUTINE mpp_lnk_2d_icb( cdname, pt2d, cd_type, psgn, kexti, kextj ) |
---|
[4990] | 1276 | !!---------------------------------------------------------------------- |
---|
| 1277 | !! *** routine mpp_lnk_2d_icb *** |
---|
| 1278 | !! |
---|
[9019] | 1279 | !! ** Purpose : Message passing management for 2d array (with extra halo for icebergs) |
---|
| 1280 | !! This routine receives a (1-kexti:jpi+kexti,1-kexti:jpj+kextj) |
---|
| 1281 | !! array (usually (0:jpi+1, 0:jpj+1)) from lbc_lnk_icb calls. |
---|
[4990] | 1282 | !! |
---|
| 1283 | !! ** Method : Use mppsend and mpprecv function for passing mask |
---|
| 1284 | !! between processors following neighboring subdomains. |
---|
| 1285 | !! domain parameters |
---|
[9019] | 1286 | !! jpi : first dimension of the local subdomain |
---|
| 1287 | !! jpj : second dimension of the local subdomain |
---|
| 1288 | !! kexti : number of columns for extra outer halo |
---|
| 1289 | !! kextj : number of rows for extra outer halo |
---|
[4990] | 1290 | !! nbondi : mark for "east-west local boundary" |
---|
| 1291 | !! nbondj : mark for "north-south local boundary" |
---|
| 1292 | !! noea : number for local neighboring processors |
---|
| 1293 | !! nowe : number for local neighboring processors |
---|
| 1294 | !! noso : number for local neighboring processors |
---|
| 1295 | !! nono : number for local neighboring processors |
---|
| 1296 | !!---------------------------------------------------------------------- |
---|
[10425] | 1297 | CHARACTER(len=*) , INTENT(in ) :: cdname ! name of the calling subroutine |
---|
[9019] | 1298 | REAL(wp), DIMENSION(1-kexti:jpi+kexti,1-kextj:jpj+kextj), INTENT(inout) :: pt2d ! 2D array with extra halo |
---|
| 1299 | CHARACTER(len=1) , INTENT(in ) :: cd_type ! nature of ptab array grid-points |
---|
| 1300 | REAL(wp) , INTENT(in ) :: psgn ! sign used across the north fold |
---|
| 1301 | INTEGER , INTENT(in ) :: kexti ! extra i-halo width |
---|
| 1302 | INTEGER , INTENT(in ) :: kextj ! extra j-halo width |
---|
| 1303 | ! |
---|
[4990] | 1304 | INTEGER :: jl ! dummy loop indices |
---|
[9019] | 1305 | INTEGER :: imigr, iihom, ijhom ! local integers |
---|
| 1306 | INTEGER :: ipreci, iprecj ! - - |
---|
[4990] | 1307 | INTEGER :: ml_req1, ml_req2, ml_err ! for key_mpi_isend |
---|
| 1308 | INTEGER, DIMENSION(MPI_STATUS_SIZE) :: ml_stat ! for key_mpi_isend |
---|
| 1309 | !! |
---|
[9019] | 1310 | REAL(wp), DIMENSION(1-kexti:jpi+kexti,nn_hls+kextj,2) :: r2dns, r2dsn |
---|
| 1311 | REAL(wp), DIMENSION(1-kextj:jpj+kextj,nn_hls+kexti,2) :: r2dwe, r2dew |
---|
[4990] | 1312 | !!---------------------------------------------------------------------- |
---|
| 1313 | |
---|
[9019] | 1314 | ipreci = nn_hls + kexti ! take into account outer extra 2D overlap area |
---|
| 1315 | iprecj = nn_hls + kextj |
---|
[4990] | 1316 | |
---|
[10425] | 1317 | IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, 1, 1, 1, ld_lbc = .TRUE. ) |
---|
[4990] | 1318 | |
---|
| 1319 | ! 1. standard boundary treatment |
---|
| 1320 | ! ------------------------------ |
---|
| 1321 | ! Order matters Here !!!! |
---|
| 1322 | ! |
---|
| 1323 | ! ! East-West boundaries |
---|
| 1324 | ! !* Cyclic east-west |
---|
[9667] | 1325 | IF( l_Iperio ) THEN |
---|
| 1326 | pt2d(1-kexti: 1 ,:) = pt2d(jpim1-kexti: jpim1 ,:) ! east |
---|
| 1327 | pt2d( jpi :jpi+kexti,:) = pt2d( 2 :2+kexti,:) ! west |
---|
[4990] | 1328 | ! |
---|
| 1329 | ELSE !* closed |
---|
[9667] | 1330 | IF( .NOT. cd_type == 'F' ) pt2d( 1-kexti :nn_hls ,:) = 0._wp ! east except at F-point |
---|
| 1331 | pt2d(jpi-nn_hls+1:jpi+kexti,:) = 0._wp ! west |
---|
[4990] | 1332 | ENDIF |
---|
[9667] | 1333 | ! ! North-South boundaries |
---|
| 1334 | IF( l_Jperio ) THEN !* cyclic (only with no mpp j-split) |
---|
| 1335 | pt2d(:,1-kextj: 1 ) = pt2d(:,jpjm1-kextj: jpjm1) ! north |
---|
| 1336 | pt2d(:, jpj :jpj+kextj) = pt2d(:, 2 :2+kextj) ! south |
---|
| 1337 | ELSE !* closed |
---|
| 1338 | IF( .NOT. cd_type == 'F' ) pt2d(:, 1-kextj :nn_hls ) = 0._wp ! north except at F-point |
---|
| 1339 | pt2d(:,jpj-nn_hls+1:jpj+kextj) = 0._wp ! south |
---|
| 1340 | ENDIF |
---|
[4990] | 1341 | ! |
---|
| 1342 | |
---|
| 1343 | ! north fold treatment |
---|
| 1344 | ! ----------------------- |
---|
| 1345 | IF( npolj /= 0 ) THEN |
---|
| 1346 | ! |
---|
| 1347 | SELECT CASE ( jpni ) |
---|
[9019] | 1348 | CASE ( 1 ) ; CALL lbc_nfd ( pt2d(1:jpi,1:jpj+kextj), cd_type, psgn, kextj ) |
---|
| 1349 | CASE DEFAULT ; CALL mpp_lbc_north_icb( pt2d(1:jpi,1:jpj+kextj), cd_type, psgn, kextj ) |
---|
[4990] | 1350 | END SELECT |
---|
| 1351 | ! |
---|
| 1352 | ENDIF |
---|
| 1353 | |
---|
| 1354 | ! 2. East and west directions exchange |
---|
| 1355 | ! ------------------------------------ |
---|
| 1356 | ! we play with the neigbours AND the row number because of the periodicity |
---|
| 1357 | ! |
---|
| 1358 | SELECT CASE ( nbondi ) ! Read Dirichlet lateral conditions |
---|
| 1359 | CASE ( -1, 0, 1 ) ! all exept 2 (i.e. close case) |
---|
[9019] | 1360 | iihom = jpi-nreci-kexti |
---|
[4990] | 1361 | DO jl = 1, ipreci |
---|
[9019] | 1362 | r2dew(:,jl,1) = pt2d(nn_hls+jl,:) |
---|
[4990] | 1363 | r2dwe(:,jl,1) = pt2d(iihom +jl,:) |
---|
| 1364 | END DO |
---|
| 1365 | END SELECT |
---|
| 1366 | ! |
---|
| 1367 | ! ! Migrations |
---|
[9019] | 1368 | imigr = ipreci * ( jpj + 2*kextj ) |
---|
[4990] | 1369 | ! |
---|
[10425] | 1370 | IF( ln_timing ) CALL tic_tac(.TRUE.) |
---|
| 1371 | ! |
---|
[4990] | 1372 | SELECT CASE ( nbondi ) |
---|
| 1373 | CASE ( -1 ) |
---|
[9019] | 1374 | CALL mppsend( 2, r2dwe(1-kextj,1,1), imigr, noea, ml_req1 ) |
---|
| 1375 | CALL mpprecv( 1, r2dew(1-kextj,1,2), imigr, noea ) |
---|
[4990] | 1376 | IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) |
---|
| 1377 | CASE ( 0 ) |
---|
[9019] | 1378 | CALL mppsend( 1, r2dew(1-kextj,1,1), imigr, nowe, ml_req1 ) |
---|
| 1379 | CALL mppsend( 2, r2dwe(1-kextj,1,1), imigr, noea, ml_req2 ) |
---|
| 1380 | CALL mpprecv( 1, r2dew(1-kextj,1,2), imigr, noea ) |
---|
| 1381 | CALL mpprecv( 2, r2dwe(1-kextj,1,2), imigr, nowe ) |
---|
[4990] | 1382 | IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) |
---|
| 1383 | IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) |
---|
| 1384 | CASE ( 1 ) |
---|
[9019] | 1385 | CALL mppsend( 1, r2dew(1-kextj,1,1), imigr, nowe, ml_req1 ) |
---|
| 1386 | CALL mpprecv( 2, r2dwe(1-kextj,1,2), imigr, nowe ) |
---|
[4990] | 1387 | IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) |
---|
| 1388 | END SELECT |
---|
| 1389 | ! |
---|
[10425] | 1390 | IF( ln_timing ) CALL tic_tac(.FALSE.) |
---|
| 1391 | ! |
---|
[4990] | 1392 | ! ! Write Dirichlet lateral conditions |
---|
[9019] | 1393 | iihom = jpi - nn_hls |
---|
[4990] | 1394 | ! |
---|
| 1395 | SELECT CASE ( nbondi ) |
---|
| 1396 | CASE ( -1 ) |
---|
| 1397 | DO jl = 1, ipreci |
---|
| 1398 | pt2d(iihom+jl,:) = r2dew(:,jl,2) |
---|
| 1399 | END DO |
---|
| 1400 | CASE ( 0 ) |
---|
| 1401 | DO jl = 1, ipreci |
---|
[9019] | 1402 | pt2d(jl-kexti,:) = r2dwe(:,jl,2) |
---|
| 1403 | pt2d(iihom+jl,:) = r2dew(:,jl,2) |
---|
[4990] | 1404 | END DO |
---|
| 1405 | CASE ( 1 ) |
---|
| 1406 | DO jl = 1, ipreci |
---|
[9019] | 1407 | pt2d(jl-kexti,:) = r2dwe(:,jl,2) |
---|
[4990] | 1408 | END DO |
---|
| 1409 | END SELECT |
---|
| 1410 | |
---|
| 1411 | |
---|
| 1412 | ! 3. North and south directions |
---|
| 1413 | ! ----------------------------- |
---|
| 1414 | ! always closed : we play only with the neigbours |
---|
| 1415 | ! |
---|
| 1416 | IF( nbondj /= 2 ) THEN ! Read Dirichlet lateral conditions |
---|
[9019] | 1417 | ijhom = jpj-nrecj-kextj |
---|
[4990] | 1418 | DO jl = 1, iprecj |
---|
| 1419 | r2dsn(:,jl,1) = pt2d(:,ijhom +jl) |
---|
[9019] | 1420 | r2dns(:,jl,1) = pt2d(:,nn_hls+jl) |
---|
[4990] | 1421 | END DO |
---|
| 1422 | ENDIF |
---|
| 1423 | ! |
---|
| 1424 | ! ! Migrations |
---|
[9019] | 1425 | imigr = iprecj * ( jpi + 2*kexti ) |
---|
[4990] | 1426 | ! |
---|
[10425] | 1427 | IF( ln_timing ) CALL tic_tac(.TRUE.) |
---|
| 1428 | ! |
---|
[4990] | 1429 | SELECT CASE ( nbondj ) |
---|
| 1430 | CASE ( -1 ) |
---|
[9019] | 1431 | CALL mppsend( 4, r2dsn(1-kexti,1,1), imigr, nono, ml_req1 ) |
---|
| 1432 | CALL mpprecv( 3, r2dns(1-kexti,1,2), imigr, nono ) |
---|
[4990] | 1433 | IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) |
---|
| 1434 | CASE ( 0 ) |
---|
[9019] | 1435 | CALL mppsend( 3, r2dns(1-kexti,1,1), imigr, noso, ml_req1 ) |
---|
| 1436 | CALL mppsend( 4, r2dsn(1-kexti,1,1), imigr, nono, ml_req2 ) |
---|
| 1437 | CALL mpprecv( 3, r2dns(1-kexti,1,2), imigr, nono ) |
---|
| 1438 | CALL mpprecv( 4, r2dsn(1-kexti,1,2), imigr, noso ) |
---|
[4990] | 1439 | IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) |
---|
| 1440 | IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) |
---|
| 1441 | CASE ( 1 ) |
---|
[9019] | 1442 | CALL mppsend( 3, r2dns(1-kexti,1,1), imigr, noso, ml_req1 ) |
---|
| 1443 | CALL mpprecv( 4, r2dsn(1-kexti,1,2), imigr, noso ) |
---|
[4990] | 1444 | IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) |
---|
| 1445 | END SELECT |
---|
| 1446 | ! |
---|
[10425] | 1447 | IF( ln_timing ) CALL tic_tac(.FALSE.) |
---|
| 1448 | ! |
---|
[4990] | 1449 | ! ! Write Dirichlet lateral conditions |
---|
[9019] | 1450 | ijhom = jpj - nn_hls |
---|
[4990] | 1451 | ! |
---|
| 1452 | SELECT CASE ( nbondj ) |
---|
| 1453 | CASE ( -1 ) |
---|
| 1454 | DO jl = 1, iprecj |
---|
| 1455 | pt2d(:,ijhom+jl) = r2dns(:,jl,2) |
---|
| 1456 | END DO |
---|
| 1457 | CASE ( 0 ) |
---|
| 1458 | DO jl = 1, iprecj |
---|
[9019] | 1459 | pt2d(:,jl-kextj) = r2dsn(:,jl,2) |
---|
| 1460 | pt2d(:,ijhom+jl) = r2dns(:,jl,2) |
---|
[4990] | 1461 | END DO |
---|
| 1462 | CASE ( 1 ) |
---|
| 1463 | DO jl = 1, iprecj |
---|
[9019] | 1464 | pt2d(:,jl-kextj) = r2dsn(:,jl,2) |
---|
[4990] | 1465 | END DO |
---|
| 1466 | END SELECT |
---|
[9019] | 1467 | ! |
---|
[4990] | 1468 | END SUBROUTINE mpp_lnk_2d_icb |
---|
[10425] | 1469 | |
---|
| 1470 | |
---|
[10437] | 1471 | SUBROUTINE mpp_report( cdname, kpk, kpl, kpf, ld_lbc, ld_glb, ld_dlg ) |
---|
[10425] | 1472 | !!---------------------------------------------------------------------- |
---|
| 1473 | !! *** routine mpp_report *** |
---|
| 1474 | !! |
---|
| 1475 | !! ** Purpose : report use of mpp routines per time-setp |
---|
| 1476 | !! |
---|
| 1477 | !!---------------------------------------------------------------------- |
---|
| 1478 | CHARACTER(len=*), INTENT(in ) :: cdname ! name of the calling subroutine |
---|
| 1479 | INTEGER , OPTIONAL, INTENT(in ) :: kpk, kpl, kpf |
---|
[10437] | 1480 | LOGICAL , OPTIONAL, INTENT(in ) :: ld_lbc, ld_glb, ld_dlg |
---|
[10425] | 1481 | !! |
---|
[11081] | 1482 | CHARACTER(len=128) :: ccountname ! name of a subroutine to count communications |
---|
[10437] | 1483 | LOGICAL :: ll_lbc, ll_glb, ll_dlg |
---|
[11081] | 1484 | INTEGER :: ji, jj, jk, jh, jf, jcount ! dummy loop indices |
---|
[10425] | 1485 | !!---------------------------------------------------------------------- |
---|
| 1486 | ! |
---|
| 1487 | ll_lbc = .FALSE. |
---|
| 1488 | IF( PRESENT(ld_lbc) ) ll_lbc = ld_lbc |
---|
| 1489 | ll_glb = .FALSE. |
---|
| 1490 | IF( PRESENT(ld_glb) ) ll_glb = ld_glb |
---|
[10437] | 1491 | ll_dlg = .FALSE. |
---|
| 1492 | IF( PRESENT(ld_dlg) ) ll_dlg = ld_dlg |
---|
[10425] | 1493 | ! |
---|
| 1494 | ! find the smallest common frequency: default = frequency product, if multiple, choose the larger of the 2 frequency |
---|
| 1495 | IF( ncom_dttrc /= 1 ) CALL ctl_stop( 'STOP', 'mpp_report, ncom_dttrc /= 1 not coded...' ) |
---|
| 1496 | ncom_freq = ncom_fsbc |
---|
| 1497 | ! |
---|
| 1498 | IF ( ncom_stp == nit000+ncom_freq ) THEN ! avoid to count extra communications in potential initializations at nit000 |
---|
| 1499 | IF( ll_lbc ) THEN |
---|
| 1500 | IF( .NOT. ALLOCATED(ncomm_sequence) ) ALLOCATE( ncomm_sequence(ncom_rec_max,2) ) |
---|
| 1501 | IF( .NOT. ALLOCATED( crname_lbc) ) ALLOCATE( crname_lbc(ncom_rec_max ) ) |
---|
| 1502 | n_sequence_lbc = n_sequence_lbc + 1 |
---|
| 1503 | IF( n_sequence_lbc > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' ) ! deadlock |
---|
| 1504 | crname_lbc(n_sequence_lbc) = cdname ! keep the name of the calling routine |
---|
| 1505 | ncomm_sequence(n_sequence_lbc,1) = kpk*kpl ! size of 3rd and 4th dimensions |
---|
| 1506 | ncomm_sequence(n_sequence_lbc,2) = kpf ! number of arrays to be treated (multi) |
---|
| 1507 | ENDIF |
---|
| 1508 | IF( ll_glb ) THEN |
---|
| 1509 | IF( .NOT. ALLOCATED(crname_glb) ) ALLOCATE( crname_glb(ncom_rec_max) ) |
---|
| 1510 | n_sequence_glb = n_sequence_glb + 1 |
---|
| 1511 | IF( n_sequence_glb > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' ) ! deadlock |
---|
| 1512 | crname_glb(n_sequence_glb) = cdname ! keep the name of the calling routine |
---|
| 1513 | ENDIF |
---|
[10437] | 1514 | IF( ll_dlg ) THEN |
---|
| 1515 | IF( .NOT. ALLOCATED(crname_dlg) ) ALLOCATE( crname_dlg(ncom_rec_max) ) |
---|
| 1516 | n_sequence_dlg = n_sequence_dlg + 1 |
---|
| 1517 | IF( n_sequence_dlg > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' ) ! deadlock |
---|
| 1518 | crname_dlg(n_sequence_dlg) = cdname ! keep the name of the calling routine |
---|
| 1519 | ENDIF |
---|
[10425] | 1520 | ELSE IF ( ncom_stp == nit000+2*ncom_freq ) THEN |
---|
| 1521 | CALL ctl_opn( numcom, 'communication_report.txt', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE., narea ) |
---|
| 1522 | WRITE(numcom,*) ' ' |
---|
| 1523 | WRITE(numcom,*) ' ------------------------------------------------------------' |
---|
| 1524 | WRITE(numcom,*) ' Communication pattern report (second oce+sbc+top time step):' |
---|
| 1525 | WRITE(numcom,*) ' ------------------------------------------------------------' |
---|
| 1526 | WRITE(numcom,*) ' ' |
---|
| 1527 | WRITE(numcom,'(A,I4)') ' Exchanged halos : ', n_sequence_lbc |
---|
| 1528 | jj = 0; jk = 0; jf = 0; jh = 0 |
---|
| 1529 | DO ji = 1, n_sequence_lbc |
---|
| 1530 | IF ( ncomm_sequence(ji,1) .GT. 1 ) jk = jk + 1 |
---|
| 1531 | IF ( ncomm_sequence(ji,2) .GT. 1 ) jf = jf + 1 |
---|
| 1532 | IF ( ncomm_sequence(ji,1) .GT. 1 .AND. ncomm_sequence(ji,2) .GT. 1 ) jj = jj + 1 |
---|
| 1533 | jh = MAX (jh, ncomm_sequence(ji,1)*ncomm_sequence(ji,2)) |
---|
| 1534 | END DO |
---|
| 1535 | WRITE(numcom,'(A,I3)') ' 3D Exchanged halos : ', jk |
---|
| 1536 | WRITE(numcom,'(A,I3)') ' Multi arrays exchanged halos : ', jf |
---|
| 1537 | WRITE(numcom,'(A,I3)') ' from which 3D : ', jj |
---|
| 1538 | WRITE(numcom,'(A,I10)') ' Array max size : ', jh*jpi*jpj |
---|
| 1539 | WRITE(numcom,*) ' ' |
---|
| 1540 | WRITE(numcom,*) ' lbc_lnk called' |
---|
[11081] | 1541 | DO ji = 1, n_sequence_lbc - 1 |
---|
| 1542 | IF ( crname_lbc(ji) /= 'already counted' ) THEN |
---|
| 1543 | ccountname = crname_lbc(ji) |
---|
| 1544 | crname_lbc(ji) = 'already counted' |
---|
| 1545 | jcount = 1 |
---|
| 1546 | DO jj = ji + 1, n_sequence_lbc |
---|
| 1547 | IF ( ccountname == crname_lbc(jj) ) THEN |
---|
| 1548 | jcount = jcount + 1 |
---|
| 1549 | crname_lbc(jj) = 'already counted' |
---|
| 1550 | END IF |
---|
| 1551 | END DO |
---|
| 1552 | WRITE(numcom,'(A, I4, A, A)') ' - ', jcount,' times by subroutine ', TRIM(ccountname) |
---|
[10425] | 1553 | END IF |
---|
| 1554 | END DO |
---|
[11081] | 1555 | IF ( crname_lbc(n_sequence_lbc) /= 'already counted' ) THEN |
---|
| 1556 | WRITE(numcom,'(A, I4, A, A)') ' - ', 1,' times by subroutine ', TRIM(crname_lbc(ncom_rec_max)) |
---|
| 1557 | END IF |
---|
[10425] | 1558 | WRITE(numcom,*) ' ' |
---|
| 1559 | IF ( n_sequence_glb > 0 ) THEN |
---|
| 1560 | WRITE(numcom,'(A,I4)') ' Global communications : ', n_sequence_glb |
---|
| 1561 | jj = 1 |
---|
| 1562 | DO ji = 2, n_sequence_glb |
---|
| 1563 | IF( crname_glb(ji-1) /= crname_glb(ji) ) THEN |
---|
| 1564 | WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_glb(ji-1)) |
---|
| 1565 | jj = 0 |
---|
| 1566 | END IF |
---|
| 1567 | jj = jj + 1 |
---|
| 1568 | END DO |
---|
| 1569 | WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_glb(n_sequence_glb)) |
---|
| 1570 | DEALLOCATE(crname_glb) |
---|
| 1571 | ELSE |
---|
| 1572 | WRITE(numcom,*) ' No MPI global communication ' |
---|
| 1573 | ENDIF |
---|
| 1574 | WRITE(numcom,*) ' ' |
---|
[10437] | 1575 | IF ( n_sequence_dlg > 0 ) THEN |
---|
| 1576 | WRITE(numcom,'(A,I4)') ' Delayed global communications : ', n_sequence_dlg |
---|
| 1577 | jj = 1 |
---|
| 1578 | DO ji = 2, n_sequence_dlg |
---|
| 1579 | IF( crname_dlg(ji-1) /= crname_dlg(ji) ) THEN |
---|
| 1580 | WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_dlg(ji-1)) |
---|
| 1581 | jj = 0 |
---|
| 1582 | END IF |
---|
| 1583 | jj = jj + 1 |
---|
| 1584 | END DO |
---|
| 1585 | WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_dlg(n_sequence_dlg)) |
---|
| 1586 | DEALLOCATE(crname_dlg) |
---|
| 1587 | ELSE |
---|
| 1588 | WRITE(numcom,*) ' No MPI delayed global communication ' |
---|
| 1589 | ENDIF |
---|
| 1590 | WRITE(numcom,*) ' ' |
---|
[10425] | 1591 | WRITE(numcom,*) ' -----------------------------------------------' |
---|
| 1592 | WRITE(numcom,*) ' ' |
---|
| 1593 | DEALLOCATE(ncomm_sequence) |
---|
| 1594 | DEALLOCATE(crname_lbc) |
---|
| 1595 | ENDIF |
---|
| 1596 | END SUBROUTINE mpp_report |
---|
| 1597 | |
---|
[6140] | 1598 | |
---|
[10425] | 1599 | SUBROUTINE tic_tac (ld_tic, ld_global) |
---|
| 1600 | |
---|
| 1601 | LOGICAL, INTENT(IN) :: ld_tic |
---|
| 1602 | LOGICAL, OPTIONAL, INTENT(IN) :: ld_global |
---|
| 1603 | REAL(wp), DIMENSION(2), SAVE :: tic_wt |
---|
| 1604 | REAL(wp), SAVE :: tic_ct = 0._wp |
---|
| 1605 | INTEGER :: ii |
---|
| 1606 | |
---|
| 1607 | IF( ncom_stp <= nit000 ) RETURN |
---|
| 1608 | IF( ncom_stp == nitend ) RETURN |
---|
| 1609 | ii = 1 |
---|
| 1610 | IF( PRESENT( ld_global ) ) THEN |
---|
| 1611 | IF( ld_global ) ii = 2 |
---|
| 1612 | END IF |
---|
| 1613 | |
---|
| 1614 | IF ( ld_tic ) THEN |
---|
| 1615 | tic_wt(ii) = MPI_Wtime() ! start count tic->tac (waiting time) |
---|
| 1616 | IF ( tic_ct > 0.0_wp ) compute_time = compute_time + MPI_Wtime() - tic_ct ! cumulate count tac->tic |
---|
| 1617 | ELSE |
---|
| 1618 | waiting_time(ii) = waiting_time(ii) + MPI_Wtime() - tic_wt(ii) ! cumulate count tic->tac |
---|
| 1619 | tic_ct = MPI_Wtime() ! start count tac->tic (waiting time) |
---|
| 1620 | ENDIF |
---|
| 1621 | |
---|
| 1622 | END SUBROUTINE tic_tac |
---|
| 1623 | |
---|
| 1624 | |
---|
[13] | 1625 | #else |
---|
| 1626 | !!---------------------------------------------------------------------- |
---|
| 1627 | !! Default case: Dummy module share memory computing |
---|
| 1628 | !!---------------------------------------------------------------------- |
---|
[2715] | 1629 | USE in_out_manager |
---|
[1976] | 1630 | |
---|
[13] | 1631 | INTERFACE mpp_sum |
---|
[10425] | 1632 | MODULE PROCEDURE mppsum_int, mppsum_a_int, mppsum_real, mppsum_a_real, mppsum_realdd, mppsum_a_realdd |
---|
[13] | 1633 | END INTERFACE |
---|
| 1634 | INTERFACE mpp_max |
---|
[681] | 1635 | MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real |
---|
[13] | 1636 | END INTERFACE |
---|
| 1637 | INTERFACE mpp_min |
---|
| 1638 | MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real |
---|
| 1639 | END INTERFACE |
---|
[1344] | 1640 | INTERFACE mpp_minloc |
---|
| 1641 | MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d |
---|
| 1642 | END INTERFACE |
---|
| 1643 | INTERFACE mpp_maxloc |
---|
| 1644 | MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d |
---|
| 1645 | END INTERFACE |
---|
[3] | 1646 | |
---|
[13] | 1647 | LOGICAL, PUBLIC, PARAMETER :: lk_mpp = .FALSE. !: mpp flag |
---|
[4147] | 1648 | LOGICAL, PUBLIC :: ln_nnogather !: namelist control of northfold comms (needed here in case "key_mpp_mpi" is not used) |
---|
[9570] | 1649 | INTEGER, PUBLIC :: mpi_comm_oce ! opa local communicator |
---|
[10425] | 1650 | |
---|
| 1651 | INTEGER, PARAMETER, PUBLIC :: nbdelay = 0 ! make sure we don't enter loops: DO ji = 1, nbdelay |
---|
| 1652 | CHARACTER(len=32), DIMENSION(1), PUBLIC :: c_delaylist = 'empty' |
---|
| 1653 | CHARACTER(len=32), DIMENSION(1), PUBLIC :: c_delaycpnt = 'empty' |
---|
[10528] | 1654 | LOGICAL, PUBLIC :: l_full_nf_update = .TRUE. |
---|
[10425] | 1655 | TYPE :: DELAYARR |
---|
| 1656 | REAL( wp), POINTER, DIMENSION(:) :: z1d => NULL() |
---|
| 1657 | COMPLEX(wp), POINTER, DIMENSION(:) :: y1d => NULL() |
---|
| 1658 | END TYPE DELAYARR |
---|
| 1659 | TYPE( DELAYARR ), DIMENSION(1), PUBLIC :: todelay |
---|
| 1660 | INTEGER, PUBLIC, DIMENSION(1) :: ndelayid = -1 |
---|
[2715] | 1661 | !!---------------------------------------------------------------------- |
---|
[13] | 1662 | CONTAINS |
---|
[3] | 1663 | |
---|
[2715] | 1664 | INTEGER FUNCTION lib_mpp_alloc(kumout) ! Dummy function |
---|
| 1665 | INTEGER, INTENT(in) :: kumout |
---|
| 1666 | lib_mpp_alloc = 0 |
---|
| 1667 | END FUNCTION lib_mpp_alloc |
---|
| 1668 | |
---|
[5407] | 1669 | FUNCTION mynode( ldtxt, ldname, kumnam_ref, knumnam_cfg, kumond , kstop, localComm ) RESULT (function_value) |
---|
[1579] | 1670 | INTEGER, OPTIONAL , INTENT(in ) :: localComm |
---|
[3764] | 1671 | CHARACTER(len=*),DIMENSION(:) :: ldtxt |
---|
[5407] | 1672 | CHARACTER(len=*) :: ldname |
---|
[4314] | 1673 | INTEGER :: kumnam_ref, knumnam_cfg , kumond , kstop |
---|
[9570] | 1674 | IF( PRESENT( localComm ) ) mpi_comm_oce = localComm |
---|
[5412] | 1675 | function_value = 0 |
---|
[1579] | 1676 | IF( .FALSE. ) ldtxt(:) = 'never done' |
---|
[5407] | 1677 | CALL ctl_opn( kumond, TRIM(ldname), 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. , 1 ) |
---|
[13] | 1678 | END FUNCTION mynode |
---|
[3] | 1679 | |
---|
[13] | 1680 | SUBROUTINE mppsync ! Dummy routine |
---|
| 1681 | END SUBROUTINE mppsync |
---|
[3] | 1682 | |
---|
[10425] | 1683 | !!---------------------------------------------------------------------- |
---|
| 1684 | !! *** mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real *** |
---|
| 1685 | !! |
---|
| 1686 | !!---------------------------------------------------------------------- |
---|
| 1687 | !! |
---|
| 1688 | # define OPERATION_MAX |
---|
| 1689 | # define INTEGER_TYPE |
---|
| 1690 | # define DIM_0d |
---|
| 1691 | # define ROUTINE_ALLREDUCE mppmax_int |
---|
| 1692 | # include "mpp_allreduce_generic.h90" |
---|
| 1693 | # undef ROUTINE_ALLREDUCE |
---|
| 1694 | # undef DIM_0d |
---|
| 1695 | # define DIM_1d |
---|
| 1696 | # define ROUTINE_ALLREDUCE mppmax_a_int |
---|
| 1697 | # include "mpp_allreduce_generic.h90" |
---|
| 1698 | # undef ROUTINE_ALLREDUCE |
---|
| 1699 | # undef DIM_1d |
---|
| 1700 | # undef INTEGER_TYPE |
---|
| 1701 | ! |
---|
| 1702 | # define REAL_TYPE |
---|
| 1703 | # define DIM_0d |
---|
| 1704 | # define ROUTINE_ALLREDUCE mppmax_real |
---|
| 1705 | # include "mpp_allreduce_generic.h90" |
---|
| 1706 | # undef ROUTINE_ALLREDUCE |
---|
| 1707 | # undef DIM_0d |
---|
| 1708 | # define DIM_1d |
---|
| 1709 | # define ROUTINE_ALLREDUCE mppmax_a_real |
---|
| 1710 | # include "mpp_allreduce_generic.h90" |
---|
| 1711 | # undef ROUTINE_ALLREDUCE |
---|
| 1712 | # undef DIM_1d |
---|
| 1713 | # undef REAL_TYPE |
---|
| 1714 | # undef OPERATION_MAX |
---|
| 1715 | !!---------------------------------------------------------------------- |
---|
| 1716 | !! *** mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real *** |
---|
| 1717 | !! |
---|
| 1718 | !!---------------------------------------------------------------------- |
---|
| 1719 | !! |
---|
| 1720 | # define OPERATION_MIN |
---|
| 1721 | # define INTEGER_TYPE |
---|
| 1722 | # define DIM_0d |
---|
| 1723 | # define ROUTINE_ALLREDUCE mppmin_int |
---|
| 1724 | # include "mpp_allreduce_generic.h90" |
---|
| 1725 | # undef ROUTINE_ALLREDUCE |
---|
| 1726 | # undef DIM_0d |
---|
| 1727 | # define DIM_1d |
---|
| 1728 | # define ROUTINE_ALLREDUCE mppmin_a_int |
---|
| 1729 | # include "mpp_allreduce_generic.h90" |
---|
| 1730 | # undef ROUTINE_ALLREDUCE |
---|
| 1731 | # undef DIM_1d |
---|
| 1732 | # undef INTEGER_TYPE |
---|
| 1733 | ! |
---|
| 1734 | # define REAL_TYPE |
---|
| 1735 | # define DIM_0d |
---|
| 1736 | # define ROUTINE_ALLREDUCE mppmin_real |
---|
| 1737 | # include "mpp_allreduce_generic.h90" |
---|
| 1738 | # undef ROUTINE_ALLREDUCE |
---|
| 1739 | # undef DIM_0d |
---|
| 1740 | # define DIM_1d |
---|
| 1741 | # define ROUTINE_ALLREDUCE mppmin_a_real |
---|
| 1742 | # include "mpp_allreduce_generic.h90" |
---|
| 1743 | # undef ROUTINE_ALLREDUCE |
---|
| 1744 | # undef DIM_1d |
---|
| 1745 | # undef REAL_TYPE |
---|
| 1746 | # undef OPERATION_MIN |
---|
[3] | 1747 | |
---|
[10425] | 1748 | !!---------------------------------------------------------------------- |
---|
| 1749 | !! *** mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real *** |
---|
| 1750 | !! |
---|
| 1751 | !! Global sum of 1D array or a variable (integer, real or complex) |
---|
| 1752 | !!---------------------------------------------------------------------- |
---|
| 1753 | !! |
---|
| 1754 | # define OPERATION_SUM |
---|
| 1755 | # define INTEGER_TYPE |
---|
| 1756 | # define DIM_0d |
---|
| 1757 | # define ROUTINE_ALLREDUCE mppsum_int |
---|
| 1758 | # include "mpp_allreduce_generic.h90" |
---|
| 1759 | # undef ROUTINE_ALLREDUCE |
---|
| 1760 | # undef DIM_0d |
---|
| 1761 | # define DIM_1d |
---|
| 1762 | # define ROUTINE_ALLREDUCE mppsum_a_int |
---|
| 1763 | # include "mpp_allreduce_generic.h90" |
---|
| 1764 | # undef ROUTINE_ALLREDUCE |
---|
| 1765 | # undef DIM_1d |
---|
| 1766 | # undef INTEGER_TYPE |
---|
| 1767 | ! |
---|
| 1768 | # define REAL_TYPE |
---|
| 1769 | # define DIM_0d |
---|
| 1770 | # define ROUTINE_ALLREDUCE mppsum_real |
---|
| 1771 | # include "mpp_allreduce_generic.h90" |
---|
| 1772 | # undef ROUTINE_ALLREDUCE |
---|
| 1773 | # undef DIM_0d |
---|
| 1774 | # define DIM_1d |
---|
| 1775 | # define ROUTINE_ALLREDUCE mppsum_a_real |
---|
| 1776 | # include "mpp_allreduce_generic.h90" |
---|
| 1777 | # undef ROUTINE_ALLREDUCE |
---|
| 1778 | # undef DIM_1d |
---|
| 1779 | # undef REAL_TYPE |
---|
| 1780 | # undef OPERATION_SUM |
---|
[3] | 1781 | |
---|
[10425] | 1782 | # define OPERATION_SUM_DD |
---|
| 1783 | # define COMPLEX_TYPE |
---|
| 1784 | # define DIM_0d |
---|
| 1785 | # define ROUTINE_ALLREDUCE mppsum_realdd |
---|
| 1786 | # include "mpp_allreduce_generic.h90" |
---|
| 1787 | # undef ROUTINE_ALLREDUCE |
---|
| 1788 | # undef DIM_0d |
---|
| 1789 | # define DIM_1d |
---|
| 1790 | # define ROUTINE_ALLREDUCE mppsum_a_realdd |
---|
| 1791 | # include "mpp_allreduce_generic.h90" |
---|
| 1792 | # undef ROUTINE_ALLREDUCE |
---|
| 1793 | # undef DIM_1d |
---|
| 1794 | # undef COMPLEX_TYPE |
---|
| 1795 | # undef OPERATION_SUM_DD |
---|
[3] | 1796 | |
---|
[10425] | 1797 | !!---------------------------------------------------------------------- |
---|
| 1798 | !! *** mpp_minloc2d, mpp_minloc3d, mpp_maxloc2d, mpp_maxloc3d |
---|
| 1799 | !! |
---|
| 1800 | !!---------------------------------------------------------------------- |
---|
| 1801 | !! |
---|
| 1802 | # define OPERATION_MINLOC |
---|
| 1803 | # define DIM_2d |
---|
| 1804 | # define ROUTINE_LOC mpp_minloc2d |
---|
| 1805 | # include "mpp_loc_generic.h90" |
---|
| 1806 | # undef ROUTINE_LOC |
---|
| 1807 | # undef DIM_2d |
---|
| 1808 | # define DIM_3d |
---|
| 1809 | # define ROUTINE_LOC mpp_minloc3d |
---|
| 1810 | # include "mpp_loc_generic.h90" |
---|
| 1811 | # undef ROUTINE_LOC |
---|
| 1812 | # undef DIM_3d |
---|
| 1813 | # undef OPERATION_MINLOC |
---|
[2480] | 1814 | |
---|
[10425] | 1815 | # define OPERATION_MAXLOC |
---|
| 1816 | # define DIM_2d |
---|
| 1817 | # define ROUTINE_LOC mpp_maxloc2d |
---|
| 1818 | # include "mpp_loc_generic.h90" |
---|
| 1819 | # undef ROUTINE_LOC |
---|
| 1820 | # undef DIM_2d |
---|
| 1821 | # define DIM_3d |
---|
| 1822 | # define ROUTINE_LOC mpp_maxloc3d |
---|
| 1823 | # include "mpp_loc_generic.h90" |
---|
| 1824 | # undef ROUTINE_LOC |
---|
| 1825 | # undef DIM_3d |
---|
| 1826 | # undef OPERATION_MAXLOC |
---|
[13] | 1827 | |
---|
[10425] | 1828 | SUBROUTINE mpp_delay_sum( cdname, cdelay, y_in, pout, ldlast, kcom ) |
---|
| 1829 | CHARACTER(len=*), INTENT(in ) :: cdname ! name of the calling subroutine |
---|
| 1830 | CHARACTER(len=*), INTENT(in ) :: cdelay ! name (used as id) of the delayed operation |
---|
| 1831 | COMPLEX(wp), INTENT(in ), DIMENSION(:) :: y_in |
---|
| 1832 | REAL(wp), INTENT( out), DIMENSION(:) :: pout |
---|
| 1833 | LOGICAL, INTENT(in ) :: ldlast ! true if this is the last time we call this routine |
---|
| 1834 | INTEGER, INTENT(in ), OPTIONAL :: kcom |
---|
| 1835 | ! |
---|
| 1836 | pout(:) = REAL(y_in(:), wp) |
---|
| 1837 | END SUBROUTINE mpp_delay_sum |
---|
[3764] | 1838 | |
---|
[10425] | 1839 | SUBROUTINE mpp_delay_max( cdname, cdelay, p_in, pout, ldlast, kcom ) |
---|
| 1840 | CHARACTER(len=*), INTENT(in ) :: cdname ! name of the calling subroutine |
---|
| 1841 | CHARACTER(len=*), INTENT(in ) :: cdelay ! name (used as id) of the delayed operation |
---|
| 1842 | REAL(wp), INTENT(in ), DIMENSION(:) :: p_in |
---|
| 1843 | REAL(wp), INTENT( out), DIMENSION(:) :: pout |
---|
| 1844 | LOGICAL, INTENT(in ) :: ldlast ! true if this is the last time we call this routine |
---|
| 1845 | INTEGER, INTENT(in ), OPTIONAL :: kcom |
---|
| 1846 | ! |
---|
| 1847 | pout(:) = p_in(:) |
---|
| 1848 | END SUBROUTINE mpp_delay_max |
---|
[3294] | 1849 | |
---|
[10425] | 1850 | SUBROUTINE mpp_delay_rcv( kid ) |
---|
| 1851 | INTEGER,INTENT(in ) :: kid |
---|
| 1852 | WRITE(*,*) 'mpp_delay_rcv: You should not have seen this print! error?', kid |
---|
| 1853 | END SUBROUTINE mpp_delay_rcv |
---|
| 1854 | |
---|
| 1855 | SUBROUTINE mppstop( ldfinal, ld_force_abort ) |
---|
| 1856 | LOGICAL, OPTIONAL, INTENT(in) :: ldfinal ! source process number |
---|
| 1857 | LOGICAL, OPTIONAL, INTENT(in) :: ld_force_abort ! source process number |
---|
[3799] | 1858 | STOP ! non MPP case, just stop the run |
---|
[51] | 1859 | END SUBROUTINE mppstop |
---|
| 1860 | |
---|
[2715] | 1861 | SUBROUTINE mpp_ini_znl( knum ) |
---|
| 1862 | INTEGER :: knum |
---|
| 1863 | WRITE(*,*) 'mpp_ini_znl: You should not have seen this print! error?', knum |
---|
[1345] | 1864 | END SUBROUTINE mpp_ini_znl |
---|
| 1865 | |
---|
[1344] | 1866 | SUBROUTINE mpp_comm_free( kcom ) |
---|
[869] | 1867 | INTEGER :: kcom |
---|
[1344] | 1868 | WRITE(*,*) 'mpp_comm_free: You should not have seen this print! error?', kcom |
---|
[869] | 1869 | END SUBROUTINE mpp_comm_free |
---|
[9019] | 1870 | |
---|
[3] | 1871 | #endif |
---|
[2715] | 1872 | |
---|
[13] | 1873 | !!---------------------------------------------------------------------- |
---|
[4147] | 1874 | !! All cases: ctl_stop, ctl_warn, get_unit, ctl_opn, ctl_nam routines |
---|
[2715] | 1875 | !!---------------------------------------------------------------------- |
---|
| 1876 | |
---|
| 1877 | SUBROUTINE ctl_stop( cd1, cd2, cd3, cd4, cd5 , & |
---|
| 1878 | & cd6, cd7, cd8, cd9, cd10 ) |
---|
| 1879 | !!---------------------------------------------------------------------- |
---|
| 1880 | !! *** ROUTINE stop_opa *** |
---|
| 1881 | !! |
---|
[3764] | 1882 | !! ** Purpose : print in ocean.outpput file a error message and |
---|
[2715] | 1883 | !! increment the error number (nstop) by one. |
---|
| 1884 | !!---------------------------------------------------------------------- |
---|
| 1885 | CHARACTER(len=*), INTENT(in), OPTIONAL :: cd1, cd2, cd3, cd4, cd5 |
---|
| 1886 | CHARACTER(len=*), INTENT(in), OPTIONAL :: cd6, cd7, cd8, cd9, cd10 |
---|
| 1887 | !!---------------------------------------------------------------------- |
---|
| 1888 | ! |
---|
[3764] | 1889 | nstop = nstop + 1 |
---|
[10425] | 1890 | |
---|
| 1891 | ! force to open ocean.output file |
---|
| 1892 | IF( numout == 6 ) CALL ctl_opn( numout, 'ocean.output', 'APPEND', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. ) |
---|
| 1893 | |
---|
| 1894 | WRITE(numout,cform_err) |
---|
| 1895 | IF( PRESENT(cd1 ) ) WRITE(numout,*) TRIM(cd1) |
---|
| 1896 | IF( PRESENT(cd2 ) ) WRITE(numout,*) TRIM(cd2) |
---|
| 1897 | IF( PRESENT(cd3 ) ) WRITE(numout,*) TRIM(cd3) |
---|
| 1898 | IF( PRESENT(cd4 ) ) WRITE(numout,*) TRIM(cd4) |
---|
| 1899 | IF( PRESENT(cd5 ) ) WRITE(numout,*) TRIM(cd5) |
---|
| 1900 | IF( PRESENT(cd6 ) ) WRITE(numout,*) TRIM(cd6) |
---|
| 1901 | IF( PRESENT(cd7 ) ) WRITE(numout,*) TRIM(cd7) |
---|
| 1902 | IF( PRESENT(cd8 ) ) WRITE(numout,*) TRIM(cd8) |
---|
| 1903 | IF( PRESENT(cd9 ) ) WRITE(numout,*) TRIM(cd9) |
---|
| 1904 | IF( PRESENT(cd10) ) WRITE(numout,*) TRIM(cd10) |
---|
| 1905 | |
---|
[2715] | 1906 | CALL FLUSH(numout ) |
---|
| 1907 | IF( numstp /= -1 ) CALL FLUSH(numstp ) |
---|
[9019] | 1908 | IF( numrun /= -1 ) CALL FLUSH(numrun ) |
---|
[2715] | 1909 | IF( numevo_ice /= -1 ) CALL FLUSH(numevo_ice) |
---|
| 1910 | ! |
---|
| 1911 | IF( cd1 == 'STOP' ) THEN |
---|
[10425] | 1912 | WRITE(numout,*) 'huge E-R-R-O-R : immediate stop' |
---|
| 1913 | CALL mppstop(ld_force_abort = .true.) |
---|
[2715] | 1914 | ENDIF |
---|
| 1915 | ! |
---|
| 1916 | END SUBROUTINE ctl_stop |
---|
| 1917 | |
---|
| 1918 | |
---|
| 1919 | SUBROUTINE ctl_warn( cd1, cd2, cd3, cd4, cd5, & |
---|
| 1920 | & cd6, cd7, cd8, cd9, cd10 ) |
---|
| 1921 | !!---------------------------------------------------------------------- |
---|
| 1922 | !! *** ROUTINE stop_warn *** |
---|
| 1923 | !! |
---|
[3764] | 1924 | !! ** Purpose : print in ocean.outpput file a error message and |
---|
[2715] | 1925 | !! increment the warning number (nwarn) by one. |
---|
| 1926 | !!---------------------------------------------------------------------- |
---|
| 1927 | CHARACTER(len=*), INTENT(in), OPTIONAL :: cd1, cd2, cd3, cd4, cd5 |
---|
| 1928 | CHARACTER(len=*), INTENT(in), OPTIONAL :: cd6, cd7, cd8, cd9, cd10 |
---|
| 1929 | !!---------------------------------------------------------------------- |
---|
[3764] | 1930 | ! |
---|
| 1931 | nwarn = nwarn + 1 |
---|
[2715] | 1932 | IF(lwp) THEN |
---|
| 1933 | WRITE(numout,cform_war) |
---|
[10425] | 1934 | IF( PRESENT(cd1 ) ) WRITE(numout,*) TRIM(cd1) |
---|
| 1935 | IF( PRESENT(cd2 ) ) WRITE(numout,*) TRIM(cd2) |
---|
| 1936 | IF( PRESENT(cd3 ) ) WRITE(numout,*) TRIM(cd3) |
---|
| 1937 | IF( PRESENT(cd4 ) ) WRITE(numout,*) TRIM(cd4) |
---|
| 1938 | IF( PRESENT(cd5 ) ) WRITE(numout,*) TRIM(cd5) |
---|
| 1939 | IF( PRESENT(cd6 ) ) WRITE(numout,*) TRIM(cd6) |
---|
| 1940 | IF( PRESENT(cd7 ) ) WRITE(numout,*) TRIM(cd7) |
---|
| 1941 | IF( PRESENT(cd8 ) ) WRITE(numout,*) TRIM(cd8) |
---|
| 1942 | IF( PRESENT(cd9 ) ) WRITE(numout,*) TRIM(cd9) |
---|
| 1943 | IF( PRESENT(cd10) ) WRITE(numout,*) TRIM(cd10) |
---|
[2715] | 1944 | ENDIF |
---|
| 1945 | CALL FLUSH(numout) |
---|
| 1946 | ! |
---|
| 1947 | END SUBROUTINE ctl_warn |
---|
| 1948 | |
---|
| 1949 | |
---|
| 1950 | SUBROUTINE ctl_opn( knum, cdfile, cdstat, cdform, cdacce, klengh, kout, ldwp, karea ) |
---|
| 1951 | !!---------------------------------------------------------------------- |
---|
| 1952 | !! *** ROUTINE ctl_opn *** |
---|
| 1953 | !! |
---|
| 1954 | !! ** Purpose : Open file and check if required file is available. |
---|
| 1955 | !! |
---|
| 1956 | !! ** Method : Fortan open |
---|
| 1957 | !!---------------------------------------------------------------------- |
---|
| 1958 | INTEGER , INTENT( out) :: knum ! logical unit to open |
---|
| 1959 | CHARACTER(len=*) , INTENT(in ) :: cdfile ! file name to open |
---|
| 1960 | CHARACTER(len=*) , INTENT(in ) :: cdstat ! disposition specifier |
---|
| 1961 | CHARACTER(len=*) , INTENT(in ) :: cdform ! formatting specifier |
---|
| 1962 | CHARACTER(len=*) , INTENT(in ) :: cdacce ! access specifier |
---|
| 1963 | INTEGER , INTENT(in ) :: klengh ! record length |
---|
| 1964 | INTEGER , INTENT(in ) :: kout ! number of logical units for write |
---|
| 1965 | LOGICAL , INTENT(in ) :: ldwp ! boolean term for print |
---|
| 1966 | INTEGER, OPTIONAL, INTENT(in ) :: karea ! proc number |
---|
[5836] | 1967 | ! |
---|
[2715] | 1968 | CHARACTER(len=80) :: clfile |
---|
| 1969 | INTEGER :: iost |
---|
| 1970 | !!---------------------------------------------------------------------- |
---|
[5836] | 1971 | ! |
---|
[2715] | 1972 | ! adapt filename |
---|
| 1973 | ! ---------------- |
---|
| 1974 | clfile = TRIM(cdfile) |
---|
| 1975 | IF( PRESENT( karea ) ) THEN |
---|
| 1976 | IF( karea > 1 ) WRITE(clfile, "(a,'_',i4.4)") TRIM(clfile), karea-1 |
---|
| 1977 | ENDIF |
---|
| 1978 | #if defined key_agrif |
---|
| 1979 | IF( .NOT. Agrif_Root() ) clfile = TRIM(Agrif_CFixed())//'_'//TRIM(clfile) |
---|
| 1980 | knum=Agrif_Get_Unit() |
---|
| 1981 | #else |
---|
| 1982 | knum=get_unit() |
---|
| 1983 | #endif |
---|
[10425] | 1984 | IF( TRIM(cdfile) == '/dev/null' ) clfile = TRIM(cdfile) ! force the use of /dev/null |
---|
[5836] | 1985 | ! |
---|
[2715] | 1986 | iost=0 |
---|
[10425] | 1987 | IF( cdacce(1:6) == 'DIRECT' ) THEN ! cdacce has always more than 6 characters |
---|
| 1988 | OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat, RECL=klengh , ERR=100, IOSTAT=iost ) |
---|
| 1989 | ELSE IF( TRIM(cdstat) == 'APPEND' ) THEN ! cdstat can have less than 6 characters |
---|
| 1990 | OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS='UNKNOWN', POSITION='APPEND', ERR=100, IOSTAT=iost ) |
---|
[2715] | 1991 | ELSE |
---|
[10425] | 1992 | OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat , ERR=100, IOSTAT=iost ) |
---|
[2715] | 1993 | ENDIF |
---|
[10425] | 1994 | IF( iost /= 0 .AND. TRIM(clfile) == '/dev/null' ) & ! for windows |
---|
| 1995 | & OPEN(UNIT=knum,FILE='NUL', FORM=cdform, ACCESS=cdacce, STATUS=cdstat , ERR=100, IOSTAT=iost ) |
---|
[2715] | 1996 | IF( iost == 0 ) THEN |
---|
| 1997 | IF(ldwp) THEN |
---|
[10425] | 1998 | WRITE(kout,*) ' file : ', TRIM(clfile),' open ok' |
---|
[2715] | 1999 | WRITE(kout,*) ' unit = ', knum |
---|
| 2000 | WRITE(kout,*) ' status = ', cdstat |
---|
| 2001 | WRITE(kout,*) ' form = ', cdform |
---|
| 2002 | WRITE(kout,*) ' access = ', cdacce |
---|
| 2003 | WRITE(kout,*) |
---|
| 2004 | ENDIF |
---|
| 2005 | ENDIF |
---|
| 2006 | 100 CONTINUE |
---|
| 2007 | IF( iost /= 0 ) THEN |
---|
| 2008 | IF(ldwp) THEN |
---|
| 2009 | WRITE(kout,*) |
---|
[10425] | 2010 | WRITE(kout,*) ' ===>>>> : bad opening file: ', TRIM(clfile) |
---|
[2715] | 2011 | WRITE(kout,*) ' ======= === ' |
---|
| 2012 | WRITE(kout,*) ' unit = ', knum |
---|
| 2013 | WRITE(kout,*) ' status = ', cdstat |
---|
| 2014 | WRITE(kout,*) ' form = ', cdform |
---|
| 2015 | WRITE(kout,*) ' access = ', cdacce |
---|
| 2016 | WRITE(kout,*) ' iostat = ', iost |
---|
| 2017 | WRITE(kout,*) ' we stop. verify the file ' |
---|
| 2018 | WRITE(kout,*) |
---|
[9438] | 2019 | ELSE !!! Force writing to make sure we get the information - at least once - in this violent STOP!! |
---|
| 2020 | WRITE(*,*) |
---|
[10425] | 2021 | WRITE(*,*) ' ===>>>> : bad opening file: ', TRIM(clfile) |
---|
[9438] | 2022 | WRITE(*,*) ' ======= === ' |
---|
| 2023 | WRITE(*,*) ' unit = ', knum |
---|
| 2024 | WRITE(*,*) ' status = ', cdstat |
---|
| 2025 | WRITE(*,*) ' form = ', cdform |
---|
| 2026 | WRITE(*,*) ' access = ', cdacce |
---|
| 2027 | WRITE(*,*) ' iostat = ', iost |
---|
| 2028 | WRITE(*,*) ' we stop. verify the file ' |
---|
| 2029 | WRITE(*,*) |
---|
[2715] | 2030 | ENDIF |
---|
[9019] | 2031 | CALL FLUSH( kout ) |
---|
[2715] | 2032 | STOP 'ctl_opn bad opening' |
---|
| 2033 | ENDIF |
---|
[5836] | 2034 | ! |
---|
[2715] | 2035 | END SUBROUTINE ctl_opn |
---|
| 2036 | |
---|
[5836] | 2037 | |
---|
[4147] | 2038 | SUBROUTINE ctl_nam ( kios, cdnam, ldwp ) |
---|
| 2039 | !!---------------------------------------------------------------------- |
---|
| 2040 | !! *** ROUTINE ctl_nam *** |
---|
| 2041 | !! |
---|
| 2042 | !! ** Purpose : Informations when error while reading a namelist |
---|
| 2043 | !! |
---|
| 2044 | !! ** Method : Fortan open |
---|
| 2045 | !!---------------------------------------------------------------------- |
---|
[5836] | 2046 | INTEGER , INTENT(inout) :: kios ! IO status after reading the namelist |
---|
| 2047 | CHARACTER(len=*), INTENT(in ) :: cdnam ! group name of namelist for which error occurs |
---|
[7646] | 2048 | CHARACTER(len=5) :: clios ! string to convert iostat in character for print |
---|
[5836] | 2049 | LOGICAL , INTENT(in ) :: ldwp ! boolean term for print |
---|
[4147] | 2050 | !!---------------------------------------------------------------------- |
---|
[5836] | 2051 | ! |
---|
[7646] | 2052 | WRITE (clios, '(I5.0)') kios |
---|
[4147] | 2053 | IF( kios < 0 ) THEN |
---|
[5836] | 2054 | CALL ctl_warn( 'end of record or file while reading namelist ' & |
---|
| 2055 | & // TRIM(cdnam) // ' iostat = ' // TRIM(clios) ) |
---|
[4147] | 2056 | ENDIF |
---|
[5836] | 2057 | ! |
---|
[4147] | 2058 | IF( kios > 0 ) THEN |
---|
[5836] | 2059 | CALL ctl_stop( 'misspelled variable in namelist ' & |
---|
| 2060 | & // TRIM(cdnam) // ' iostat = ' // TRIM(clios) ) |
---|
[4147] | 2061 | ENDIF |
---|
| 2062 | kios = 0 |
---|
| 2063 | RETURN |
---|
[5836] | 2064 | ! |
---|
[4147] | 2065 | END SUBROUTINE ctl_nam |
---|
| 2066 | |
---|
[5836] | 2067 | |
---|
[2715] | 2068 | INTEGER FUNCTION get_unit() |
---|
| 2069 | !!---------------------------------------------------------------------- |
---|
| 2070 | !! *** FUNCTION get_unit *** |
---|
| 2071 | !! |
---|
| 2072 | !! ** Purpose : return the index of an unused logical unit |
---|
| 2073 | !!---------------------------------------------------------------------- |
---|
[3764] | 2074 | LOGICAL :: llopn |
---|
[2715] | 2075 | !!---------------------------------------------------------------------- |
---|
| 2076 | ! |
---|
| 2077 | get_unit = 15 ! choose a unit that is big enough then it is not already used in NEMO |
---|
| 2078 | llopn = .TRUE. |
---|
| 2079 | DO WHILE( (get_unit < 998) .AND. llopn ) |
---|
| 2080 | get_unit = get_unit + 1 |
---|
| 2081 | INQUIRE( unit = get_unit, opened = llopn ) |
---|
| 2082 | END DO |
---|
| 2083 | IF( (get_unit == 999) .AND. llopn ) THEN |
---|
| 2084 | CALL ctl_stop( 'get_unit: All logical units until 999 are used...' ) |
---|
| 2085 | get_unit = -1 |
---|
| 2086 | ENDIF |
---|
| 2087 | ! |
---|
| 2088 | END FUNCTION get_unit |
---|
| 2089 | |
---|
| 2090 | !!---------------------------------------------------------------------- |
---|
[3] | 2091 | END MODULE lib_mpp |
---|