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

Generated with Doxygen.