00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 subroutine prism_put_inquire ( field_id, date, date_bounds, info, ierror )
00014
00015
00016
00017 use PRISM, dummy_interface => prism_put_inquire
00018 use PRISM_calendar
00019 use PSMILe
00020 use PSMILe_SMIOC, only : sga_smioc_comp, transient
00021
00022 implicit none
00023
00024
00025
00026 Integer, Intent (In) :: field_id
00027
00028
00029
00030 Type(PRISM_Time_Struct), Intent (In) :: date
00031
00032
00033
00034 Type(PRISM_Time_Struct), Intent (In) :: date_bounds(2)
00035
00036
00037
00038
00039
00040
00041
00042
00043 Integer, Intent (Out) :: info
00044
00045
00046
00047 Integer, Intent (Out) :: ierror
00048
00049
00050
00051
00052
00053
00054
00055 Type (GridFunction), Pointer :: fp
00056 Type (transient), Pointer :: sga_smioc_transi(:)
00057
00058 Double Precision :: julian_day, julian_dayb(2)
00059 Double Precision :: julian_sec, julian_secb(2)
00060
00061 Double Precision :: delta_sec
00062
00063 Logical :: action(3)
00064 Logical :: flag(3)
00065
00066
00067 Logical :: precise = .false.
00068
00069 Integer :: add_days
00070 Integer :: lag
00071
00072 Integer :: i, j
00073 Integer :: il_nb_transi_out
00074
00075 Integer, Parameter :: nerrp = 2
00076 Integer :: ierrp (nerrp)
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097 Character(len=len_cvs_string), save :: mycvs =
00098 '$Id: prism_put_inquire.F90 2325 2010-04-21 15:00:07Z valcke $'
00099
00100
00101
00102 #ifdef VERBOSE
00103 print *, trim(ch_id), ': prism_put_inquire: field_id', field_id
00104 call psmile_flushstd
00105 #endif /* VERBOSE */
00106
00107
00108
00109
00110
00111 ierror = 0
00112 info = PRISM_NOACTION
00113 flag = .false.
00114 action = .false.
00115
00116
00117
00118
00119
00120 if ( Fields(field_id)%status /= PSMILe_status_defined ) then
00121
00122 ierror = PRISM_Error_Arg
00123
00124 print *, trim(ch_id), ': prism_put_inquire: eof field_id not defined'
00125 call psmile_flushstd
00126
00127 return
00128
00129 endif
00130
00131 fp => Fields(field_id)
00132 sga_smioc_transi => sga_smioc_comp(Fields(field_id)%comp_id)%sga_smioc_transi
00133 il_nb_transi_out = sga_smioc_transi(fp%smioc_loc)%ig_nb_transi_out
00134
00135
00136
00137
00138
00139 if ( il_nb_transi_out < 1 ) then
00140 #ifdef VERBOSE
00141 print *, trim(ch_id), ': prism_put_inquire: eof ig_nb_transi_out ', il_nb_transi_out
00142 call psmile_flushstd
00143 #endif /* VERBOSE */
00144 return
00145 endif
00146
00147 if ( fp%smioc_loc == PRISM_UNDEFINED ) then
00148
00149 ierror = PRISM_Error_Arg
00150
00151 print *, trim(ch_id), ': prism_put_inquire: WARNING: smioc_loc undefined'
00152 call psmile_flushstd
00153
00154 return
00155
00156 endif
00157
00158
00159
00160
00161
00162
00163
00164 call psmile_date2ju ( date, julian_day, julian_sec )
00165 call psmile_date2ju ( date_bounds(1), julian_dayb(1), julian_secb(1))
00166 call psmile_date2ju ( date_bounds(2), julian_dayb(2), julian_secb(2))
00167
00168
00169
00170 if ( julian_dayb(2) < julian_dayb(1) .or. &
00171 ( julian_dayb(1) == julian_dayb(2) .and. &
00172 julian_secb(2) < julian_secb(1) ) ) then
00173
00174 ierror = PRISM_Error_Date
00175
00176 print *, trim(ch_id), ': prism_put_inquire: WARNING: upper bound < lower bound'
00177 call psmile_flushstd
00178
00179 return
00180
00181 endif
00182
00183
00184
00185 if ( ( julian_dayb(1) > julian_day .or. &
00186 julian_day > julian_dayb(2) ) .or. &
00187 ( julian_dayb(1) == julian_day .and. &
00188 julian_sec < julian_secb(1) ) .or. &
00189 ( julian_dayb(2) == julian_day .and. &
00190 julian_sec > julian_secb(2) ) ) then
00191
00192 ierrp (1) = field_id
00193
00194 ierror = PRISM_Error_Date
00195
00196 print *, trim(ch_id), ': prism_put_inquire: WARNING: date out of bounds'
00197 call psmile_flushstd
00198
00199 return
00200
00201 endif
00202
00203
00204
00205
00206
00207 do i = 1, il_nb_transi_out
00208 if ( sga_smioc_transi(fp%smioc_loc)%sga_transi_out(i)%ig_src_timeop /= PSMILe_undef ) then
00209 info = info + 1
00210 #ifdef VERBOSE
00211 print *, trim(ch_id), ': prism_put_inquire: eof summation required'
00212 call psmile_flushstd
00213 #endif /* VERBOSE */
00214 exit
00215 endif
00216 enddo
00217
00218
00219
00220
00221
00222 do i = 1, il_nb_transi_out
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232 lag = sga_smioc_transi(fp%smioc_loc)%sga_transi_out(i)%ig_lag
00233
00234 if ( lag /= PSMILe_undef ) then
00235
00236 #ifdef DEBUG
00237 print *, trim(ch_id), ': prism_put_inquire: lag is active and set to ', lag
00238 #endif
00239 delta_sec = (julian_dayb(2) - julian_day) * 86400.0 &
00240 + julian_secb(2) - julian_secb(1)
00241
00242 julian_sec = julian_sec + lag * delta_sec
00243 add_days = floor(julian_sec / 86400.0)
00244 julian_day = julian_day + add_days
00245 julian_sec = julian_sec - float(add_days) * 86400.0
00246
00247 julian_secb(1) = julian_secb(1) + lag * delta_sec
00248 add_days = floor(julian_secb(1) / 86400.0)
00249 julian_dayb(1) = julian_dayb(1) + add_days
00250 julian_secb(1) = julian_secb(1) - float(add_days) * 86400.0
00251
00252 julian_secb(2) = julian_secb(2) + lag * delta_sec
00253 add_days = floor(julian_secb(2) / 86400.0)
00254 julian_dayb(2) = julian_dayb(2) + add_days
00255 julian_secb(2) = julian_secb(2) - float(add_days) * 86400.0
00256
00257 endif
00258
00259
00260
00261
00262
00263
00264 call psmile_check_action ( field_id, i, precise, &
00265 julian_day, julian_dayb(1), &
00266 julian_sec, julian_secb(1), &
00267 action )
00268
00269
00270 do j = 1, 3
00271 if (action(j)) flag(j) = .true.
00272 enddo
00273
00274 enddo
00275
00276 if ( flag(1) ) info = info + 1000
00277 if ( flag(2) ) info = info + 100
00278 if ( flag(3) ) info = info + 10
00279
00280 #ifdef VERBOSE
00281 print *, trim(ch_id), ': prism_put_inquire: eof ierror ', ierror
00282 call psmile_flushstd
00283 #endif /* VERBOSE */
00284
00285 end subroutine prism_put_inquire