1 | ! scale_LUT_io: Contains subroutines to load and save scaling Look Up Tables (LUTs) to a file |
---|
2 | |
---|
3 | ! June 2010 Written by Roj Marchand |
---|
4 | |
---|
5 | module scale_LUTs_io |
---|
6 | implicit none |
---|
7 | |
---|
8 | contains |
---|
9 | |
---|
10 | subroutine load_scale_LUTs(hp) |
---|
11 | |
---|
12 | use radar_simulator_types |
---|
13 | |
---|
14 | type(class_param), intent(inout) :: hp |
---|
15 | |
---|
16 | logical :: LUT_file_exists |
---|
17 | integer :: i,j,k,ind |
---|
18 | |
---|
19 | ! load scale LUT from file |
---|
20 | |
---|
21 | inquire(file=trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat', & |
---|
22 | exist=LUT_file_exists) |
---|
23 | |
---|
24 | if(.not.LUT_file_exists) then |
---|
25 | |
---|
26 | write(*,*) '*************************************************' |
---|
27 | write(*,*) 'Warning: Could NOT FIND radar LUT file: ', & |
---|
28 | trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat' |
---|
29 | write(*,*) 'Will calculated LUT values as needed' |
---|
30 | write(*,*) '*************************************************' |
---|
31 | |
---|
32 | return |
---|
33 | else |
---|
34 | |
---|
35 | OPEN(unit=12,file=trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat',& |
---|
36 | form='unformatted', & |
---|
37 | err= 89, & |
---|
38 | access='DIRECT',& |
---|
39 | recl=28) |
---|
40 | |
---|
41 | write(*,*) 'Loading radar LUT file: ', & |
---|
42 | trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat' |
---|
43 | |
---|
44 | do i=1,maxhclass |
---|
45 | do j=1,mt_ntt |
---|
46 | do k=1,nRe_types |
---|
47 | |
---|
48 | ind = i+(j-1)*maxhclass+(k-1)*(nRe_types*mt_ntt) |
---|
49 | |
---|
50 | read(12,rec=ind) hp%Z_scale_flag(i,j,k), & |
---|
51 | hp%Ze_scaled(i,j,k), & |
---|
52 | hp%Zr_scaled(i,j,k), & |
---|
53 | hp%kr_scaled(i,j,k) |
---|
54 | |
---|
55 | ! if(ind==1482329) then |
---|
56 | ! write (*,*) ind, hp%Z_scale_flag(i,j,k), & |
---|
57 | ! hp%Ze_scaled(i,j,k), & |
---|
58 | ! hp%Zr_scaled(i,j,k), & |
---|
59 | ! hp%kr_scaled(i,j,k) |
---|
60 | !endif |
---|
61 | |
---|
62 | enddo |
---|
63 | enddo |
---|
64 | enddo |
---|
65 | |
---|
66 | close(unit=12) |
---|
67 | return |
---|
68 | endif |
---|
69 | |
---|
70 | 89 write(*,*) 'Error: Found but could NOT READ radar LUT file: ', & |
---|
71 | trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat' |
---|
72 | stop |
---|
73 | |
---|
74 | end subroutine load_scale_LUTs |
---|
75 | |
---|
76 | subroutine save_scale_LUTs(hp) |
---|
77 | |
---|
78 | use radar_simulator_types |
---|
79 | |
---|
80 | type(class_param), intent(inout) :: hp |
---|
81 | |
---|
82 | logical :: LUT_file_exists |
---|
83 | integer :: i,j,k,ind |
---|
84 | |
---|
85 | inquire(file=trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat', & |
---|
86 | exist=LUT_file_exists) |
---|
87 | |
---|
88 | OPEN(unit=12,file=trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat',& |
---|
89 | form='unformatted', & |
---|
90 | err= 99, & |
---|
91 | access='DIRECT',& |
---|
92 | recl=28) |
---|
93 | |
---|
94 | write(*,*) 'Creating or Updating radar LUT file: ', & |
---|
95 | trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat' |
---|
96 | |
---|
97 | do i=1,maxhclass |
---|
98 | do j=1,mt_ntt |
---|
99 | do k=1,nRe_types |
---|
100 | |
---|
101 | ind = i+(j-1)*maxhclass+(k-1)*(nRe_types*mt_ntt) |
---|
102 | |
---|
103 | if(.not.LUT_file_exists .or. hp%Z_scale_added_flag(i,j,k)) then |
---|
104 | |
---|
105 | hp%Z_scale_added_flag(i,j,k)=.false. |
---|
106 | |
---|
107 | write(12,rec=ind) hp%Z_scale_flag(i,j,k), & |
---|
108 | hp%Ze_scaled(i,j,k), & |
---|
109 | hp%Zr_scaled(i,j,k), & |
---|
110 | hp%kr_scaled(i,j,k) |
---|
111 | |
---|
112 | ! 1482329 T 0.170626345026495 0.00000000000000 1.827402935860823E-003 |
---|
113 | |
---|
114 | !if(ind==1482329) then |
---|
115 | ! write (*,*) ind, hp%Z_scale_flag(i,j,k), & |
---|
116 | ! hp%Ze_scaled(i,j,k), & |
---|
117 | ! hp%Zr_scaled(i,j,k), & |
---|
118 | ! hp%kr_scaled(i,j,k) |
---|
119 | !endif |
---|
120 | endif |
---|
121 | enddo |
---|
122 | enddo |
---|
123 | enddo |
---|
124 | |
---|
125 | close(unit=12) |
---|
126 | return |
---|
127 | |
---|
128 | 99 write(*,*) 'Error: Unable to create/update radar LUT file: ', & |
---|
129 | trim(hp%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat' |
---|
130 | return |
---|
131 | |
---|
132 | end subroutine save_scale_LUTs |
---|
133 | |
---|
134 | |
---|
135 | end module scale_LUTs_io |
---|
136 | |
---|