Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add STRM and VPOT at 200mb to SFS #1072

Open
wants to merge 12 commits into
base: develop
Choose a base branch
from
16 changes: 16 additions & 0 deletions parm/post_avblflds.xml
Original file line number Diff line number Diff line change
Expand Up @@ -8479,5 +8479,21 @@
<scale>3.0</scale>
</param>

<param>
<post_avblfldidx>1021</post_avblfldidx>
<shortname>VPOT_ON_ISOBARIC_SFC_FROM_WIND_FLD</shortname>
<pname>VPOT</pname>
<fixed_sfc1_type>isobaric_sfc</fixed_sfc1_type>
<scale>3.0</scale>
</param>

<param>
<post_avblfldidx>1022</post_avblfldidx>
<shortname>STRM_ON_ISOBARIC_SFC_FROM_WIND_FLD</shortname>
<pname>STRM</pname>
<fixed_sfc1_type>isobaric_sfc</fixed_sfc1_type>
<scale>3.0</scale>
</param>

</post_avblflds>
</postxml>
12 changes: 12 additions & 0 deletions parm/sfs/postcntrl_sfs.xml
Original file line number Diff line number Diff line change
Expand Up @@ -68,6 +68,18 @@
<scale>5.0</scale>
</param>

<param>
<shortname>VPOT_ON_ISOBARIC_SFC_FROM_WIND_FLD</shortname>
<level>20000.</level>
<scale>3.0</scale>
</param>

<param>
<shortname>STRM_ON_ISOBARIC_SFC_FROM_WIND_FLD</shortname>
<level>20000.</level>
<scale>3.0</scale>
</param>

<param>
<shortname>MSLET_ON_MEAN_SEA_LVL</shortname>
<scale>6.0</scale>
Expand Down
86 changes: 85 additions & 1 deletion parm/sfs/postxconfig-NT-sfs.txt
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
1
114
116
GFSPRS
0
ncep_nco
Expand Down Expand Up @@ -352,6 +352,90 @@ isobaric_sfc
?
?
?
1021
VPOT_ON_ISOBARIC_SFC_FROM_WIND_FLD
?
1
tmpl4_0
VPOT
?
?
isobaric_sfc
0
?
1
20000.
?
0
?
0
?
?
?
?
0
0.0
0
0.0
?
0
0.0
0
0.0
0
0.0
0
0.0
1
3.0
0
0
0
?
?
?
1022
STRM_ON_ISOBARIC_SFC_FROM_WIND_FLD
?
1
tmpl4_0
STRM
?
?
isobaric_sfc
0
?
1
20000.
?
0
?
0
?
?
?
?
0
0.0
0
0.0
?
0
0.0
0
0.0
0
0.0
0
0.0
1
3.0
0
0
0
?
?
?
23
MSLET_ON_MEAN_SEA_LVL
?
Expand Down
13 changes: 12 additions & 1 deletion sorc/ncep_post.fd/COLLECT_LOC.f
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,8 @@
!> Date | Programmer | Comments
!> -----------|---------------------|----------
!> 2000-01-06 | Jim Tuccillo | Initial
!> 2021-06-01 | George Vandenberghe | 2D Decomposition
!> 2021-06-01 | George Vandenberghe | 2D Decomposition
!> 2024-11-19 | George Vandenberghe | Add timers
!>
!> @author Jim Tuccillo IBM @date 2000-01-06
!--------------------------------------------------------------------------------
Expand All @@ -33,10 +34,12 @@ SUBROUTINE COLLECT_LOC ( A, B )
real, dimension(ista_2l:iend_2u,jsta_2l:jend_2u), intent(in) :: a
real, dimension(im,jm), intent(out) :: b
integer ierr,n
real*8 ta,tb,tc,td,te
real, allocatable :: rbufs(:)
allocate(buff(im*jm))
jj=( jexa(me)-jsxa(me)+1) * (iexa(me)-isxa(me)+1)
allocate( rbufs(( jexa(me)-jsxa(me)+1) * (iexa(me)-isxa(me)+1)) )
ta=mpi_wtime()
!
if ( num_procs <= 1 ) then
b = a
Expand Down Expand Up @@ -81,6 +84,9 @@ SUBROUTINE COLLECT_LOC ( A, B )
deallocate(buff)
deallocate(rbufs)

