UPP (develop)
Loading...
Searching...
No Matches
getVariable.f
1!!!@PROCESS NOEXTCHK
2subroutine getvariable(fileName,DateStr,dh,VarName,VarBuff,IM,JSTA_2L,JEND_2U,LM,IM1,JS,JE,LM1)
3
4
5! SUBPROGRAM DOCUMENTATION BLOCK
6! . . .
7! SUBPROGRAM: getVariable Read data from WRF output
8! PRGRMMR: MIKE BALDWIN ORG: NSSL/SPC DATE: 2002-04-08
9!
10! ABSTRACT: THIS ROUTINE READS DATA FROM A WRF OUTPUT FILE
11! USING WRF I/O API.
12! .
13!
14! PROGRAM HISTORY LOG:
15!
16! Date | Programmer | Comments
17! -----------|---------------|---------
18! 2024-08-06 | Jaymes Kenyon | Read-in netCDF fill values for MPAS applications
19! 2024-09-05 | Jaymes Kenyon | Limiting write statements to process 0 only
20!
21! USAGE: CALL getVariable(fileName,DateStr,dh,VarName,VarBuff,IM,JSTA_2L,JEND_2U,LM,IM1,JS,JE,LM1)
22!
23! INPUT ARGUMENT LIST:
24! fileName : Character(len=256) : name of WRF output file
25! DateStr : Character(len=19) : date/time of requested variable
26! dh : integer : data handle
27! VarName : Character(len=31) : variable name
28! IM : integer : X dimension of data array
29! JSTA_2L : integer : start Y dimension of data array
30! JEND_2U : integer : end Y dimension of data array
31! LM : integer : Z dimension of data array
32! IM1 : integer : amount of data pulled in X dimension
33! JS : integer : start Y dimension of amount of data array pulled
34! JE : integer : end Y dimension of amount of data array pulled
35! LM1 : integer : amount of data pulled in Z dimension
36!
37! data is flipped in the Z dimension from what is originally given
38! the code requires the Z dimension to increase with pressure
39!
40! OUTPUT ARGUMENT LIST:
41! VarBuff : real(IM,JSTA_2L:JEND_2U,LM) : requested data array
42!
43! OUTPUT FILES:
44! NONE
45!
46! SUBPROGRAMS CALLED:
47! UTILITIES:
48! NONE
49! LIBRARY:
50! WRF I/O API
51! NETCDF
52
53 ! This subroutine reads the values of the variable named VarName into the buffer
54 ! VarBuff. VarBuff is filled with data only for I=1,IM1 and for J=JS,JE
55 ! and for L=1,Lm1, presumably this will be
56 ! the portion of VarBuff that is needed for this task.
57 ! use mpi
58 use wrf_io_flags_mod, only: wrf_real, wrf_real8
59 use ctlblk_mod, only: me, spval, submodelname
60!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61 implicit none
62
63 include "mpif.h"
64
65!
66 character(len=256) ,intent(in) :: fileName
67 character(len=19) ,intent(in) :: DateStr
68 integer ,intent(in) :: dh
69 character(*) ,intent(in) :: VarName
70 real,intent(out) :: VarBuff(IM,JSTA_2L:JEND_2U,LM)
71 integer,intent(in) :: IM,LM,JSTA_2L,JEND_2U
72 integer,intent(in) :: IM1,LM1,JS,JE
73 integer :: ndim
74 integer :: WrfType,i,j,l,ll
75 integer, dimension(4) :: start_index, end_index
76 character (len= 4) :: staggering
77 character (len= 3) :: ordering
78 character (len=80), dimension(3) :: dimnames
79 real, allocatable, dimension(:,:,:,:) :: data
80 integer :: ierr,size,mype,idsize,ier
81 character(len=132) :: Stagger
82 real :: FillValue
83 integer :: OutCount
84
85! call set_wrf_debug_level ( 1 )
86 call mpi_comm_rank(mpi_comm_world,mype,ier)
87 start_index = 1
88 end_index = 1
89! print*,'SPVAL in getVariable = ',SPVAL
90 call ext_ncd_get_var_info(dh,trim(varname),ndim,ordering,stagger,start_index,end_index,wrftype,ierr)
91 IF ( ierr /= 0 ) THEN
92 if (me==0) write(*,*)'Error: ',ierr,trim(varname),' not found in ',filename
93 varbuff=0.
94 return
95 ENDIF
96 allocate(data (end_index(1), end_index(2), end_index(3), 1))
97 if( wrftype /= wrf_real .AND. wrftype /= wrf_real8 ) then !Ignore if not a real variable
98 if (me==0) write(*,*) 'Error: Not a real variable',wrftype
99 return
100 endif
101! write(*,'(A9,1x,I1,3(1x,I3),1x,A,1x,A)')&
102! trim(VarName), ndim, end_index(1), end_index(2), end_index(3), &
103! trim(ordering), trim(DateStr)
104! allocate(data (end_index(1), end_index(2), end_index(3), 1))
105! call ext_ncd_read_field(dh,DateStr,TRIM(VarName),data,WrfType,0,0,0,ordering,&
106! CHANGE WrfType to WRF_REAL BECAUSE THIS TELLS WRF IO API TO CONVERT TO REAL
107! print *,' GWVX XT_NCD GET FIELD',size(data), size(varbuff),mype
108 idsize=size(data)
109 if(mype == 0) then
110 call ext_ncd_read_field(dh,datestr,trim(varname),data,wrftype,0,0,0,ordering,&
111 staggering, dimnames , &
112 start_index,end_index, & !dom
113 start_index,end_index, & !mem
114 start_index,end_index, & !pat
115 ierr)
116 endif
117 call mpi_bcast(data,idsize,mpi_real,0,mpi_comm_world,ierr)
118 IF ( ierr /= 0 ) THEN
119 if (me==0) then
120 write(*,*)'Error reading ',varname,' from ',filename
121 write(*,*)' ndim = ', ndim
122 write(*,*)' end_index(1) ',end_index(1)
123 write(*,*)' end_index(2) ',end_index(2)
124 write(*,*)' end_index(3) ',end_index(3)
125 endif
126 varbuff = 0.0
127 return
128 ENDIF
129
130 if (me==0) then
131 if (im1>end_index(1)) write(*,*) 'Err:',varname,' IM1=',im1,&
132 ' but data dim=',end_index(1)
133 if (je>end_index(2)) write(*,*) 'Err:',varname,' JE=',je,&
134 ' but data dim=',end_index(2)
135 if (lm1>end_index(3)) write(*,*) 'Err:',varname,' LM1=',lm1,&
136 ' but data dim=',end_index(3)
137 if (ndim>3) then
138 write(*,*) 'Error: ndim = ',ndim
139 endif
140 endif
141
142 if (submodelname=='MPAS') then
143 ! For MPAS: determine the fill value associated with the variable
144 call ext_ncd_get_var_ti_real(dh,"_FillValue",trim(varname),fillvalue,1,outcount,ierr)
145 do l=1,lm1
146 ll=lm1-l+1 ! flip the z axis not sure about soil
147 do i=1,im1
148 do j=js,je
149 if (data(i,j,ll,1) /= fillvalue) then
150 varbuff(i,j,l)=data(i,j,ll,1)
151 else ! For MPAS: assign SPVAL where FillValue is present
152 varbuff(i,j,l)=spval
153 endif
154 enddo
155 enddo
156 enddo
157 else
158 do l=1,lm1
159 ll=lm1-l+1 ! flip the z axis not sure about soil
160 do i=1,im1
161 do j=js,je
162 varbuff(i,j,l)=data(i,j,ll,1)
163 enddo
164 enddo
165! write(*,*) Varname,' L ',l,': = ',data(1,1,ll,1)
166 enddo
167 endif
168
169 deallocate(data)
170 return
171
172end subroutine getvariable