UPP  V11.0.0
 All Data Structures Files Functions Pages
getIVariableN.f
1 !!!@PROCESS NOEXTCHK
2 subroutine 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 ! USAGE: CALL getVariable(fileName,DateStr,dh,VarName,VarBuff,IM,JSTA_2L,JEND_2U,LM,IM1,JS,JE,LM1)
17 !
18 ! INPUT ARGUMENT LIST:
19 ! fileName : Character(len=256) : name of WRF output file
20 ! DateStr : Character(len=19) : date/time of requested variable
21 ! dh : integer : data handle
22 ! VarName : Character(len=31) : variable name
23 ! IM : integer : X dimension of data array
24 ! JSTA_2L : integer : start Y dimension of data array
25 ! JEND_2U : integer : end Y dimension of data array
26 ! LM : integer : Z dimension of data array
27 ! IM1 : integer : amount of data pulled in X dimension
28 ! JS : integer : start Y dimension of amount of data array pulled
29 ! JE : integer : end Y dimension of amount of data array pulled
30 ! LM1 : integer : amount of data pulled in Z dimension
31 !
32 ! data is flipped in the Z dimension from what is originally given
33 ! the code requires the Z dimension to increase with pressure
34 !
35 ! OUTPUT ARGUMENT LIST:
36 ! VarBuff : real(IM,JSTA_2L:JEND_2U,LM) : requested data array
37 !
38 ! OUTPUT FILES:
39 ! NONE
40 !
41 ! SUBPROGRAMS CALLED:
42 ! UTILITIES:
43 ! NONE
44 ! LIBRARY:
45 ! WRF I/O API
46 ! NETCDF
47 
48  ! This subroutine reads the values of the variable named VarName into the buffer
49  ! VarBuff. VarBuff is filled with data only for I=1,IM1 and for J=JS,JE
50  ! and for L=1,Lm1, presumably this will be
51  ! the portion of VarBuff that is needed for this task.
52 
54 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55  implicit none
56 !
57  character(len=256) ,intent(in) :: filename
58  character(len=19) ,intent(in) :: datestr
59  integer ,intent(in) :: dh
60  character(*) ,intent(in) :: varname
61  integer,intent(out) :: varbuff(im,jsta_2l:jend_2u,lm)
62  integer,intent(in) :: im,lm,jsta_2l,jend_2u
63  integer,intent(in) :: im1,lm1,js,je
64  integer :: ndim
65  integer :: wrftype,i,j,l,ll
66  integer, dimension(4) :: start_index, end_index
67  character (len= 4) :: staggering
68  character (len= 3) :: ordering
69  character (len=80), dimension(3) :: dimnames
70  integer, allocatable, dimension(:,:,:,:) :: data
71 ! real, allocatable, dimension(:,:,:,:) :: data
72  integer :: ierr
73  character(len=132) :: stagger
74 
75 ! call set_wrf_debug_level ( 1 )
76 
77  start_index = 1
78  end_index = 1
79  call ext_ncd_get_var_info(dh,trim(varname),ndim,ordering,stagger,start_index,end_index,wrftype,ierr)
80  IF ( ierr /= 0 ) THEN
81  write(*,*)'Error: ',ierr,trim(varname),' not found in ',filename
82  varbuff=0.
83  return
84  ENDIF
85  allocate(data (end_index(1), end_index(2), end_index(3), 1))
86  write(*,*)'WrfType in getIVariable= ',wrftype
87 ! if( WrfType /= WRF_REAL .AND. WrfType /= WRF_REAL8 ) then !Ignore if not a real variable
88 ! write(*,*) 'Error: Not a real variable',WrfType
89 ! return
90 ! endif
91 ! write(*,'(A9,1x,I1,3(1x,I3),1x,A,1x,A)')&
92 ! trim(VarName), ndim, end_index(1), end_index(2), end_index(3), &
93 ! trim(ordering), trim(DateStr)
94 ! allocate(data (end_index(1), end_index(2), end_index(3), 1))
95 ! call ext_ncd_read_field(dh,DateStr,TRIM(VarName),data,WrfType,0,0,0,ordering,&
96 ! CHANGE WrfType to WRF_REAL BECAUSE THIS TELLS WRF IO API TO CONVERT TO REAL
97  call ext_ncd_read_field(dh,datestr,trim(varname),data,wrftype,0,0,0,ordering,&
98  staggering, dimnames , &
99  start_index,end_index, & !dom
100  start_index,end_index, & !mem
101  start_index,end_index, & !pat
102  ierr)
103  IF ( ierr /= 0 ) THEN
104  write(*,*)'Error reading ',varname,' from ',filename
105  write(*,*)' ndim = ', ndim
106  write(*,*)' end_index(1) ',end_index(1)
107  write(*,*)' end_index(2) ',end_index(2)
108  write(*,*)' end_index(3) ',end_index(3)
109  varbuff = 0.0
110  return
111  ENDIF
112  if (im1>end_index(1)) write(*,*) 'Err:',varname,' IM1=',im1,&
113  ' but data dim=',end_index(1)
114  if (je>end_index(2)) write(*,*) 'Err:',varname,' JE=',je,&
115  ' but data dim=',end_index(2)
116  if (lm1>end_index(3)) write(*,*) 'Err:',varname,' LM1=',lm1,&
117  ' but data dim=',end_index(3)
118  if (ndim>3) then
119  write(*,*) 'Error: ndim = ',ndim
120  endif
121  do l=1,lm1
122  ll=lm1-l+1 ! flip the z axis not sure about soil
123  do i=1,im1
124  do j=js,je
125  varbuff(i,j,l)=data(i,j,ll,1)
126  enddo
127  enddo
128 ! write(*,*) Varname,' L ',l,': = ',data(1,1,ll,1)
129  enddo
130  deallocate(data)
131  return
132 
133 end subroutine getivariablen