00001
00002
00003
00004
00005
00006 subroutine psmile_store_source_locs_1d (found, loc, ibeg, len, &
00007 send_info, nloc, opt, &
00008 ialloc, ipart, nprev, nadd, ierror)
00009
00010
00011
00012 use PRISM_constants
00013
00014 use PSMILe, dummy_interface => PSMILe_store_source_locs_1d
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 (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 Type(Send_information), Intent(InOut) :: send_info
00065
00066
00067
00068 Integer, Intent (Out) :: nadd
00069
00070
00071
00072 integer, Intent (Out) :: ierror
00073
00074
00075
00076
00077
00078
00079
00080 Integer :: i, n
00081
00082 Integer, Pointer :: srcloc (:)
00083
00084
00085
00086 Integer, Parameter :: nerrp = 2
00087 Integer :: ierrp (nerrp)
00088
00089
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_source_locs_1d.F90 2325 2010-04-21 15:00:07Z valcke $'
00112
00113
00114
00115
00116
00117 #ifdef VERBOSE
00118 print 9990, trim(ch_id), opt, ibeg, len
00119
00120 call psmile_flushstd
00121 #endif /* VERBOSE */
00122
00123 ierror = 0
00124
00125 if (nprev == 0 .and. &
00126 .not. Associated (send_info%srclocs(ialloc, ipart)%vector) ) then
00127
00128 Allocate (send_info%srclocs(ialloc, ipart)%vector(nloc), STAT = ierror)
00129 if ( ierror > 0 ) then
00130 ierrp (1) = ierror
00131 ierrp (2) = nloc
00132 call psmile_error ( PRISM_Error_Alloc, 'send_info%srclocs()%vector', &
00133 ierrp, 2, __FILE__, __LINE__ )
00134 return
00135 endif
00136
00137 endif
00138
00139
00140
00141
00142
00143
00144 srcloc => send_info%srclocs(ialloc, ipart)%vector
00145
00146 n = nprev
00147
00148 if (opt == 0) then
00149
00150 do i = ibeg, len
00151 if (abs(found (i)) == 1) then
00152 srcloc (n+1) = loc (i)
00153 n = n + 1
00154 endif
00155 end do
00156
00157 else
00158
00159
00160 do i = ibeg, len
00161 if (found (i) == opt) then
00162 srcloc (n+1) = loc (i)
00163 n = n + 1
00164 endif
00165 end do
00166 endif
00167
00168 #ifdef PRISM_ASSERTION
00169 if (nloc < n) then
00170 print *, "nloc, n", nloc, n
00171 call psmile_assert ( __FILE__, __LINE__, "nloc < n")
00172 endif
00173 #endif
00174
00175
00176
00177 nadd = n - nprev
00178
00179 send_info%npoints(ialloc, ipart) = send_info%npoints(ialloc, ipart) + nadd
00180
00181
00182 #ifdef DEBUGX
00183 n = nprev
00184 if (opt == 0) then
00185 do i = ibeg, len
00186 if (abs(found (i)) == 1) then
00187 n = n + 1
00188 print *, i, loc(i), n, srcloc(n)
00189 endif
00190 end do
00191 else
00192 do i = ibeg, len
00193 if (found (i) == opt) then
00194 n = n + 1
00195 print *, i, loc(i), n, srcloc(n)
00196 endif
00197 end do
00198 endif
00199 #endif
00200
00201 #ifdef VERBOSE
00202 print 9980, trim(ch_id), ierror
00203
00204 call psmile_flushstd
00205 #endif /* VERBOSE */
00206
00207
00208
00209 #ifdef VERBOSE
00210
00211 9990 format (1x, a, ': psmile_store_source_locs_1d: opt ', i2, &
00212 '; ibeg, len ', 2i6)
00213 9980 format (1x, a, ': psmile_store_source_locs_1d: eof ierror =', i3)
00214
00215 #endif /* VERBOSE */
00216
00217 end subroutine PSMILe_store_source_locs_1d