libsim  Versione7.2.6
dballe_test2.F03
1 #include "config.h"
2 
3 program dballe_test
4 
5 use dballe_class
7 use log4fortran
9 implicit none
10 
11 type(dbasession) :: session, sessionto
12 type(dbaconnection) :: connection
13 integer :: ier, category
14 CHARACTER(len=512) :: a_name
15 
16 !questa chiamata prende dal launcher il nome univoco
17 call l4f_launcher(a_name,a_name_force="dballe_test")
18 !init di log4fortran
19 ier=l4f_init()
20 !imposta a_name
21 category=l4f_category_get(a_name//".main")
22 
23 ! connect to dsn type DBA
24 connection=dbaconnection(dsn="sqlite:dballe_test.sqlite")
25 session=dbasession(connection,wipe=.true.,write=.true.)
26 call write1() ! write etherogeneous ensamble of data with attributes and constant data using macro object
27 call delete3() ! delete some attributes from one var only where are some defined metadati
28 call writeattronly() ! append some attributes to station 1
29 
30 
31 ! connect to dsn type BUFR file for write
32 sessionto=dbasession(filename="dballe_test2.bufr",wipe=.true.,write=.true.,memdb=.false.,template="generic")
33 call export2bufr()
34 
35 # ifndef F2003_FULL_FEATURES
36 call sessionto%delete()
37 #endif
38 
39 ! connect to dsn type BUFR file for write with memdb
40 sessionto=dbasession(filename="dballe_test2_memdb.bufr",wipe=.true.,write=.true.,memdb=.true.,template="generic")
41 call export2bufr()
42 
43 
44 # ifndef F2003_FULL_FEATURES
45 !close everythings
46 call sessionto%delete()
47 call session%delete()
48 call connection%delete()
49 #endif
50 
51 !chiudo il logger
52 CALL l4f_category_delete(category)
53 ier=l4f_fini()
54 
55 contains
56 
57 subroutine write1()
58 
59 type(dbametaanddata),allocatable :: metaanddata(:)
60 type(dbadcv) :: attrv
61 integer :: i
62 
63 print *,"----------------------------------------------"
64 print *,"--------------- write1 ------------------------"
65 
66 allocate(metaanddata(5)) ! one metadata for data and one for constant data
67 
68 metaanddata(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 
75 metaanddata(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 
82 metaanddata(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 
89 metaanddata(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 
96 metaanddata(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
105 allocate (metaanddata(1)%dataattrv%dataattr(1))
106 
107 ! first data
108 allocate (metaanddata(1)%dataattrv%dataattr(1)%dat,source=dbadatai("B13003",85))
109 
110 ! create an etherogeneous ensamble of attr
111 allocate (attrv%dcv(3))
112 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",30.))
113 allocate (attrv%dcv(2)%dat,source=dbadatai("*B33193",50))
114 allocate (attrv%dcv(3)%dat,source=dbadatar("*B33194",70.))
115 !assemble data and attribute
116 metaanddata(1)%dataattrv%dataattr(1)%attrv=attrv
117 
118 ! second station
119 allocate (metaanddata(2)%dataattrv%dataattr(1))
120 
121 allocate (metaanddata(2)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",27315))
122 ! create an etherogeneous ensamble of attr
123 deallocate(attrv%dcv)
124 allocate (attrv%dcv(2))
125 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",30.))
126 allocate (attrv%dcv(2)%dat,source=dbadatai("*B33193",50))
127 !assemble data and attribute
128 metaanddata(2)%dataattrv%dataattr(1)%attrv=attrv
129 
130 
131 ! 3d station with all data and attributes missing
132 allocate (metaanddata(3)%dataattrv%dataattr(1))
133 
134 allocate (metaanddata(3)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",imiss))
135 ! create an etherogeneous ensamble of attr
136 deallocate(attrv%dcv)
137 allocate (attrv%dcv(1))
138 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",rmiss))
139 !assemble data and attribute
140 metaanddata(3)%dataattrv%dataattr(1)%attrv=attrv
141 
142 
143 ! 4d station with data present and attributes missing
144 allocate (metaanddata(4)%dataattrv%dataattr(1))
145 
146 allocate (metaanddata(4)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",27315))
147 ! create an etherogeneous ensamble of attr
148 deallocate(attrv%dcv)
149 allocate (attrv%dcv(1))
150 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",rmiss))
151 !assemble data and attribute
152 metaanddata(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
157 allocate (metaanddata(5)%dataattrv%dataattr(1))
158 
159 ! first data
160 allocate (metaanddata(5)%dataattrv%dataattr(1)%dat,source=dbadatai("B13003",-5))
161 
162 ! create an attr
163 deallocate (attrv%dcv)
164 allocate (attrv%dcv(1))
165 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33192",30.))
166 !assemble data and attribute
167 metaanddata(5)%dataattrv%dataattr(1)%attrv=attrv
168 
169 
170 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
171 
172 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
173 ! display and save everythings
174 do i=1,size(metaanddata)
175  call metaanddata(i)%display()
176  !call session%extrude(metaanddata=metaanddata(i))
177  call metaanddata(i)%extrude(session)
178 end do
179 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
180 
181 end subroutine write1
182 
183 
184 
185 subroutine delete3()
186 
187 type(dbametadata),allocatable :: metadata(:)
188 type(dbafilter) :: filter
189 
190 print *,"----------------------------------------------"
191 print *,"--------------- delete ----------------------"
192 
193 
194 allocate(metadata(1))
195 
196 metadata(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 
203 filter=dbafilter(var="B13003",starvarlist="*B33193,*B33194")
204 call session%set(filter=filter)
205 call session%dissolveattr(metadata)
206 
207 deallocate(metadata)
208 
209 end subroutine delete3
210 
211 subroutine writeattronly()
212 
213 type(dbametaanddata),allocatable :: metaanddata(:)
214 type(dbadcv) :: attrv
215 integer :: i
216 
217 print *,"----------------------------------------------"
218 print *,"--------------- writeattronly ------------------------"
219 
220 allocate(metaanddata(1)) ! one metadata for data and one for constant data
221 
222 metaanddata(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
231 allocate (metaanddata(1)%dataattrv%dataattr(1))
232 
233 ! first data
234 allocate (metaanddata(1)%dataattrv%dataattr(1)%dat,source=dbadatai("B12101",28315))
235 
236 ! create an etherogeneous ensamble of attr
237 allocate (attrv%dcv(2))
238 allocate (attrv%dcv(1)%dat,source=dbadatar("*B33193"))
239 allocate (attrv%dcv(2)%dat,source=dbadatar("*B33194",40.))
240 !assemble data and attribute
241 metaanddata(1)%dataattrv%dataattr(1)%attrv=attrv
242 
243 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
244 
245 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
246 ! display and save everythings
247 do i=1,size(metaanddata)
248  call metaanddata(i)%display()
249  !call session%extrude(metaanddata=metaanddata(i))
250  call metaanddata(i)%extrude(session,attronly=.true.)
251 end do
252 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
253 
254 end subroutine writeattronly
255 
256 subroutine export2bufr()
257 type(dbametaanddata), allocatable:: metaanddatav(:)
258 integer :: i
259 
260 print *,"----------------------------------------------"
261 print *,"--------------- export2bufr ------------------------"
262 
263 call session%ingest(metaanddatav)
264 
265 do 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
274 end do
275 
276 end subroutine export2bufr
277 
278 
279 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.