da_setup_pseudo_obs.inc

References to this file elsewhere.
1 subroutine da_setup_pseudo_obs(iv, ob)
2 
3    !-------------------------------------------------------------------------
4    ! Purpose: Sets up pseudo ob part of observation structure.
5    !-------------------------------------------------------------------------
6 
7    implicit none
8 
9    type(iv_type),    intent(inout) :: iv   ! Obs and header structure.
10    type(y_type),     intent(inout) :: ob   ! (Smaller) observation structure.
11 
12    integer                       :: n    ! Loop counters.
13 
14    if (trace_use_dull) call da_trace_entry("da_setup_pseudo_obs")
15 
16    allocate (iv%pseudo(1:iv%info(pseudo)%nlocal))
17 
18    do n=1, iv%info(pseudo)%nlocal
19       iv%pseudo(n) % u % inv = missing_r
20       iv%pseudo(n) % v % inv = missing_r
21       iv%pseudo(n) % t % inv = missing_r
22       iv%pseudo(n) % p % inv = missing_r
23       iv%pseudo(n) % q % inv = missing_r
24 
25       iv%pseudo(n) % u % error = missing_r
26       iv%pseudo(n) % v % error = missing_r
27       iv%pseudo(n) % t % error = missing_r
28       iv%pseudo(n) % p % error = missing_r
29       iv%pseudo(n) % q % error = missing_r
30 
31       iv%pseudo(n) % u % qc  = missing_data
32       iv%pseudo(n) % v % qc  = missing_data
33       iv%pseudo(n) % t % qc  = missing_data
34       iv%pseudo(n) % p % qc  = missing_data
35       iv%pseudo(n) % q % qc  = missing_data
36 
37       ob%pseudo(n) % u = missing_r
38       ob%pseudo(n) % v = missing_r
39       ob%pseudo(n) % t = missing_r
40       ob%pseudo(n) % p = missing_r
41       ob%pseudo(n) % q = missing_r
42 
43       !---------------------------------------------------------------
44       ! [1.0] Initialise components of innovation vector:
45       !---------------------------------------------------------------
46 
47       iv%info(pseudo)%x(:,n)  = pseudo_x
48       iv%info(pseudo)%y(:,n)  = pseudo_y
49       iv%info(pseudo)%zk(:,n) = pseudo_z
50 
51       iv%info(pseudo)%i(:,n) = int(pseudo_x)
52       iv%info(pseudo)%j(:,n) = int(pseudo_y)
53 
54       iv%info(pseudo)%dx(:,n) = pseudo_x-real(iv%info(pseudo)%i(:,n))
55       iv%info(pseudo)%dy(:,n) = pseudo_y-real(iv%info(pseudo)%j(:,n))
56       iv%info(pseudo)%dxm(:,n)=1.0-iv%info(pseudo)%dx(:,n)
57       iv%info(pseudo)%dym(:,n)=1.0-iv%info(pseudo)%dy(:,n)
58 
59       if (pseudo_var(1:1) == 'u' .or. pseudo_var(1:1) == 'U') then
60          iv%pseudo(n) % u % inv = pseudo_val
61          iv%pseudo(n) % u % error = pseudo_err
62          iv%pseudo(n) % u % qc = 0
63       else if (pseudo_var(1:1) == 'v' .or. pseudo_var(1:1) == 'V') then
64          iv%pseudo(n) % v % inv = pseudo_val
65          iv%pseudo(n) % v % error = pseudo_err
66          iv%pseudo(n) % v % qc = 0
67       else if (pseudo_var(1:1) == 't' .or. pseudo_var(1:1) == 'T') then
68          iv%pseudo(n) % t % inv = pseudo_val
69          iv%pseudo(n) % t % error = pseudo_err
70          iv%pseudo(n) % t % qc = 0
71       else if (pseudo_var(1:1) == 'p' .or. pseudo_var(1:1) == 'P') then
72          iv%pseudo(n) % p % inv = pseudo_val
73          iv%pseudo(n) % p % error = pseudo_err
74          iv%pseudo(n) % p % qc = 0
75       else if (pseudo_var(1:1) == 'q' .or. pseudo_var(1:1) == 'Q') then
76          iv%pseudo(n) % q % inv = pseudo_val
77          iv%pseudo(n) % q % error = pseudo_err
78          iv%pseudo(n) % q % qc = 0
79       end if 
80       
81       write(unit=stdout,fmt='(a4,2f15.5)')pseudo_var, pseudo_val, pseudo_err
82       write(unit=stdout,fmt='(3f15.5)')pseudo_x, pseudo_y, pseudo_z
83    end do
84 
85    if (trace_use_dull) call da_trace_exit("da_setup_pseudo_obs")
86 
87 end subroutine da_setup_pseudo_obs
88 
89