00001 subroutine psmile_store_mask_locs_21d ( range, control, found, &
00002 rangez, controlz, foundz, &
00003 send_info, ipart, ncpl, ncplz, ierror )
00004
00005 use PRISM_constants
00006 use PSMILe
00007
00008 Implicit none
00009
00010 Integer, Intent(In) :: range (2,ndim_2d)
00011 Integer, Intent(In) :: control(2,ndim_2d)
00012 Integer, Intent(In) :: found ( range(1,1):range(2,1),
00013 range(1,2):range(2,2))
00014 Integer, Intent(In) :: rangez(2)
00015 Integer, Intent(In) :: controlz(2)
00016 Integer, Intent(In) :: foundz ( rangez(1):rangez(2) )
00017 Type(Send_information), Intent(InOut) :: send_info
00018 Integer, Intent(In) :: ipart
00019 Integer, Intent(In) :: ncpl
00020 Integer, Intent(In) :: ncplz
00021
00022 Integer, Intent(Out) :: ierror
00023
00024
00025
00026 Integer, Parameter :: indl = 1
00027 Integer, Parameter :: indz = 2
00028
00029 Integer :: i, ii, jj, n
00030
00031
00032
00033 Integer, Parameter :: nerrp = 2
00034 Integer :: ierrp (nerrp)
00035
00036
00037
00038 #ifdef VERBOSE
00039 print 9990, trim(ch_id), ipart
00040
00041 call psmile_flushstd
00042 #endif /* VERBOSE */
00043
00044 ierror = 0
00045
00046
00047
00048 i = indl
00049
00050 if ( .not. Associated (send_info%msklocs(i,ipart)%vector) ) then
00051 Allocate (send_info%msklocs(i,ipart)%vector(ncpl), &
00052 STAT = ierror)
00053
00054 if ( ierror > 0 ) then
00055 ierrp (1) = ierror
00056 ierrp (2) = ncpl
00057 ierror = PRISM_Error_Alloc
00058 call psmile_error ( ierror, 'send_info%msklocs%vector', &
00059 ierrp, 2, __FILE__, __LINE__ )
00060 return
00061 endif
00062 endif
00063
00064 n = 0
00065 do jj = range(1,2), range(2,2)
00066 do ii = range(1,1), range(2,1)
00067 if (abs(found(ii,jj)) == 1 .or. found(ii,jj) == 0 ) then
00068 n = n + 1
00069 if ( ii >= control(1,1) .and. ii <= control(2,1) .and. &
00070 jj >= control(1,2) .and. jj <= control(2,2) ) then
00071 send_info%msklocs(i,ipart)%vector(n) = .true.
00072 else
00073 send_info%msklocs(i,ipart)%vector(n) = .false.
00074 endif
00075 endif
00076 end do
00077 end do
00078
00079
00080
00081
00082
00083 i = indz
00084
00085 if ( .not. Associated (send_info%msklocs(i, ipart)%vector) ) then
00086 Allocate (send_info%msklocs(i, ipart)%vector(ncplz), &
00087 STAT = ierror)
00088
00089 if ( ierror > 0 ) then
00090 ierrp (1) = ierror
00091 ierrp (2) = ncplz
00092 ierror = PRISM_Error_Alloc
00093 call psmile_error ( ierror, 'send_info%msklocs(indz,:)%vector', &
00094 ierrp, 2, __FILE__, __LINE__ )
00095 return
00096 endif
00097
00098 endif
00099
00100 n = 0
00101 do jj = rangez(1), rangez(2)
00102 if (abs(foundz(jj)) == 1 .or. foundz(jj) == 0 ) then
00103 n = n + 1
00104 if ( jj >= controlz(1) .and. jj <= controlz(2) ) then
00105 send_info%msklocs(i,ipart)%vector(n) = .true.
00106 else
00107 send_info%msklocs(i,ipart)%vector(n) = .false.
00108 endif
00109 endif
00110 end do
00111
00112
00113
00114
00115 #ifdef VERBOSE
00116 print 9980, trim(ch_id), ierror
00117
00118 call psmile_flushstd
00119 #endif /* VERBOSE */
00120
00121
00122
00123 #ifdef VERBOSE
00124
00125 9990 format (1x, a, ': psmile_store_mask_locs_21d: ipart', i3)
00126
00127 9980 format (1x, a, ': psmile_store_mask_locs_21d: eof ierror =', i3)
00128
00129 #endif /* VERBOSE */
00130
00131 end subroutine psmile_store_mask_locs_21d