00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_get_dble ( field_id, julian_day, julian_sec, &
00012 julian_dayb, julian_secb, data_array, cpl_action, io_action, &
00013 info, ierror )
00014
00015
00016
00017 use PRISM_constants
00018
00019 use PSMILe
00020 use PSMILe_SMIOC, only : sga_smioc_comp, transient, transient_in
00021
00022 implicit none
00023
00024
00025
00026 Integer, Intent (In) :: field_id
00027
00028
00029
00030 Real (PSMILe_float_kind), Intent (In) :: julian_day, julian_dayb(2)
00031 Real (PSMILe_float_kind), Intent (In) :: julian_sec, julian_secb(2)
00032
00033
00034
00035
00036
00037 Logical, Intent (In) :: cpl_action
00038
00039 Logical, Intent (In) :: io_action
00040
00041
00042
00043
00044
00045 Double Precision, Intent (InOut) :: data_array(*)
00046
00047
00048
00049 Integer, Intent (InOut) :: info
00050
00051
00052
00053 Integer, Intent (Out) :: ierror
00054
00055
00056
00057
00058
00059
00060
00061 Double Precision, Allocatable :: data_scattered(:)
00062
00063 Integer :: len_gathered, len_scattered, i, j
00064
00065 Integer :: len, len_3d
00066
00067 Integer :: nbr_fields
00068
00069 Logical :: local_timeop
00070 Logical :: local_add
00071 Logical :: local_multiply
00072 Logical :: local_gather
00073
00074 Integer :: task_id
00075
00076 Integer :: n_grid_dim
00077
00078 Type (GridFunction), Pointer :: fp
00079
00080
00081
00082 Type (transient), Pointer :: sga_smioc_transi(:)
00083 Type (transient_in), Pointer :: sg_transi_in
00084
00085
00086
00087 #ifdef __PSMILE_WITH_IO
00088 Logical :: debug_action
00089 Integer :: il_taskid,il_transiouts,il_smioc_loc
00090 #endif
00091
00092
00093
00094 Integer :: stat_nsum
00095 Integer :: nsum
00096
00097 Integer :: shape_in(2,6)
00098 Integer :: shape_out(2,6)
00099
00100 Logical :: mask_set
00101 Logical :: mask_needed
00102 Logical, Pointer :: mask_array(:)
00103 Logical, Pointer :: mask_aux(:)
00104
00105 Double Precision, Allocatable :: recv_buf(:)
00106 Double Precision, Allocatable :: stat_max(:)
00107 Double Precision, Allocatable :: stat_min(:)
00108 Double Precision, Allocatable :: stat_sum(:)
00109 Double Precision, Allocatable :: stat_mean(:)
00110
00111 Character(len=6), Save :: cstats (3)
00112
00113
00114
00115 Integer, Parameter :: nerrp = 2
00116 Integer :: ierr, ierrp (nerrp)
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141 Character(len=len_cvs_string), save :: mycvs =
00142 '$Id: psmile_get_dble.F90 2687 2010-10-28 15:15:52Z coquart $'
00143
00144 Data cstats /'masked', 'valid', 'all'/
00145
00146
00147
00148 #ifdef VERBOSE
00149 print *, trim(ch_id), ': psmile_get_dble: field_id', field_id
00150
00151 call psmile_flushstd
00152 #endif /* VERBOSE */
00153
00154
00155
00156
00157
00158 ierror = 0
00159 task_id = 1
00160 local_timeop = .false.
00161
00162 fp => Fields(field_id)
00163 sga_smioc_transi => sga_smioc_comp(Fields(field_id)%comp_id)%sga_smioc_transi
00164 sg_transi_in => sga_smioc_transi(fp%smioc_loc)%sg_transi_in
00165
00166 n_grid_dim = Grids(Methods(fp%method_id)%grid_id)%n_dim
00167
00168
00169
00170
00171
00172 len = fp%size
00173
00174 if ( fp%transi_type == PSMILe_bundle ) then
00175 nbr_fields = fp%var_shape(2,n_grid_dim+1)
00176 else
00177 nbr_fields = 1
00178 endif
00179
00180 len_3d = len / nbr_fields
00181
00182
00183
00184
00185
00186 local_add = sg_transi_in%sg_tgt_local_trans%dg_add_scalar /= PSMILe_dundef
00187 local_multiply = sg_transi_in%sg_tgt_local_trans%dg_mult_scalar /= PSMILe_dundef
00188 local_gather = sg_transi_in%sg_tgt_local_trans%ig_gather == PSMILe_gath
00189
00190 mask_set = fp%mask_id /= PRISM_UNDEFINED
00191
00192 if (mask_set) mask_set = Masks(fp%mask_id)%status /= PSMILe_status_free
00193
00194 #ifdef __PSMILE_WITH_IO
00195
00196
00197
00198
00199
00200 if ( io_action ) then
00201
00202 if ( local_gather ) then
00203
00204 call psmile_read_byid_dble ( field_id, task_id, data_scattered, &
00205 julian_day, julian_sec, &
00206 julian_dayb,julian_secb,local_timeop,ierror )
00207 else
00208
00209 call psmile_read_byid_dble ( field_id, task_id, data_array, &
00210 julian_day, julian_sec, &
00211 julian_dayb,julian_secb,local_timeop,ierror )
00212 endif
00213
00214 if ( local_timeop ) info = info + 1
00215
00216 endif
00217
00218 #endif
00219
00220
00221
00222
00223
00224 if ( cpl_action ) then
00225
00226 if ( local_gather ) then
00227
00228 call psmile_get_field_dble (field_id, data_scattered, len_3d, &
00229 nbr_fields, ierror)
00230
00231 else
00232
00233 call psmile_get_field_dble (field_id, data_array, len_3d, &
00234 nbr_fields, ierror)
00235
00236 endif
00237
00238 endif
00239
00240
00241
00242
00243
00244
00245 IF ( local_multiply ) THEN
00246
00247 IF ( local_gather ) THEN
00248
00249 data_scattered(1:len) = data_scattered(1:len) * &
00250 sg_transi_in%sg_tgt_local_trans%dg_mult_scalar
00251
00252 ELSE
00253
00254 data_array(1:len) = data_array(1:len) * &
00255 sg_transi_in%sg_tgt_local_trans%dg_mult_scalar
00256
00257 ENDIF
00258
00259 ENDIF
00260
00261
00262
00263 IF ( local_add ) THEN
00264
00265 IF ( local_gather ) THEN
00266
00267 data_scattered(1:len) = data_scattered(1:len) + &
00268 sg_transi_in%sg_tgt_local_trans%dg_add_scalar
00269 ELSE
00270
00271 data_array(1:len) = data_array(1:len) + &
00272 sg_transi_in%sg_tgt_local_trans%dg_add_scalar
00273 ENDIF
00274
00275 ENDIF
00276
00277
00278
00279 #ifdef __PSMILE_WITH_IO
00280
00281
00282
00283 il_smioc_loc=fp%smioc_loc
00284 debug_action= sg_transi_in%ig_debugmode .eq. PSMILe_true
00285 if ( debug_action ) then
00286 il_transiouts=0
00287 if(associated(sga_smioc_transi(il_smioc_loc)%sga_transi_out)) &
00288 il_transiouts=size(sga_smioc_transi(il_smioc_loc)%sga_transi_out)
00289
00290
00291
00292 il_taskid=2*il_transiouts+task_id
00293
00294
00295
00296 if(il_taskid.le.size( fp%io_task_lookup)) &
00297 debug_action= fp%io_task_lookup(il_taskid).gt.0
00298 endif
00299
00300 #ifdef VERBOSE
00301 print *, trim(ch_id), ': psmile_get_dble: debug_action', debug_action
00302
00303 call psmile_flushstd
00304 #endif /* VERBOSE */
00305
00306
00307
00308
00309
00310
00311 if ( debug_action ) then
00312
00313 if ( local_gather ) then
00314
00315 call psmile_write_byid_dble ( field_id, il_taskid, data_scattered, &
00316 julian_day, julian_sec, ierror )
00317 else
00318
00319
00320 call psmile_write_byid_dble ( field_id, il_taskid, data_array, &
00321 julian_day, julian_sec, ierror )
00322
00323 endif
00324
00325 endif
00326
00327 #endif
00328
00329
00330
00331
00332
00333
00334 if ( local_gather ) then
00335
00336
00337
00338
00339
00340 if ( .not. mask_set ) then
00341
00342 len_gathered = len_3d
00343 len_scattered = len_gathered
00344
00345 else
00346
00347 len_gathered = 0
00348 len_scattered = 1
00349
00350 do i = 1, n_grid_dim
00351
00352 len_scattered = len_scattered &
00353 * ( Masks(fp%mask_id)%mask_shape(2,i) - &
00354 Masks(fp%mask_id)%mask_shape(1,i) + 1 )
00355
00356 if ( len_scattered /= len_3d ) then
00357
00358 ierror = PRISM_Error_Parameter
00359
00360 ierrp(1) = len_scattered
00361 ierrp(2) = len_3d
00362
00363 call psmile_error ( PRISM_Error_Parameter, 'Size of Mask and Field', &
00364 ierrp, nerrp, __FILE__, __LINE__ )
00365
00366 return
00367 endif
00368
00369 enddo
00370
00371 len_gathered = count ( Masks(fp%mask_id)%mask_array )
00372
00373 endif
00374
00375
00376
00377 call psmile_loc_trans_dble ( PSMILe_gath, nbr_fields, &
00378 len_scattered, data_scattered, &
00379 len_gathered, data_array, field_id )
00380 endif
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391 mask_needed = sg_transi_in%iga_stats(1) == PSMILe_true .or. &
00392 sg_transi_in%iga_stats(2) == PSMILe_true .or. &
00393 sg_transi_in%iga_stats(3) == PSMILe_true
00394
00395 if ( mask_needed ) then
00396
00397
00398
00399 shape_out = 1
00400 shape_out(2,6) = nbr_fields
00401
00402 shape_in = 1
00403
00404 do i = 1, n_grid_dim
00405 shape_in(1,i) = fp%var_shape(1,i)
00406 shape_in(2,i) = fp%var_shape(2,i)
00407 enddo
00408
00409 shape_in(2,6)=nbr_fields
00410
00411
00412
00413 allocate (mask_aux(len/nbr_fields), STAT=ierr)
00414 if ( ierr /= 0 ) then
00415 ierrp (1) = 1
00416 ierror = PRISM_Error_Alloc
00417 call psmile_error ( ierror, 'mask_aux', ierrp(1), 1, &
00418 __FILE__, __LINE__ )
00419 return
00420 endif
00421
00422 allocate (recv_buf(1:shape_out(2,6)), &
00423 stat_max(1:shape_out(2,6)), &
00424 stat_min(1:shape_out(2,6)), &
00425 stat_sum(1:shape_out(2,6)), &
00426 stat_mean(1:shape_out(2,6)), STAT=ierror)
00427 if ( ierror /= 0 ) then
00428 ierrp (1) = nbr_fields * 5
00429 ierror = PRISM_Error_Alloc
00430 call psmile_error ( ierror, 'recv_buf, stat_{max,min,sum,mean}', &
00431 ierrp(1), 1, __FILE__, __LINE__ )
00432 return
00433 endif
00434
00435
00436
00437 do j = 1, 3
00438
00439 if ( sg_transi_in%iga_stats(j) == PSMILe_true ) then
00440
00441
00442
00443
00444
00445 if ( mask_set ) then
00446
00447 select case (j)
00448 case (1)
00449 mask_aux = .not. Masks(fp%mask_id)%mask_array
00450 mask_array => mask_aux
00451 case (2)
00452 mask_array => Masks(fp%mask_id)%mask_array
00453 case (3)
00454 mask_aux = .true.
00455 mask_array => mask_aux
00456 end select
00457
00458 else
00459
00460 mask_aux = j > 1
00461 mask_array => mask_aux
00462 endif
00463
00464 if ( local_gather ) then
00465
00466 call psmile_multi_reduce_dble ( PSMILe_max, shape_in, data_scattered, &
00467 shape_out, stat_max, mask_array, ierror )
00468
00469 call psmile_multi_reduce_dble ( PSMILe_min, shape_in, data_scattered, &
00470 shape_out, stat_min, mask_array, ierror )
00471
00472 call psmile_multi_reduce_dble ( PSMILe_integral, shape_in, data_scattered, &
00473 shape_out, stat_sum, mask_array, ierror )
00474 else
00475
00476 call psmile_multi_reduce_dble ( PSMILe_max, shape_in, data_array, &
00477 shape_out, stat_max, mask_array, ierror )
00478
00479 call psmile_multi_reduce_dble ( PSMILe_min, shape_in, data_array, &
00480 shape_out, stat_min, mask_array, ierror )
00481
00482 call psmile_multi_reduce_dble ( PSMILe_integral, shape_in, data_array, &
00483 shape_out, stat_sum, mask_array, ierror )
00484
00485 endif
00486
00487 if ( Comps(fp%comp_id)%act_comm /= MPI_COMM_NULL ) then
00488
00489 call MPI_Allreduce ( stat_max, recv_buf, shape_out(2,6), MPI_DOUBLE_PRECISION, &
00490 MPI_MAX, Comps(fp%comp_id)%act_comm, ierror )
00491 stat_max (:) = recv_buf (:)
00492
00493 call MPI_Allreduce ( stat_min, recv_buf, shape_out(2,6), MPI_DOUBLE_PRECISION, &
00494 MPI_MIN, Comps(fp%comp_id)%act_comm, ierror )
00495 stat_min (:) = recv_buf (:)
00496
00497 call MPI_Allreduce ( stat_sum, recv_buf, shape_out(2,6), MPI_DOUBLE_PRECISION, &
00498 MPI_SUM, Comps(fp%comp_id)%act_comm, ierror )
00499
00500 stat_sum (:) = recv_buf (:)
00501
00502
00503 if (j == 3) then
00504 stat_nsum = len_3d
00505 else
00506 stat_nsum = count(mask_array)
00507 endif
00508
00509 call MPI_Allreduce ( stat_nsum, nsum, 1, MPI_INTEGER, &
00510 MPI_SUM, Comps(fp%comp_id)%act_comm, ierror )
00511
00512 if (nsum > 0) then
00513 stat_mean(:) = stat_sum(:) / nsum
00514 else
00515 stat_mean = 0.0
00516 endif
00517
00518 endif
00519
00520
00521
00522 write (*, 9990) trim(ch_id)
00523
00524 write (*, 9980) trim(ch_id), trim(cstats (j))
00525
00526 write (*,*) trim(ch_id), &
00527 ': ... for field ', trim(fp%local_name)
00528
00529 write (*,'(1x,a,a)') trim(ch_id), &
00530 ': BundleNr. Min Max Sum Mean'
00531
00532 write (*, 9950) trim(ch_id)
00533
00534 do i = 1, shape_out(2,6)
00535 write(*,'(1x,a,a2,i3,6x,4(1x,e14.6))') trim(ch_id), &
00536 ': ', i, stat_min(i), stat_max(i), stat_sum(i), stat_mean(i)
00537 enddo
00538
00539 write (*, 9990) trim(ch_id)
00540
00541 call psmile_flushstd
00542
00543 endif
00544
00545 enddo
00546
00547
00548
00549 deallocate (recv_buf, stat_min, &
00550 stat_max, stat_sum, stat_mean, STAT = ierror)
00551 #if defined ( VERBOSE)
00552 if ( ierror /= 0 ) then
00553 ierrp (1) = nbr_fields
00554 ierror = PRISM_Error_Dealloc
00555 call psmile_error ( ierror, 'recv_buf, stat_{min,max,sum,mean}', &
00556 ierrp(1), 1, __FILE__, __LINE__ )
00557 return
00558 endif
00559 #endif
00560
00561 deallocate (mask_aux, STAT = ierror)
00562 #if defined ( VERBOSE)
00563 if ( ierror /= 0 ) then
00564 ierrp (1) = 1
00565 ierror = PRISM_Error_Dealloc
00566 call psmile_error ( ierror, 'mask_aux', &
00567 ierrp(1), 1, __FILE__, __LINE__ )
00568 return
00569 endif
00570 #endif
00571
00572 endif
00573
00574
00575 #ifdef VERBOSE
00576 print *, trim(ch_id), ': psmile_get_dble: eof ierror ', ierror
00577
00578 call psmile_flushstd
00579 #endif /* VERBOSE */
00580
00581
00582
00583 9990 format (1x, a, ': ', 65('='))
00584 9980 format (1x, a, ': Statistics over ', a, ' points')
00585 9950 format (1x, a, ': ', 65('-'))
00586
00587 end subroutine psmile_get_dble