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

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 474 for trunk/NEMO/OPA_SRC/DYN – NEMO

Changeset 474 for trunk/NEMO/OPA_SRC/DYN


Ignore:
Timestamp:
2006-05-11T17:24:19+02:00 (18 years ago)
Author:
opalod
Message:

nemo_v1_update_061: SM: end of ctl_stop + mpi optimization in _bilap

Location:
trunk/NEMO/OPA_SRC/DYN
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/DYN/dynhpg.F90

    r455 r474  
    203203      IF( ln_hpg_djc )   ioptio = ioptio + 1 
    204204      IF( ln_hpg_rot )   ioptio = ioptio + 1 
    205       IF ( ioptio > 1 ) THEN 
    206           IF(lwp) WRITE(numout,cform_err) 
    207           IF(lwp) WRITE(numout,*) ' several hydrostatic pressure gradient options used' 
    208           nstop = nstop + 1 
    209       ENDIF 
     205      IF ( ioptio > 1 )   & 
     206           &   CALL ctl_stop( ' several hydrostatic pressure gradient options used' ) 
    210207 
    211208      IF( lk_dynhpg_jki ) THEN 
  • trunk/NEMO/OPA_SRC/DYN/dynldf.F90

    r456 r474  
    152152      IF( ln_dynldf_lap   )   ioptio = ioptio + 1 
    153153      IF( ln_dynldf_bilap )   ioptio = ioptio + 1 
    154       IF( ioptio /= 1 )   THEN 
    155           IF(lwp) WRITE(numout,cform_err) 
    156           IF(lwp) WRITE(numout,*) '          use ONE of the 2 lap/bilap operator type on dynamics' 
    157           nstop = nstop + 1 
    158       ENDIF 
     154      IF( ioptio /= 1 ) CALL ctl_stop( '          use ONE of the 2 lap/bilap operator type on dynamics' ) 
    159155      ioptio = 0 
    160156      IF( ln_dynldf_level )   ioptio = ioptio + 1 
    161157      IF( ln_dynldf_hor   )   ioptio = ioptio + 1 
    162158      IF( ln_dynldf_iso   )   ioptio = ioptio + 1 
    163       IF( ioptio /= 1 ) THEN 
    164          IF(lwp) WRITE(numout,cform_err) 
    165          IF(lwp) WRITE(numout,*) '          use only ONE direction (level/hor/iso)' 
    166          nstop = nstop + 1 
    167       ENDIF 
     159      IF( ioptio /= 1 ) CALL ctl_stop( '          use only ONE direction (level/hor/iso)' ) 
    168160 
    169161      ! defined the type of lateral diffusion from ln_dynldf_... logicals 
     
    205197      ENDIF 
    206198 
    207       IF( ierr == 1 ) THEN 
    208          IF(lwp) WRITE(numout,cform_err) 
    209          IF(lwp) WRITE(numout,*) ' iso-level in z-coordinate - partial step, not allowed' 
    210          nstop = nstop + 1 
    211       ENDIF 
    212       IF( ierr == 2 ) THEN 
    213          IF(lwp) WRITE(numout,cform_err) 
    214          IF(lwp) WRITE(numout,*) ' isoneutral bilaplacian operator does not exist' 
    215          nstop = nstop + 1 
    216       ENDIF 
     199      IF( ierr == 1 )   & 
     200           &   CALL ctl_stop( ' iso-level in z-coordinate - partial step, not allowed' ) 
     201      IF( ierr == 2 )   & 
     202           &   CALL ctl_stop( ' isoneutral bilaplacian operator does not exist' ) 
    217203      IF( nldf == 1 .OR. nldf == 3 ) THEN      ! rotation 
    218          IF( .NOT.lk_ldfslp ) THEN 
    219             IF(lwp) WRITE(numout,cform_err) 
    220             IF(lwp) WRITE(numout,*) '          the rotation of the diffusive tensor require key_ldfslp' 
    221             nstop = nstop + 1 
    222          ENDIF 
     204         IF( .NOT.lk_ldfslp )   & 
     205           &   CALL ctl_stop(  '          the rotation of the diffusive tensor require key_ldfslp' ) 
    223206      ENDIF 
    224207 
  • trunk/NEMO/OPA_SRC/DYN/dynldf_bilap.F90

    r455 r474  
    8686      REAL(wp) ::   zua, zva, zbt, ze2u, ze2v ! temporary scalar 
    8787      REAL(wp), DIMENSION(jpi,jpj) ::   & 
    88          zuf, zut, zlu, zlv, zcu, zcv         ! temporary workspace 
     88         zcu, zcv                             ! temporary workspace 
     89      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   & 
     90         zuf, zut, zlu, zlv                   ! temporary workspace 
    8991      !!---------------------------------------------------------------------- 
    9092      !!  OPA 8.5, LODYC-IPSL (2002) 
     
    9698         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~' 
    9799      ENDIF 
    98       zuf(:,:) = 0.e0 
    99       zut(:,:) = 0.e0 
    100       zlu(:,:) = 0.e0 
    101       zlv(:,:) = 0.e0 
     100 
     101!!bug gm this should be enough 
     102!!$      zuf(:,:,jpk) = 0.e0 
     103!!$      zut(:,:,jpk) = 0.e0 
     104!!$      zlu(:,:,jpk) = 0.e0 
     105!!$      zlv(:,:,jpk) = 0.e0 
     106      zuf(:,:,:) = 0.e0 
     107      zut(:,:,:) = 0.e0 
     108      zlu(:,:,:) = 0.e0 
     109      zlv(:,:,:) = 0.e0 
    102110 
    103111      !                                                ! =============== 
     
    108116 
    109117         IF( ln_sco .OR. ln_zps ) THEN   ! s-coordinate or z-coordinate with partial steps 
    110             zuf(:,:) = rotb(:,:,jk) * fse3f(:,:,jk) 
     118            zuf(:,:,jk) = rotb(:,:,jk) * fse3f(:,:,jk) 
    111119            DO jj = 2, jpjm1 
    112120               DO ji = fs_2, fs_jpim1   ! vector opt. 
    113                   zlu(ji,jj) = - ( zuf(ji,jj) - zuf(ji,jj-1) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
     121                  zlu(ji,jj,jk) = - ( zuf(ji,jj,jk) - zuf(ji,jj-1,jk) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
    114122                     &         + ( hdivb(ji+1,jj,jk) - hdivb(ji,jj,jk) ) / e1u(ji,jj) 
    115123    
    116                   zlv(ji,jj) = + ( zuf(ji,jj) - zuf(ji-1,jj) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
     124                  zlv(ji,jj,jk) = + ( zuf(ji,jj,jk) - zuf(ji-1,jj,jk) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
    117125                     &         + ( hdivb(ji,jj+1,jk) - hdivb(ji,jj,jk) ) / e2v(ji,jj) 
    118126               END DO 
     
    121129            DO jj = 2, jpjm1 
    122130               DO ji = fs_2, fs_jpim1   ! vector opt. 
    123                   zlu(ji,jj) = - ( rotb (ji  ,jj,jk) - rotb (ji,jj-1,jk) ) / e2u(ji,jj)   & 
     131                  zlu(ji,jj,jk) = - ( rotb (ji  ,jj,jk) - rotb (ji,jj-1,jk) ) / e2u(ji,jj)   & 
    124132                     &         + ( hdivb(ji+1,jj,jk) - hdivb(ji,jj  ,jk) ) / e1u(ji,jj) 
    125133    
    126                   zlv(ji,jj) = + ( rotb (ji,jj  ,jk) - rotb (ji-1,jj,jk) ) / e1v(ji,jj)   & 
     134                  zlv(ji,jj,jk) = + ( rotb (ji,jj  ,jk) - rotb (ji-1,jj,jk) ) / e1v(ji,jj)   & 
    127135                     &         + ( hdivb(ji,jj+1,jk) - hdivb(ji  ,jj,jk) ) / e2v(ji,jj) 
    128136               END DO   
    129137            END DO   
    130138         ENDIF 
    131  
    132          ! Boundary conditions on the laplacian  (zlu,zlv) 
    133          CALL lbc_lnk( zlu, 'U', -1. ) 
    134          CALL lbc_lnk( zlv, 'V', -1. ) 
    135           
    136           
     139      ENDDO 
     140 
     141      ! Boundary conditions on the laplacian  (zlu,zlv) 
     142      CALL lbc_lnk( zlu, 'U', -1. ) 
     143      CALL lbc_lnk( zlv, 'V', -1. ) 
     144          
     145          
     146      DO jk = 1, jpkm1 
     147    
    137148         ! Third derivative 
    138149         ! ---------------- 
    139150          
    140151         ! Multiply by the eddy viscosity coef. (at u- and v-points) 
    141          zlu(:,:) = zlu(:,:) * fsahmu(:,:,jk) 
    142          zlv(:,:) = zlv(:,:) * fsahmv(:,:,jk) 
     152         zlu(:,:,jk) = zlu(:,:,jk) * fsahmu(:,:,jk) 
     153         zlv(:,:,jk) = zlv(:,:,jk) * fsahmv(:,:,jk) 
    143154          
    144155         ! Contravariant "laplacian" 
    145          zcu(:,:) = e1u(:,:) * zlu(:,:) 
    146          zcv(:,:) = e2v(:,:) * zlv(:,:) 
     156         zcu(:,:) = e1u(:,:) * zlu(:,:,jk) 
     157         zcv(:,:) = e2v(:,:) * zlv(:,:,jk) 
    147158          
    148159         ! Laplacian curl ( * e3f if s-coordinates or z-coordinate with partial steps) 
    149160         DO jj = 1, jpjm1 
    150161            DO ji = 1, fs_jpim1   ! vector opt. 
    151                zuf(ji,jj) = fmask(ji,jj,jk) * (  zcv(ji+1,jj  ) - zcv(ji,jj)      & 
     162               zuf(ji,jj,jk) = fmask(ji,jj,jk) * (  zcv(ji+1,jj  ) - zcv(ji,jj)      & 
    152163                  &                            - zcu(ji  ,jj+1) + zcu(ji,jj)  )   & 
    153164#if defined key_zco 
     
    163174            DO ji = 1, fs_jpim1   ! vector opt. 
    164175#if defined key_zco 
    165                zlu(ji,jj) = e2u(ji,jj) * zlu(ji,jj) 
    166                zlv(ji,jj) = e1v(ji,jj) * zlv(ji,jj) 
    167 #else 
    168                zlu(ji,jj) = e2u(ji,jj) * fse3u(ji,jj,jk) * zlu(ji,jj) 
    169                zlv(ji,jj) = e1v(ji,jj) * fse3v(ji,jj,jk) * zlv(ji,jj) 
     176               zlu(ji,jj,jk) = e2u(ji,jj) * zlu(ji,jj,jk) 
     177               zlv(ji,jj,jk) = e1v(ji,jj) * zlv(ji,jj,jk) 
     178#else 
     179               zlu(ji,jj,jk) = e2u(ji,jj) * fse3u(ji,jj,jk) * zlu(ji,jj,jk) 
     180               zlv(ji,jj,jk) = e1v(ji,jj) * fse3v(ji,jj,jk) * zlv(ji,jj,jk) 
    170181#endif 
    171182            END DO 
     
    180191               zbt = e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) 
    181192#endif 
    182                zut(ji,jj) = (  zlu(ji,jj) - zlu(ji-1,jj  )   & 
    183                   &          + zlv(ji,jj) - zlv(ji  ,jj-1) ) / zbt 
     193               zut(ji,jj,jk) = (  zlu(ji,jj,jk) - zlu(ji-1,jj  ,jk)   & 
     194                  &             + zlv(ji,jj,jk) - zlv(ji  ,jj-1,jk) ) / zbt 
    184195            END DO 
    185196         END DO 
     197      END DO 
    186198 
    187199 
    188200      ! boundary conditions on the laplacian curl and div (zuf,zut) 
     201!!bug gm no need to do this 2 following lbc... 
    189202      CALL lbc_lnk( zuf, 'F', 1. ) 
    190203      CALL lbc_lnk( zut, 'T', 1. ) 
    191204 
    192           
     205      DO jk = 1, jpkm1       
     206    
    193207         ! Bilaplacian 
    194208         ! ----------- 
     
    204218#endif 
    205219               ! horizontal biharmonic diffusive trends 
    206                zua = - ( zuf(ji  ,jj) - zuf(ji,jj-1) ) / ze2u   & 
    207                   &  + ( zut(ji+1,jj) - zut(ji,jj  ) ) / e1u(ji,jj) 
    208  
    209                zva = + ( zuf(ji,jj  ) - zuf(ji-1,jj) ) / ze2v   & 
    210                   &  + ( zut(ji,jj+1) - zut(ji  ,jj) ) / e2v(ji,jj) 
     220               zua = - ( zuf(ji  ,jj,jk) - zuf(ji,jj-1,jk) ) / ze2u   & 
     221                  &  + ( zut(ji+1,jj,jk) - zut(ji,jj  ,jk) ) / e1u(ji,jj) 
     222 
     223               zva = + ( zuf(ji,jj  ,jk) - zuf(ji-1,jj,jk) ) / ze2v   & 
     224                  &  + ( zut(ji,jj+1,jk) - zut(ji  ,jj,jk) ) / e2v(ji,jj) 
    211225               ! add it to the general momentum trends 
    212226               ua(ji,jj,jk) = ua(ji,jj,jk) + zua 
  • trunk/NEMO/OPA_SRC/DYN/dynspg.F90

    r372 r474  
    168168      IF(lk_dynspg_rl )   ioptio = ioptio + 1 
    169169 
    170       IF( ( ioptio > 1 .AND. .NOT. lk_esopa ) .OR. ioptio == 0 ) THEN 
    171          IF(lwp) WRITE(numout,cform_err) 
    172          IF(lwp) WRITE(numout,*) ' Choose only one surface pressure gradient scheme with a key cpp' 
    173          nstop = nstop + 1 
    174       ENDIF 
     170      IF( ( ioptio > 1 .AND. .NOT. lk_esopa ) .OR. ioptio == 0 )   & 
     171           &   CALL ctl_stop( ' Choose only one surface pressure gradient scheme with a key cpp' ) 
    175172 
    176173      IF( lk_esopa     )   nspg = -1 
     
    199196      ! -------------------------- 
    200197      IF( lk_dynspg_ts ) THEN 
    201          IF( MOD( rdt , rdtbt ) /= 0. ) THEN 
    202             IF(lwp) WRITE(numout,cform_err) 
    203             IF(lwp) WRITE(numout,*) ' The barotropic timestep must be an integer divisor of the baroclinic timestep' 
    204             nstop = nstop + 1 
    205          ENDIF 
     198         IF( MOD( rdt , rdtbt ) /= 0. )   & 
     199           &   CALL ctl_stop( ' The barotropic timestep must be an integer divisor of the baroclinic timestep' ) 
    206200      ENDIF 
    207201 
  • trunk/NEMO/OPA_SRC/DYN/dynspg_flt.F90

    r455 r474  
    286286            CALL sol_sor_e( kindic ) 
    287287         ELSE                          ! e r r o r in nsolv namelist parameter 
    288             IF(lwp) WRITE(numout,cform_err) 
    289             IF(lwp) WRITE(numout,*) ' dyn_spg_flt : e r r o r, nsolv = 1, 2 ,3 or 4' 
    290             IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~                not = ', nsolv 
    291             nstop = nstop + 1 
     288            WRITE(ctmp1,*) ' ~~~~~~~~~~~                not = ', nsolv 
     289            CALL ctl_stop( ' dyn_spg_flt : e r r o r, nsolv = 1, 2 ,3 or 4', ctmp1 ) 
    292290         ENDIF 
    293291      ENDIF 
  • trunk/NEMO/OPA_SRC/DYN/dynspg_flt_jki.F90

    r455 r474  
    300300            CALL sol_sor_e( kindic ) 
    301301         ELSE                          ! e r r o r in nsolv namelist parameter 
    302             IF(lwp) WRITE(numout,cform_err) 
    303             IF(lwp) WRITE(numout,*) ' dyn_spg_flt_jki : e r r o r, nsolv = 1, 2, 3 or 4' 
    304             IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~~~~~~                not = ', nsolv 
    305             nstop = nstop + 1 
    306          ENDIF 
     302            WRITE(ctmp1,*) ' ~~~~~~~~~~~~~~~~                not = ', nsolv 
     303            CALL ctl_stop( ' dyn_spg_flt_jki : e r r o r, nsolv = 1, 2, 3 or 4', ctmp1 ) 
    307304      ENDIF 
    308305 
  • trunk/NEMO/OPA_SRC/DYN/dynspg_rl.F90

    r359 r474  
    229229            CALL sol_sor_e( kindic ) 
    230230         CASE DEFAULT                  ! e r r o r in nsolv namelist parameter 
    231             IF(lwp) WRITE(numout,cform_err) 
    232             IF(lwp) WRITE(numout,*) ' dyn_spg_rl : e r r o r, nsolv = 1, 2 ,3 or 4' 
    233             IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~                not = ', nsolv 
    234             nstop = nstop + 1 
     231            WRITE(ctmp1,*) ' ~~~~~~~~~~                not = ', nsolv 
     232            CALL ctl_stop( ' dyn_spg_rl : e r r o r, nsolv = 1, 2 ,3 or 4', ctmp1 ) 
    235233         END SELECT 
    236234      ENDIF 
  • trunk/NEMO/OPA_SRC/DYN/dynvor.F90

    r455 r474  
    677677         ioptio = ioptio + 1 
    678678      ENDIF 
    679       IF( ioptio /= 1 .AND. .NOT. lk_esopa ) THEN 
    680          WRITE(numout,cform_err) 
    681          IF(lwp) WRITE(numout,*) ' use ONE and ONLY one vorticity scheme' 
    682          nstop = nstop + 1 
    683       ENDIF 
     679      IF( ioptio /= 1 .AND. .NOT. lk_esopa ) CALL ctl_stop( ' use ONE and ONLY one vorticity scheme' ) 
    684680      IF( lk_esopa ) THEN 
    685681         nvor = -1 
Note: See TracChangeset for help on using the changeset viewer.