Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_define_structures / da_allocate_y.inc
blob13935e1a5213bc219903b04ac831f66e60b55ba7
1 subroutine da_allocate_y (iv, y)
3    !---------------------------------------------------------------------------
4    ! Purpose: Allocate arrays used in y and residual obs structures.
5    !---------------------------------------------------------------------------
7    implicit none
8    
9    type (iv_type), intent(in)    :: iv      ! Ob type input.
10    type (y_type),  intent(inout) :: y       ! Residual type structure.
12    integer :: n, i    ! Loop counter.
13    integer :: nlevels ! Number of levels.
15    !---------------------------------------------------------------------------
16    !  [1.0] Copy number of observations:
17    !---------------------------------------------------------------------------
19    if (trace_use) call da_trace_entry("da_allocate_y")
21    y % nlocal(:) = iv%info(:)%nlocal
22    y % ntotal(:) = iv%info(:)%ntotal
24    y % num_inst     = iv % num_inst
26   !---------------------------------------------------------------------------
27   ! [2.0] Allocate:
28   !---------------------------------------------------------------------------
30    if (y % nlocal(synop) > 0) then
31       allocate (y % synop(1:y % nlocal(synop)))
32       y % synop(1:y % nlocal(synop)) % u = 0.0
33       y % synop(1:y % nlocal(synop)) % v = 0.0
34       y % synop(1:y % nlocal(synop)) % t = 0.0
35       y % synop(1:y % nlocal(synop)) % p = 0.0
36       y % synop(1:y % nlocal(synop)) % q = 0.0
37    end if
39    if (y % nlocal(ships) > 0) then
40       allocate (y % ships(1:y % nlocal(ships)))
41       y % ships(1:y % nlocal(ships)) % u = 0.0
42       y % ships(1:y % nlocal(ships)) % v = 0.0
43       y % ships(1:y % nlocal(ships)) % t = 0.0
44       y % ships(1:y % nlocal(ships)) % p = 0.0
45       y % ships(1:y % nlocal(ships)) % q = 0.0
46    end if
48    if (y % nlocal(metar) > 0) then
49       allocate (y % metar(1:y % nlocal(metar)))
50       y % metar(1:y % nlocal(metar)) % u = 0.0
51       y % metar(1:y % nlocal(metar)) % v = 0.0
52       y % metar(1:y % nlocal(metar)) % t = 0.0
53       y % metar(1:y % nlocal(metar)) % p = 0.0
54       y % metar(1:y % nlocal(metar)) % q = 0.0
55    end if
57    if (y % nlocal(geoamv) > 0) then
58       allocate (y % geoamv(1:y % nlocal(geoamv)))
59       do n = 1, y % nlocal(geoamv)
60          nlevels = iv%info(geoamv)%levels(n)
61          allocate (y % geoamv(n)%u(1:nlevels))
62          allocate (y % geoamv(n)%v(1:nlevels))
63          y % geoamv(n) % u(1:nlevels) = 0.0
64          y % geoamv(n) % v(1:nlevels) = 0.0
65       end do
66    end if
68    if (y % nlocal(polaramv) > 0) then
69       allocate (y % polaramv(1:y % nlocal(polaramv)))
70       do n = 1, y % nlocal(polaramv)
71          nlevels = iv%info(polaramv)%levels(n)
72          allocate (y % polaramv(n)%u(1:nlevels))
73          allocate (y % polaramv(n)%v(1:nlevels))
74          y % polaramv(n) % u(1:nlevels) = 0.0
75          y % polaramv(n) % v(1:nlevels) = 0.0
76       end do
77    end if
79    if (y % nlocal(gpspw) > 0) then
80       allocate (y % gpspw(1:y % nlocal(gpspw)))
81       y % gpspw(1:y % nlocal(gpspw)) % tpw = 0.0
82    end if
84    if (y % nlocal(gpsref) > 0) then
85       allocate (y % gpsref(1:y % nlocal(gpsref)))
86       do n = 1, y % nlocal(gpsref)
87          nlevels = iv%info(gpsref)%levels(n)
88          allocate (y % gpsref(n)%ref(1:nlevels))
89          allocate (y % gpsref(n)%  p(1:nlevels))
90          allocate (y % gpsref(n)%  t(1:nlevels))
91          allocate (y % gpsref(n)%  q(1:nlevels))
93          y % gpsref(n) % ref(1:nlevels) = 0.0
94          y % gpsref(n) %   p(1:nlevels) = 0.0
95          y % gpsref(n) %   t(1:nlevels) = 0.0
96          y % gpsref(n) %   q(1:nlevels) = 0.0
97       end do
98    end if
100    if (y % nlocal(gpseph) > 0) then
101       allocate (y % gpseph(1:y % nlocal(gpseph)))
102       do n = 1, y % nlocal(gpseph)
103          nlevels = iv%info(gpseph)%levels(n)
104          if ( nlevels > 0 ) then
105             allocate (y % gpseph(n)%eph(1:nlevels))
106             y % gpseph(n) % eph(1:nlevels) = 0.0
107          end if
108       end do
109    end if
111    if (y % nlocal(sound) > 0) then
112       allocate (y % sound(1:y % nlocal(sound)))
113       do n = 1, y % nlocal(sound)
114          nlevels = max(1,iv%info(sound)%levels(n))
115          allocate (y % sound(n)%u(1:nlevels))
116          allocate (y % sound(n)%v(1:nlevels))
117          allocate (y % sound(n)%t(1:nlevels))
118          allocate (y % sound(n)%q(1:nlevels))
119          y % sound(n) % u(1:nlevels) = 0.0
120          y % sound(n) % v(1:nlevels) = 0.0
121          y % sound(n) % t(1:nlevels) = 0.0
122          y % sound(n) % q(1:nlevels) = 0.0
123       end do
124    end if
126    if (y % nlocal(sonde_sfc) > 0) then
127       allocate (y % sonde_sfc(1:y % nlocal(sonde_sfc)))
129       y % sonde_sfc(1:y % nlocal(sonde_sfc)) % u = 0.0
130       y % sonde_sfc(1:y % nlocal(sonde_sfc)) % v = 0.0
131       y % sonde_sfc(1:y % nlocal(sonde_sfc)) % t = 0.0
132       y % sonde_sfc(1:y % nlocal(sonde_sfc)) % p = 0.0
133       y % sonde_sfc(1:y % nlocal(sonde_sfc)) % q = 0.0
134    end if
135      
136    if (y % nlocal(mtgirs) > 0) then
137       allocate (y % mtgirs(1:y % nlocal(mtgirs)))
138       do n = 1, y % nlocal(mtgirs)
139          nlevels = max(1,iv%info(mtgirs)%levels(n))
140          allocate (y % mtgirs(n)%u(1:nlevels))
141          allocate (y % mtgirs(n)%v(1:nlevels))
142          allocate (y % mtgirs(n)%t(1:nlevels))
143          allocate (y % mtgirs(n)%q(1:nlevels))
144          y % mtgirs(n) % u(1:nlevels) = 0.0
145          y % mtgirs(n) % v(1:nlevels) = 0.0
146          y % mtgirs(n) % t(1:nlevels) = 0.0
147          y % mtgirs(n) % q(1:nlevels) = 0.0
148       end do
149    end if
151    if (y % nlocal(tamdar) > 0) then
152       allocate (y % tamdar(1:y % nlocal(tamdar)))
153       do n = 1, y % nlocal(tamdar)
154          nlevels = max(1,iv%info(tamdar)%levels(n))
155          allocate (y % tamdar(n)%u(1:nlevels))
156          allocate (y % tamdar(n)%v(1:nlevels))
157          allocate (y % tamdar(n)%t(1:nlevels))
158          allocate (y % tamdar(n)%q(1:nlevels))
159          y % tamdar(n) % u(1:nlevels) = 0.0
160          y % tamdar(n) % v(1:nlevels) = 0.0
161          y % tamdar(n) % t(1:nlevels) = 0.0
162          y % tamdar(n) % q(1:nlevels) = 0.0
163       end do
164    end if
165    if (y % nlocal(tamdar_sfc) > 0) then
166       allocate (y % tamdar_sfc(1:y % nlocal(tamdar_sfc)))
168       y % tamdar_sfc(1:y % nlocal(tamdar_sfc)) % u = 0.0
169       y % tamdar_sfc(1:y % nlocal(tamdar_sfc)) % v = 0.0
170       y % tamdar_sfc(1:y % nlocal(tamdar_sfc)) % t = 0.0
171       y % tamdar_sfc(1:y % nlocal(tamdar_sfc)) % p = 0.0
172       y % tamdar_sfc(1:y % nlocal(tamdar_sfc)) % q = 0.
173    end if
175    if (y % nlocal(pilot) > 0) then
176       allocate (y % pilot(1:y % nlocal(pilot)))
177       do n = 1, y % nlocal(pilot)
178          nlevels = iv%info(pilot)%levels(n)
179          allocate (y % pilot(n)%u(1:nlevels))
180          allocate (y % pilot(n)%v(1:nlevels))
181          y % pilot(n) % u(1:nlevels) = 0.0
182          y % pilot(n) % v(1:nlevels) = 0.0
183       end do
184    end if
186    if (y % nlocal(radar) > 0) then
187       allocate (y % radar(1:y % nlocal(radar)))
188       do n = 1, y % nlocal(radar)
189          nlevels = iv%info(radar)%levels(n)
190          allocate (y % radar(n)%rv(1:nlevels))
191          allocate (y % radar(n)%rf(1:nlevels))
192          if ( use_radar_rhv ) then
193             allocate (y % radar(n)%rrn(1:nlevels))
194             allocate (y % radar(n)%rsn(1:nlevels))
195             allocate (y % radar(n)%rgr(1:nlevels))
196          end if
197          if ( use_radar_rqv ) then
198             allocate (y % radar(n)%rqv(1:nlevels))
199          end if
201          y % radar(n) % rv(1:nlevels)  = 0.0
202          y % radar(n) % rf(1:nlevels)  = 0.0
203          if ( use_radar_rhv ) then
204             y % radar(n) % rrn(1:nlevels) = 0.0
205             y % radar(n) % rsn(1:nlevels) = 0.0
206             y % radar(n) % rgr(1:nlevels) = 0.0
207          end if
208          if ( use_radar_rqv ) then
209             y % radar(n) % rqv(1:nlevels) = 0.0
210          end if
211       end do
212    end if
214    if (y % nlocal(airep) > 0) then
215       allocate (y % airep(1:y % nlocal(airep)))
216       do n = 1, y % nlocal(airep)
217          nlevels = iv%info(airep)%levels(n)
218          allocate (y % airep(n)%u(1:nlevels))
219          allocate (y % airep(n)%v(1:nlevels))
220          allocate (y % airep(n)%t(1:nlevels))
221          allocate (y % airep(n)%q(1:nlevels))
222          y % airep(n) % u(1:nlevels) = 0.0
223          y % airep(n) % v(1:nlevels) = 0.0
224          y % airep(n) % t(1:nlevels) = 0.0
225          y % airep(n) % q(1:nlevels) = 0.0
226       end do
227    end if
229    if (y % nlocal(bogus) > 0) then
230       allocate (y % bogus(1:y % nlocal(bogus)))
231       do n = 1, y % nlocal(bogus)
232          nlevels = iv%info(bogus)%levels(n)
233          allocate (y % bogus(n)%u(1:nlevels))
234          allocate (y % bogus(n)%v(1:nlevels))
235          allocate (y % bogus(n)%t(1:nlevels))
236          allocate (y % bogus(n)%q(1:nlevels))
237          y % bogus(n) % u(1:nlevels) = 0.0
238          y % bogus(n) % v(1:nlevels) = 0.0
239          y % bogus(n) % t(1:nlevels) = 0.0
240          y % bogus(n) % q(1:nlevels) = 0.0
241       end do
243       y % bogus(1:y % nlocal(bogus)) % slp = 0.0
244    end if
246    if (y % nlocal(satem) > 0) then
247       allocate (y % satem(1:y % nlocal(satem)))
248       do n = 1, y % nlocal(satem)
249          nlevels = iv%info(satem)%levels(n)
250          allocate (y % satem(n) % thickness(1:nlevels))
251          y % satem(n) % thickness(1:nlevels) = 0.0
252       end do
253    end if
255    if (y % nlocal(ssmi_tb) > 0) then
256       allocate (y % ssmi_tb(1:y % nlocal(ssmi_tb)))
257       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb19v = 0.0
258       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb19h = 0.0
259       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb22v = 0.0
260       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb37v = 0.0
261       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb37h = 0.0
262       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb85v = 0.0
263       y % ssmi_tb(1:y % nlocal(ssmi_tb)) % tb85h = 0.0
264    end if
266    if (y % nlocal(ssmi_rv) > 0) then
267         allocate (y % ssmi_rv(1:y % nlocal(ssmi_rv)))
268         y % ssmi_rv(1:y % nlocal(ssmi_rv)) % tpw = 0.0
269         y % ssmi_rv(1:y % nlocal(ssmi_rv)) % Speed = 0.0
270    end if
271    
272    if (y % nlocal(ssmt1) > 0) then
273       allocate (y % ssmt1(1:y % nlocal(ssmt1)))
274       do n = 1, y % nlocal(ssmt1)
275          nlevels = iv%info(ssmt1)%levels(n)
276          allocate (y % ssmt1(n) % t(1:nlevels))
277          y % ssmt1(n) % t(1:nlevels) = 0.0
278       end do
279    end if
280    
281    if (y % nlocal(ssmt2) > 0) then
282       allocate (y % ssmt2(1:y % nlocal(ssmt2)))
283       do n = 1, y % nlocal(ssmt2)
284          nlevels=iv%info(ssmt2)%levels(n)
285          allocate (y % ssmt2(n) % rh(1:nlevels))
286          y % ssmt2(n) % rh(1:nlevels) = 0.0
287       end do
288    end if
289    
290    if (y % nlocal(pseudo) > 0) then
291         allocate (y % pseudo(1:y % nlocal(pseudo)))
292         y % pseudo(1:y % nlocal(pseudo)) % u = 0.0
293         y % pseudo(1:y % nlocal(pseudo)) % v = 0.0
294         y % pseudo(1:y % nlocal(pseudo)) % t = 0.0
295         y % pseudo(1:y % nlocal(pseudo)) % p = 0.0
296         y % pseudo(1:y % nlocal(pseudo)) % q = 0.0
297    end if
299    if (y % nlocal(qscat) > 0) then
300       allocate (y % qscat(1:y % nlocal(qscat)))
301       y % qscat(1:y % nlocal(qscat)) % u = 0.0
302       y % qscat(1:y % nlocal(qscat)) % v = 0.0
303    end if
304       
305    if (y % nlocal(profiler) > 0) then
306       allocate (y % profiler(1:y % nlocal(profiler)))
307       do n = 1, y % nlocal(profiler)
308          nlevels = iv%info(profiler)%levels(n)
309          allocate (y % profiler(n)%u(1:nlevels))
310          allocate (y % profiler(n)%v(1:nlevels))
311          y % profiler(n) % u(1:nlevels) = 0.0
312          y % profiler(n) % v(1:nlevels) = 0.0
313       end do
314    end if
316    if (y % nlocal(buoy) > 0) then
317       allocate (y % buoy(1:y % nlocal(buoy)))
318       y % buoy(1:y % nlocal(buoy)) % u = 0.0
319       y % buoy(1:y % nlocal(buoy)) % v = 0.0
320       y % buoy(1:y % nlocal(buoy)) % t = 0.0
321       y % buoy(1:y % nlocal(buoy)) % p = 0.0
322       y % buoy(1:y % nlocal(buoy)) % q = 0.0
323    end if
325    if (y % nlocal(rain) > 0) then
326       allocate (y % rain(1:y % nlocal(rain)))
327       y % rain(1:y % nlocal(rain)) % rain = 0.0
328    end if
330    if (y % num_inst > 0) then
331       allocate (y % instid(1:y % num_inst))
332       do i = 1,  y % num_inst
333          y % instid(i) % num_rad = iv % instid(i) % num_rad
334          y % instid(i) % nchan   = iv % instid(i) % nchan
335          ! allocate (y % instid(i) % ichan(1:y % instid(i) % nchan))
336          ! do n = 1, y % instid(i) % nchan
337          !     y % instid(i) % ichan(n) = n
338          ! end do
339          if (y % instid(i) % num_rad < 1)  then
340             nullify (y % instid(i) % tb)
341             cycle
342          end if
343          allocate (y % instid(i) % tb(1:y % instid(i) % nchan, y % instid(i) % num_rad))
344          y % instid(i) % tb(:,:) = 0.0
345       end do
346    end if
348    if (y % nlocal(airsr) > 0) then
349       allocate (y % airsr(1:y % nlocal(airsr)))
350       do n = 1, y % nlocal(airsr)
351          nlevels = iv%info(airsr)%levels(n)
352          allocate (y % airsr(n)%t(1:nlevels))
353          allocate (y % airsr(n)%q(1:nlevels))
354          y % airsr(n) % t(1:nlevels) = 0.0
355          y % airsr(n) % q(1:nlevels) = 0.0
356       end do
357    end if
359    if (trace_use) call da_trace_exit("da_allocate_y")
361 end subroutine da_allocate_y