00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 subroutine prism_put ( field_id, date, date_bounds, data_array, info, ierror )
00014
00015
00016
00017 use PRISM_constants
00018 use PRISM_calendar
00019
00020 use PSMILe
00021 use PSMILe_SMIOC, only : sga_smioc_comp, transient
00022
00023 implicit none
00024
00025
00026
00027 Integer, Intent (InOut) :: field_id
00028
00029
00030
00031 Type(PRISM_Time_Struct), Intent (In) :: date
00032
00033
00034
00035 Type(PRISM_Time_Struct), Intent (In) :: date_bounds(2)
00036
00037
00038
00039
00040 Double Precision, Intent (In) :: data_array(*)
00041
00042
00043
00044
00045
00046
00047 integer, Intent (Out) :: info
00048
00049
00050
00051 integer, Intent (Out) :: ierror
00052
00053
00054
00055
00056
00057
00058
00059 DOUBLE PRECISION :: my_dble
00060 Double Precision, Parameter :: dble_huge = huge(my_dble)
00061
00062 Real (PSMILe_float_kind) :: julian_day, julian_dayb(2)
00063 Real (PSMILe_float_kind) :: julian_sec, julian_secb(2)
00064
00065 Real (PSMILe_float_kind) :: julian_day_save, julian_dayb_save(2)
00066 Real (PSMILe_float_kind) :: julian_sec_save, julian_secb_save(2)
00067
00068 Real (PSMILe_float_kind) :: delta_sec
00069
00070 Logical :: action(3)
00071 Logical :: flag(3)
00072
00073 Logical :: precise = .false.
00074
00075 Integer :: i, j
00076 Integer :: add_days
00077 Integer :: lag
00078 Integer, Parameter :: nerrp = 2
00079 Integer :: ierrp (nerrp)
00080 Integer :: il_o, il_omax, nbr_out
00081
00082 Type (transient), pointer :: sga_smioc_transi(:)
00083
00084
00085
00086 Type (Userdef), pointer :: ug
00087 Type (GridFunction), pointer :: fp
00088
00089 Integer :: il_udef_id
00090 Integer :: il_userdef
00091 Integer :: il_dim1, il_dim2, il_nbfld
00092 Integer :: field_id_1, il_side
00093 Integer :: field_id_2, il_fsize, il_gsppp
00094 Integer :: il_size1
00095
00096 Integer, Allocatable :: ila_ch_act(:,:)
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121 Character(len=len_cvs_string), save :: mycvs =
00122 '$Id: prism_put.F90 2805 2010-12-07 10:16:28Z coquart $'
00123
00124
00125 #ifdef VERBOSE
00126 print *, trim(ch_id), ': prism_put: field_id', field_id
00127
00128 call psmile_flushstd
00129 #endif /* VERBOSE */
00130
00131
00132
00133
00134
00135 ierror = 0
00136 info = PRISM_NOACTION
00137 action = .false.
00138 flag = .false.
00139
00140
00141
00142
00143
00144 if ( field_id == PRISM_UNDEFINED ) then
00145 #ifdef VERBOSE
00146 PRINT*, TRIM(ch_id), ': prism_put we return because field_id undefined:', field_id
00147 print *, trim(ch_id), ': prism_put: eof ierror', ierror
00148 call psmile_flushstd
00149 #endif /* VERBOSE */
00150 return
00151 endif
00152
00153 if ( field_id > size (Fields) .or. field_id < 1) then
00154
00155 ierrp (1) = field_id
00156 ierror = PRISM_Error_Arg
00157 call psmile_error ( ierror, 'field_id', &
00158 ierrp, 1, __FILE__, __LINE__ )
00159 return
00160 endif
00161
00162
00163
00164
00165
00166 if ( Fields(field_id)%status /= PSMILe_status_defined ) then
00167 ierrp (1) = field_id
00168
00169 ierror = PRISM_Error_Arg
00170
00171 call psmile_error ( ierror, 'field_id', &
00172 ierrp, 1, __FILE__, __LINE__ )
00173 return
00174 endif
00175
00176
00177
00178
00179
00180 if ( Fields(field_id)%smioc_loc == PRISM_UNDEFINED .or. &
00181 (.not. Fields(field_id)%used_for_coupling .and. &
00182 .not. Fields(field_id)%used_for_io )) then
00183 #ifdef VERBOSE
00184 print *, trim(ch_id), ': prism_put: nothing to do'
00185 print *, trim(ch_id), ': prism_put: eof ierror ', ierror
00186 call psmile_flushstd
00187 #endif /* VERBOSE */
00188 return
00189 endif
00190
00191
00192
00193
00194
00195
00196
00197 call psmile_date2ju ( date, julian_day, julian_sec )
00198 call psmile_date2ju ( date_bounds(1), julian_dayb(1), julian_secb(1))
00199 call psmile_date2ju ( date_bounds(2), julian_dayb(2), julian_secb(2))
00200
00201
00202
00203 if ( julian_dayb(2) < julian_dayb(1) .or. &
00204 ( julian_dayb(1) == julian_dayb(2) .and. &
00205 julian_secb(2) < julian_secb(1) ) ) then
00206
00207 ierrp (1) = field_id
00208
00209 ierror = PRISM_Error_Date
00210
00211 call psmile_error ( ierror, 'upper bound < lower bound', &
00212 ierrp, 1, __FILE__, __LINE__ )
00213 return
00214
00215 endif
00216
00217
00218
00219 if ( ( julian_dayb(1) > julian_day .or. &
00220 julian_day > julian_dayb(2) ) .or. &
00221 ( julian_dayb(1) == julian_day .and. &
00222 julian_sec < julian_secb(1) ) .or. &
00223 ( julian_dayb(2) == julian_day .and. &
00224 julian_sec > julian_secb(2) ) ) then
00225
00226 ierrp (1) = field_id
00227
00228 ierror = PRISM_Error_Date
00229
00230 call psmile_error ( ierror, 'date out of bounds', &
00231 ierrp, 1, __FILE__, __LINE__ )
00232 return
00233
00234 endif
00235
00236
00237
00238
00239
00240 field_id_1 = field_id
00241 fp => Fields(field_id)
00242
00243
00244 nbr_out = 0
00245 if ( Associated(fp%Taskout) ) then
00246 nbr_out = size (fp%Taskout)
00247 endif
00248 #ifdef DEBUG
00249 print *, trim(ch_id),': prism_put: Field_id = ',field_id, nbr_out
00250 #endif /* DEBUG */
00251 #ifdef PRISM_ASSERTION
00252 if (nbr_out == 0) then
00253 print *, trim(ch_id), ": prism_put: nbr_out", nbr_out
00254 call PSMILE_Assert (__FILE__, __LINE__, &
00255 "nbr_out == 0")
00256 endif
00257 #endif
00258
00259
00260
00261 Allocate ( ila_ch_act(nbr_out,4), STAT=ierror )
00262 if ( ierror > 0 ) then
00263 ierrp (1) = ierror
00264 ierrp (2) = nbr_out
00265 call PSMILe_Error ( PRISM_Error_Alloc, 'ila_ch_act', &
00266 ierrp, 2, __FILE__, __LINE__ )
00267 return
00268 endif
00269
00270
00271
00272 do il_o = 1, nbr_out
00273
00274 ila_ch_act(il_o,1) = field_id
00275 ila_ch_act(il_o,2) = il_o
00276 ila_ch_act(il_o,3) = PSMILe_false
00277
00278 il_udef_id = fp%Taskout(il_o)%userdef_id
00279 ila_ch_act(il_o,4) = il_udef_id
00280
00281 if ( il_udef_id /= PSMILe_undef ) then
00282
00283 ug => Userdefs(il_udef_id)
00284 il_side = ug%ig_transi_side
00285 il_dim1 = size ( fp%var_shape(:,:), dim=1 )
00286 il_dim2 = size ( fp%var_shape(:,:), dim=2 )
00287 field_id_2 = fp%Taskout(1)%assoc_var_id
00288
00289 ila_ch_act(il_o,1) = field_id_2
00290 ila_ch_act(il_o,2) = 1
00291 ila_ch_act(il_o,3) = PSMILe_true
00292 ila_ch_act(il_o,4) = il_udef_id
00293 #ifdef DEBUG
00294 PRINT *, ' Field_id, Userdef_id field_id_2 = ', field_id, il_udef_id, field_id_2
00295 PRINT *, ' il_side, number of dimensions of data array = ',il_side, il_dim1, il_dim2
00296 PRINT *, ' content of var_shape array (geogr.) = ', fp%var_shape(:,:)
00297 #endif
00298
00299
00300
00301 il_fsize = Fields(field_id_2)%size
00302 il_gsppp = ug%ig_nb_ppp
00303 il_nbfld = ug%ig_nbr_fields
00304
00305 il_size1 = il_fsize / il_nbfld
00306 #ifdef DEBUG
00307 PRINT *, ' prism_put : il_nbfld il_fsize = ', il_nbfld, il_fsize
00308 PRINT *, ' prism_put : il_size1 il_gsppp = ', il_size1, il_gsppp
00309 #endif
00310 #ifdef PRISM_ASSERTION
00311 if ( il_size1 /= il_gsppp ) then
00312 call PSMILe_Assert (__FILE__, __LINE__, "Incorrect data size")
00313 endif
00314 #endif
00315
00316
00317
00318 if ( Fields(field_id)%dataType == PRISM_Real ) then
00319 Allocate ( ug%real_gridless(1:il_gsppp,1,1,il_nbfld), STAT=ierror )
00320 elseif ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00321 Allocate ( ug%dble_gridless(1:il_gsppp,1,1,il_nbfld), STAT=ierror )
00322 endif
00323
00324
00325 if ( fp%dataType == PRISM_Real ) then
00326 call psmile_gridless_func_real ( field_id, il_udef_id, il_side, &
00327 data_array, ierror )
00328 else if ( fp%dataType == PRISM_Double_Precision ) then
00329 call psmile_gridless_func_dble ( field_id, il_udef_id, il_side, &
00330 data_array, ierror)
00331 endif
00332
00333 Nullify (ug)
00334
00335 endif
00336
00337 enddo
00338
00339 sga_smioc_transi => sga_smioc_comp(Fields(field_id)%comp_id)%sga_smioc_transi
00340
00341
00342
00343
00344
00345 julian_day_save = julian_day
00346 julian_dayb_save = julian_dayb
00347 julian_sec_save = julian_sec
00348 julian_secb_save = julian_secb
00349
00350 il_omax = sga_smioc_transi(Fields(field_id)%smioc_loc)%ig_nb_transi_out
00351 do il_o = 1, il_omax
00352
00353 #ifdef DEBUG
00354 print *, " prism_put : il_omax nbr_out = ",il_omax, nbr_out
00355 #endif
00356 julian_day = julian_day_save
00357 julian_dayb = julian_dayb_save
00358 julian_sec = julian_sec_save
00359 julian_secb = julian_secb_save
00360
00361
00362
00363 field_id = ila_ch_act(il_o,1)
00364 i = ila_ch_act(il_o,2)
00365 il_userdef = ila_ch_act(il_o,3)
00366 il_udef_id = ila_ch_act(il_o,4)
00367
00368
00369
00370
00371
00372
00373 if ( Fields(field_id)%Taskout(i)%Judate_Ubnd%days /= dble_huge ) then
00374 if ( julian_dayb(1) /= Fields(field_id)%Taskout(i)%Judate_Ubnd%days .or. &
00375 julian_secb(1) /= Fields(field_id)%Taskout(i)%Judate_Ubnd%secs ) then
00376
00377 print *, trim(ch_id), ': prism_put upper bound from previous put:'
00378 print *, trim(ch_id), ': - days: ', Fields(field_id)%Taskout(i)%Judate_Ubnd%days
00379 print *, trim(ch_id), ': - secs: ', Fields(field_id)%Taskout(i)%Judate_Ubnd%secs
00380
00381 print *, trim(ch_id), ': prism_put lower bound from this call:'
00382 print *, trim(ch_id), ': - days: ', julian_dayb(1)
00383 print *, trim(ch_id), ': - secs: ', julian_secb(1)
00384
00385 ierrp (1) = field_id
00386 ierror = PRISM_Error_Date
00387 call psmile_error ( ierror, 'inconsistent date bounds for prism_put', &
00388 ierrp, 1, __FILE__, __LINE__ )
00389 return
00390 endif
00391 endif
00392
00393
00394
00395 Fields(field_id)%Taskout(i)%Judate_Lbnd%days=julian_dayb(1)
00396 Fields(field_id)%Taskout(i)%Judate_Lbnd%secs=julian_secb(1)
00397
00398 Fields(field_id)%Taskout(i)%Judate_Ubnd%days=julian_dayb(2)
00399 Fields(field_id)%Taskout(i)%Judate_Ubnd%secs=julian_secb(2)
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409 lag = sga_smioc_transi(Fields(field_id)%smioc_loc)%sga_transi_out(i)%ig_lag
00410
00411 if ( lag /= PSMILe_undef ) then
00412
00413 #ifdef DEBUG
00414 print *, trim(ch_id), ': prism_put: lag is active and set to ', lag
00415 #endif
00416 delta_sec = (julian_dayb(2) - julian_dayb(1)) * 86400.0 &
00417 + julian_secb(2) - julian_secb(1)
00418
00419 julian_sec = julian_sec + lag * delta_sec
00420 add_days = floor(julian_sec / 86400.0)
00421 julian_day = julian_day + add_days
00422 julian_sec = julian_sec - float(add_days) * 86400.0
00423
00424 julian_secb(1) = julian_secb(1) + lag * delta_sec
00425 add_days = floor(julian_secb(1) / 86400.0)
00426 julian_dayb(1) = julian_dayb(1) + add_days
00427 julian_secb(1) = julian_secb(1) - float(add_days) * 86400.0
00428
00429 julian_secb(2) = julian_secb(2) + lag * delta_sec
00430 add_days = floor(julian_secb(2) / 86400.0)
00431 julian_dayb(2) = julian_dayb(2) + add_days
00432 julian_secb(2) = julian_secb(2) - float(add_days) * 86400.0
00433
00434 endif
00435
00436
00437
00438
00439
00440
00441 call psmile_check_action ( field_id, i, precise, &
00442 julian_day, julian_dayb(1), &
00443 julian_sec, julian_secb(1), &
00444 action )
00445
00446
00447 do j = 1, 3
00448 if (action(j)) flag(j) = .true.
00449 enddo
00450
00451
00452
00453
00454
00455
00456 if ( il_userdef .EQ. PSMILe_true ) then
00457
00458 ug => Userdefs(il_udef_id)
00459
00460
00461 if ( Fields(field_id)%dataType == PRISM_Real ) then
00462 call psmile_put_real ( field_id, i, julian_day, julian_sec, &
00463 julian_dayb, julian_secb, ug%real_gridless, action, info, ierror )
00464 else if ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00465 call psmile_put_dble ( field_id, i, julian_day, julian_sec, &
00466 julian_dayb, julian_secb, ug%dble_gridless, action, info, ierror )
00467 endif
00468
00469
00470
00471
00472 if ( Fields(field_id)%dataType == PRISM_Real ) then
00473 Deallocate (ug%real_gridless, STAT=ierror)
00474 else if ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00475 Deallocate (ug%dble_gridless, STAT=ierror)
00476 endif
00477
00478 if (ierror > 0) then
00479 ierrp (1) = ierror
00480 ierrp (2) = il_udef_id
00481 ierror = PRISM_Error_Dealloc
00482 call psmile_error ( ierror, 'Userdefs', &
00483 ierrp, 2, __FILE__, __LINE__ )
00484 return
00485 endif
00486
00487 Nullify (ug%real_gridless)
00488 Nullify (ug%dble_gridless)
00489
00490 else
00491
00492 if ( Fields(field_id)%dataType == PRISM_Integer ) then
00493 call psmile_put_int ( field_id, i, julian_day, julian_sec, &
00494 julian_dayb, julian_secb, data_array, action, info, ierror )
00495
00496 else if ( Fields(field_id)%dataType == PRISM_Real ) then
00497 call psmile_put_real ( field_id, i, julian_day, julian_sec, &
00498 julian_dayb, julian_secb, data_array, action, info, ierror )
00499
00500 else if ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00501 call psmile_put_dble ( field_id, i, julian_day, julian_sec, &
00502 julian_dayb, julian_secb, data_array, action, info, ierror )
00503 endif
00504
00505 endif
00506
00507 enddo
00508
00509 info = min(info,1)
00510
00511 if ( flag(1) ) info = info + 1000
00512 if ( flag(2) ) info = info + 100
00513 if ( flag(3) ) info = info + 10
00514
00515
00516
00517 Deallocate ( ila_ch_act, STAT=ierror )
00518 if ( ierror > 0 ) then
00519 ierrp (1) = ierror
00520 call PSMILe_Error ( PRISM_Error_Dealloc, 'ila_ch_act', &
00521 ierrp, 1, __FILE__, __LINE__ )
00522 return
00523 endif
00524
00525
00526
00527 Nullify (fp)
00528 field_id = field_id_1
00529
00530
00531 #ifdef VERBOSE
00532 print *, trim(ch_id), ': prism_put: eof ierror ', ierror
00533 call psmile_flushstd
00534 #endif /* VERBOSE */
00535
00536 end subroutine prism_put