updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / hydro / nudging / module_nudging_utils.F
blob884a7234de84c136f0d2ea5704e14c458488bfa3
1 !  Program Name:
2 !  Author(s)/Contact(s):
3 !  Abstract:
4 !  History Log:
5
6 !  Usage:
7 !  Parameters: <Specify typical arguments passed>
8 !  Input Files:
9 !        <list file names and briefly describe the data they include>
10 !  Output Files:
11 !        <list file names and briefly describe the information they include>
12
13 !  Condition codes:
14 !        <list exit condition or error codes returned >
15 !        If appropriate, descriptive troubleshooting instructions or
16 !        likely causes for failures could be mentioned here with the
17 !        appropriate error code
18
19 !  User controllable options: <if applicable>
21 module module_nudging_utils
23 real :: totalNudgeTime
24 integer :: sysClockCountRate, sysClockCountMax
25 character(len=4) :: clockType
27 contains
29 !===================================================================================================
30 ! NOTE for whichUtilites
31 ! whUniLoop was fastest for single index searches.
32 ! I still havent tested multiple index searches (which and whichLoop)
35 !===================================================================================================
36 ! Program Names: 
37 !   functions: whichPack and whichLoop 
38 ! Author(s)/Contact(s): 
39 !   James L McCreight <jamesmcc><ucar><edu>
40 ! Abstract: 
41 !   Identify indices in a vector which are TRUE, reutrns zero length vector
42 !   if there are no matches.
43 ! History Log: 
44 !   6/04/15 -Created, JLM.
45 ! Usage: 
46 ! Parameters:  
47 ! Input Files:
48 ! Output Files: 
49 ! Condition codes: 
50 ! User controllable options: None. 
51 ! Notes: 
52 !   JLM: Recent catastrophic failure reported for pack on ifort, with work arround.
53 !   JLM: https://software.intel.com/en-us/forums/topic/559308#comments
55 subroutine whichPack(theMask, which, nWhich)
56 implicit none
57 logical, allocatable, intent(in),  dimension(:) :: theMask
58 integer,              intent(out), dimension(:) :: which
59 integer,              intent(out)               :: nwhich
61 integer :: ii
62 which = -9999
63 nWhich = sum( (/ (1, ii=1,size(theMask)) /), mask=theMask)
64 if(nWhich .gt. size(which)) then
65    which = -9999
66    return
67 end if
68 which(1:nWhich) = pack( (/ (ii, ii=1,size(theMask)) /), mask=theMask)
69 end subroutine whichPack
71 subroutine whichLoop(theMask, which, nWhich)
72 implicit none
73 logical, intent(in),  dimension(:) :: theMask
74 integer, intent(out), dimension(:) :: which
75 integer, intent(out)               :: nwhich
77 integer :: ii
78 which = -9999
79 nWhich = 1
80 do ii=1,size(theMask)
81    if(nWhich .gt. size(which)) then
82       which = -9999
83       return
84    end if
85    if(theMask(ii)) then 
86       which(nWhich)=ii
87       nWhich = nWhich + 1
88    endif
89 end do
90 nWhich = nWhich-1
91 end subroutine whichLoop
93 !===================================================================================================
94 ! Program Names: 
95 !   function: whUnique
96 ! Author(s)/Contact(s): 
97 !   James L McCreight <jamesmcc><ucar><edu>
98 ! Abstract: 
99 !   Identify THE index in a logical vector which is TRUE. Returns
100 !   -1 if not unique or none are true.
101 ! History Log: 
102 !   6/04/15 -Created, JLM.
103 ! Usage: 
104 ! Parameters:  
105 ! Input Files:
106 ! Output Files: 
107 ! Condition codes: 
108 ! User controllable options: None. 
109 ! Notes: 
111 function whUnique(theMask, unsafe)
112   implicit none
113   integer                             :: whUnique !! return value
114   logical, allocatable, dimension(:), intent(in)  :: theMask
115   logical, optional, intent(in)  :: unsafe
116   integer, allocatable, dimension(:) :: whUniques
117   integer :: i, nMatches
118   if(present(unsafe)) then
119      !whUniques=pack( (/ (i, i=1,size(theMask)) /), mask= theMask)     
120      !whUnique = whUniques(1)
121      whUnique=sum( (/ (i, i=1,size(theMask)) /), mask= theMask)     
122   else 
123      nMatches = sum( (/ (1, i=1,size(theMask)) /), mask= theMask )
124      if (nMatches .gt. 1 .OR. nMatches .eq. 0) then
125         whUnique=-1
126      else 
127         whUnique=sum( (/ (i, i=1,size(theMask)) /), mask= theMask)     
128      end if
129   end if
130 end function whUnique
133 !===================================================================================================
134 ! Program Names: 
135 !   function: whUnique
136 ! Author(s)/Contact(s): 
137 !   James L McCreight <jamesmcc><ucar><edu>
138 ! Abstract: 
139 !   Simply returns the first match, no check for uniques. On gfortran this
140 !    was the fastest of the bunch even/especially for max indices on huge arrays.
141 ! History Log: 
142 !   6/04/15 -Created, JLM.
143 ! Usage: 
144 ! Parameters:  
145 ! Input Files:
146 ! Output Files: 
147 ! Condition codes: 
148 ! User controllable options: None. 
149 ! Notes: 
151 function whUniLoop(theMask)
152   implicit none
153   integer                                         :: whUniLoop !! return value
154   logical, allocatable, dimension(:), intent(in)  :: theMask
155   integer :: ii
156   whUniLoop = -9999
157   do ii=1,size(theMask)
158      if(theMask(ii)) then
159         whUniLoop = ii
160         return
161      end if
162   end do
163 end function whUniLoop
165 !===================================================================================================
166 ! Program Names: 
167 !   function: whInLoop
168 ! Author(s)/Contact(s): 
169 !   James L McCreight <jamesmcc><ucar><edu>
170 ! Abstract: 
171 !   Identify the indices of elements in a first vector which are present in the 
172 !   second vector, returns 0 for no matches. This can be slow, it's a double do/for loop.
173 ! History Log: 
174 !   6/04/15 -Created, JLM.
175 ! Usage: 
176 ! Parameters:  
177 ! Input Files:
178 ! Output Files: 
179 ! Condition codes: 
180 ! User controllable options: None. 
181 ! Notes: Can be slow, use with caution.
183 ! parallelize this? ||||||||||||||||||||||||||||||||||
184 subroutine whichInLoop(vecToSearch, vecToMatch, which, nWhich)
185 implicit none
186 character(len=15), intent(in),  dimension(:) :: vecToSearch
187 character(len=15), intent(in),  dimension(:) :: vecToMatch
188 integer, intent(out), dimension(:) :: which
189 integer, intent(out)               :: nWhich
190 integer :: ii, jj
191 which = -9999
192 nWhich = 0
193 do ii=1,size(vecToSearch)
194    do jj=1,size(vecToMatch)
195       if(trim(adjustl(vecToSearch(ii))) .eq. trim(adjustl(vecToMatch(jj)))) then
196          which(ii)=ii
197          nWhich = nWhich + 1
198          exit
199       end if
200    end do
201 end do
202 end subroutine whichInLoop
205 ! parallelize this? ||||||||||||||||||||||||||||||||||
206 subroutine whichInLoop2(vecToSearch, vecToMatch, which, nWhich)
207 implicit none
208 character(len=15), intent(in),  dimension(:) :: vecToSearch
209 character(len=15), intent(in),  dimension(:) :: vecToMatch
210 integer, intent(out), dimension(:) :: which
211 integer, intent(out)               :: nWhich
212 integer :: ii, jj
213 which = -9999
214 nWhich = 0
215 do ii=1,size(vecToSearch)
216    if(any(vecToMatch .eq. vecToSearch(ii))) then
217       which(ii)=ii
218       nWhich = nWhich + 1
219    end if
220 end do
221 end subroutine whichInLoop2
224 !===================================================================================================
225 ! Program Names: 
226 !   accum_nudging_time
227 ! Author(s)/Contact(s): 
228 !   James L McCreight <jamesmcc><ucar><edu>
229 ! Abstract: 
230 !   Tally up the total cpu or wall time used by nudging.
231 ! History Log: 
232 !   8/20/15 -Created, JLM.
233 ! Usage: 
234 ! Parameters:  
235 !   start, end: real times for end-diff timing & accumulation
236 !   sectionLabel: prints a message with the timing for the section
237 !      print*, 'Ndg: ' // sectionLabel // '(seconds ' // trim(clockType) // ' time):', diff
238 !   optional - accum: accumulate this towards the overall time or simply print the above 
239 !      message? Do not accum for nested sections of code, but still give the diagnostic.
240 ! Input Files:
241 ! Output Files: 
242 ! Condition codes: 
243 ! User controllable options: None. 
244 ! Notes: 
245 subroutine accum_nudging_time(start, end, sectionLabel, accum)
246 implicit none
247 real,              intent(in) :: start, end
248 character(len=*),  intent(in) :: sectionLabel
249 logical, optional, intent(in):: accum
250 logical :: accumLocal
251 real :: diff
252 accumLocal=.TRUE.
253 if(present(accum)) accumLocal = accum
254 diff=end-start
255 if(clockType.eq.'wall') then
256    if(diff .lt. 0) diff = diff + sysClockCountMax
257    diff=diff/sysClockCountRate
258 endif
259 if (accumLocal) totalNudgeTime = totalNudgeTime + diff
261 print*,'Ndg: Timing: ' // sectionLabel // ' (' // trim(clockType) // ' time, seconds):', diff
263 if(accumLocal) print*,'Ndg: Timing: accum totalNudgeTime: ',totalNudgeTime
264 end subroutine accum_nudging_time
267 !===================================================================================================
268 ! Program Names: 
269 !   nudging_timer
270 ! Author(s)/Contact(s): 
271 !   James L McCreight <jamesmcc><ucar><edu>
272 ! Abstract: 
273 !   Return your choice of cpu time or wall time
274 ! History Log: 
275 !   8/20/15 -Created, JLM.
276 ! Usage: 
277 ! Parameters:  
278 ! Input Files:
279 ! Output Files: 
280 ! Condition codes: 
281 ! User controllable options: None. 
282 ! Notes: 
283 subroutine nudging_timer(time)
284 implicit none
285 real, intent(inout) :: time
286 integer :: count 
287 if(clockType.eq.'cpu') call cpu_time(time)
288 if(clockType.eq.'wall') then
289    call system_clock(count=count) 
290    time=real(count)
291 end if
292 end subroutine nudging_timer
295 !===================================================================================================
296 end module module_nudging_utils