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 8885 for branches/2017/dev_CNRS_2017/NEMOGCM/NEMO – NEMO

Ignore:
Timestamp:
2017-12-04T10:41:40+01:00 (6 years ago)
Author:
clem
Message:

remove useless references to clem's comments

Location:
branches/2017/dev_CNRS_2017/NEMOGCM/NEMO
Files:
15 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/LIM_SRC_3/icedyn_adv_umx.F90

    r8882 r8885  
    2424   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero) 
    2525   USE lbclnk         ! lateral boundary conditions (or mpp links) 
    26    USE timing         ! Timing 
    2726 
    2827   IMPLICIT NONE 
     
    261260      !!---------------------------------------------------------------------- 
    262261      ! 
    263 !!      IF( nn_timing == 1 )  CALL timing_start('macho') 
    264       ! 
    265262      IF( MOD( (kt - 1) / nn_fsbc , 2 ) == 0 ) THEN         !==  odd ice time step:  adv_x then adv_y  ==! 
    266263         ! 
     
    301298      ENDIF       
    302299      ! 
    303 !!      IF( nn_timing == 1 )  CALL timing_stop('macho') 
    304       ! 
    305300   END SUBROUTINE macho 
    306301 
     
    327322      REAL(wp), DIMENSION(jpi,jpj) :: ztu1, ztu2, ztu3, ztu4 
    328323      !!---------------------------------------------------------------------- 
    329       ! 
    330 !!      IF( nn_timing == 1 )  CALL timing_start('ultimate_x') 
    331324      ! 
    332325      !                                                     !--  Laplacian in i-direction  --! 
     
    423416      END SELECT 
    424417      ! 
    425 !!      IF( nn_timing == 1 )  CALL timing_stop('ultimate_x') 
    426       ! 
    427418   END SUBROUTINE ultimate_x 
    428419    
     
    449440      REAL(wp), DIMENSION(jpi,jpj) :: ztv1, ztv2, ztv3, ztv4 
    450441      !!---------------------------------------------------------------------- 
    451       ! 
    452 !!      IF( nn_timing == 1 )  CALL timing_start('ultimate_y') 
    453442      ! 
    454443      !                                                     !--  Laplacian in j-direction  --! 
     
    543532      END SELECT 
    544533      ! 
    545 !!      IF( nn_timing == 1 )  CALL timing_stop('ultimate_y') 
    546       ! 
    547534   END SUBROUTINE ultimate_y 
    548535    
     
    572559      !!---------------------------------------------------------------------- 
    573560      ! 
    574 !!      IF( nn_timing == 1 )  CALL timing_start('nonosc_2d') 
    575       ! 
    576561      zbig = 1.e+40_wp 
    577562      zsml = 1.e-15_wp 
    578563 
    579       ! clem test 
     564      ! test on divergence 
    580565      DO jj = 2, jpjm1 
    581566         DO ji = fs_2, fs_jpim1   ! vector opt.   
     
    646631      END DO 
    647632      ! 
    648 !!      IF( nn_timing == 1 )  CALL timing_stop('nonosc_2d') 
    649       ! 
    650633   END SUBROUTINE nonosc_2d 
    651634 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/LIM_SRC_3/iceistate.F90

    r8882 r8885  
    226226                        zh_i_ini(ji,jj,i_fill) = ( zvt_i_ini(ji,jj) - zV ) / MAX( za_i_ini(ji,jj,i_fill), epsi10 )  
    227227 
    228                         ! clem: correction if concentration of upper cat is greater than lower cat 
    229                         !       (it should be a gaussian around jl0 but sometimes it is not) 
     228                        ! correction if concentration of upper cat is greater than lower cat 
     229                        !   (it should be a gaussian around jl0 but sometimes it is not) 
    230230                        IF ( jl0 /= jpl ) THEN 
    231231                           DO jl = jpl, jl0+1, -1 
     
    471471      v_ice_b(:,:)     = v_ice(:,:) 
    472472 
    473 !!!clem 
    474 !!      ! Output the initial state and forcings 
     473!!clem: output of initial state should be written here but it is impossible because 
     474!!      the ocean and ice are in the same file 
    475475!!      CALL dia_wri_state( 'output.init', nit000 ) 
    476 !!!       
    477476      ! 
    478477   END SUBROUTINE ice_istate 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/LIM_SRC_3/iceitd.F90

    r8882 r8885  
    231231                           h_i_1d(ji) = h_i_1d(ji) * a_i_1d(ji) / ( a_i_1d(ji) - zda0 ) 
    232232                           a_i_1d(ji)  = a_i_1d(ji) - zda0 
    233                            v_i_1d(ji)  = a_i_1d(ji) * h_i_1d(ji) ! clem-useless ? 
     233                           v_i_1d(ji)  = a_i_1d(ji) * h_i_1d(ji) ! useless ? 
    234234                        ENDIF 
    235235                        ! 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/LIM_SRC_3/icethd_dh.F90

    r8882 r8885  
    270270      !------------------------------ 
    271271      ! qla_ice is always >=0 (upwards), heat goes to the atmosphere, therefore snow sublimates 
    272       ! clem comment: not counted in mass/heat exchange in iceupdate.F90 since this is an exchange with atm. (not ocean) 
     272      !   comment: not counted in mass/heat exchange in iceupdate.F90 since this is an exchange with atm. (not ocean) 
    273273      zdeltah(1:npti,:) = 0._wp 
    274274      DO ji = 1, npti 
     
    332332               hfx_res_1d(ji) = hfx_res_1d(ji) + zfmdt * a_i_1d(ji) * zEi * r1_rdtice 
    333333                
    334                ! Contribution to salt flux (clem: using s_i_1d and not sz_i_1d(jk) is ok) 
     334               ! Contribution to salt flux (using s_i_1d and not sz_i_1d(jk) is ok) 
    335335               sfx_res_1d(ji) = sfx_res_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * s_i_1d(ji) * r1_rdtice 
    336336                
     
    358358               zQm            = zfmdt * zEw                           ! Energy of the melt water sent to the ocean [J/m2, <0] 
    359359                
    360                ! Contribution to salt flux >0 (clem: using s_i_1d and not sz_i_1d(jk) is ok) 
     360               ! Contribution to salt flux >0 (using s_i_1d and not sz_i_1d(jk) is ok) 
    361361               sfx_sum_1d(ji) = sfx_sum_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * s_i_1d(ji) * r1_rdtice 
    362362                
     
    377377            zdeltah(ji,jk) = zdeltah(ji,jk) + zdum 
    378378            dh_i_sub(ji)  = dh_i_sub(ji) + zdum 
    379             ! Salt flux > 0 (clem2016: flux is sent to the ocean for simplicity but salt should remain in the ice except if all ice is melted. 
     379            ! Salt flux > 0 (clem: flux is sent to the ocean for simplicity but salt should remain in the ice except if all ice is melted. 
    380380            !                          It must be corrected at some point) 
    381381            sfx_sub_1d(ji) = sfx_sub_1d(ji) - rhoic * a_i_1d(ji) * zdum * s_i_1d(ji) * r1_rdtice 
     
    526526                  hfx_res_1d(ji) = hfx_res_1d(ji) + zfmdt * a_i_1d(ji) * zEi * r1_rdtice 
    527527 
    528                   ! Contribution to salt flux (clem: using s_i_1d and not sz_i_1d(jk) is ok) 
     528                  ! Contribution to salt flux (using s_i_1d and not sz_i_1d(jk) is ok) 
    529529                  sfx_res_1d(ji) = sfx_res_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * s_i_1d(ji) * r1_rdtice 
    530530                                     
     
    559559                  hfx_thd_1d(ji)  = hfx_thd_1d(ji) + zfmdt * a_i_1d(ji) * zEw * r1_rdtice 
    560560 
    561                   ! Contribution to salt flux (clem: using s_i_1d and not sz_i_1d(jk) is ok) 
     561                  ! Contribution to salt flux (using s_i_1d and not sz_i_1d(jk) is ok) 
    562562                  sfx_bom_1d(ji)  = sfx_bom_1d(ji) - rhoic *  a_i_1d(ji) * zdeltah(ji,jk) * s_i_1d(ji) * r1_rdtice 
    563563                   
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/LIM_SRC_3/icevar.F90

    r8882 r8885  
    615615                  zh_i(ji,i_fill) = ( zhti(ji) * zati(ji) - zV ) / MAX( za_i(ji,i_fill), epsi10 )  
    616616                  ! 
    617                   ! clem: correction if concentration of upper cat is greater than lower cat 
    618                   !       (it should be a gaussian around jl0 but sometimes it is not) 
     617                  ! correction if concentration of upper cat is greater than lower cat 
     618                  !    (it should be a gaussian around jl0 but sometimes it is not) 
    619619                  IF ( jl0 /= jpl ) THEN 
    620620                     DO jl = jpl, jl0+1, -1 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/NST_SRC/agrif_lim3_interp.F90

    r8882 r8885  
    5252      !!----------------------------------------------------------------------- 
    5353      ! 
    54       IF( Agrif_Root() .OR. nn_ice==0 )  RETURN   ! clem2017: do not interpolate if inside Parent domain or if child domain does not have ice 
     54      IF( Agrif_Root() .OR. nn_ice==0 )  RETURN   ! do not interpolate if inside Parent domain or if child domain does not have ice 
    5555      ! 
    5656      SELECT CASE( cd_type ) 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/NST_SRC/agrif_lim3_update.F90

    r8882 r8885  
    5454      IF( ( MOD( (kt-nit000)/nn_fsbc + 1, Agrif_irhot() * Agrif_Parent(nn_fsbc) / nn_fsbc ) /=0 ) .AND. (kt /= 0) ) RETURN ! do not update if nb of child time steps differ from time refinement 
    5555                                                                                                                           ! i.e. update only at the parent time step 
    56       IF( nn_ice == 0 ) RETURN   ! clem2017: do not update if child domain does not have ice 
     56      IF( nn_ice == 0 ) RETURN   ! do not update if child domain does not have ice 
    5757      ! 
    5858      Agrif_SpecialValueFineGrid    = -9999. 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/BDY/bdyini.F90

    r8882 r8885  
    502502            DO igrd = 1, jpbgrd 
    503503               id_dummy = iom_varid( inum, 'nbi'//cgrid(igrd), kdimsz=kdimsz )   
    504                !clem nblendta(igrd,ib_bdy) = kdimsz(1) 
    505                !clem jpbdtau = MAX(jpbdtau, kdimsz(1)) 
    506504               nblendta(igrd,ib_bdy) = MAXVAL(kdimsz) 
    507505               jpbdtau = MAX(jpbdtau, MAXVAL(kdimsz)) 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r8882 r8885  
    249249      ENDIF 
    250250          
    251       ! clem: heat and salt content 
     251      ! heat and salt contents 
    252252      IF( iom_use("heatc") ) THEN 
    253253         z2d(:,:)  = 0._wp  
     
    967967 
    968968#if defined key_lim3 
    969       IF( nn_ice == 2 ) THEN   ! clem2017: condition in case agrif + lim but no-ice in child grid 
     969      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid 
    970970         CALL ice_wri_state( kt, id_i, nh_i ) 
    971971      ENDIF 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk.F90

    r8884 r8885  
    111111   LOGICAL  ::   ln_taudif      ! logical flag to use the "mean of stress module - module of mean stress" data 
    112112   REAL(wp) ::   rn_pfac        ! multiplication factor for precipitation 
    113    REAL(wp) ::   rn_efac        ! multiplication factor for evaporation (clem) 
    114    REAL(wp) ::   rn_vfac        ! multiplication factor for ice/ocean velocity in the calculation of wind stress (clem) 
     113   REAL(wp) ::   rn_efac        ! multiplication factor for evaporation 
     114   REAL(wp) ::   rn_vfac        ! multiplication factor for ice/ocean velocity in the calculation of wind stress 
    115115   REAL(wp) ::   rn_zqt         ! z(q,t) : height of humidity and temperature measurements 
    116116   REAL(wp) ::   rn_zu          ! z(u)   : height of wind measurements 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/SBC/sbcmod.F90

    r8882 r8885  
    312312      ! 
    313313#if defined key_lim3 
    314       IF( lk_agrif .AND. nn_ice == 0 ) THEN 
    315                           IF( sbc_ice_alloc() /= 0 )   CALL ctl_stop('STOP', 'sbc_ice_alloc : unable to allocate arrays' )  ! clem2017: allocate ice arrays in case agrif + lim + no-ice in child grid 
     314      IF( lk_agrif .AND. nn_ice == 0 ) THEN            ! allocate ice arrays in case agrif + ice-model + no-ice in child grid 
     315                          IF( sbc_ice_alloc() /= 0 )   CALL ctl_stop('STOP', 'sbc_ice_alloc : unable to allocate arrays' ) 
    316316      ELSEIF( nn_ice == 2 ) THEN 
    317                           CALL ice_init                ! LIM3 initialization 
     317                          CALL ice_init                ! ICE initialization 
    318318      ENDIF 
    319319#endif 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/nemogcm.F90

    r8882 r8885  
    203203#if defined key_agrif 
    204204      CALL Agrif_ParentGrid_To_ChildGrid() 
    205       IF( ln_diaobs ) CALL dia_obs_wri 
     205      IF( ln_diaobs )        CALL dia_obs_wri 
    206206      IF( nn_timing == 1 )   CALL timing_finalize 
    207207      CALL Agrif_ChildGrid_To_ParentGrid() 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/SAS_SRC/diawri.F90

    r8882 r8885  
    396396 
    397397#if defined key_lim3 
    398       IF( nn_ice == 2 ) THEN   ! clem2017: condition in case agrif + lim but no-ice in child grid 
     398      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + lim but no-ice in child grid 
    399399         CALL ice_wri_state( kt, id_i, nh_i ) 
    400400      ENDIF 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/SAS_SRC/nemogcm.F90

    r8882 r8885  
    3131   USE icbini         ! handle bergs, initialisation 
    3232   USE icbstp         ! handle bergs, calving, themodynamics and transport 
    33    USE bdyini         ! open boundary cond. setting       (bdy_init routine). clem: mandatory for LIM3 
    34    USE bdydta         ! open boundary cond. setting   (bdy_dta_init routine). clem: mandatory for LIM3 
     33   USE bdyini         ! open boundary cond. setting       (bdy_init routine). mandatory for sea-ice 
     34   USE bdydta         ! open boundary cond. setting   (bdy_dta_init routine). mandatory for sea-ice 
    3535   ! 
    3636   USE lib_mpp        ! distributed memory computing 
     
    132132      ! 
    133133#if defined key_agrif 
    134 !!clem2017      IF( .NOT. Agrif_Root() ) THEN 
    135134         CALL Agrif_ParentGrid_To_ChildGrid() 
    136135         IF( ln_timing )   CALL timing_finalize 
    137136         CALL Agrif_ChildGrid_To_ParentGrid() 
    138 !!clem2017      ENDIF 
    139137#endif 
    140138      IF( ln_timing )   CALL timing_finalize 
     
    366364                           CALL sbc_init   ! Forcings : surface module  
    367365 
    368       ! ==> clem: open boundaries init. is mandatory for LIM3 because ice BDY is not decoupled from   
     366      ! ==> clem: open boundaries init. is mandatory for sea-ice because ice BDY is not decoupled from   
    369367      !           the environment of ocean BDY. Therefore bdy is called in both OPA and SAS modules.  
    370368      !           This is not clean and should be changed in the future.  
     
    511509      USE dom_oce   , ONLY: dom_oce_alloc 
    512510      USE bdy_oce   , ONLY: ln_bdy, bdy_oce_alloc 
    513       USE oce         ! clem: mandatory for LIM3 because needed for bdy arrays 
     511      USE oce       ! mandatory for sea-ice because needed for bdy arrays 
    514512      ! 
    515513      INTEGER :: ierr 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/SAS_SRC/step.F90

    r8882 r8885  
    2424   USE diawri           ! Standard run outputs             (dia_wri routine) 
    2525   USE bdy_oce   , ONLY: ln_bdy 
    26    USE bdydta           ! clem: mandatory for ESIM 
     26   USE bdydta           ! mandatory for sea-ice 
    2727   USE stpctl           ! time stepping control            (stp_ctl routine) 
    2828   ! 
     
    3737 
    3838#if defined key_agrif 
    39    USE agrif_oce, ONLY: lk_agrif_debug  !clem 
     39   USE agrif_oce, ONLY: lk_agrif_debug 
    4040#endif 
    4141    
     
    8989                             CALL iom_setkt( kstp - nit000 + 1, cxios_context )   ! tell iom we are at time step kstp 
    9090 
    91       ! ==> clem: open boundaries is mandatory for ESIM because ice BDY is not decoupled from   
     91      ! ==> clem: open boundaries is mandatory for sea-ice because ice BDY is not decoupled from   
    9292      !           the environment of ocean BDY. Therefore bdy is called in both OPA and SAS modules. 
    9393      !           From SAS: ocean bdy data are wrong  (but we do not care) and ice bdy data are OK.   
     
    114114                             CALL dia_wri_state( 'output.abort', kstp ) 
    115115      ENDIF 
    116       IF( kstp == nit000   ) CALL iom_close( numror )     ! close input  ocean restart file (clem: not sure...) 
     116      IF( kstp == nit000   ) CALL iom_close( numror )     ! close input  ocean restart file 
    117117       
    118118      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
Note: See TracChangeset for help on using the changeset viewer.