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
00026
00027
00028 subroutine prism_def_partition (grid_id, nbr_blocks, &
00029 partition_array, extent_array, ierror )
00030
00031 use prism_constants, only : prism_error_arg, prism_gaussreduced_regvrt
00032 use psmile_common, only : ch_id
00033 use psmile_user_data, only : psmile_store_data_partition, test_user_grid_id, &
00034 grid_data, get_grid_type
00035 use psmile_geometry, only : intersect
00036 use psmile_grid, only : get_size_of_shape
00037
00038 implicit none
00039
00040 integer, intent (in) :: grid_id
00041 integer, intent (in) :: nbr_blocks
00042 integer, intent (in) :: partition_array(1:nbr_blocks,*)
00043 integer, intent (in) :: extent_array(1:nbr_blocks,*)
00044 integer, intent (out) :: ierror
00045
00046 integer :: shape_size(2)
00047 integer :: i, j
00048 integer, allocatable :: extent(:), range_a(:,:), range_b(:,:)
00049 logical, allocatable :: varying_dimension(:)
00050
00051 #ifdef VERBOSE
00052 print 9990, trim(ch_id)
00053 call psmile_flushstd
00054 #endif /* VERBOSE */
00055
00056 ierror = 0
00057
00058
00059 if ( nbr_blocks < 1 ) then
00060
00061 ierror = PRISM_Error_Arg
00062 call psmile_error (ierror, 'nbr_blocks', (/grid_id, nbr_blocks/), 2, &
00063 __FILE__, __LINE__ )
00064 return
00065 endif
00066
00067
00068 call test_user_grid_id(grid_id, ierror)
00069
00070 shape_size = get_size_of_shape(get_grid_type(grid_id))
00071
00072
00073
00074
00075 if (nbr_blocks > 1 .and. &
00076 has_normal_extent_array(get_grid_type(grid_id))) then
00077
00078 allocate (extent(shape_size(2)), varying_dimension(shape_size(2)))
00079
00080
00081 extent = extent_array(1, 1:shape_size(2))
00082 varying_dimension = .false.
00083
00084 do i = 2, nbr_blocks
00085 varying_dimension = varying_dimension .or. &
00086 extent /= extent_array(i, 1:shape_size(2))
00087 enddo
00088
00089
00090 if (count (varying_dimension) > 1) then
00091
00092 ierror = prism_error_arg
00093 call psmile_error (ierror, "'More than one varying dimension specified'", &
00094 (/grid_id, count (varying_dimension)/), 2, &
00095 __FILE__, __LINE__ )
00096 endif
00097
00098 deallocate (extent, varying_dimension)
00099
00100 allocate (range_a(2, shape_size(2)), range_b(2, shape_size(2)))
00101
00102 do i = 1, nbr_blocks - 1
00103 range_a (1,:) = partition_array(i,1:shape_size(2)) + 1
00104 range_a (2,:) = partition_array(i,1:shape_size(2)) + extent_array(i,1:shape_size(2))
00105 do j = i + 1, nbr_blocks
00106 range_b (1,:) = partition_array(j,1:shape_size(2)) + 1
00107 range_b (2,:) = partition_array(j,1:shape_size(2)) + extent_array(j,1:shape_size(2))
00108
00109 if (intersect(range_a, range_b)) then
00110
00111 ierror = prism_error_arg
00112 print *, "block a", range_a
00113 print *, "block b", range_b
00114 call psmile_flushstd
00115 call psmile_assert (__FILE__, __LINE__, &
00116 "Blocks intersect in global index space!")
00117 endif
00118 enddo
00119 enddo
00120
00121 deallocate (range_a, range_b)
00122 endif
00123
00124 call psmile_store_data_partition (grid_id, nbr_blocks, &
00125 partition_array, extent_array, ierror)
00126
00127 #ifdef VERBOSE
00128 print 9980, trim(ch_id), ierror
00129 call psmile_flushstd
00130 #endif /* VERBOSE */
00131
00132 9990 format (1x, a, ': prism_def_partition: ')
00133 9980 format (1x, a, ': prism_def_partition: eof ierror =', i5)
00134
00135 contains
00136
00137 function has_normal_extent_array(grid_type)
00138 integer, intent(in) :: grid_type
00139 logical :: has_normal_extent_array
00140
00141 has_normal_extent_array = grid_type /= prism_gaussreduced_regvrt
00142
00143 end function has_normal_extent_array
00144
00145 end subroutine prism_def_partition