00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_put_dble ( field_id, task_id, julian_day, julian_sec, &
00012 julian_dayb, julian_secb, data_array, action, info, ierror )
00013
00014
00015
00016 use PRISM_constants
00017 use PRISM_calendar
00018
00019 use PSMILe
00020 use PSMILe_SMIOC, only : sga_smioc_comp, transient, transient_out
00021
00022 implicit none
00023
00024
00025
00026 Integer, Intent (In) :: field_id
00027
00028
00029
00030 Integer, Intent (In) :: task_id
00031
00032
00033
00034 Real (PSMILe_float_kind), Intent (In) :: julian_day, julian_dayb(2)
00035 Real (PSMILe_float_kind), Intent (In) :: julian_sec, julian_secb(2)
00036
00037
00038
00039
00040
00041 Double Precision, Intent (In) :: data_array(*)
00042
00043
00044
00045 Logical, Intent (In) :: action(3)
00046
00047
00048
00049
00050
00051
00052
00053 Integer, Intent (InOut) :: info
00054
00055
00056
00057 Integer, Intent (Out) :: ierror
00058
00059
00060
00061
00062
00063
00064
00065 Double Precision :: add_scalar
00066 Double Precision :: mul_scalar
00067
00068 Double Precision, Allocatable :: data_scattered(:)
00069 Double Precision, Allocatable :: data_reduced(:)
00070
00071 Integer :: shape_in(2,6)
00072
00073 Integer :: i, j
00074 Integer :: len_gathered, len_scattered
00075 Integer :: len, len_3d, len_reduced
00076 Integer :: nbr_fields
00077
00078 Integer :: rdim(6)
00079 Integer :: nbr_reductions
00080
00081 Integer :: local_reduce
00082 Integer :: local_timeop
00083
00084 Logical :: local_scatter
00085 Logical :: local_multiply
00086 Logical :: local_add
00087
00088 Logical :: local_operation = .false.
00089 Logical :: mask_set
00090 Logical :: mask_needed
00091 Logical, Pointer :: mask_array(:)
00092 Logical, Pointer :: mask_aux (:)
00093
00094 Type (GridFunction), Pointer :: fp
00095 Type (Taskout_type), Pointer :: Taskout
00096
00097
00098
00099 Type (transient), Pointer :: sga_smioc_transi(:)
00100 Type (transient_out), Pointer :: sga_transi_out
00101
00102
00103
00104 #ifdef __PSMILE_WITH_IO
00105 Real (PSMILe_float_kind) :: jd_end,js_end,delta_sec
00106 Real (PSMILe_float_kind) :: jd_cur,js_cur
00107
00108 Integer :: il_taskid,il_transiouts,il_smioc_loc
00109 Integer :: il_taskid_restr
00110 Integer :: lag,add_days
00111
00112 Logical :: debug_action, restart_action
00113 #endif
00114
00115
00116
00117 Integer :: stat_nsum
00118 Integer :: nsum
00119
00120 Integer :: shape_out(2,6)
00121
00122 Double Precision, Allocatable :: recv_buf(:)
00123 Double Precision, Allocatable :: stat_max(:)
00124 Double Precision, Allocatable :: stat_min(:)
00125 Double Precision, Allocatable :: stat_sum(:)
00126 Double Precision, Allocatable :: stat_mean(:)
00127
00128 Logical :: stats
00129
00130 Character(len=6), Save :: cstats (3)
00131
00132 Integer :: n_grid_dim
00133
00134
00135
00136 Integer :: ierr
00137 Integer, Parameter :: nerrp = 3
00138 Integer :: ierrp (nerrp)
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170 Data cstats /'masked', 'valid', 'all'/
00171
00172 #ifdef VERBOSE
00173 print *, trim(ch_id), ': psmile_put_dble: field_id', field_id
00174
00175 call psmile_flushstd
00176 #endif /* VERBOSE */
00177
00178
00179
00180
00181
00182 ierror = 0
00183 local_reduce = PSMILe_undef
00184
00185 local_operation = .false.
00186
00187 Nullify(mask_array)
00188 Nullify(mask_aux)
00189
00190 fp => Fields(field_id)
00191 Taskout => Fields(field_id)%Taskout(task_id)
00192 sga_smioc_transi => sga_smioc_comp(Fields(field_id)%comp_id)%sga_smioc_transi
00193 sga_transi_out => sga_smioc_transi(fp%smioc_loc)%sga_transi_out(task_id)
00194 n_grid_dim = Grids(Methods(fp%method_id)%grid_id)%n_dim
00195
00196
00197
00198
00199
00200 len = fp%size
00201
00202 if ( fp%transi_type == PSMILe_bundle ) then
00203 nbr_fields = fp%var_shape(2,n_grid_dim+1)
00204 else
00205 nbr_fields = 1
00206 endif
00207
00208 len_3d = len / nbr_fields
00209
00210
00211
00212
00213
00214 local_scatter = sga_transi_out%sg_src_local_trans%ig_scatter == PSMILe_true
00215 local_timeop = sga_transi_out%ig_src_timeop
00216
00217
00218
00219 mul_scalar = sga_transi_out%sg_src_local_trans%dg_mult_scalar
00220 add_scalar = sga_transi_out%sg_src_local_trans%dg_add_scalar
00221
00222 local_add = add_scalar /= PSMILe_dundef
00223 local_multiply = mul_scalar /= PSMILe_dundef
00224
00225
00226
00227 mask_set = fp%mask_id /= PRISM_UNDEFINED
00228
00229 if (mask_set) then
00230 mask_set = Masks(fp%mask_id)%status /= PSMILe_status_free
00231 endif
00232
00233
00234
00235
00236
00237 stats = sga_transi_out%iga_stats(1) == PSMILe_true .or. &
00238 sga_transi_out%iga_stats(2) == PSMILe_true .or. &
00239 sga_transi_out%iga_stats(3) == PSMILe_true
00240
00241 mask_needed = stats .or. local_scatter .or. &
00242 local_reduce == PSMILe_max .or. local_reduce == PSMILe_min
00243
00244 if ( mask_set ) mask_array => Masks(fp%mask_id)%mask_array
00245
00246 if ( mask_needed ) then
00247 allocate (mask_aux(len_3d), STAT=ierr)
00248 if ( ierr /= 0 ) then
00249 ierrp (1) = 1
00250 ierror = PRISM_Error_Alloc
00251 call psmile_error ( ierror, 'mask_aux', ierrp(1), 1, &
00252 __FILE__, __LINE__ )
00253 return
00254 endif
00255
00256 mask_aux = .true.
00257
00258 if (.not. mask_set) mask_array => mask_aux
00259 endif
00260
00261
00262
00263
00264
00265 if ( Taskout%nsum == 0 ) then
00266 Taskout%start_day = julian_dayb(1)
00267 Taskout%start_sec = julian_secb(1)
00268 endif
00269
00270
00271
00272
00273
00274
00275 if ( local_timeop /= PSMILe_undef .or. local_add .or. local_multiply ) then
00276
00277 local_operation = .true.
00278
00279 if ( .not. associated(Taskout%buffer_dble) ) then
00280 Allocate(Taskout%buffer_dble(1:len), STAT = ierr )
00281 if ( ierr /= 0 ) then
00282 ierrp (1) = 1
00283 ierror = PRISM_Error_Alloc
00284 call psmile_error ( ierror, 'Taskout%buffer_dble', &
00285 ierrp(1), 1, __FILE__, __LINE__ )
00286 return
00287 endif
00288
00289 Taskout%buffer_dble = 0.0
00290 Taskout%nsum = 0
00291 endif
00292
00293 endif
00294
00295
00296
00297
00298
00299
00300 if ( local_timeop /= PSMILe_undef ) then
00301
00302 Taskout%nsum=Taskout%nsum+1
00303
00304 Taskout%buffer_dble(1:len) = &
00305 Taskout%buffer_dble(1:len) + data_array(1:len)
00306
00307 info = info + 1
00308
00309 endif
00310
00311
00312
00313
00314
00315 Taskout%end_day = julian_dayb(2)
00316 Taskout%end_sec = julian_secb(2)
00317
00318
00319
00320
00321
00322 if ( .not. action(1) .and. .not. action(2) .and. .not. action(3) ) then
00323 #ifdef VERBOSE
00324 print *, trim(ch_id), ': psmile_put_dble: eof nothing to do! ierror ', ierror
00325 #endif /* VERBOSE */
00326 return
00327 endif
00328
00329
00330
00331
00332
00333
00334
00335 IF ( local_multiply ) THEN
00336
00337
00338
00339 if (( local_timeop /= PSMILe_undef) .OR. ( local_add ) ) then
00340
00341 Taskout%buffer_dble(1:len) = Taskout%buffer_dble(1:len) * mul_scalar
00342
00343 else
00344
00345
00346
00347 Taskout%buffer_dble(1:len) = data_array(1:len) * mul_scalar
00348
00349 endif
00350
00351 ENDIF
00352
00353
00354
00355
00356 if ( local_add ) then
00357
00358
00359
00360 if ( local_timeop /= PSMILe_undef ) then
00361
00362 Taskout%buffer_dble(1:len) = Taskout%buffer_dble(1:len) &
00363 + (add_scalar*Taskout%nsum)
00364
00365 else
00366
00367
00368
00369 Taskout%buffer_dble(1:len) = data_array(1:len) + add_scalar
00370
00371 endif
00372
00373 endif
00374
00375
00376
00377
00378 if ( local_timeop == PSMILe_tave .and. Taskout%nsum > 0 ) then
00379
00380 Taskout%buffer_dble(1:len) = Taskout%buffer_dble(1:len) &
00381 / dble(Taskout%nsum)
00382 endif
00383
00384
00385
00386
00387
00388 if ( local_scatter ) then
00389
00390
00391
00392
00393
00394
00395
00396
00397 if ( .not. mask_set ) then
00398
00399 len_gathered = len_3d
00400 len_scattered = len_gathered
00401
00402 else
00403 len = len_3d
00404 len_gathered = 0
00405 len_scattered = 1
00406
00407 do i = 1, n_grid_dim
00408
00409 len_scattered = len_scattered &
00410 * ( Masks(fp%mask_id)%mask_shape(2,i) - &
00411 Masks(fp%mask_id)%mask_shape(1,i) + 1 )
00412 enddo
00413
00414 if ( len_scattered /= len ) then
00415 ierrp (1) = field_id
00416 ierrp (2) = len
00417 ierrp (3) = len_scattered
00418 ierror = PRISM_Error_Size
00419
00420 call psmile_error ( ierror, fp%local_name, ierrp, 3, &
00421 __FILE__, __LINE__ )
00422 return
00423 endif
00424
00425 len_gathered = count(mask_array)
00426
00427 endif
00428
00429
00430
00431
00432
00433 Allocate (data_scattered(len_scattered), STAT = ierr )
00434
00435 if ( ierr /= 0 ) then
00436 ierrp (1) = len_scattered
00437 ierror = PRISM_Error_Alloc
00438 call psmile_error ( ierror, 'data_scattered', &
00439 ierrp(1), 1, __FILE__, __LINE__ )
00440 return
00441 endif
00442
00443 if ( local_operation ) then
00444
00445 call psmile_loc_trans_dble ( PSMILe_scat, nbr_fields, &
00446 len_gathered, Taskout%buffer_dble, &
00447 len_scattered, data_scattered, field_id )
00448 else
00449
00450 call psmile_loc_trans_dble ( PSMILe_scat, nbr_fields, &
00451 len_gathered, data_array, &
00452 len_scattered, data_scattered, field_id )
00453
00454 endif
00455
00456 endif
00457
00458
00459
00460
00461
00462
00463
00464 shape_in = 1
00465
00466 do i = 1, n_grid_dim
00467 shape_in(1,i) = fp%var_shape(1,i)
00468 shape_in(2,i) = fp%var_shape(2,i)
00469 enddo
00470
00471 if ( Fields(field_id)%transi_type == PSMILe_bundle ) &
00472 shape_in(2,6) = nbr_fields
00473
00474
00475
00476
00477
00478 if ( local_reduce == PSMILe_max .or. local_reduce == PSMILe_min ) then
00479
00480
00481
00482 shape_out = shape_in
00483 len_reduced = 1
00484
00485 nbr_reductions = 0
00486 rdim = 6
00487
00488 do i = 1, nbr_reductions
00489 shape_out(1,rdim(i)) = 1
00490 shape_out(2,rdim(i)) = 1
00491 enddo
00492
00493 do i = 1, 6
00494 len_reduced = len_reduced * ( shape_out(2,i)-shape_out(1,i) + 1 )
00495 enddo
00496
00497 allocate (data_reduced(len_reduced), STAT=ierr)
00498 if ( ierr /= 0 ) then
00499 ierrp (1) = len_reduced
00500 ierror = PRISM_Error_Alloc
00501 call psmile_error ( ierror, 'data_reduced', ierrp(1), 1, &
00502 __FILE__, __LINE__ )
00503 return
00504 endif
00505
00506 if ( local_scatter ) then
00507
00508
00509
00510 call psmile_multi_reduce_dble ( local_reduce, shape_in, &
00511 data_scattered, shape_out, data_reduced, mask_array, ierror )
00512 else
00513
00514 if ( local_operation ) then
00515
00516
00517
00518 call psmile_multi_reduce_dble ( local_reduce, shape_in, &
00519 Taskout%buffer_dble, &
00520 shape_out, data_reduced, mask_array, ierror )
00521 else
00522
00523
00524
00525 call psmile_multi_reduce_dble ( local_reduce,shape_in, &
00526 data_array, shape_out, data_reduced, mask_array, ierror )
00527
00528 endif
00529
00530 endif
00531
00532 endif
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543 if (stats) then
00544
00545 shape_out = 1
00546 shape_out(2,6) = nbr_fields
00547
00548
00549
00550 allocate (recv_buf(1:shape_out(2,6)), &
00551 stat_max(1:shape_out(2,6)), &
00552 stat_min(1:shape_out(2,6)), &
00553 stat_sum(1:shape_out(2,6)), &
00554 stat_mean(1:shape_out(2,6)), STAT=ierr)
00555 if ( ierr /= 0 ) then
00556 ierrp (1) = nbr_fields * 5
00557 ierror = PRISM_Error_Alloc
00558 call psmile_error ( ierror, 'recv_buf', ierrp(1), 1, __FILE__, __LINE__ )
00559 return
00560 endif
00561
00562
00563
00564 do j = 1, 3
00565
00566 if ( sga_transi_out%iga_stats(j) == PSMILe_true ) then
00567
00568
00569
00570
00571
00572 if ( mask_set ) then
00573
00574 select case (j)
00575 case (1)
00576 mask_aux = .not. Masks(fp%mask_id)%mask_array
00577 mask_array => mask_aux
00578 case (2)
00579 mask_array => Masks(fp%mask_id)%mask_array
00580 case (3)
00581 mask_aux = .true.
00582 mask_array => mask_aux
00583 end select
00584
00585 else
00586
00587 mask_aux = j > 1
00588 mask_array => mask_aux
00589
00590 endif
00591
00592 if ( local_scatter ) then
00593
00594
00595
00596 call psmile_multi_reduce_dble ( PSMILe_max, shape_in, &
00597 data_scattered, shape_out, stat_max, mask_array, ierror )
00598
00599 call psmile_multi_reduce_dble ( PSMILe_min, shape_in, &
00600 data_scattered, shape_out, stat_min, mask_array, ierror )
00601
00602 call psmile_multi_reduce_dble ( PSMILe_integral, shape_in, &
00603 data_scattered, shape_out, stat_sum, mask_array, ierror )
00604 else
00605
00606 if ( local_operation ) then
00607
00608
00609
00610 call psmile_multi_reduce_dble ( PSMILe_max, shape_in, &
00611 Taskout%buffer_dble, &
00612 shape_out, stat_max, mask_array, ierror )
00613
00614 call psmile_multi_reduce_dble ( PSMILe_min, shape_in, &
00615 Taskout%buffer_dble, &
00616 shape_out, stat_min, mask_array, ierror )
00617
00618 call psmile_multi_reduce_dble ( PSMILe_integral, shape_in, &
00619 Taskout%buffer_dble, &
00620 shape_out, stat_sum, mask_array, ierror )
00621 else
00622
00623
00624
00625 call psmile_multi_reduce_dble ( PSMILe_max, &
00626 shape_in, data_array, shape_out, stat_max, &
00627 mask_array, ierror )
00628
00629 call psmile_multi_reduce_dble ( PSMILe_min, &
00630 shape_in, data_array, shape_out, stat_min, &
00631 mask_array, ierror )
00632
00633 call psmile_multi_reduce_dble ( PSMILe_integral, &
00634 shape_in, data_array, shape_out, stat_sum, &
00635 mask_array, ierror )
00636
00637 endif
00638
00639 endif
00640
00641 if ( Comps(fp%comp_id)%act_comm /= MPI_COMM_NULL ) then
00642
00643 call MPI_Allreduce ( stat_max, recv_buf, shape_out(2,6), MPI_DOUBLE_PRECISION, &
00644 MPI_MAX, Comps(fp%comp_id)%act_comm, ierror )
00645
00646 stat_max (:) = recv_buf(:)
00647
00648 call MPI_Allreduce ( stat_min, recv_buf, shape_out(2,6), MPI_DOUBLE_PRECISION, &
00649 MPI_MIN, Comps(fp%comp_id)%act_comm, ierror )
00650
00651 stat_min (:) = recv_buf(:)
00652
00653 call MPI_Allreduce ( stat_sum, recv_buf, shape_out(2,6), MPI_DOUBLE_PRECISION, &
00654 MPI_SUM, Comps(fp%comp_id)%act_comm, ierror )
00655
00656 stat_sum (:) = recv_buf(:)
00657
00658
00659 if (j == 3) then
00660 stat_nsum = len_3d
00661 else
00662 stat_nsum = count(mask_array)
00663 endif
00664
00665 call MPI_Allreduce ( stat_nsum, nsum, 1, MPI_INTEGER, &
00666 MPI_SUM, Comps(fp%comp_id)%act_comm, ierror )
00667 if (nsum > 0) then
00668 stat_mean(:) = stat_sum(:) / nsum
00669 else
00670 stat_mean = 0.0
00671 endif
00672
00673 endif
00674
00675 write (*, 9990) trim(ch_id)
00676 write (*, 9980) trim(ch_id), trim(cstats (j))
00677
00678 write (*,*) trim(ch_id), &
00679 ': ... for field ', trim(fp%local_name)
00680
00681 write (*,'(1x,a,a)') trim(ch_id), &
00682 ': BundleNr. Min Max Sum Mean'
00683
00684 write (*, 9950)
00685
00686 do i = 1, shape_out(2,6)
00687 write(*,'(1x,a,a2,i3,6x,4(1x,e14.6))') trim(ch_id), &
00688 ': ', i, stat_min(i), stat_max(i), stat_sum(i), stat_mean(i)
00689 enddo
00690
00691 write (*, 9990) trim(ch_id)
00692
00693 call psmile_flushstd
00694
00695 endif
00696
00697 enddo
00698
00699
00700
00701 Deallocate (recv_buf, stat_min, stat_max, stat_sum, stat_mean, STAT=ierror)
00702 #if defined ( VERBOSE)
00703 if ( ierror /= 0 ) then
00704 ierrp (1) = nbr_fields
00705 ierror = PRISM_Error_Dealloc
00706 call psmile_error ( ierror, 'recv_buf, stat_{min,max,sum,mean}', &
00707 ierrp(1), 1, __FILE__, __LINE__ )
00708 return
00709 endif
00710 #endif
00711 end if
00712
00713
00714 #ifdef __PSMILE_WITH_IO
00715
00716
00717
00718
00719
00720
00721 #ifdef VERBOSE
00722 print *, trim(ch_id), ': psmile_put_dble: io_action', action(2)
00723
00724 call psmile_flushstd
00725 #endif /* VERBOSE */
00726
00727 lag = sga_transi_out%ig_lag
00728
00729 if ( action(2) ) then
00730
00731
00732
00733
00734 if(lag /= PSMILe_undef) then
00735
00736 delta_sec = (julian_dayb(2) - julian_day) * 86400.0 &
00737 + julian_secb(2) - julian_secb(1)
00738
00739 js_cur = julian_sec - lag * delta_sec
00740 add_days = floor(js_end / 86400.0)
00741 jd_cur = julian_day + add_days
00742 js_cur = js_cur - dble(add_days) * 86400.0
00743
00744 else
00745
00746 js_cur=julian_sec
00747 jd_cur=julian_day
00748
00749 endif
00750
00751 if ( local_reduce /= PSMILe_undef ) then
00752
00753
00754
00755
00756 print *, 'psmile_put_dble: output of reduced field not supported.'
00757
00758 else if ( local_scatter ) then
00759
00760 call psmile_write_byid_dble ( field_id, task_id, data_scattered, &
00761 jd_cur, js_cur, ierror )
00762
00763 else
00764
00765 if ( local_operation ) then
00766
00767 call psmile_write_byid_dble ( field_id, task_id, &
00768 Taskout%buffer_dble, &
00769 jd_cur, js_cur, ierror )
00770
00771 Taskout%buffer_dble = 0.0
00772 Taskout%nsum = 0
00773 else
00774
00775 call psmile_write_byid_dble ( field_id, task_id, data_array, &
00776 jd_cur, js_cur, ierror )
00777
00778 endif
00779
00780 endif
00781
00782 endif
00783
00784
00785
00786
00787
00788 il_smioc_loc=fp%smioc_loc
00789 debug_action=sga_transi_out%ig_debugmode.eq.PSMILe_true
00790 if ( debug_action ) then
00791 il_transiouts=0
00792 if(associated(sga_smioc_transi(il_smioc_loc)%sga_transi_out)) &
00793 il_transiouts=size(sga_smioc_transi(il_smioc_loc)%sga_transi_out)
00794
00795
00796
00797 il_taskid=il_transiouts+task_id
00798
00799
00800
00801 if(il_taskid.le.size( fp%io_task_lookup)) &
00802 debug_action= fp%io_task_lookup(il_taskid).gt.0
00803 endif
00804
00805 #ifdef VERBOSE
00806 print *, trim(ch_id), ': psmile_put_dble: debug_action', debug_action
00807
00808 call psmile_flushstd
00809 #endif /* VERBOSE */
00810
00811 if ( debug_action ) then
00812
00813 if ( local_reduce /= PSMILe_undef ) then
00814
00815
00816
00817
00818 print *, 'psmile_put_dble: output of reduced field not supported.'
00819
00820 else if ( local_scatter ) then
00821
00822 call psmile_write_byid_dble ( field_id, il_taskid, data_scattered, &
00823 julian_day, julian_sec, ierror )
00824
00825 else
00826
00827 if ( local_operation ) then
00828
00829 call psmile_write_byid_dble ( field_id, il_taskid, &
00830 Taskout%buffer_dble, &
00831 julian_day, julian_sec, ierror )
00832 else
00833
00834 call psmile_write_byid_dble ( field_id, il_taskid, data_array, &
00835 julian_day, julian_sec, ierror )
00836
00837 endif
00838
00839 endif
00840
00841 endif
00842
00843 #endif
00844
00845
00846
00847
00848
00849 if ( action(1) ) then
00850
00851 if ( local_reduce /= PSMILe_undef ) then
00852
00853
00854
00855
00856 print *, 'psmile_put_dble: coupling of reduced field not supported.'
00857
00858 else if ( local_scatter ) then
00859
00860 call psmile_put_field_dble (field_id, task_id, data_scattered, &
00861 len_3d, nbr_fields, ierror)
00862
00863 else
00864
00865 if ( local_operation ) then
00866
00867 call psmile_put_field_dble (field_id, task_id, &
00868 Taskout%buffer_dble, &
00869 len_3d, nbr_fields, ierror)
00870
00871
00872 if ( .not. action(3) ) then
00873 Taskout%buffer_dble = 0.0
00874 Taskout%nsum = 0
00875 endif
00876
00877 else
00878
00879 call psmile_put_field_dble (field_id, task_id, data_array, &
00880 len_3d, nbr_fields, ierror)
00881
00882 endif
00883
00884 endif
00885
00886 endif
00887
00888 #ifdef __PSMILE_WITH_IO
00889
00890
00891
00892
00893
00894 if ( action(3) ) then
00895
00896 il_transiouts=0
00897 if(associated(sga_smioc_transi(il_smioc_loc)%sga_transi_out)) &
00898 il_transiouts=size(sga_smioc_transi(il_smioc_loc)%sga_transi_out)
00899
00900
00901
00902 il_taskid_restr=3*il_transiouts+task_id+1
00903
00904
00905
00906 restart_action = .false.
00907 if(il_taskid_restr.le.size( fp%io_task_lookup)) &
00908 restart_action = fp%io_task_lookup(il_taskid_restr).gt.0
00909
00910
00911
00912
00913
00914 if(restart_action) &
00915 call psmile_date2ju ( PRISM_Jobend_date, jd_end, js_end )
00916
00917 if ( local_reduce /= PSMILe_undef ) then
00918
00919
00920
00921
00922 print *, 'psmile_put_dble: restart of reduced field not supported.'
00923
00924 else if ( local_scatter ) then
00925
00926 if ( restart_action ) &
00927 call psmile_write_byid_dble ( field_id, il_taskid_restr, &
00928 data_scattered, &
00929 jd_end, js_end, ierror )
00930
00931 else
00932
00933 if ( local_operation ) then
00934
00935 if ( restart_action ) &
00936 call psmile_write_byid_dble ( field_id, il_taskid_restr, &
00937 Taskout%buffer_dble, &
00938 jd_end, js_end, ierror )
00939
00940 Taskout%buffer_dble = 0.0
00941 Taskout%nsum = 0
00942
00943 else
00944
00945 if ( restart_action ) &
00946 call psmile_write_byid_dble ( field_id, il_taskid_restr, &
00947 data_array, &
00948 jd_end, js_end, ierror )
00949
00950 endif
00951
00952 endif
00953
00954 endif
00955
00956 #endif
00957
00958
00959
00960
00961
00962 if ( allocated(data_reduced) ) deallocate(data_reduced)
00963 if ( allocated(data_scattered) ) deallocate(data_scattered)
00964 if ( associated(mask_aux) ) deallocate(mask_aux)
00965
00966 #ifdef VERBOSE
00967 print *, trim(ch_id), ': psmile_put_dble: eof ierror ', ierror
00968
00969 call psmile_flushstd
00970 #endif /* VERBOSE */
00971
00972
00973
00974 9990 format (1x, a, ': ', 65('='))
00975 9980 format (1x, a, ': Statistics over ', a, ' points')
00976 9950 format (1x, a, ': ', 65('-'))
00977
00978
00979 end subroutine psmile_put_dble