00001 
00002 
00003 
00004 
00005 
00006       subroutine psmile_store_dest_locs_3d (found, loc, range, control,  &
00007                                             send_info, nloc, opt,        &
00008                                             nprev, nadd, ierror)
00009 
00010 
00011 
00012       use PRISM_constants
00013 
00014       use PSMILe, dummy_interface => PSMILe_Store_dest_locs_3d
00015 #ifdef DEBUG_TRACE
00016       use psmile_debug_trace
00017 #endif
00018 
00019       implicit none
00020 
00021 
00022 
00023       Integer, Intent (In)            :: range (2, ndim_3d)
00024 
00025 
00026 
00027       Integer, Intent (In)            :: control (2, ndim_3d)
00028 
00029 
00030 
00031       Integer, Intent (In)            :: found (range(1,1):range(2,1), 
00032                                                 range(1,2):range(2,2), 
00033                                                 range(1,3):range(2,3))
00034 
00035 
00036 
00037 
00038 
00039       Integer, Intent (In)            :: loc   (ndim_3d,               
00040                                                 range(1,1):range(2,1), 
00041                                                 range(1,2):range(2,2), 
00042                                                 range(1,3):range(2,3))
00043 
00044 
00045 
00046       Integer, Intent (In)            :: nloc
00047 
00048 
00049 
00050       Integer, Intent (In)            :: opt
00051 
00052 
00053 
00054 
00055 
00056 
00057 
00058 
00059       Integer, Intent (In)            :: nprev
00060 
00061 
00062 
00063 
00064 
00065       Type(Send_information), Intent(Inout) :: send_info
00066 
00067 
00068 
00069       Integer, Intent (Out)           :: nadd
00070 
00071 
00072 
00073       integer, Intent (Out)           :: ierror
00074 
00075 
00076 
00077 
00078 
00079 
00080 
00081       Integer, Parameter              :: val_both    =  0
00082       Integer, Parameter              :: val_abs     =  1
00083 
00084       Integer                         :: i, j, k, n
00085 
00086 
00087 
00088       Integer, Parameter              :: nerrp = 2
00089       Integer                         :: ierrp (nerrp)
00090 
00091 
00092 
00093 
00094 
00095 
00096 
00097 
00098 
00099 
00100 
00101 
00102 
00103 
00104 
00105 
00106 
00107 
00108 
00109 
00110    Character(len=len_cvs_string), save :: mycvs = 
00111        '$Id: psmile_store_dest_locs_3d.F90 2927 2011-01-28 14:04:12Z hanke $'
00112 
00113 
00114 
00115 
00116 
00117 #ifdef VERBOSE
00118       print 9990, trim(ch_id), opt, control
00119 
00120       call psmile_flushstd
00121 #endif /* VERBOSE */
00122 
00123 #ifdef DEBUG_TRACE
00124       if (range(1,1) <= ictl_ind(1) .and. ictl_ind(1) <= range(2,1) .and. &
00125           range(1,2) <= ictl_ind(2) .and. ictl_ind(2) <= range(2,2) .and. &
00126           range(1,3) <= ictl_ind(3) .and. ictl_ind(3) <= range(2,3)) then
00127          print *, '### prepared to store dest: ictl_ind found, loc', &
00128                   found(   ictl_ind(1),ictl_ind(2),ictl_ind(3)), &
00129                     loc(:, ictl_ind(1),ictl_ind(2),ictl_ind(3))
00130       endif
00131 #endif
00132 
00133       ierror = 0
00134 
00135       if (nprev == 0 .and. .not. Associated (send_info%dstijk) ) then
00136          Allocate (send_info%dstijk(1:ndim_3d, 1:nloc), STAT = ierror)
00137 
00138 
00139          if ( ierror > 0 ) then
00140             ierrp (1) = ierror
00141             ierrp (2) = nloc * ndim_3d
00142 
00143             ierror = PRISM_Error_Alloc
00144 
00145             call psmile_error ( ierror, 'send_info%dstijk', &
00146                                 ierrp, 2, __FILE__, __LINE__ )
00147             return
00148          endif
00149 
00150 #ifdef DEBUG
00151 
00152          send_info%dstijk = PSMILe_Undef
00153 #endif
00154       endif
00155 
00156 
00157 
00158 
00159       n = nprev
00160 
00161          if (opt == val_both) then
00162             do k = control (1, 3), control (2, 3)
00163                do j = control (1, 2), control (2, 2)
00164 
00165                   do i = control (1, 1), control (2, 1)
00166                      if (abs(found (i,j,k)) == val_abs) then
00167                         n = n + 1
00168 #ifdef DEBUG_TRACE
00169                         if ( i == ictl_ind(1) .and. j == ictl_ind(2) .and. k == ictl_ind(3) ) then
00170                            print *, '### store dest a): ictl_ind found, loc', n, &
00171                                 found( ictl_ind(1),ictl_ind(2),ictl_ind(3)), &
00172                                 loc(:, ictl_ind(1),ictl_ind(2),ictl_ind(3))
00173                         endif
00174 #endif
00175                         send_info%dstijk (1, n) = i
00176                         send_info%dstijk (2, n) = j
00177                         send_info%dstijk (3, n) = k
00178                      endif
00179                   end do
00180                end do
00181             end do
00182 
00183          else
00184 
00185             do k = control (1, 3), control (2, 3)
00186                do j = control (1, 2), control (2, 2)
00187 
00188                   do i = control (1, 1), control (2, 1)
00189                      if (found (i,j,k) == opt) then
00190                         n = n + 1
00191 #ifdef DEBUG_TRACE
00192                         if ( i == ictl_ind(1) .and. j == ictl_ind(2) .and. k == ictl_ind(3) ) then
00193                            print *, '### store dest b): ictl_ind found, loc', n, &
00194                                 found( ictl_ind(1),ictl_ind(2),ictl_ind(3)),  &
00195                                 loc(:, ictl_ind(1),ictl_ind(2),ictl_ind(3))
00196                         endif
00197 #endif
00198                         send_info%dstijk (1, n) = i
00199                         send_info%dstijk (2, n) = j
00200                         send_info%dstijk (3, n) = k
00201                      endif
00202                   end do
00203                end do
00204             end do
00205          endif
00206 
00207 #ifdef PRISM_ASSERTION
00208       if (nloc < n) then
00209          print *, "nloc, n", nloc, n
00210          call psmile_assert ( __FILE__, __LINE__, "nloc < n")
00211       endif
00212 #endif
00213 
00214 
00215 
00216       nadd = n - nprev
00217 
00218 #ifdef VERBOSE
00219       print 9980, trim(ch_id), ierror, nadd
00220 
00221       call psmile_flushstd
00222 #endif /* VERBOSE */
00223 
00224 
00225 
00226 #ifdef VERBOSE
00227 
00228 9990 format (1x, a, ': psmile_store_dest_locs_3d: opt ', i2, &
00229              '; control ', 6i6)
00230 9980 format (1x, a, ': psmile_store_dest_locs_3d: eof ierror =', i3, &
00231                     '; nadd =', i8)
00232 
00233 #endif /* VERBOSE */
00234 
00235       end subroutine PSMILe_Store_dest_locs_3d