1 subroutine da_res_generic_set_info( res_generic, proc_domain, &
4 !---------------------------------------------------------------------------
5 ! Purpose: Eliminate redundant code for many obs types.
7 ! Method: Set common fields in res_generic.
8 !---------------------------------------------------------------------------
12 type(residual_generic_type), intent(inout) :: res_generic ! generic type
13 logical, intent(in ) :: proc_domain
14 integer, intent(in ) :: obs_global_index
16 res_generic%proc_domain = proc_domain
17 res_generic%obs_global_index = obs_global_index
18 end subroutine da_res_generic_set_info
21 subroutine da_res_sound_create_template( template)
23 !---------------------------------------------------------------------------
24 ! Purpose: Return storage template for specific type stored as
25 ! residual_generic_type.
26 !---------------------------------------------------------------------------
30 type(residual_template_type), intent(inout) :: template
32 ! only vector data for this type
36 end subroutine da_res_sound_create_template
39 subroutine da_res_sound_to_generic( res_specific, iv_num_levels, &
42 !---------------------------------------------------------------------------
43 ! Purpose: Eliminate redundant code for many obs types.
45 ! Method: Specific type res_specific is translated to generic type
46 ! res_generic. Pointer manipulation is used for vector data, no
47 ! vector data is copied. Scalar data is copied. This routine
48 ! allocates memory for res_generic. The caller must ensure that
49 ! memory is deallocated later.
50 ! iv_num_levels is used as a sanity check and is ignored for
51 ! generic types that do not contain vector data.
52 !---------------------------------------------------------------------------
56 type(residual_sound_type), intent(in ) :: res_specific ! specific type
57 integer, intent(in ) :: iv_num_levels ! levels
58 type(residual_generic_type), intent(inout) :: res_generic ! generic type
60 type(residual_template_type) :: template
63 call da_res_sound_create_template( template)
64 allocate( res_generic%values(template%lbnd:template%ubnd))
65 ! only vector data for this type
66 ! store references to vector data
67 res_generic%values(1)%ptr => res_specific%u
68 res_generic%values(2)%ptr => res_specific%v
69 res_generic%values(3)%ptr => res_specific%t
70 res_generic%values(4)%ptr => res_specific%q
73 !TBH: NOTE: We could handle iv_num_levels < size(res_generic%values(n)%ptr)
74 !TBH: HOWEVER, we would have to add "num_levels" state to
75 !TBH: residual_generic_type AND send this around. Better to fix
76 !TBH: allocation in specific types to avoid wasting memory!
77 if (size(res_generic%values(n)%ptr) /= iv_num_levels) then
78 call da_error(__FILE__,__LINE__, &
79 (/'residual_sound_to_generic: mismatched levels'/))
83 end subroutine da_res_sound_to_generic
86 subroutine da_res_sound_from_generic( res_generic, res_specific)
88 !---------------------------------------------------------------------------
89 ! Purpose: Eliminate redundant code for many obs types.
91 ! Method: Generic type res_generic is translated to specific type
92 ! res_specific. Pointer manipulation is used for vector data, no
93 ! vector data is copied. Scalar data is copied.
94 !---------------------------------------------------------------------------
98 type(residual_generic_type), intent(in ) :: res_generic ! generic type
99 type(residual_sound_type), intent(inout) :: res_specific ! specific type
101 ! only vector data for this type
102 ! store references to vector data
103 res_specific%u => res_generic%values(1)%ptr
104 res_specific%v => res_generic%values(2)%ptr
105 res_specific%t => res_generic%values(3)%ptr
106 res_specific%q => res_generic%values(4)%ptr
108 end subroutine da_res_sound_from_generic
111 subroutine da_res_synop_create_template( template)
113 !---------------------------------------------------------------------------
114 ! Purpose: Return storage template for specific type stored as
115 ! residual_generic_type.
116 !---------------------------------------------------------------------------
120 type(residual_template_type), intent(inout) :: template
122 ! only scalar data for this type
126 end subroutine da_res_synop_create_template
129 subroutine da_res_synop_to_generic( res_specific, iv_num_levels, &
132 !---------------------------------------------------------------------------
133 ! Purpose: Eliminate redundant code for many obs types.
135 ! Method: Specific type res_specific is translated to generic type
136 ! res_generic. Pointer manipulation is used for vector data, no
137 ! vector data is copied. Scalar data is copied. This routine
138 ! allocates memory for res_generic. The caller must ensure that
139 ! memory is deallocated later.
140 ! iv_num_levels is used as a sanity check and is ignored for
141 ! generic types that do not contain vector data.
142 !---------------------------------------------------------------------------
146 type(residual_synop_type), intent(in ) :: res_specific ! specific type
147 integer, intent(in ) :: iv_num_levels ! levels
148 type(residual_generic_type), intent(inout) :: res_generic ! generic type
150 type(residual_template_type) :: template
152 ! use to avoid compiler warning
153 if (iv_num_levels==0) then
156 call da_res_synop_create_template( template)
157 allocate( res_generic%values(template%lbnd:template%ubnd))
158 ! only scalar data for this type
159 allocate( res_generic%values(0)%ptr(5))
161 res_generic%values(0)%ptr(1) = res_specific%u
162 res_generic%values(0)%ptr(2) = res_specific%v
163 res_generic%values(0)%ptr(3) = res_specific%t
164 res_generic%values(0)%ptr(4) = res_specific%p
165 res_generic%values(0)%ptr(5) = res_specific%q
167 end subroutine da_res_synop_to_generic
170 subroutine da_res_synop_from_generic( res_generic, res_specific)
172 !---------------------------------------------------------------------------
173 ! Purpose: Eliminate redundant code for many obs types.
175 ! Method: Generic type res_generic is translated to specific type
176 ! res_specific. Pointer manipulation is used for vector data, no
177 ! vector data is copied. Scalar data is copied.
178 !---------------------------------------------------------------------------
182 type(residual_generic_type), intent(in ) :: res_generic ! generic type
183 type(residual_synop_type), intent(inout) :: res_specific ! specific type
185 ! only scalar data for this type
187 res_specific%u = res_generic%values(0)%ptr(1)
188 res_specific%v = res_generic%values(0)%ptr(2)
189 res_specific%t = res_generic%values(0)%ptr(3)
190 res_specific%p = res_generic%values(0)%ptr(4)
191 res_specific%q = res_generic%values(0)%ptr(5)
193 end subroutine da_res_synop_from_generic
196 subroutine da_y_facade_create( slice, num_obs, num_obs_glo, template)
198 !---------------------------------------------------------------------------
199 ! Purpose: Create a y_facade_type containing specified number of
200 ! residual_generic_type objects.
202 ! Method: Allocate memory and call residual_generic_create.
203 ! Call y_facade_free() to deallocate memory allocated here.
204 ! If template is not present, residual_generic_type objects will be
205 ! left uninitialized. If template is present, then
206 ! residual_generic_type objects will be allocated but no values
207 ! will be copied into them.
209 !---------------------------------------------------------------------------
213 type(y_facade_type), intent(inout) :: slice ! selected residual obs
214 integer, intent(in ) :: num_obs
215 integer, intent(in ) :: num_obs_glo
216 type(residual_template_type), optional, intent(in ) :: template
220 slice%num_obs = num_obs
221 slice%num_obs_glo = num_obs_glo
222 allocate( slice%obs(slice%num_obs))
223 do n=1, slice%num_obs
224 call da_res_generic_create( slice%obs(n), template=template)
227 end subroutine da_y_facade_create
230 subroutine da_res_generic_create( res_generic, template)
232 !---------------------------------------------------------------------------
233 ! Purpose: Create a residual_generic_type object. This must be called
234 ! before any other operation is performed on a
235 ! residual_generic_type object. Do not pass an already-allocated
236 ! object into this routine as res_generic or you will cause a
239 ! Method: If template is not present, create an empty object.
240 ! If template is present, create an uninitialized object with
241 ! space allocated to match template. Caller is responsible
242 ! for deallocation via call to residual_generic_free(). Note
243 ! that this is *NOT* a copy-constructor because values are not
244 ! copied. Also, proc_domain and obs_global_index fields are
245 ! not copied from the template. Finally, memory is not allocated
246 ! for vector or scalar data, these pointers
247 ! (res_generic%values(:)%ptr) are nullified.
248 !---------------------------------------------------------------------------
252 type(residual_generic_type), intent(inout) :: res_generic ! generic type
253 type(residual_template_type), optional, intent(in ) :: template
257 nullify( res_generic%values)
258 if (present( template)) then
259 allocate( res_generic%values(template%lbnd:template%ubnd))
260 do i=template%lbnd, template%ubnd
261 nullify( res_generic%values(i)%ptr)
264 end subroutine da_res_generic_create
267 subroutine da_res_generic_alloc_and_set( res_generic, val_index, values)
269 ! TOdo: replace this with a full-blown copy-constructor!
270 !---------------------------------------------------------------------------
271 ! Purpose: Allocates and initializes one of the values(either scalar or
272 ! vector) in an already-created residual_generic_type object.
274 ! Method: Allocate pointer and copy data from values. Note that a call
275 ! to residual_generic_free() will deallocate scalar data but leave
276 ! vector data alone. It is still the callers responsibility to
277 ! deallocate pointers to vector data. In this particular case,
278 ! any vector data allocated here is later passed by reference to
279 ! a global specific object via call to residual_*_from_generic()
280 ! (called from y_type_insert_sound_global()) and later explictly
281 ! deallocated via call to free_global_*().
282 !---------------------------------------------------------------------------
286 type(residual_generic_type), intent(inout) :: res_generic ! generic type
287 integer, intent(in ) :: val_index ! which value
288 real, intent(in ) :: values(:) ! values
290 if (( val_index < LBOUND(res_generic%values,1)) .OR. &
291 ( val_index > UBOUND(res_generic%values,1))) then
292 call da_error(__FILE__,__LINE__, &
293 (/'residual_generic_alloc_and_set: bad val_index'/))
295 allocate( res_generic%values(val_index)%ptr(size(values)))
296 res_generic%values(val_index)%ptr(:) = values(:)
298 end subroutine da_res_generic_alloc_and_set
301 logical function da_res_generic_has_vector( res_generic)
303 !---------------------------------------------------------------------------
304 ! Purpose: Returns .true. iff res_generic stores vector values
305 !---------------------------------------------------------------------------
309 type(residual_generic_type), intent(in) :: res_generic ! generic type
310 da_res_generic_has_vector =( UBOUND(res_generic%values,1) > 0)
312 end function da_res_generic_has_vector
315 logical function da_res_generic_has_scalar( res_generic)
317 !---------------------------------------------------------------------------
318 ! Purpose: Returns .true. iff res_generic stores scalar values.
319 !---------------------------------------------------------------------------
323 type(residual_generic_type), intent(in) :: res_generic ! generic type
324 da_res_generic_has_scalar =( LBOUND(res_generic%values,1) == 0)
326 end function da_res_generic_has_scalar
329 subroutine da_res_generic_free( res_generic)
331 !---------------------------------------------------------------------------
332 ! Purpose: Free memory for residual_generic_type object.
334 ! Method: pointers to vector values are assumed to point to memory allocated
335 ! by others and are not deallocated.
336 ! This will fail if called for a residual_generic_type object that
337 ! has never been created.
338 !---------------------------------------------------------------------------
342 type(residual_generic_type), intent(inout) :: res_generic ! generic type
344 logical :: oktofreevalues
346 oktofreevalues = .false.
347 if (da_res_generic_has_scalar( res_generic)) then
349 deallocate( res_generic%values(0)%ptr)
350 oktofreevalues = .true.
352 if (da_res_generic_has_vector( res_generic)) then
353 oktofreevalues = .true.
355 if (oktofreevalues) then
356 deallocate( res_generic%values)
359 end subroutine da_res_generic_free
362 subroutine da_y_facade_free( slice)
364 !---------------------------------------------------------------------------
365 ! Purpose: Free memory for y_facade_type object.
367 ! Method: Brute force. May want to make this smarter...
368 !---------------------------------------------------------------------------
372 type(y_facade_type), intent(inout) :: slice
375 if (associated( slice%obs)) then
376 do n=1, size(slice%obs)
377 call da_res_generic_free( slice%obs(n))
379 deallocate( slice%obs)
381 end subroutine da_y_facade_free