Opened 7 years ago
Closed 7 years ago
#751 closed upgrade (Done)
Important: Added new vertical stretching functional
Reported by: | arango | Owned by: | |
---|---|---|---|
Priority: | major | Milestone: | Release ROMS/TOMS 3.7 |
Component: | Nonlinear | Version: | 3.7 |
Keywords: | Cc: |
Description
Added new vertical stretching option, Vstretching=5, based on a quadratic Legendre polynomial function that allows higher resolution near the surface developed by Souza et al. (2015). Many thanks to Joao Sousa and Brian Powell for coding and testing this new vertical stretching. For more documentation, please check WikiROMS or the publication reference.
- The new vertical stretching functional is coded in set_scoord.F:
! !---------------------------------------------------------------------- ! Stretching 5 case using a quadratic Legendre polynomial function ! aproach for the s-coordinate to enhance the surface exchange layer. ! ! J. Souza, B.S. Powell, A.C. Castillo-Trujillo, and P. Flament, 2015: ! The Vorticity Balance of the Ocean Surface in Hawaii from a ! Regional Reanalysis.'' J. Phys. Oceanogr., 45, 424-440. ! ! Added by Joao Marcos Souza - SOEST - 05/07/2012. !----------------------------------------------------------------------- ! ELSE IF (Vstretching(ng).eq.5) THEN SCALARS(ng)%sc_w(N(ng))=0.0_r8 SCALARS(ng)%Cs_w(N(ng))=0.0_r8 DO k=N(ng)-1,1,-1 rk=REAL(k,r8) rN=REAL(N(ng),r8) sc_w=-(rk*rk - 2.0_r8*rk*rN + rk + rN*rN - rN)/(rN*rN - rN)- & 0.01_r8*(rk*rk - rk*rN)/(1.0_r8 - rN) SCALARS(ng)%sc_w(k)=sc_w IF (theta_s(ng).gt.0.0_r8) THEN Csur=(1.0_r8-COSH(theta_s(ng)*sc_w))/ & & (COSH(theta_s(ng))-1.0_r8) ELSE Csur=-sc_w**2 END IF IF (theta_b(ng).gt.0.0_r8) THEN Cbot=(EXP(theta_b(ng)*Csur)-1.0_r8)/ & & (1.0_r8-EXP(-theta_b(ng))) SCALARS(ng)%Cs_w(k)=Cbot ELSE SCALARS(ng)%Cs_w(k)=Csur END IF END DO SCALARS(ng)%sc_w(0)=-1.0_r8 SCALARS(ng)%Cs_w(0)=-1.0_r8 ! DO k=1,N(ng) rk=REAL(k,r8)-0.5_r8 rN=REAL(N(ng),r8) sc_r=-(rk*rk - 2.0_r8*rk*rN + rk + rN*rN - rN)/(rN*rN - rN)- & 0.01_r8*(rk*rk - rk*rN)/(1.0_r8 - rN) SCALARS(ng)%sc_r(k)=sc_r IF (theta_s(ng).gt.0.0_r8) THEN Csur=(1.0_r8-COSH(theta_s(ng)*sc_r))/ & & (COSH(theta_s(ng))-1.0_r8) ELSE Csur=-sc_r**2 END IF IF (theta_b(ng).gt.0.0_r8) THEN Cbot=(EXP(theta_b(ng)*Csur)-1.0_r8)/ & & (1.0_r8-EXP(-theta_b(ng))) SCALARS(ng)%Cs_r(k)=Cbot ELSE SCALARS(ng)%Cs_r(k)=Csur END IF END DO END IF
The new functional is similar to Vstretching=4, but with a quadratic instead or linearly varying fractional vertical coordinate, s.
- The Matlab scripts stretching.m, set_depth.m, and scoord.m were updated with the new vertical stretching option. The difference between Vstretchin=4 and Vstretching=5 is illustrated below.
For Vtransform=2 and Vstretching=4, we get almost 4 vertical levels in the top 50 m:
Contrarily, with Vtransform=2 and Vstretching=5, we double the number of vertical levels to 8 in the top 50 m:
![]()
Notice that the increased vertical resolution at the near-surface causes a loss of resolution at the bottom.
- All the input ocean_*.in scripts were modified to update the documentation about Vstretching.
Reference:
Souza, J. M., B. Powell, A. C. Castillo-Trujillo, and P. Flament, 2015: The Vorticity Balance of the Ocean Surface in Hawaii from a Regional Reanalysis, J. Phys. Oceanogr., 45, 424-440, doi:10.1175/JPO-D-14-0074.1.