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 8197 for branches/UKMO/GO6_dyn_vrt_diag/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_flt.F90 – NEMO

Ignore:
Timestamp:
2017-06-21T11:39:54+02:00 (7 years ago)
Author:
glong
Message:

Changed id's to be chars e.g. hpg to more easily identify output (and updated field_def.xml accordingly). Also rearranged scaling factors in dyn_vrt_dia subroutines in divcur.F90

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/GO6_dyn_vrt_diag/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_flt.F90

    r6486 r8197  
    4848   USE lib_fortran 
    4949   USE timing          ! Timing 
     50   USE divcur          ! for dyn_vrt_dia_3d 
    5051#if defined key_agrif 
    5152   USE agrif_opa_interp 
     
    108109      INTEGER, INTENT(  out) ::   kindic   ! solver convergence flag (<0 if not converge) 
    109110      ! 
     111      CHARACTER(len=3) ::  id_vrt_dia_spg = "spg" ! TODO remove once flags set properly 
    110112      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    111113      REAL(wp) ::   z2dt, z2dtg, zgcb, zbtd, ztdgu, ztdgv   ! local scalars 
     
    130132         !                                                        ! gcx, gcxb 
    131133      ENDIF 
     134 
     135      IF ( l_trddyn .OR. ( id_vrt_dia_spg == "spg" ) )   THEN 
     136         CALL wrk_alloc( jpi, jpj, jpk, ztrdu, ztrdv ) 
     137      END IF 
    132138 
    133139      ! Local constant initialization 
     
    184190         END DO 
    185191         ! 
    186          IF( l_trddyn )   THEN                      ! temporary save of spg trends 
    187             CALL wrk_alloc( jpi, jpj, jpk, ztrdu, ztrdv ) 
     192                                          ! temporary save of spg trends 
     193         IF ( l_trddyn .OR. ( id_vrt_dia_spg == "spg" ) )   THEN 
    188194            DO jk = 1, jpkm1              ! unweighted time stepping  
    189195               DO jj = 2, jpjm1 
     
    328334#endif       
    329335 
    330       IF( l_trddyn )   THEN                      
     336      IF ( l_trddyn .OR. ( id_vrt_dia_spg == "spg" ) )   THEN                      
    331337         ztrdu(:,:,:) = ua(:,:,:)                 ! save the after velocity before the filtered SPG 
    332338         ztrdv(:,:,:) = va(:,:,:) 
    333339         ! 
    334          CALL wrk_alloc( jpi, jpj, zpw ) 
    335          ! 
    336          zpw(:,:) = - z2dt * gcx(:,:) 
    337          CALL iom_put( "ssh_flt" , zpw )          ! output equivalent ssh modification due to implicit filter 
    338          ! 
    339          !                                        ! save surface pressure flux: -pw at z=0 
    340          zpw(:,:) = - rau0 * grav * sshn(:,:) * wn(:,:,1) * tmask(:,:,1) 
    341          CALL iom_put( "pw0_exp" , zpw ) 
    342          zpw(:,:) = wn(:,:,1) 
    343          CALL iom_put( "w0" , zpw ) 
    344          zpw(:,:) =  rau0 * z2dtg * gcx(:,:) * wn(:,:,1) * tmask(:,:,1) 
    345          CALL iom_put( "pw0_flt" , zpw ) 
    346          ! 
    347          CALL wrk_dealloc( jpi, jpj, zpw )  
    348          !                                    
     340         IF ( l_trddyn ) THEN 
     341            CALL wrk_alloc( jpi, jpj, zpw ) 
     342            ! 
     343            zpw(:,:) = - z2dt * gcx(:,:) 
     344            CALL iom_put( "ssh_flt" , zpw )          ! output equivalent ssh modification due to implicit filter 
     345            ! 
     346            !                                        ! save surface pressure flux: -pw at z=0 
     347            zpw(:,:) = - rau0 * grav * sshn(:,:) * wn(:,:,1) * tmask(:,:,1) 
     348            CALL iom_put( "pw0_exp" , zpw ) 
     349            zpw(:,:) = wn(:,:,1) 
     350            CALL iom_put( "w0" , zpw ) 
     351            zpw(:,:) =  rau0 * z2dtg * gcx(:,:) * wn(:,:,1) * tmask(:,:,1) 
     352            CALL iom_put( "pw0_flt" , zpw ) 
     353            ! 
     354            CALL wrk_dealloc( jpi, jpj, zpw )  
     355            !                                    
     356         ENDIF 
    349357      ENDIF 
    350358       
     
    363371      END DO 
    364372 
    365       IF( l_trddyn )   THEN                      ! save the explicit SPG trends for further diagnostics 
     373      IF ( l_trddyn .OR. ( id_vrt_dia_spg == "spg" ) )   THEN                      ! save the explicit SPG trends for further diagnostics 
    366374         ztrdu(:,:,:) = ( ua(:,:,:) - ztrdu(:,:,:) ) / z2dt 
    367375         ztrdv(:,:,:) = ( va(:,:,:) - ztrdv(:,:,:) ) / z2dt 
    368          CALL trd_dyn( ztrdu, ztrdv, jpdyn_spgflt, kt ) 
     376         ! 
     377         IF ( l_trddyn ) THEN 
     378            CALL trd_dyn( ztrdu, ztrdv, jpdyn_spgflt, kt ) 
     379         END IF 
     380         ! 
     381         IF ( id_vrt_dia_spg == "spg" ) THEN 
     382            ! TODO remove kt after validation 
     383            CALL dyn_vrt_dia_3d( ztrdu, ztrdv, id_vrt_dia_spg, kt ) 
     384         END IF 
    369385         ! 
    370386         CALL wrk_dealloc( jpi, jpj, jpk, ztrdu, ztrdv )  
Note: See TracChangeset for help on using the changeset viewer.