UPP (develop)
Loading...
Searching...
No Matches
getIVariableN.f
1!!!@PROCESS NOEXTCHK
2subroutine getivariablen(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
58 use wrf_io_flags_mod
59 use ctlblk_mod, only: me, spval, submodelname
60!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61 implicit none
62!
63 character(len=256) ,intent(in) :: fileName
64 character(len=19) ,intent(in) :: DateStr
65 integer ,intent(in) :: dh
66 character(*) ,intent(in) :: VarName
67 integer,intent(out) :: VarBuff(IM,JSTA_2L:JEND_2U,LM)
68 integer,intent(in) :: IM,LM,JSTA_2L,JEND_2U
69 integer,intent(in) :: IM1,LM1,JS,JE
70 integer :: ndim
71 integer :: WrfType,i,j,l,ll
72 integer, dimension(4) :: start_index, end_index
73 character (len= 4) :: staggering
74 character (len= 3) :: ordering
75 character (len=80), dimension(3) :: dimnames
76 integer, allocatable, dimension(:,:,:,:) :: data
77! real, allocatable, dimension(:,:,:,:) :: data
78 integer :: ierr
79 character(len=132) :: Stagger
80 real :: FillValue
81 integer :: OutCount
82
83! call set_wrf_debug_level ( 1 )
84
85 start_index = 1
86 end_index = 1
87 call ext_ncd_get_var_info(dh,trim(varname),ndim,ordering,stagger,start_index,end_index,wrftype,ierr)
88 IF ( ierr /= 0 ) THEN
89 if (me==0) write(*,*)'Error: ',ierr,trim(varname),' not found in ',filename
90 varbuff=0.
91 return
92 ENDIF
93 allocate(data (end_index(1), end_index(2), end_index(3), 1))
94 if (me==0) write(*,*)'WrfType in getIVariable= ',wrftype
95! if( WrfType /= WRF_REAL .AND. WrfType /= WRF_REAL8 ) then !Ignore if not a real variable
96! write(*,*) 'Error: Not a real variable',WrfType
97! return
98! endif
99! write(*,'(A9,1x,I1,3(1x,I3),1x,A,1x,A)')&
100! trim(VarName), ndim, end_index(1), end_index(2), end_index(3), &
101! trim(ordering), trim(DateStr)
102! allocate(data (end_index(1), end_index(2), end_index(3), 1))
103! call ext_ncd_read_field(dh,DateStr,TRIM(VarName),data,WrfType,0,0,0,ordering,&
104! CHANGE WrfType to WRF_REAL BECAUSE THIS TELLS WRF IO API TO CONVERT TO REAL
105 call ext_ncd_read_field(dh,datestr,trim(varname),data,wrftype,0,0,0,ordering,&
106 staggering, dimnames , &
107 start_index,end_index, & !dom
108 start_index,end_index, & !mem
109 start_index,end_index, & !pat
110 ierr)
111 IF ( ierr /= 0 ) THEN
112 if (me==0) then
113 write(*,*)'Error reading ',varname,' from ',filename
114 write(*,*)' ndim = ', ndim
115 write(*,*)' end_index(1) ',end_index(1)
116 write(*,*)' end_index(2) ',end_index(2)
117 write(*,*)' end_index(3) ',end_index(3)
118 write(*,*)'Error reading ',varname,' from ',filename
119 endif
120 varbuff = 0.0
121 return
122 ENDIF
123
124 if (me==0) then
125 if (im1>end_index(1)) write(*,*) 'Err:',varname,' IM1=',im1,&
126 ' but data dim=',end_index(1)
127 if (je>end_index(2)) write(*,*) 'Err:',varname,' JE=',je,&
128 ' but data dim=',end_index(2)
129 if (lm1>end_index(3)) write(*,*) 'Err:',varname,' LM1=',lm1,&
130 ' but data dim=',end_index(3)
131 if (ndim>3) then
132 write(*,*) 'Error: ndim = ',ndim
133 endif
134 endif
135
136 if (submodelname=='MPAS') then
137 ! For MPAS: determine the fill value associated with the variable
138 call ext_ncd_get_var_ti_real(dh,"_FillValue",trim(varname),fillvalue,1,outcount,ierr)
139 do l=1,lm1
140 ll=lm1-l+1 ! flip the z axis not sure about soil
141 do i=1,im1
142 do j=js,je
143 if (data(i,j,ll,1) /= fillvalue) then
144 varbuff(i,j,l)=data(i,j,ll,1)
145 else ! For MPAS: assign SPVAL where FillValue is present
146 varbuff(i,j,l)=spval
147 endif
148 enddo
149 enddo
150 enddo
151 else
152 do l=1,lm1
153 ll=lm1-l+1 ! flip the z axis not sure about soil
154 do i=1,im1
155 do j=js,je
156 varbuff(i,j,l)=data(i,j,ll,1)
157 enddo
158 enddo
159! write(*,*) Varname,' L ',l,': = ',data(1,1,ll,1)
160 enddo
161 endif
162
163 deallocate(data)
164 return
165
166end subroutine getivariablen