1module sqliteFiles_mod
2 ! MODULE sqliteFiles_mod (prefix='sqlf' category='3. Observation input/output')
3 !
4 !:Purpose: To store the filenames of the sqlite observation files and call
5 ! subroutines in readSqlite to read and update sqlite files.
6 !
7 use mathPhysConstants_mod
8 use sqliteRead_mod
9 use obsSpaceData_mod
10 use fSQLite
11 use utilities_mod
12 use codePrecision_mod
13 use obsUtil_mod
14 use obsVariableTransforms_mod
15 use timeCoord_mod
16
17 implicit none
18 save
19 private
20 public :: sqlf_getDateStamp, sqlf_updateFile, sqlf_readFile, sqlf_cleanFile
21 public :: sqlf_addCloudParametersandEmissivity
22
23 type(fSQL_DATABASE) :: db ! type for SQLIte file handle
24 type(FSQL_STATUS) :: statusSqlite
25
26 contains
27
28 !--------------------------------------------------------------------------
29 ! sqlf_getDateStamp
30 !--------------------------------------------------------------------------
31 subroutine sqlf_getDateStamp(dateStamp, sqliteFileName)
32 !
33 ! Purpose: get dateStamp from an SQLite file
34 !
35
36 implicit none
37
38 ! Arguments:
39 integer , intent(out) :: dateStamp
40 character(len=*), intent(in) :: sqliteFileName
41
42 ! Locals:
43 logical :: fileExists
44 integer :: ier, imode, validTime, validDate, validDateRecv, validTimeRecv
45 integer :: newdate
46 integer, allocatable :: headDateValues(:), headTimeValues(:)
47
48 validDate = MPC_missingValue_INT
49 validTime = MPC_missingValue_INT
50
51 inquire(file = trim(sqliteFileName), exist = fileExists)
52
53 if (fileExists) then
54 call sqlr_getColumnValuesDate(headDateValues, headTimeValues, fileName=trim(sqliteFileName))
55 call tim_getValidDateTimeFromList(headDateValues, headTimeValues, validDate, validTime)
56 end if
57
58 ! Make sure all mpi tasks have a valid date (important for split sqlite files)
59 call rpn_comm_allreduce(validDate, validDateRecv, 1, "MPI_INTEGER", "MPI_MAX", "GRID", ier)
60 call rpn_comm_allreduce(validTime, validTimeRecv, 1, "MPI_INTEGER", "MPI_MAX", "GRID", ier)
61
62 if (validDateRecv == MPC_missingValue_INT .or. validTimeRecv == MPC_missingValue_INT) then
63 write(*,*) 'sqlf_getDateStamp: WARNING: Error in getting valid date and time!'
64 dateStamp = 0
65 else
66 ! printable to stamp, validTime must be multiplied with 1e6 to make newdate work
67 imode = 3
68 ier = newdate(dateStamp, validDateRecv, validTimeRecv * 1000000, imode)
69 write(*,*)'sqlf_getDateStamp: SQLite files valid date (YYYYMMDD): ', validDateRecv
70 write(*,*)'sqlf_getDateStamp: SQLite files valid time (HH): ', validTimeRecv
71 write(*,*)'sqlf_getDateStamp: SQLite files dateStamp : ', datestamp
72 end if
73
74 end subroutine sqlf_getDateStamp
75
76 !--------------------------------------------------------------------------
77 ! sqlf_readFile
78 !--------------------------------------------------------------------------
79 subroutine sqlf_readFile(obsdat, fileName, familyType, fileIndex)
80 implicit none
81
82 ! Arguments:
83 type (struct_obs), intent(inout) :: obsdat
84 character(len=*), intent(in) :: fileName
85 character(len=*), intent(in) :: familyType
86 integer , intent(in) :: fileIndex
87
88 ! Locals:
89 integer :: headerIndexBegin, headerIndexEnd, headerIndex
90 integer :: numBody, numHeader
91 character(len=*), parameter :: my_name = 'sqlf_readFile'
92 character(len=*), parameter :: my_warning = '****** '// my_name //' WARNING: '
93 character(len=*), parameter :: my_error = '******** '// my_name //' ERROR: '
94 real(pre_obsReal) :: missingValue
95
96 write(*,*)' '
97 write(*,*)' '//trim(my_name)//': Starting '
98 write(*,*)' '
99 missingValue = real(MPC_missingValue_R8,pre_obsReal)
100 write(*,*) my_name//': FileName : ', trim(FileName)
101 write(*,*) my_name//': FamilyType : ', FamilyType
102
103 headerIndexBegin = obs_numheader(obsdat) + 1
104 call sqlr_readSqlite(obsdat, trim(familyType), trim(fileName))
105 headerIndexEnd = obs_numheader(obsdat)
106 if (trim(familyType) == 'TO') then
107 call sqlr_readSqlite_avhrr(obsdat, trim(fileName), headerIndexBegin, headerIndexEnd)
108 end if
109
110 if (trim(familyType) /= 'TO') then
111 call ovt_transformObsValues (obsdat, headerIndexBegin, headerIndexEnd)
112 call ovt_adjustHumGZ (obsdat, headerIndexBegin, headerIndexEnd)
113 call obsu_computeVertCoordSurfObs(obsdat, headerIndexBegin, headerIndexEnd)
114 end if
115
116 do headerIndex = headerIndexBegin, headerIndexEnd
117 call obs_headSet_i(obsdat, OBS_IDF, headerIndex, fileIndex)
118 call obs_setFamily(obsdat, trim(familyType), headerIndex)
119 end do
120
121 ! For GP family, initialize OBS_OER to element 15032 (ZTD formal error)
122 ! for all ZTD data (element 15031)
123 if (trim(familyType) == 'GP') then
124 write(*,*)' Initializing OBS_OER for GB-GPS ZTD to formal error (ele 15032)'
125 call obsu_setGbgpsError(obsdat, headerIndexBegin, headerIndexEnd)
126 end if
127
128 numHeader = obs_numHeader(obsdat)
129 numBody = obs_numBody(obsdat)
130 write(*,*) my_name//': obs_numheader', trim(familyType), numHeader
131 write(*,*) my_name//': obs_numbody ', trim(familyType), numBody
132 write(*,*)' '
133 write(*,*)' '//trim(my_name)//' END '
134 write(*,*)' '
135
136 end subroutine sqlf_readFile
137
138 !--------------------------------------------------------------------------
139 ! sqlf_updateFile
140 !--------------------------------------------------------------------------
141 subroutine sqlf_updateFile(obsSpaceData, fileName, familyType, fileIndex)
142 implicit none
143
144 ! Arguments:
145 type (struct_obs), intent(inout) :: obsSpaceData
146 character(len=*), intent(in) :: fileName
147 character(len=*), intent(in) :: familyType
148 integer , intent(in) :: fileIndex
149
150 ! Locals:
151 character(len=*), parameter :: myName = 'sqlf_updateFile'
152 character(len=*), parameter :: myWarning = '****** '// myName //' WARNING: '
153 character(len=*), parameter :: myError = '******** '// myName //' ERROR: '
154
155 call utl_tmg_start(13,'----UpdateSqliteFile')
156 write(*,*) myName//' Starting'
157 write(*,*) myName//': FileName : ',trim(fileName)
158 write(*,*) myName//': FamilyType : ',FamilyType
159
160 call fSQL_open(db, fileName, statusSqlite)
161 if (fSQL_error(statusSqlite) /= FSQL_OK) then
162 write(*,*) 'fSQL_open: ', fSQL_errmsg(statusSqlite)
163 write(*,*) myError, fSQL_errmsg(statusSqlite)
164 end if
165
166 call sqlr_updateSqlite(db, obsSpaceData, familyType, fileName, fileIndex)
167 call sqlr_insertSqlite(db, obsSpaceData, familyType, fileName, fileIndex)
168
169 write(*,*)' closed database -->', trim(FileName)
170 call fSQL_close(db, statusSqlite)
171 write(*,*)' '
172 write(*,*)'================================================='
173 write(*,*)' '//trim(myName)//' END '
174 write(*,*)'================================================='
175 write(*,*)' '
176 call utl_tmg_stop(13)
177
178 end subroutine sqlf_updateFile
179
180 !--------------------------------------------------------------------------
181 ! sqlf_cleanFile
182 !--------------------------------------------------------------------------
183 subroutine sqlf_cleanFile(fileName, familyType)
184 !
185 !:Purpose: to reduce the number of observation data in an SQLite file
186 !
187 implicit none
188
189 ! Arguments:
190 character(len=*), intent(in) :: fileName
191 character(len=*), intent(in) :: familyType
192
193 ! Locals:
194 character(len=*), parameter :: myName = 'sqlf_cleanFile'
195
196 write(*,*) myName//': Starting'
197 write(*,*) myName//': FileName : ',trim(fileName)
198 write(*,*) myName//': FamilyType : ',FamilyType
199
200 call sqlr_cleanSqlite(db, fileName)
201
202 write(*,*)myName//': Finished'
203
204 end subroutine sqlf_cleanFile
205
206 !--------------------------------------------------------------------------
207 ! sqlf_addCloudParametersandEmissivity
208 !--------------------------------------------------------------------------
209 subroutine sqlf_addCloudParametersandEmissivity(obsSpaceData, fileIndex, fileName)
210 !
211 !:Purpose: To insert cloud parameters in obsspace data into sqlite file
212 !
213 implicit none
214
215 ! Arguments:
216 type (struct_obs), intent(inout) :: obsSpaceData
217 character(len=*), intent(in) :: fileName
218 integer, intent(in) :: fileIndex
219
220 ! Locals:
221 character(len=*), parameter :: myName = 'sqlf_addCloudParametersandEmissivity'
222 character(len=*), parameter :: myError = '******** '// myName //' ERROR: '
223
224 call fSQL_open(db, fileName, statusSqlite)
225 if (fSQL_error(statusSqlite) /= FSQL_OK) then
226 write(*,*) 'fSQL_open: ', fSQL_errmsg(statusSqlite)
227 write(*,*) myError, fSQL_errmsg(statusSqlite)
228 end if
229 call sqlr_addCloudParametersandEmissivity(db, obsSpaceData,fileIndex)
230 call fSQL_close(db, statusSqlite)
231 end subroutine sqlf_addCloudParametersandEmissivity
232
233end module sqliteFiles_mod