tb=mpi_wtime()
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@KarinaAsmar-NOAA Clean up the debugging code.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@GeorgeVandenberghe-NOAA Would you please clean up the debugging part of COLLECT_LOC.f? Let me know when done and I'll push it to this branch.

if(me .eq. 0) print 109,' GWVX COLLECT TIME ',im,jm,tb-ta
109 format(a,2i10,f20.10)
end
!
!-----------------------------------------------------------------------
Expand All @@ -104,6 +110,8 @@ SUBROUTINE COLLECT_ALL ( A, B )
real, dimension(im,jm), intent(out) :: b
integer ierr,n
real, allocatable :: rbufs(:)
real*8 tb,ta
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@KarinaAsmar-NOAA Clean up the debugging code in this routine.

ta=mpi_wtime()
allocate(buff(im*jm))
jj=( jexa(me)-jsxa(me)+1) * (iexa(me)-isxa(me)+1)
allocate( rbufs(( jexa(me)-jsxa(me)+1) * (iexa(me)-isxa(me)+1)) )
Expand Down Expand Up @@ -146,6 +154,9 @@ SUBROUTINE COLLECT_ALL ( A, B )

deallocate(buff)
deallocate(rbufs)
tb=mpi_wtime()
if(me .eq. 0) print 109,' GWVX COLLECT_ALL',tb-ta
109 format(a,f20.10)

end

112 changes: 111 additions & 1 deletion sorc/ncep_post.fd/MDL2P.f
Original file line number Diff line number Diff line change
Expand Up @@ -39,6 +39,7 @@
!> 2023-08-24 | Y Mao | Add gtg_on option for GTG interpolation
!> 2023-09-12 | J Kenyon | Prevent spurious supercooled rain and cloud water
!> 2024-04-23 | E James | Adding smoke emissions (ebb) from RRFS
!> 2024-09-23 | K Asmar | Add velocity potential and streamfunction from wind vectors
!>
!> @author T Black W/NP2 @date 1999-09-23
!--------------------------------------------------------------------------------------
Expand Down Expand Up @@ -75,7 +76,7 @@ SUBROUTINE MDL2P(iostatusD3D)
IEND_2U, slrutah_on, gtg_on
use rqstfld_mod, only: IGET, LVLS, ID, IAVBLFLD, LVLSXML
use gridspec_mod, only: GRIDTYPE, MAPTYPE, DXVAL
use upp_physics, only: FPVSNEW, CALRH, CALVOR, CALSLR_ROEBBER, CALSLR_UUTAH
use upp_physics, only: FPVSNEW, CALRH, CALVOR, CALSLR_ROEBBER, CALSLR_UUTAH, CALCHIPSI
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!
Expand Down Expand Up @@ -107,6 +108,8 @@ SUBROUTINE MDL2P(iostatusD3D)
INTEGER, dimension(ista_2l:iend_2u,jsta_2l:jend_2u) :: NL1X, NL1XF
real, dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U,LSM) :: TPRS, QPRS, FPRS
real, dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U,LSM) :: RHPRS
real, dimension(ISTA:IEND,JSTA:JEND,LSM) :: USLP, VSLP
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@KarinaAsmar-NOAA Clean up USLP and VSLP.

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@KarinaAsmar-NOAA Clean up USLP and VSLP from this routine.

