00001
00002
00003
00004
00005
00006 subroutine psmile_store_source_locs_2d (found, loc, ibeg, len, &
00007 send_info, nloc, opt, ialloc, ipart, &
00008 nprev, nadd, ierror)
00009
00010
00011
00012 use PRISM_constants
00013
00014 use PSMILe, dummy_interface => PSMILe_store_source_locs_2d
00015
00016 Implicit none
00017
00018
00019
00020 Integer, Intent (In) :: ibeg
00021
00022
00023
00024 Integer, Intent (In) :: len
00025
00026
00027
00028 Integer, Intent (In) :: found (len)
00029
00030
00031
00032
00033 Integer, Intent (In) :: loc (ndim_2d, len)
00034
00035
00036
00037 Integer, Intent (In) :: nloc
00038
00039
00040
00041 Integer, Intent (In) :: opt
00042
00043
00044
00045
00046
00047
00048
00049
00050 Integer, Intent (In) :: ialloc
00051
00052
00053
00054 Integer, Intent (In) :: ipart
00055
00056
00057
00058 Integer, Intent (In) :: nprev
00059
00060
00061
00062
00063
00064
00065 Type(Send_information), Intent(InOut) :: send_info
00066
00067
00068
00069
00070
00071 Integer, Intent (Out) :: nadd
00072
00073
00074
00075 Integer, Intent (Out) :: ierror
00076
00077
00078
00079
00080
00081
00082
00083 Integer :: i, n
00084
00085 Integer, Pointer :: srcloc (:)
00086
00087
00088
00089 Integer, Parameter :: nerrp = 2
00090 Integer :: ierrp (nerrp)
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113 Character(len=len_cvs_string), save :: mycvs =
00114 '$Id: psmile_store_source_locs_2d.F90 2933 2011-01-31 16:42:27Z hanke $'
00115
00116
00117
00118
00119
00120 #ifdef VERBOSE
00121 print 9990, trim(ch_id), opt, ibeg, len
00122
00123 call psmile_flushstd
00124 #endif /* VERBOSE */
00125
00126 ierror = 0
00127
00128 if (nprev == 0 .and. &
00129 .not. Associated (send_info%srclocs(ialloc, ipart)%vector) ) then
00130
00131
00132
00133 Allocate (send_info%srclocs(ialloc, ipart)%vector(ndim_2d*nloc), &
00134 STAT = ierror)
00135 if ( ierror > 0 ) then
00136 ierrp (1) = ierror
00137 ierrp (2) = nloc * ndim_2d
00138
00139 ierror = PRISM_Error_Alloc
00140 call psmile_error ( ierror, 'send_info%srclocs()%vector', &
00141 ierrp, 2, __FILE__, __LINE__ )
00142 return
00143 endif
00144 endif
00145
00146
00147
00148
00149
00150
00151 srcloc => send_info%srclocs(ialloc, ipart)%vector
00152
00153 n = nprev * ndim_2d
00154 if (opt == 0) then
00155
00156 do i = ibeg, len
00157 if (abs(found (i)) == 1) then
00158 srcloc (n+1) = loc (1, i)
00159 srcloc (n+2) = loc (2, i)
00160
00161 n = n + ndim_2d
00162 endif
00163 end do
00164
00165 else
00166
00167
00168 do i = ibeg, len
00169 if (found (i) == opt) then
00170 srcloc (n+1) = loc (1, i)
00171 srcloc (n+2) = loc (2, i)
00172
00173 n = n + ndim_2d
00174 endif
00175 end do
00176 endif
00177
00178 #ifdef PRISM_ASSERTION
00179 if (nloc*ndim_2d < n) then
00180 print *, "nloc, n", nloc, n/ndim_2d
00181 call psmile_assert ( __FILE__, __LINE__, "nloc < n")
00182 endif
00183 #endif
00184
00185
00186
00187 nadd = n/ndim_2d - nprev
00188
00189 send_info%npoints(ialloc, ipart) = send_info%npoints(ialloc, ipart) + nadd
00190
00191
00192 #ifdef VERBOSE
00193 print 9980, trim(ch_id), ierror
00194
00195 call psmile_flushstd
00196 #endif /* VERBOSE */
00197
00198
00199
00200 #ifdef VERBOSE
00201
00202 9990 format (1x, a, ': psmile_store_source_locs_2d: opt ', i2, &
00203 '; ibeg, len ', 2i6)
00204 9980 format (1x, a, ': psmile_store_source_locs_2d: eof ierror =', i3)
00205
00206 #endif /* VERBOSE */
00207
00208 end subroutine PSMILe_store_source_locs_2d