psmile_store_data_partition.F90
Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 subroutine psmile_store_data_partition (grid_id, nbr_blocks, &
00026 partition_array, &
00027 extent_array, ierror)
00028
00029 use psmile_user_data, only : grid_data
00030 use psmile_grid, only : get_size_of_shape
00031 use prism_constants, only : PRISM_Error_Arg
00032 use psmile_common, only : ch_id
00033
00034 implicit none
00035
00036 integer, intent (in) :: grid_id
00037 integer, intent (in) :: nbr_blocks
00038 integer, intent (in) :: partition_array(1:nbr_blocks,*)
00039 integer, intent (in) :: extent_array(1:nbr_blocks,*)
00040 integer, intent (out) :: ierror
00041
00042 integer :: shape_size(2)
00043 #ifdef DEBUG
00044 integer :: i
00045 #endif
00046
00047 #ifdef VERBOSE
00048 print 9990, trim(ch_id)
00049 call psmile_flushstd
00050 #endif /* VERBOSE */
00051
00052 ierror = 0
00053
00054 if (associated (grid_data(grid_id)%partition_array)) then
00055
00056 ierror = PRISM_Error_Arg
00057 call psmile_error ( PRISM_Error_Arg, 'partition data already existed', &
00058 (/grid_id/), 1, __FILE__, __LINE__ )
00059 return
00060 endif
00061
00062 shape_size = get_size_of_shape(grid_data(grid_id)%grid_type)
00063
00064 allocate (grid_data(grid_id)%partition_array(nbr_blocks,1:shape_size(2)), &
00065 grid_data(grid_id)%extent_array(nbr_blocks,1:shape_size(2)))
00066
00067 grid_data(grid_id)%partition_array = partition_array(:,1:shape_size(2))
00068 grid_data(grid_id)%extent_array = extent_array(:,1:shape_size(2))
00069
00070 #ifdef DEBUG
00071 DO i=1,nbr_blocks
00072 PRINT*, 'grid_data(grid_id)%partition_array defined by the user:',partition_array(i,1:shape_size(2))
00073 PRINT*, 'grid_data(grid_id)%extent_array defined by the user:',extent_array(i,1:shape_size(2))
00074 ENDDO
00075 #endif
00076
00077 #ifdef VERBOSE
00078 print 9980, trim(ch_id), ierror
00079 call psmile_flushstd
00080 #endif /* VERBOSE */
00081
00082 9990 format (1x, a, ': psmile_store_data_partition: ')
00083 9980 format (1x, a, ': psmile_store_data_partition: eof ierror =', i5)
00084 end subroutine psmile_store_data_partition