real, dimension(ista_2l:iend_2u,jsta_2l:jend_2u) :: CHI, PSI
!
INTEGER K, NSMOOTH
!
Expand Down Expand Up @@ -228,6 +231,7 @@ SUBROUTINE MDL2P(iostatusD3D)
(IGET(257) > 0) .OR. (IGET(258) > 0) .OR. &
(IGET(294) > 0) .OR. (IGET(268) > 0) .OR. &
(IGET(331) > 0) .OR. (IGET(326) > 0) .OR. &
(IGET(1021) > 0) .OR. (IGET(1022) > 0) .OR. &
! add D3D fields
(IGET(354) > 0) .OR. (IGET(355) > 0) .OR. &
(IGET(356) > 0) .OR. (IGET(357) > 0) .OR. &
Expand Down Expand Up @@ -1769,6 +1773,17 @@ SUBROUTINE MDL2P(iostatusD3D)
endif
ENDIF
ENDIF
!
! *** K. ASMAR - SAVE ALL P-LEVELS OF U/V WINDS FOR VELOCITY POTENTIAL AND STREAMFUNCTION
!
IF (IGET(1021)>0 .OR. IGET(1022)>0) THEN
DO J=JSTA,JEND
DO I=ISTA,IEND
USLP(I,J,LP)=USL(I,J)
VSLP(I,J,LP)=VSL(I,J)
ENDDO
ENDDO
ENDIF
!
!*** ABSOLUTE VORTICITY
!
Expand Down Expand Up @@ -1816,6 +1831,101 @@ SUBROUTINE MDL2P(iostatusD3D)
ENDIF
ENDIF
!
!*** CHIPSI
!
IF (IGET(1021) > 0 .or. IGET(1022) > 0) THEN
IF (LVLS(LP,IGET(1021)) > 0 .or. LVLS(LP,IGET(1022)) > 0) THEN
CALL CALCHIPSI(USL,VSL,CHI,PSI)
! print *,'me=',me,'EGRID1=',EGRID1(1:10,JSTA)
!
!*** CHI
!
IF (LVLS(LP,IGET(1021)) > 0) THEN
!$omp parallel do private(i,j)
DO J=JSTA,JEND
DO I=ISTA,IEND
GRID1(I,J) = CHI(I,J)
ENDDO
ENDDO
IF (SMFLAG .or. ioform == 'binarympiio' ) THEN
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@KarinaAsmar-NOAA @JesseMeng-NOAA May I know the reason of applying the smoothc function for velocity potential calculation?

call AllGETHERV(GRID1)
if (ioform == 'binarympiio') then
! nsmooth = max(2, min(30,nint(jm/94.0)))
! do k=1,5
CALL SMOOTHC(GRID1,SDUMMY,IM,JM,0.5)
CALL SMOOTHC(GRID1,SDUMMY,IM,JM,-0.5)
! enddo
else
NSMOOTH = nint(4.*(13500./dxm))
! endif
do k=1,NSMOOTH
CALL SMOOTH(GRID1,SDUMMY,IM,JM,0.5)
end do
endif
ENDIF
if(grib == 'grib2')then
cfld = cfld + 1
fld_info(cfld)%ifld=IAVBLFLD(IGET(1021))
fld_info(cfld)%lvl=LVLSXML(LP,IGET(1021))
!$omp parallel do private(i,j,ii,jj)
do j=1,jend-jsta+1
jj = jsta+j-1
do i=1,iend-ista+1
ii=ista+i-1
datapd(i,j,cfld) = GRID1(ii,jj)
enddo
enddo
endif
ENDIF !CHI
!
!*** PSI
!
IF (LVLS(LP,IGET(1022)) > 0) THEN
!$omp parallel do private(i,j)
DO J=JSTA,JEND
DO I=ISTA,IEND
GRID1(I,J) = PSI(I,J)
ENDDO
ENDDO
IF (SMFLAG .or. ioform == 'binarympiio' ) THEN
call AllGETHERV(GRID1)
if (ioform == 'binarympiio') then
! nsmooth = max(2, min(30,nint(jm/94.0)))
! do k=1,5
CALL SMOOTHC(GRID1,SDUMMY,IM,JM,0.5)
CALL SMOOTHC(GRID1,SDUMMY,IM,JM,-0.5)
! enddo
else
NSMOOTH = nint(4.*(13500./dxm))
! endif
do k=1,NSMOOTH
CALL SMOOTH(GRID1,SDUMMY,IM,JM,0.5)
end do
endif
ENDIF
if(grib == 'grib2')then
cfld = cfld + 1
fld_info(cfld)%ifld=IAVBLFLD(IGET(1022))
fld_info(cfld)%lvl=LVLSXML(LP,IGET(1022))
!$omp parallel do private(i,j,ii,jj)
do j=1,jend-jsta+1
jj = jsta+j-1
do i=1,iend-ista+1
ii=ista+i-1
datapd(i,j,cfld) = GRID1(ii,jj)
enddo
enddo
endif
ENDIF !PSI
ENDIF !LVLS(CHIPSI)
ENDIF !CHIPSI
!
!
! GEOSTROPHIC STREAMFUNCTION.
IF (IGET(086) > 0) THEN
IF (LVLS(LP,IGET(086)) > 0) THEN
Expand Down
Loading