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