Update version info for release v4.6.1 (#2122)
[WRF.git] / var / da / da_gen_be / da_print_be_stats_v.inc
blob4dc01d831cad146ba36894ec5aa367d677abf490
1 subroutine da_print_be_stats_v(outunit, variable, nk, num_bins2d, &
2    e_vec, e_val, e_vec_loc, e_val_loc)
4    !----------------------------------------------------------------------------
5    ! Purpose: To print out the first global and local eigenvector and 
6    !           eigenvalues for 'variable'.
7    !
8    ! Input   : outunit    --- output unit
9    !           variable   --- variable name: psi, chi_u, t_u, rh
10    !           nk         --- the number of vertical modes or levels
11    !           num_bins2d --- the number of the 2d bins
12    !           e_vec, e_val         --- global eigenvector and eigenvalues
13    !           e_vec_loc, e_val_loc --- local eigenvectors and eigenvalues
14    !
15    ! Output  : fort.174,178,182,186,(190) --- The first 5 global eigenvectors
16    !           fort.175,179,183,187,(191) --- The global eigenvalues
17    !           fort.176,180,184,188,(192) --- The first 5 local eigenvectors
18    !           fort.177,181,185,189,(193) --- The first 5 local eigenvalues
19    !      
20    !          * in parenthisis, the units for 2d fields ps_u (ps).
21    !----------------------------------------------------------------------------
23    implicit none
25    integer, intent(inout)   :: outunit                    ! Output file unit.
26    character*10, intent(in) :: variable                   ! Variable name
27    integer, intent(in)      :: nk                         ! Vertical dimension.
28    integer, intent(in)      :: num_bins2d                 ! # bins for 2D fields.
29    real, intent(in)         :: e_vec(1:nk,1:nk)           ! Domain-averaged eigenvectors.
30    real, intent(in)         :: e_val(1:nk)                ! Domain-averaged eigenvalues.
31    real, intent(in) :: e_vec_loc(1:nk,1:nk,1:num_bins2d)  ! Latitudinally varying eigenvectors.
32    real, intent(in) :: e_val_loc(1:nk,1:num_bins2d)       ! Latitudinally varying eigenvalues.
34    integer                  :: k, m, b, mn                ! Loop counters.
36    if (trace_use) call da_trace_entry("da_print_be_stats_v")
38    if (nk > 5) then
39       mn = 5
40    else
41       mn = nk
42    end if
44    ! 1, Global vectors:
45    write(unit=stdout,fmt='(3a,i5)')' First 5 Global eigenvectors for variable ', trim(variable), &
46                      ' in unit ', outunit
48    open(unit=outunit)
49    do k = 1, nk
50       write(unit=outunit,fmt='(i4,5f15.8)') k, (e_vec(k,m), m = 1, mn)
51    end do
52    close(unit=outunit)
54    ! 2, Global values:
55    outunit = outunit + 1
57    write(unit=stdout,fmt='(3a,i5)')' Global eigenvalues for variable ', trim(variable), &
58                      ' in unit ', outunit
60    open(unit=outunit)
61    do k = 1, nk
62       write(unit=outunit,fmt='(i4,1pe18.5)') k, e_val(k)
63    end do
64    close(unit=outunit)
66    ! 3, Local vectors:
68    outunit = outunit + 1
70    write(unit=stdout,fmt='(3a,i5)')' First 5 local eigenvectors for variable ', trim(variable), &
71                      ' in unit ', outunit
73    open(unit=outunit)
74    do b = 1, num_bins2d
75      write(unit=outunit,fmt='(/"bin =",i6)') b
76      do k = 1, nk
77        write(unit=outunit,fmt='(i4,5f15.8)') k, (e_vec_loc(k,m,b), m = 1, mn)
78      end do
79    end do
80    close(unit=outunit)
82    ! 4. Local values:
84    outunit = outunit + 1
86    write(unit=stdout,fmt='(3a,i5)')' First 5 local eigenvalues for variable ', trim(variable), &
87                      ' in unit ', outunit
89    open(unit=outunit)
90    do b = 1, num_bins2d
91       write(unit=outunit,fmt='(i4,5(1pe18.5))') b, (e_val_loc(m,b), m = 1, mn)
92    end do
93    close(unit=outunit)
95    outunit = outunit + 1 
96    write(unit=stdout,fmt=*) ' '
98    if (trace_use) call da_trace_exit("da_print_be_stats_v")
100 end subroutine da_print_be_stats_v