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

Feature/rtodling/cascase #226

Merged
merged 4 commits into from
Nov 3, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 13 additions & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,19 @@ and this project adheres to [Semantic Versioning](https://semver.org/spec/v2.0.0

### Removed

### Added

### Fixed

### Changed

- add Cascade knob to g5fcst_stats.pl and regrid.pl
- revised dyn_blob: more general on the blobs
- make sure echorc.x exits w/ success code when applicable

### Removed


## [1.4.10] - 2021-10-08

### Added
Expand Down
11 changes: 8 additions & 3 deletions GEOS_Util/post/g5fcst_stats.pl
Original file line number Diff line number Diff line change
Expand Up @@ -688,13 +688,18 @@ sub submit_calcjob {
$ntspn = 36;
$qos = "#SBATCH --qos=dastest"; # wired for now since only way to use SKY
$partition = "#SBATCH --partition=preops"; # wired for now since only way to use SKY
} elsif ( $npn == 48 ) {
$mynodes = "cas";
$ntspn = 42;
$qos = "#SBATCH --qos=dastest"; # wired for now since only way to use CAS
$partition = "#SBATCH --partition=preops"; # wired for now since only way to use CAS
} else {
$mynodes = "hasw";
$ntspn = 24;
$qos = "";
$partition = "";
# $qos = "#SBATCH --qos=dastest"; # wired for now since only way to use SKY
# $partition = "#SBATCH --partition=preops"; # wired for now since only way to use SKY
# $qos = "#SBATCH --qos=dastest"; # wired for now since only way to use HASW
# $partition = "#SBATCH --partition=preops"; # wired for now since only way to use HASW
}
if ( $usrnodes ne "null" ) { $mynodes = $usrnodes }; # overwrite with specification from command line

Expand Down Expand Up @@ -981,7 +986,7 @@ sub usage {
[dirname(\$FVHOME) or \$NOBACKUP]
-noarchive do not archive outputs [archives by default]

-nodes nodesname specify nodes (e.g., sky or hasw)
-nodes nodesname specify nodes (e.g., sky, hasw, or cas)
-das check for DAS hidden files before attempting to fetch files
and set no prompt; requires \$FVHOME environment variable;

Expand Down
2 changes: 2 additions & 0 deletions GEOS_Util/post/regrid.pl
Original file line number Diff line number Diff line change
Expand Up @@ -2410,6 +2410,8 @@ sub regrid_upperair_rsts_CS {
my $npn = `facter processorcount`; chomp($npn);
if ( $npn == 40 ) {
$mynodes = "sky";
} elsif ( $npn == 48 ) {
$mynodes = "cas";
} else {
$mynodes = "hasw";
}
Expand Down
129 changes: 108 additions & 21 deletions GMAO_hermes/dyn_blob.f90
Original file line number Diff line number Diff line change
Expand Up @@ -5,34 +5,44 @@ program dyn_blob
use m_dyn, only: dyn_vect
use m_dyn, only: dyn_clean
use m_const, only: radius_earth
use m_die, only: die

implicit none

character(len=*), parameter :: fname = 'bkg.eta.nc4'
character(len=*), parameter :: pname = 'fsens.eta.nc4'
character(len=*), parameter :: bname = 'blob.eta.nc4'
character(len=*), parameter :: dname = 'delta.eta.nc4'
real, parameter :: PPMV2GpG = 1.6571E-6 ! ((47.9982 g/mol)/((28.9644 g/mol))*1e-6(mol/mol)-> g/g
integer, parameter :: dyntype=5
integer nymd, nhms, ier, im, jm, km, freq, nstep
integer ii, nymdw, nhmsw
integer, parameter :: npnts=4
real blocs(2,npnts)
real,allocatable:: blocs(:,:)
logical :: advect=.false.

real,parameter:: adrate = 0.75
real,parameter:: lenfcst = 24
real,parameter:: deltlen = 0.001 ! make this really small
real,parameter:: corrlen = 800. ! to be read from file
real,parameter:: corrlength = corrlen*1000/radius_earth
real :: corrlength

integer, parameter :: zlevout = 1 ! -1 won''t do anything
integer, parameter :: zlevout = -1! won''t do anything

type(dyn_vect) :: wind
type(dyn_vect) :: pert

logical :: delta=.true. !.false.
integer :: ntimes = 12
integer :: lu=10
real,parameter :: dlon=60.
real,allocatable :: covloc(:,:,:,:)
real,allocatable :: adcovloc(:,:,:,:)
real,allocatable :: ploc(:,:,:)
integer :: npnts,nblobs_perlon

nblobs_perlon=360./dlon
npnts = 4*nblobs_perlon
allocate(blocs(2,npnts))

! read in perturbation fields
call dyn_get ( trim(pname), nymd, nhms, pert, ier, timidx=1, freq=freq, nstep=nstep, vectype=dyntype )
Expand All @@ -53,14 +63,21 @@ program dyn_blob
pert%delp = 0.0
pert%q = 0.0
pert%pt = 0.0
pert%ts = 0.0

call set_blobs_
if ( delta ) then
corrlength = deltlen*1000/radius_earth
else
corrlength = corrlen*1000/radius_earth
endif
call make_blobs_
!call make_sfc_blobs_
call dyn_put ( trim(bname), nymd, nhms, 0, pert, ier, freq=freq, nstep=nstep, vectype=dyntype )

! read in wind fields
call dyn_get ( trim(fname), nymdw, nhmsw, wind, ier, timidx=1, freq=freq, nstep=nstep, vectype=dyntype )
if (advect) then
call dyn_get ( trim(fname), nymdw, nhmsw, wind, ier, timidx=1, freq=freq, nstep=nstep, vectype=dyntype )
call advect_blobs_
endif
if (zlevout>0) then
Expand All @@ -76,6 +93,7 @@ program dyn_blob
call dyn_put ( trim(pname), nymd, nhms, 0, pert, ier, freq=freq, nstep=nstep, vectype=dyntype )
endif


contains

subroutine wrtout_(lu,fld)
Expand Down Expand Up @@ -112,19 +130,30 @@ subroutine readin_(lu,fld)
end subroutine readin_

subroutine set_blobs_
blocs(1,1)=10 ! :,1=lats, :,2=lons
blocs(2,1)= 0

! blocs(1,2)=90
blocs(1,2)=30
blocs(2,2)=180

blocs(1,3)=45
blocs(2,3)=50
integer :: np,n

! blocs(1,4)=-90
blocs(1,4)=-30
blocs(2,4)=-60
np=1
do n=1,nblobs_perlon
blocs(1,np)=45 ! 1,:=lats, 2,:=lons
blocs(2,np)= -180.+(n-1)*dlon
np=np+1
enddo
do n=1,nblobs_perlon
blocs(1,np)=10 ! 1,:=lats, 2,:=lons
blocs(2,np)= -180.+(n-1)*dlon
np=np+1
enddo
do n=1,nblobs_perlon
blocs(1,np)=-30 ! 1,:=lats, 2,:=lons
blocs(2,np)= -180.+(n-1)*dlon
np=np+1
enddo
do n=1,nblobs_perlon
blocs(1,np)=-60 ! 1,:=lats, 2,:=lons
blocs(2,np)= -180.+(n-1)*dlon
np=np+1
enddo

allocate(ploc(1,1,3))
allocate(covloc(im,jm,3,npnts))
Expand All @@ -134,29 +163,87 @@ subroutine clean_blobs_
if(advect) deallocate(adcovloc)
deallocate(covloc)
deallocate(ploc)
deallocate(blocs)
end subroutine clean_blobs_

subroutine make_blobs_
integer nn,mm,ii,jj,iii,jjj
integer nn,ma,mb,ii,jj,kk,iii,jjj,jlat,jlon
real pi,cs,sn,dist

mm=km
do nn=1,npnts
! ma=km-9 ! 850 mb
! mb=km-9 ! 850 mb
! ma=km-22 ! 500 mb
! mb=km-22 ! 500 mb
! ma=14 ! 1 mb
! mb=14 ! 1 mb
! ma=25 ! 10 mb
! mb=25 ! 10 mb
ma=km ! sfc
mb=km ! sfc
do jlat=1,4
do jlon=1,nblobs_perlon
nn=nn+1
if(nn>npnts) call die('make_blobs_','Trying to access more pnts than avail',99)

call globeloc ( ploc, blocs(1,nn:nn), blocs(2,nn:nn) )
call globeloc ( covloc(:,:,:,nn), pert%grid%lat, pert%grid%lon )
do jj=1,jm
do ii=1,im
dist = sqrt( (covloc(ii,jj,1,nn)-ploc(1,1,1))**2 + &
(covloc(ii,jj,2,nn)-ploc(1,1,2))**2 + &
(covloc(ii,jj,3,nn)-ploc(1,1,3))**2 )/corrlength
if (dist<10*corrlength) pert%pt(ii,jj,:mm) = gc(dist)
if (jlon==1) then
if (dist<10*corrlength) then
! do kk=1,km
! pert%delp(ii,jj,kk) = 0.01*gc(dist)*(pert%grid%bk(kk+1)-pert%grid%bk(kk))
! enddo
pert%ps(ii,jj) = 0.01*gc(dist) ! 1mb
endif
endif
if (jlon==2) then
if (dist<10*corrlength) pert%pt(ii,jj,ma:mb) = gc(dist)
endif
if (jlon==3) then
if (dist<10*corrlength) pert%u(ii,jj,ma:mb) = gc(dist)
endif
if (jlon==4) then
if (dist<10*corrlength) pert%v(ii,jj,ma:mb) = gc(dist)
endif
if (jlon==5) then
if (dist<10*corrlength) pert%q(ii,jj,ma:mb,1) = gc(dist)
endif
if (jlon==6) then
if (dist<10*corrlength) pert%q(ii,jj,ma:mb,2) = gc(dist)/PPMV2GpG ! 1 g/g
endif
enddo
enddo
mm=mm-1
enddo
enddo

end subroutine make_blobs_

subroutine make_sfc_blobs_
integer nn,ii,jj,iii,jjj
real pi,cs,sn,dist

do nn=1,npnts
call globeloc ( ploc, blocs(1,nn:nn), blocs(2,nn:nn) )
call globeloc ( covloc(:,:,:,nn), pert%grid%lat, pert%grid%lon )
do jj=1,jm
do ii=1,im
dist = sqrt( (covloc(ii,jj,1,nn)-ploc(1,1,1))**2 + &
(covloc(ii,jj,2,nn)-ploc(1,1,2))**2 + &
(covloc(ii,jj,3,nn)-ploc(1,1,3))**2 )/corrlength
if (dist<10*corrlength) then
! pert%ts(ii,jj) = gc(dist)
pert%ps(ii,jj) = 100.*gc(dist)
endif
enddo
enddo
enddo

end subroutine make_sfc_blobs_

subroutine advect_blobs_
integer nn,ii,jj,kk,iii,jjj
real fcstlen,dt,pi,cs,sn,dist
Expand Down
2 changes: 0 additions & 2 deletions GMAO_hermes/echorc.f90
Original file line number Diff line number Diff line change
Expand Up @@ -132,8 +132,6 @@ program echorc

endif

stop 0

CONTAINS

!-------------------------------------------------------------------------
Expand Down