libsim  Versione6.3.0
dballe_test2.F03
1 program dballe_test
2 
3 use dballe_class
5 use log4fortran
7 implicit none
8 
9 type(dbasession) :: session, sessionto
10 type(dbaconnection) :: connection
11 integer :: ier, category
12 CHARACTER(len=512) :: a_name
13 
14 !questa chiamata prende dal launcher il nome univoco
15 call l4f_launcher(a_name,a_name_force="dballe_test")
16 !init di log4fortran
17 ier=l4f_init()
18 !imposta a_name
19 category=l4f_category_get(a_name//".main")
20 
21 ! connect to dsn type DBA
22 connection=dbaconnection(dsn="sqlite:dballe_test.sqlite")
23 session=dbasession(connection,wipe=.true.,write=.true.)
24 call write1() ! write etherogeneous ensamble of data with attributes and constant data using macro object
25 call delete3() ! delete some attributes from one var only where are some defined metadati
26 call writeattronly() ! append some attributes to station 1
27 
28 
29 ! connect to dsn type BUFR file for write
30 sessionto=dbasession(filename="dballe_test2.bufr",wipe=.true.,write=.true.,memdb=.false.,template="generic")
31 call export2bufr()
32 
33 call sessionto%delete()
34 
35 ! connect to dsn type BUFR file for write with memdb
36 sessionto=dbasession(filename="dballe_test2_memdb.bufr",wipe=.true.,write=.true.,memdb=.true.,template="generic")
37 call export2bufr()
38 
39 
40 !close everythings
41 call sessionto%delete()
42 call session%delete()
43 call connection%delete()
44 
45 
46 !chiudo il logger
47 CALL l4f_category_delete(category)
48 ier=l4f_fini()
49 
50 contains
51 
52 subroutine write1()
53 
54 type(dbametaanddata),allocatable :: metaanddata(:)
55 type(dbadcv) :: attrv
56 integer :: i
57 
58 print *,"----------------------------------------------"
59 print *,"--------------- write1 ------------------------"
60 
61 allocate(metaanddata(5)) ! one metadata for data and one for constant data
62 
63 metaanddata(1)%metadata=dbametadata( &
64  level=dbalevel(level1=105, l1=2000) &
65  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
66  ,ana=dbaana(lon=10.d0,lat=45.d0) &
67  ,network=dbanetwork("generic") &
68  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
69 
70 metaanddata(2)%metadata=dbametadata( &
71  level=dbalevel(level1=105, l1=2000) &
72  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
73  ,ana=dbaana(lon=11.d0,lat=45.d0) &
74  ,network=dbanetwork("generic") &
75  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
76 
77 metaanddata(3)%metadata=dbametadata( &
78  level=dbalevel(level1=105, l1=2000) &
79  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
80  ,ana=dbaana(lon=12.d0,lat=45.d0) &
81  ,network=dbanetwork("generic") &
82  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
83 
84 metaanddata(4)%metadata=dbametadata( &
85  level=dbalevel(level1=105, l1=2000) &
86  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
87  ,ana=dbaana(lon=13.d0,lat=45.d0) &
88  ,network=dbanetwork("generic") &
89  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
90 
91 metaanddata(5)%metadata=dbametadata( &
92  level=dbalevel(level1=105, l1=2000) &
93  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
94  ,ana=dbaana(lon=14.d0,lat=45.d0) &
95  ,network=dbanetwork("generic") &
96  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
97 
98 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99 ! create an etherogeneous ensamble of data
100 allocate (metaanddata(1)%dataattrv%dataattr(1))
101 
102 ! first data
103 allocate (metaanddata(1)%dataattrv%dataattr(1)%dat,source=dbadatai("B13003",85))
104 
105 ! create an etherogeneous ensamble of attr
106 allocate (attrv%dcv(3))
107 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",30.))
108 allocate (attrv%dcv(2)%dat,source=dbadatai("*B33193",50))
109 allocate (attrv%dcv(3)%dat,source=dbadatar("*B33194",70.))
110 !assemble data and attribute
111 metaanddata(1)%dataattrv%dataattr(1)%attrv=attrv
112 
113 ! second station
114 allocate (metaanddata(2)%dataattrv%dataattr(1))
115 
116 allocate (metaanddata(2)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",27315))
117 ! create an etherogeneous ensamble of attr
118 deallocate(attrv%dcv)
119 allocate (attrv%dcv(2))
120 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",30.))
121 allocate (attrv%dcv(2)%dat,source=dbadatai("*B33193",50))
122 !assemble data and attribute
123 metaanddata(2)%dataattrv%dataattr(1)%attrv=attrv
124 
125 
126 ! 3d station with all data and attributes missing
127 allocate (metaanddata(3)%dataattrv%dataattr(1))
128 
129 allocate (metaanddata(3)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",imiss))
130 ! create an etherogeneous ensamble of attr
131 deallocate(attrv%dcv)
132 allocate (attrv%dcv(1))
133 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",rmiss))
134 !assemble data and attribute
135 metaanddata(3)%dataattrv%dataattr(1)%attrv=attrv
136 
137 
138 ! 4d station with data present and attributes missing
139 allocate (metaanddata(4)%dataattrv%dataattr(1))
140 
141 allocate (metaanddata(4)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",27315))
142 ! create an etherogeneous ensamble of attr
143 deallocate(attrv%dcv)
144 allocate (attrv%dcv(1))
145 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",rmiss))
146 !assemble data and attribute
147 metaanddata(4)%dataattrv%dataattr(1)%attrv=attrv
148 
149 
150 
151 ! ! 5d station with data invalid (Value -5 is outside the range [0,126] for 013003 (RELATIVE HUMIDITY) ) and attribute
152 allocate (metaanddata(5)%dataattrv%dataattr(1))
153 
154 ! first data
155 allocate (metaanddata(5)%dataattrv%dataattr(1)%dat,source=dbadatai("B13003",-5))
156 
157 ! create an attr
158 deallocate (attrv%dcv)
159 allocate (attrv%dcv(1))
160 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",30.))
161 !assemble data and attribute
162 metaanddata(5)%dataattrv%dataattr(1)%attrv=attrv
163 
164 
165 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
166 
167 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
168 ! display and save everythings
169 do i=1,size(metaanddata)
170  call metaanddata(i)%display()
171  !call session%extrude(metaanddata=metaanddata(i))
172  call metaanddata(i)%extrude(session)
173 end do
174 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
175 
176 end subroutine write1
177 
178 
179 
180 subroutine delete3()
181 
182 type(dbametadata),allocatable :: metadata(:)
183 type(dbafilter) :: filter
184 
185 print *,"----------------------------------------------"
186 print *,"--------------- delete ----------------------"
187 
188 
189 allocate(metadata(1))
190 
191 metadata(1)=dbametadata( &
192  level=dbalevel(level1=105, l1=2000) &
193  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
194  ,ana=dbaana(lon=10.d0,lat=45.d0) &
195  ,network=dbanetwork("generic") &
196  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
197 
198 filter=dbafilter(var="B13003",starvarlist="*B33193,*B33194")
199 call session%set(filter=filter)
200 call session%dissolveattr(metadata)
201 
202 deallocate(metadata)
203 
204 end subroutine delete3
205 
206 subroutine writeattronly()
207 
208 type(dbametaanddata),allocatable :: metaanddata(:)
209 type(dbadcv) :: attrv
210 integer :: i
211 
212 print *,"----------------------------------------------"
213 print *,"--------------- writeattronly ------------------------"
214 
215 allocate(metaanddata(1)) ! one metadata for data and one for constant data
216 
217 metaanddata(1)%metadata=dbametadata( &
218  level=dbalevel(level1=105, l1=2000) &
219  ,timerange=dbatimerange(timerange=4, p1=3600,p2=7200) &
220  ,ana=dbaana(lon=11.d0,lat=45.d0) &
221  ,network=dbanetwork("generic") &
222  ,datetime=dbadatetime(datetime_new(2014,01,06,18,00)))
223 
224 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
225 ! create an etherogeneous ensamble of data
226 allocate (metaanddata(1)%dataattrv%dataattr(1))
227 
228 ! first data
229 allocate (metaanddata(1)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",28315))
230 
231 ! create an etherogeneous ensamble of attr
232 allocate (attrv%dcv(2))
233 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33193"))
234 allocate (attrv%dcv(2)%dat,source=dbadatar("*B33194",40.))
235 !assemble data and attribute
236 metaanddata(1)%dataattrv%dataattr(1)%attrv=attrv
237 
238 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
239 
240 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
241 ! display and save everythings
242 do i=1,size(metaanddata)
243  call metaanddata(i)%display()
244  !call session%extrude(metaanddata=metaanddata(i))
245  call metaanddata(i)%extrude(session,attronly=.true.)
246 end do
247 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
248 
249 end subroutine writeattronly
250 
251 subroutine export2bufr()
252 type(dbametaanddata), allocatable:: metaanddatav(:)
253 integer :: i
254 
255 print *,"----------------------------------------------"
256 print *,"--------------- export2bufr ------------------------"
257 
258 call session%ingest(metaanddatav)
259 
260 do i=1,size(metaanddatav)
261  call metaanddatav(i)%display()
262  !call sessionto%extrude(metaanddata=metaanddatav(i))
263  call metaanddatav(i)%extrude(sessionto)
264  if (sessionto%memdb) then
265  call sessionto%messages_write_next("generic")
266  !clean memdb
267  call sessionto%remove_all()
268  endif
269 end do
270 
271 end subroutine export2bufr
272 
273 
274 end program dballe_test
Class for expressing an absolute time value.
datetime metadata
Classi per la gestione delle coordinate temporali.
manage session handle
summ of all metadata pieces
filter to apply before ingest data
log4fortran destructor
manage connection handle to a DSN
Definitions of constants and functions for working with missing values.
integer version for dbadata
class for import and export data from e to DB-All.e.
classe per la gestione del logging
Global log4fortran constructor.
timerange metadata
real version for dbadata

Generated with Doxygen.