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 subroutine prism_set_mask (mask_id, grid_id, mask_actual_shape, &
00028 mask_array, new_mask, ierror)
00029
00030 use prism_constants, only : prism_error_invalid_arg, prism_gaussreduced_regvrt, &
00031 prism_undefined
00032 use psmile, only : ch_id
00033 use psmile_user_data, only : psmile_store_data_mask, &
00034 test_user_grid_id, test_user_mask_id, &
00035 get_grid_valid_shape, get_grid_type
00036 use psmile_grid, only : get_size_of_shape
00037
00038 implicit none
00039
00040 integer, intent(inout) :: mask_id
00041 integer, intent(in) :: grid_id
00042 integer, intent(in) :: mask_actual_shape(2, *)
00043 logical, intent(in) :: mask_array (*)
00044 logical, intent(in) :: new_mask
00045 integer, intent(out) :: ierror
00046
00047 integer :: size_of_shape (2)
00048 integer :: i
00049 integer, allocatable :: grid_valid_shape(:,:)
00050
00051 #ifdef VERBOSE
00052 print 9990, trim(ch_id)
00053 call psmile_flushstd
00054 #endif /* VERBOSE */
00055
00056 ierror = 0
00057
00058 call test_user_grid_id (grid_id, ierror)
00059
00060 if (.not. new_mask) call test_user_mask_id (mask_id, ierror)
00061
00062 size_of_shape = get_size_of_shape(get_grid_type(grid_id))
00063
00064 allocate (grid_valid_shape(2, size_of_shape(2)))
00065
00066 grid_valid_shape = get_grid_valid_shape(grid_id, size_of_shape)
00067
00068
00069 do i = 1, size_of_shape(2)
00070 if (grid_valid_shape(1,i) /= prism_undefined .and. &
00071 (grid_valid_shape(1,i) < mask_actual_shape(1,i) .or. &
00072 grid_valid_shape(2,i) > mask_actual_shape(2,i))) then
00073
00074 ierror = prism_error_invalid_arg
00075 call psmile_error (ierror, "mask_actual_shape", (/0/), 1, &
00076 __FILE__, __LINE__ )
00077 endif
00078 enddo
00079
00080 deallocate (grid_valid_shape)
00081
00082 call psmile_store_data_mask (mask_id, grid_id, mask_actual_shape, &
00083 mask_array, new_mask, ierror)
00084
00085 #ifdef VERBOSE
00086 print 9980, trim(ch_id), ierror
00087 call psmile_flushstd
00088 #endif /* VERBOSE */
00089
00090 9990 format (1x, a, ': prism_set_mask: ')
00091 9980 format (1x, a, ': prism_set_mask: eof ierror =', i5)
00092
00093 end subroutine prism_set_mask