Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! particle_info_ops.f90 is part of AstroBEAR. 00008 ! 00009 ! AstroBEAR is free software: you can redistribute it and/or modify 00010 ! it under the terms of the GNU General Public License as published by 00011 ! the Free Software Foundation, either version 3 of the License, or 00012 ! (at your option) any later version. 00013 ! 00014 ! AstroBEAR is distributed in the hope that it will be useful, 00015 ! but WITHOUT ANY WARRANTY; without even the implied warranty of 00016 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00017 ! GNU General Public License for more details. 00018 ! 00019 ! You should have received a copy of the GNU General Public License 00020 ! along with AstroBEAR. If not, see <http://www.gnu.org/licenses/>. 00021 ! 00022 !######################################################################### 00025 00029 00032 MODULE ParticleInfoOps 00033 USE ParticleDeclarations 00034 USE DataDeclarations 00035 USE GlobalDeclarations 00036 USE CommonFunctions 00037 USE PointGravitySrc 00038 USE Bondi 00039 IMPLICIT NONE 00040 PUBLIC CalcMoment, CheckForNewParticle, DoAccretion, GetGasForce, SinkSetErrFlag, BondiAccretionRate 00041 00042 PRIVATE 00043 CONTAINS 00044 00047 00048 00051 SUBROUTINE SinkSetErrFlag(Info) 00052 TYPE(InfoDef) :: Info 00053 TYPE(ParticleListDef), POINTER :: ParticleList 00054 TYPE(ParticleDef), POINTER :: Particle 00055 INTEGER, DIMENSION(:,:,:), POINTER :: mSs 00056 INTEGER :: nOverlaps, nGhost 00057 INTEGER, DIMENSION(3,2) :: mS 00058 REAL(KIND=qPREC), DIMENSION(:,:), POINTER :: offsets 00059 REAL(KIND=qPREC), DIMENSION(3,2) :: xBounds 00060 INTEGER :: i 00061 00062 ParticleList=>SinkParticles 00063 DO WHILE (ASSOCIATED(ParticleList)) 00064 Particle=>ParticleList%self 00065 xBounds=spread(particle%xloc,2,2) 00066 xBounds(1:nDim,1)=particle%xloc(1:nDim)-particle%buffer(Info%level)*levels(Info%level)%dx 00067 xBounds(1:nDim,2)=particle%xloc(1:nDim)+particle%buffer(Info%level)*levels(Info%level)%dx 00068 CALL CalcPhysicalOverlaps(Info, xBounds, mSs, nOverlaps, offsets, iEVERYWHERE, lHydroPeriodic, 0) 00069 IF (nOverlaps > 0) THEN 00070 DO i=1,nOverlaps 00071 mS=mSs(i,:,:) 00072 Info%ErrFlag(ms(1,1):mS(1,2),mS(2,1):mS(2,2),mS(3,1):ms(3,2))=1 00073 END DO 00074 DEALLOCATE(mSs, offsets) 00075 NULLIFY(mSs, offsets) 00076 END IF 00077 ParticleList=>ParticleList%next 00078 END DO 00079 END SUBROUTINE SinkSetErrFlag 00080 00081 00082 00083 00084 00089 SUBROUTINE CalcMoment(Particle, Info, lGhost) 00090 TYPE(ParticleDef) :: Particle 00091 TYPE(InfoDef) :: Info 00092 LOGICAL :: lPeriodic, lGhost 00093 INTEGER, DIMENSION(:,:,:), POINTER :: mSs 00094 REAL(KIND=qPREC), DIMENSION(:,:), POINTER :: offsets 00095 INTEGER :: nOverlaps, nGhost 00096 INTEGER, DIMENSION(3,2) :: mS 00097 REAL(KIND=qPREC), DIMENSION(3) :: offset 00098 INTEGER :: i 00099 ! For new particles - include ghost cells and don't consider periodic overlaps 00100 ! For old particles - don't include ghost cells - except at physical boundaries 00101 IF (lGhost) THEN 00102 CALL CalcPhysicalOverlaps(Info, GetParticleBounds(Particle), mSs, nOverlaps, offsets, iEVERYWHERE, lHydroPeriodic) 00103 ELSE 00104 CALL CalcPhysicalOverlaps(Info, GetParticleBounds(Particle), mSs, nOverlaps, offsets, iEVERYWHERE, lHydroPeriodic, 0) 00105 END IF 00106 00107 ! IF (nOverlaps > 0) write(*,*) "GridParticleOverlaps", nOverlaps, mSs 00108 IF (nOverlaps > 0) THEN 00109 DO i=1,nOverlaps 00110 mS=mSs(i,:,:) 00111 offset=offsets(i,:) 00112 CALL CalcMomentContributions(Info, Particle, mS, offset) 00113 END DO 00114 DEALLOCATE(mSs) 00115 DEALLOCATE(offsets) 00116 NULLIFY(mSs) 00117 NULLIFY(offsets) 00118 END IF 00119 ! write(*,*) particle%moments(1) 00120 END SUBROUTINE CalcMoment 00121 00122 00125 SUBROUTINE CheckForNewParticle(Info) 00126 ! Federrath et al list the following checks: 00127 ! * Density Threshold (TrueLove Criteria -> JeansLength = sqrt(pi*cs^2/G/rho)=sqrt(pi*gamma*Press/G)/rho > 4 dx or rho < sqrt(pi*gamma*Press/G)/(4 dx) = sqrt(pi*gamma/G)/(4 dx) * sqrt(Press) 00128 ! * Refinement Check (not necessary as this routine is only called on the maxlevel) 00129 ! * Converging Flow in each direction 00130 ! * Minimum in the total Potential 00131 ! * Jeans Instability Check 00132 ! * Bound State within control volume 00133 ! * Proximity check - don't create sink particles that are close to already existing particles 00134 USE EOS 00135 INTEGER :: i,j,k,l,m,n,ll,mm,nn,ii,jj,kk 00136 TYPE(InfoDef) :: Info 00137 REAL(KIND=qPREC), POINTER, DIMENSION(:,:,:,:) :: q 00138 REAL(KIND=qPREC) :: rhojeans, dx, eth, egrav,mass,ekin,emag 00139 INTEGER, DIMENSION(3,2) :: ip,ioffset 00140 LOGICAL :: lsamex, lsamey, lConverging 00141 INTEGER, DIMENSION(3) :: ir, il, ipos 00142 REAL(KIND=qPREC), DIMENSION(3) :: pos,mom,poffset 00143 TYPE(ParticleListDef), POINTER :: ParticleList 00144 TYPE(ParticleDef), POINTER :: Particle, NewParticle 00145 LOGICAL :: lFoundNearbyParticle 00146 q=>Info%q 00147 dx=levels(maxlevel)%dx 00148 ioffset=0 00149 WHERE (lHydroPeriodic(1:nDim)) ioffset(1:nDim,2)=2 !Do two copies... 00150 ioffset(1:nDim,1)=-ioffset(1:nDim,2) 00151 00152 ! write(*,*) "JeansFact=", JeansFact,JeansFact*sqrt(press(q(10,10,10,:))) 00153 DO i=1, Info%mX(1) 00154 DO j=1, Info%mX(2) 00155 DO k=1, Info%mX(3) 00156 ! * Density Threshold 00157 rhoJeans=JeansFact*sqrt(press(q(i,j,k,:))) 00158 IF (q(i,j,k,1) > rhoJeans) THEN !TrueLove threshold reached 00159 ! Check for Proximity to existing particles 00160 ! write(*,'(A5,3I8,A,4E18.5)') "cell", Info%mGlobal(:,1)-1+(/i,j,k/), ' violates density threshhold', rhoJeans, q(i,j,k,1) 00161 ipos=(/i,j,k/) 00162 pos=Cellpos(Info,i,j,k) 00163 particlelist=>SinkParticles 00164 lFoundNearbyParticle=.false. 00165 DO WHILE (ASSOCIATED(particlelist)) 00166 particle=>particlelist%self 00167 IF (particle%iaccrete == FEDERRATH_ACCRETION) THEN 00168 ! write(*,*) "proximity check", sqrt(SUM((pos(1:nDim)-particle%xloc(1:nDim))**2)), sqrt(r_acc2) 00169 DO ii=ioffset(1,1),ioffset(1,2) 00170 DO jj=ioffset(2,1),ioffset(2,2) 00171 DO kk=ioffset(3,1),ioffset(3,2) 00172 pos=Cellpos(Info, i, j, k)+(/ii,jj,kk/)*(GxBounds(:,2)-GxBounds(:,1)) 00173 IF (SUM((pos(1:nDim)-particle%xloc(1:nDim))**2) <= r_acc2) THEN 00174 lFoundNearbyParticle=.true. 00175 END IF 00176 END DO 00177 END DO 00178 END DO 00179 IF (lFoundNearbyParticle) EXIT 00180 END IF 00181 particlelist=>particlelist%next 00182 END DO 00183 00184 IF (.NOT. lFoundNearbyParticle) THEN !Did not find an existing particle... 00185 ! write(*,*) "did not find an existing particle" 00186 ip=spread(ipos,2,2) 00187 ! Check for converging flow 00188 lConverging=.true. 00189 ir=ipos 00190 il=ipos 00191 DO l=1,nDim 00192 il(l)=ipos(l)-1 !Left indices 00193 ir(l)=ipos(l)+1 !Right indices 00194 !div(velocity) = v(ir)-v(il) < 0 00195 lConverging=q(ir(1),ir(2),ir(3),imom(l))/q(ir(1),ir(2),ir(3),irho) < & 00196 q(il(1),il(2),il(3),imom(l))/q(il(1),il(2),il(3),irho) 00197 IF (.NOT. lConverging) EXIT 00198 il(l)=ipos(l) 00199 ir(l)=ipos(l) 00200 END DO 00201 IF (lConverging) THEN 00202 ! write(*,*) "converging pass" 00203 ! Check for Gravitational Minimum 00204 ip(1:nDim,1)=ipos(1:nDim)-IR_ACC 00205 ip(1:nDim,2)=ipos(1:nDim)+IR_ACC 00206 IF (q(ipos(1),ipos(2),ipos(3),1) == MAXVAL(q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),1),lControlVolume) .AND. COUNT(q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),1)==q(ipos(1),ipos(2),ipos(3),1) .AND. lControlVolume)==1) THEN 00207 ! write(*,*) 'A', q(ipos(1),ipos(2),ipos(3),iPhi) 00208 ! write(*,*) 'B', MAXVAL(q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),iPhi),lControlVolume) 00209 ! write(*,*) 'C', lControlVolume 00210 ! write(*,*) "minimum potential pass at ", Info%mGlobal(:,1)-1+ipos! !NrSinkParticles 00211 ! write(*,'(9E13.2)') q(ipos(1)-1:ipos(1)+1,ipos(2)-1:ipos(2)+1,ipos(3)-1:ipos(3)+1,iphi) 00212 ! Jeans Instability Check 00213 eth=0 00214 egrav=0 00215 mass=0 00216 mom=0 00217 emag=0 00218 !Sum up the total energies (divided by dV) 00219 DO l=ip(1,1),ip(1,2) 00220 DO m=ip(2,1),ip(2,2) 00221 DO n=ip(3,1),ip(3,2) 00222 IF (.NOT. lControlVolume(l-ipos(1),m-ipos(2),n-ipos(3))) CYCLE 00223 mass=mass+q(l,m,n,1) 00224 mom(1:nDim)=mom(1:nDim)+q(l,m,n,imom(1:nDim)) 00225 eth=eth+1.5d0*Press(q(l,m,n,:)) 00226 IF (lMHD) emag=emag+half*SUM(q(l,m,n,iBx:iBz)**2) 00227 DO ll=ip(1,1),ip(1,2) 00228 lsamex=ll==l 00229 DO mm=ip(2,1),ip(2,2) 00230 lsamey=lsamex .AND. mm==m 00231 DO nn=ip(3,1),ip(3,2) 00232 IF (lControlVolume(ll-ipos(1),mm-ipos(2),nn-ipos(3))) THEN 00233 IF (lsamey .AND. nn==n) THEN 00234 IF (nDim == 3) THEN 00235 egrav=egrav-0.94115625*q(l,m,n,1)**2*ScaleGravdV/sink_dx ! Self-energy of a uniform cube... from arxiv.org/pdf/astro-ph/0002496 00236 ELSE 00237 egrav=egrav-2.94298464974921*q(l,m,n,1)**2*ScaleGravdV/sink_dx ! Self-energy of a uniform square.. 00238 END IF 00239 ELSE 00240 IF (nDim == 3 .OR. iCylindrical /= NOCYL) THEN 00241 egrav=egrav-half*q(l,m,n,1)*q(ll,mm,nn,1)*ScaleGravdV/sqrt(REAL((ll-l)**2+(mm-m)**2+(nn-n)**2))/sink_dx 00242 ELSE 00243 ! In 2D everything is bound - nothing can truly escape to infinity since potential goes like ln(r) so normalize potential to go to zero at r ~ box-size 00244 egrav=egrav+half*q(l,m,n,1)*q(ll,mm,nn,1)*2d0*ScaleGravdV*log(sqrt(REAL((ll-l)**2+(mm-m)**2+(nn-n)**2))*sink_dx/sqrt(product(GxBounds(1:2,2)-GxBounds(1:2,1)))) 00245 END IF 00246 END IF 00247 END IF 00248 END DO 00249 END DO 00250 END DO 00251 END DO 00252 END DO 00253 END DO 00254 IF (abs(egrav) > 2d0*eth+emag) THEN !Jeans unstable - not sure if this should be 2d0*(eth+emag) instead of 2d0*eth+emag 00255 ! write(*,'(A,6E20.12)') "Jeans unstable", egrav, eth, emag,mom/mass 00256 00257 mom=mom/mass !mom is now center of motion 00258 ekin=0 00259 DO l=ip(1,1),ip(1,2) 00260 DO m=ip(2,1),ip(2,2) 00261 DO n=ip(3,1),ip(3,2) 00262 IF (.NOT. lControlVolume(l-ipos(1),m-ipos(2),n-ipos(3))) CYCLE 00263 Ekin=Ekin+half*q(l,m,n,1)*SUM((q(l,m,n,imom(1:nDim))/q(l,m,n,1)-mom(1:nDim))**2) 00264 END DO 00265 END DO 00266 END DO 00267 IF (Egrav+Eth+Ekin+Emag < 0d0) THEN 00268 ! write(*,'(A,5E20.12)') "Bound System", Egrav, Eth, Emag, Ekin, Egrav+Eth+Ekin+Emag 00269 ! We made a new particle 00270 NULLIFY(NewParticle) 00271 CALL CreateParticle(NewParticle) 00272 NewParticle%xloc(1:nDim)=pos(1:nDim) 00273 CALL CalcMoment(NewParticle, Info,.true.) 00274 CALL AddNewSinkParticle(NewParticle) 00275 ELSE 00276 ! write(*,'(A,5E20.12)') "UnBound System", Egrav, Eth, Emag, Ekin, Egrav+Eth+Ekin+Emag 00277 END IF 00278 ELSE 00279 ! write(*,'(A,4E20.12)') "Jeans stable", egrav, eth, emag 00280 END IF 00281 END IF 00282 ELSE 00283 ! write(*,*) "converging failed" 00284 END IF 00285 ELSE 00286 ! write(*,*) "near an existing particle" 00287 END IF 00288 END IF 00289 END DO 00290 END DO 00291 END DO 00292 END SUBROUTINE CheckForNewParticle 00293 00294 00295 SUBROUTINE BondiAccretionRate(Particle, Info) 00296 TYPE(InfoDef) :: Info 00297 TYPE(ParticleDef) :: Particle 00298 INTEGER :: ipos(3), mb(3,2), i, j, k 00299 REAL(KIND=qPREC) :: v_inf, c_inf, r_bh, r_kernel, r, w, rho_inf, w_sum,dipos(3) 00300 IF (ParticleInGrid(Particle, Info, ipos)) THEN 00301 v_inf=sqrt(sum((Info%q(ipos(1),ipos(2),ipos(3),imom(1:nDim))/Info%q(ipos(1),ipos(2),ipos(3),1)-Particle%q(imom(1:nDim))/Particle%q(1))**2)) 00302 c_inf=soundspeed(Info%q(ipos(1),ipos(2),ipos(3),:)) 00303 r_bh=ScaleGrav*Particle%q(1)/(v_inf**2+c_inf**2) 00304 r_kernel=min(max(r_bh, sink_dx/4d0), IR_ACC*sink_dx/2d0) 00305 mb=1 00306 mb(1:nDim,1)=ipos(1:nDim)-IR_ACC 00307 mB(1:nDim,2)=ipos(1:nDim)+IR_ACC 00308 w_sum=0d0 00309 ! rho_inf=0d0 00310 dipos=(Particle%xloc-Info%xBounds(:,1))/sink_dx+half 00311 DO i=mb(1,1), mB(1,2) 00312 DO j=mB(2,1), mB(2,2) 00313 DO k=mb(3,1), mB(3,2) 00314 r=sqrt(sum((real((/i,j,k/),8)-dipos)**2 ))*sink_dx 00315 IF (r <= REAL(IR_ACC)*sink_dx) THEN 00316 IF (iCylindrical /= NOCYL) THEN 00317 w=exp(-(r/r_kernel)**2)*abs(Info%xBounds(1,1)+(real(i,8)-half)*levels(Info%level)%dx) 00318 ELSE 00319 w=exp(-(r/r_kernel)**2) 00320 END IF 00321 w_sum=w_sum+w 00322 ! rho_inf=rho_inf+w*Info%q(i,j,k,1) 00323 END IF 00324 END DO 00325 END DO 00326 END DO 00327 ! rho_inf=rho_inf/w_sum/BH_alpha(1.2*sink_dx/r_bh) 00328 Particle%AccretionRate=4d0*Pi*r_bh**2*sqrt(Bondi_lambda2*c_inf**2+v_inf**2)/BH_alpha(1.2*sink_dx/r_bh)/(w_sum*sink_dv) !normalized 00329 Particle%Bondi_kernel=r_kernel 00330 ! write(*,*) 'ParticleAccretionRate=', Particle%AccretionRate, 1.2*sink_dx/r_bh 00331 ! write(*,'(A,E24.15)') 'r_bh=', r_bh 00332 ! write(*,'(A,E24.15)') 'r_kernel', r_kernel 00333 ! write(*,'(A,E24.15)') 'c_inf', c_inf 00334 ! write(*,'(A,2E24.15)') 'rho_inf', rho_inf, rho_inf*BH_alpha(1.2*sink_dx/r_bh) 00335 ! write(*,'(A,E24.15)') 'v_inf', v_inf 00336 00337 ! write(*,'(A,E24.15)') 'Accretion Coeff', Particle%AccretionRate 00338 ! write(*,'(A,E24.15)') 'Accretion Rate', 4d0*Pi*r_bh**2*rho_inf*sqrt(Bondi_lambda2*c_inf**2+v_inf**2) 00339 ! write(*,*) 'rate/coeff= ', 4d0*Pi*r_bh**2*rho_inf*sqrt(Bondi_lambda2*c_inf**2+v_inf**2)/Particle%AccretionRate 00340 END IF 00341 END SUBROUTINE BondiAccretionRate 00342 00343 00346 SUBROUTINE DoAccretion(Info) 00347 TYPE(InfoDef) :: Info 00348 TYPE(ParticleListDef), POINTER :: ParticleList 00349 TYPE(ParticleDef), POINTER :: Particle 00350 ParticleList=>SinkParticles 00351 DO WHILE (ASSOCIATED(ParticleList)) 00352 Particle=>ParticleList%self 00353 SELECT CASE(Particle%iAccrete) 00354 CASE(FEDERRATH_ACCRETION) 00355 CALL SelfGravAccretion(Info, Particle) 00356 CASE(KRUMHOLZ_ACCRETION) 00357 CALL BondiAccretion(Info, Particle) 00358 END SELECT 00359 ParticleList=>ParticleList%next 00360 END DO 00361 END SUBROUTINE DoAccretion 00362 00363 SUBROUTINE SelfGravAccretion(Info, Particle) 00364 USE EOS 00365 TYPE(InfoDef) :: Info 00366 REAL(KIND=qPREC), DIMENSION(:,:,:,:), POINTER :: q 00367 INTEGER, DIMENSION(3,2) :: ip 00368 INTEGER, DIMENSION(:,:,:), POINTER :: mSs 00369 REAL(KIND=qPREC), DIMENSION(:,:), POINTER :: offsets 00370 INTEGER :: nParticles 00371 INTEGER :: nOverlaps 00372 INTEGER :: i,j,k,l,iOverlap,rmbc,ii,jj,kk,ioffset(3,2) 00373 INTEGER, DIMENSION(3) :: ipos 00374 REAL(KIND=qPREC), DIMENSION(3) :: pos, offset, poffset 00375 REAL(KIND=qPREC) :: rhoJeans, minCellBindingEnergy, minRadius, r, CellBindingEnergy,mass,rmass(3),p 00376 TYPE(ParticleDef), POINTER :: Particle, MostBoundParticle,ClosestParticle,TempParticle 00377 TYPE(ParticleListDef), POINTER :: ParticleList, TempParticleList 00378 LOGICAL :: lAccrete 00379 ioffset=0 00380 WHERE (lHydroPeriodic(1:nDim)) ioffset(1:nDim,2)=2 !Do two copies... 00381 ioffset(1:nDim,1)=-ioffset(1:nDim,2) 00382 00383 q=>Info%q 00384 ! ip(nDim+1:3,:)=1 00385 ! ip(1:nDim,1)=Info%mGlobal(1:nDim,1)-IR_ACC 00386 ! ip(1:nDim,2)=Info%mGlobal(1:nDim,2)+IR_ACC 00387 ! rmbc=levels(Info%level)%gmbc(levels(Info%level)%step) 00388 00389 ! CALL CreateOverlapList(ip, GridParticles, ips) 00390 ! nParticles=size(ips,1) 00391 00392 ! DO iParticle=1, nParticles 00393 00394 ! write(*,*) "checking on particle", loc(particle) 00395 CALL CalcPhysicalOverlaps(Info, GetParticleBounds(Particle), mSs, nOverlaps, offsets, iEVERYWHERE, lHydroPeriodic) 00396 ! write(*,*) "found overlaps", nOverlaps 00397 ! write(*,*) Info%xbounds 00398 ! write(*,*) Particle%xloc 00399 ! DO iOverlap=1,nOverlaps 00400 ! write(*,*) mSs(iOVerlap,:,:) 00401 ! END DO 00402 DO iOverlap=1,nOverlaps 00403 offset=offsets(iOverlap,:) 00404 ip=mSs(iOverlap,:,:) !overlap bounds for this particle 00405 ! CALL CalcOverlaps(Info%mGlobal, Particle%xloc, ip) 00406 00407 ! ALLOCATE(potentialoverlaps(ParticleCount(ParticleOverlaps))) 00408 ! write(*,*) "checking cells", ip 00409 00410 DO k=ip(3,1),ip(3,2) 00411 DO j=ip(2,1),ip(2,2) 00412 DO i=ip(1,1),ip(1,2) 00413 IF (q(i,j,k,1) <= MinDensity) CYCLE 00414 p=Press(q(i,j,k,:)) 00415 IF (p < q(i,j,k,1)*MinTemp) CYCLE 00416 rhoJeans=(JeansFact)**2*p/q(i,j,k,1) 00417 ! rhoJeans=JeansFact*sqrt(press(q(i,j,k,:))) 00418 IF (q(i,j,k,1) > rhoJeans) THEN !TrueLove threshold reached 00419 ipos=(/i,j,k/) 00420 ! write(*,*) "cell at",(/i,j,k/), "violates truelove", rhoJeans,press(q(i,j,k,:)),q(i,j,k,1) 00421 ! write(*,*) rhoJeans, q(i,j,k,1) 00422 pos=(REAL(ipos)-half)*sink_dx + Info%xbounds(:,1)-offset 00423 r=sqrt(sum((pos(1:nDim)-Particle%xloc(1:nDim))**2)) 00424 IF (r <= r_acc) THEN 00425 00426 !For this cell - check to see if it is most tightly bound to current particle 00427 !Need Binding Energy to be negative otherwise cell is not bound to any particle 00428 MinCellBindingEnergy=0d0 00429 NULLIFY(MostBoundParticle, ClosestParticle) 00430 minradius=huge(1d0) 00431 TempParticleList=>SinkParticles 00432 DO WHILE (ASSOCIATED(TempParticleList)) 00433 TempParticle=>TempParticleList%self 00434 IF (TempParticle%iAccrete==FEDERRATH_ACCRETION) THEN 00435 DO ii=ioffset(1,1),ioffset(1,2) 00436 DO jj=ioffset(2,1),ioffset(2,2) 00437 DO kk=ioffset(3,1),ioffset(3,2) 00438 pOffSet=(/ii,jj,kk/)*(GxBounds(:,2)-GxBounds(:,1)) 00439 r=sqrt(sum((pos(1:nDim)+pOffSet(1:nDim)-TempParticle%xloc(1:nDim))**2)) 00440 CellBindingEnergy = GetBindingEnergy(Tempparticle, pos+pOffset, q(i,j,k,:),r) 00441 IF (SUM((pos(1:nDim)+pOffset(1:nDim)-TempParticle%xloc(1:nDim))*(q(i,j,k,imom(1:nDim))/q(i,j,k,1)-TempParticle%Q(imom(1:nDim)))) <= 0) THEN 00442 IF (CellBindingEnergy < MinCellBindingEnergy) THEN 00443 MostBoundParticle=>TempParticle 00444 MinCellBindingEnergy=CellBindingEnergy 00445 END IF 00446 END IF 00447 IF (r < minradius) THEN 00448 ClosestParticle=>TempParticle 00449 minradius=r 00450 END IF 00451 END DO 00452 END DO 00453 END DO 00454 END IF 00455 TempParticleList=>TempParticleList%next 00456 END DO 00457 lAccrete=.false. 00458 ! write(*,*) minradius, r_inner_acc 00459 IF (minradius < r_inner_acc) THEN !Don't check for most bound state - since this cell is too close to a particle 00460 lAccrete = associated(ClosestParticle, target=Particle) !Automatically accrete from this particle 00461 ELSE 00462 IF (associated(MostBoundParticle, target=Particle)) THEN !this cell is most tightly bound to this particle 00463 !Check to see if it's radial velocity is toward the particle 00464 !Already did this in determining whether it was bound 00465 lAccrete = .true. !(SUM((q(i,j,k,imom(1:nDim))/q(i,j,k,1)-Particle%vel(1:nDim)) * (pos(1:nDim)-Particle%xloc(1:nDim))) < 0) 00466 END IF 00467 END IF 00468 IF (lAccrete) THEN 00469 ! write(*,*) "accreting from cell", i,j,k 00470 IF (ALL(ipos(1:nDim) >= 1) .AND. ALL(ipos(1:nDim) <= Info%mX(1:nDim))) THEN 00471 Particle%dQ(1:NrHydroVars)=Particle%dQ(1:NrHydrovars)+Info%q(i,j,k,1:NrHydroVars)*sink_dv 00472 00473 mass=sink_dv*(q(i,j,k,1)-rhoJeans) 00474 rmass(1:nDim)=mass*(pos(1:nDim)-Particle%xloc(1:nDim)) 00475 Particle%drmass=Particle%drmass+rmass 00476 IF (nDim == 2) THEN 00477 Particle%dJ=Particle%dJ+Cross2D(rmass(1:2),q(i,j,k,imom(1:2))/q(i,j,k,1)) 00478 ELSEIF (nDim == 3) THEN 00479 Particle%dJ(1:nDim)=Particle%dJ+Cross3D(rmass,q(i,j,k,imom(1:3))/q(i,j,k,1)) 00480 END IF 00481 END IF 00482 00483 IF (iE /= 0) THEN 00484 IF (lMHD) q(i,j,k,iE)=q(i,j,k,iE)-half*SUM(q(i,j,k,iBx:iBz)**2) !Subtract off magnetic energy 00485 q(i,j,k,iE)=q(i,j,k,iE)-half*SUM(q(i,j,k,m_low:m_high)**2)/q(i,j,k,1) !Subtract off kinetic energy 00486 END IF 00487 q(i,j,k,1:NrHydroVars)=q(i,j,k,1:NrHydroVars)*rhoJeans/q(i,j,k,1) !Everything else should be accreted proportionally to density... 00488 IF (lMHD) THEN 00489 q(i,j,k,iBx)=.5*(Info%aux(i,j,k,1)+Info%aux(i+1,j,k,1)) 00490 q(i,j,k,iBy)=.5*(Info%aux(i,j,k,2)+Info%aux(i,j+1,k,2)) 00491 IF (nDim == 3) q(i,j,k,iBz)=.5*(Info%aux(i,j,k,3)+Info%aux(i,j,k+1,3)) 00492 END IF 00493 IF (iE /= 0) THEN 00494 q(i,j,k,iE)=q(i,j,k,iE)+half*SUM(q(i,j,k,m_low:m_high)**2)/q(i,j,k,1) !Add in kinetic energy 00495 IF (lMHD) q(i,j,k,iE)=q(i,j,k,iE)+half*SUM(q(i,j,k,iBx:iBz)**2) !Add magnetic energy back in. 00496 END IF 00497 IF (ALL(ipos(1:nDim) >= 1) .AND. ALL(ipos(1:nDim) <= Info%mX(1:nDim))) & 00498 Particle%dQ(1:NrHydroVars)=Particle%dQ(1:NrHydrovars)-Info%q(i,j,k,1:NrHydroVars)*sink_dv 00499 00500 END IF 00501 00502 END IF 00503 END IF 00504 END DO 00505 END DO 00506 END DO 00507 END DO 00508 IF (nOverlaps > 0) THEN 00509 DEALLOCATE(mSs, offsets) 00510 NULLIFY(mSs, offsets) 00511 END IF 00512 00513 END SUBROUTINE SelfGravAccretion 00514 00515 00516 SUBROUTINE BondiAccretion(Info,Particle) 00517 TYPE(ParticleDef) :: Particle 00518 TYPE(InfoDef) :: Info 00519 REAL(KIND=qPREC) :: GM, q_fact, sKE, r, r2, esp, jsp2, sie, r_surf, dx, dz, w, drho, mass 00520 REAL(KIND=qPREC), DIMENSION(3) :: xpos, pos, v, v_rel, v_r, v_phi, sample_fact, rmass 00521 INTEGER :: i, j, k, ii, jj, kk, iOverlap, nOverlaps 00522 INTEGER, DIMENSION(3) :: ipos, sample_res, iipos 00523 INTEGER, DIMENSION(3,2) :: mB, mS 00524 INTEGER, DIMENSION(3) :: offset 00525 INTEGER, DIMENSION(:,:), POINTER :: offsets 00526 INTEGER, DIMENSION(:,:,:), POINTER :: mSs, count 00527 REAL(KIND=qPREC), DIMENSION(:,:,:,:), POINTER :: q 00528 ipos = PosCell(Info, Particle%xloc(1), Particle%xloc(2), Particle%xloc(3)) 00529 mB=spread(Info%mGlobal(:,1)-1+ipos,2,2) 00530 mB(:,1)=Info%mGlobal(:,1)-1+ipos(:)-IR_ACC 00531 mb(:,2)=Info%mGlobal(:,1)-1+ipos(:)+IR_ACC 00532 CALL CalcCellOverlaps(Info, mB, mSs, nOverlaps, offsets, iEVERYWHERE, lHydroPeriodic) 00533 IF (nOverlaps > 0) THEN 00534 q=>Info%q 00535 GM=Particle%q(1)*ScaleGrav 00536 dx=sink_dx 00537 dz=merge(dx, 0d0, nDim == 3) 00538 sample_res=1 00539 sample_fact=0d0 00540 sample_res(1:nDim)=8 00541 sample_fact(1:nDim)=1d0/REAL(sample_res(1:nDim),8) 00542 q_fact=product(sample_fact(1:nDim)) 00543 r_surf=.25*sink_dx 00544 DO iOverlap=1,nOverlaps 00545 offset=offsets(iOverlap,:) 00546 mS=mSs(iOverlap,:,:) !overlap bounds for this particle 00547 00548 ALLOCATE(count(ms(1,1):mS(1,2),mS(2,1):mS(2,2),mS(3,1):mS(3,2))) 00549 count=0 00550 00551 DO k=mS(3,1),mS(3,2) 00552 xpos(3)=Info%xBounds(3,1)+(k-1+offset(3))*dz 00553 DO j=mS(2,1),mS(2,2) 00554 xpos(2)=Info%xBounds(2,1)+(j-1+offset(2))*dx 00555 DO i=mS(1,1),mS(1,2) 00556 xpos(1)=Info%xBounds(1,1)+(i-1+offset(1))*dx 00557 v(1:nDim)=q(i,j,k,imom(1:nDim))/q(i,j,k,1) !cell velocity 00558 v_rel(1:nDim)=v(1:nDim) - Particle%q(imom(1:nDim)) !velocity w/respect to sink source 00559 sKE=half*DOT_PRODUCT(v_rel(1:nDim),v_rel(1:nDim)) 00560 IF (iCylindrical == 2) sKE=half*(q(i,j,k,iAngMom)/q(i,j,k,1))**2 !Add angular momentum to energy 00561 IF (sqrt(SUM((xpos(:)+half*(/dx,dx,dz/)-Particle%xloc(:))**2)) > IR_ACC*sink_dx) THEN 00562 count(i,j,k)=-1 00563 CYCLE 00564 END IF 00565 IF (q(i,j,k,1) <= MinDensity) CYCLE 00566 DO kk=1,sample_res(3) 00567 pos(3)=xpos(3)+(REAL(kk, 8)-half)*dz*sample_fact(3) 00568 DO jj=1,sample_res(2) 00569 pos(2)=xpos(2)+(REAL(jj, 8)-half)*dx*sample_fact(2) 00570 DO ii=1,sample_res(1) 00571 pos(1)=xpos(1)+(REAL(ii, 8)-half)*dx*sample_fact(1) 00572 r2=sum((pos(1:nDim)-Particle%xloc(1:nDim))**2) 00573 r=sqrt(r2) 00574 esp=-GM/r+sKE 00575 IF (esp < 0d0) THEN !gas parcel is bound to sink particle 00576 v_r(1:nDim)=DOT_PRODUCT(v_rel(1:nDim), pos(1:nDim)-Particle%xloc(1:nDim))*(pos(1:nDim)-Particle%xloc(1:nDim))/r2 00577 v_phi(1:nDim)=v_rel(1:nDim)-v_r(1:nDim) 00578 jsp2=SUM(v_phi(1:nDim)**2)*r2 00579 IF (iCylindrical == WITHANGMOM) jsp2=jsp2+(q(i,j,k,iAngMom)/q(i,j,k,1))**2*r2 00580 IF (jsp2 < 2d0*r_surf*(esp*r_surf+GM)) THEN !gas parcel will come within "surface" of sink particle 00581 count(i,j,k)=count(i,j,k)+1 !fraction of cell that is accumulated 00582 END IF 00583 END IF 00584 END DO 00585 END DO 00586 END DO 00587 END DO 00588 END DO 00589 END DO 00590 ! write(*,'(9I5)') count 00591 DO k=mS(3,1),mS(3,2) 00592 xpos(3)=Info%xBounds(3,1)+(real(k+offset(3),8)-half)*dz 00593 iipos(3)=k 00594 DO j=mS(2,1),mS(2,2) 00595 xpos(2)=Info%xBounds(2,1)+(real(j+offset(2),8)-half)*dx 00596 iipos(2)=j 00597 DO i=mS(1,1),mS(1,2) 00598 xpos(1)=Info%xBounds(1,1)+(real(i+offset(1),8)-half)*dx 00599 iipos(1)=i 00600 IF (count(i,j,k) <= 0) CYCLE 00601 r2=sum((xpos(1:nDim)-Particle%xloc(1:nDim))**2) 00602 w=exp(-r2/Particle%Bondi_Kernel**2) 00603 drho=min(.25d0, Particle%AccretionRate*w*REAL(count(i,j,k),KIND=qPREC)*q_fact*levels(Info%level)%dt)*Info%q(i,j,k,1) !fraction of cell to accumulate 00604 IF (drho > 0d0) THEN 00605 IF (ALL(iipos(1:nDim) >= 1) .AND. ALL(iipos(1:nDim) <= Info%mX(1:nDim))) THEN 00606 IF (iCylindrical == NOCYL) THEN 00607 Particle%dQ(1:NrHydroVars)=Particle%dQ(1:NrHydrovars)+q(i,j,k,1:NrHydroVars)*sink_dv 00608 mass=sink_dv*(drho) 00609 rmass(1:nDim)=mass*(xpos(1:nDim)-Particle%xloc(1:nDim)) 00610 Particle%drmass=Particle%drmass+rmass 00611 IF (nDim == 2) THEN 00612 Particle%dJ=Particle%dJ+Cross2D(rmass(1:2),q(i,j,k,imom(1:2))/q(i,j,k,1)) 00613 ELSEIF (nDim == 3) THEN 00614 Particle%dJ(1:nDim)=Particle%dJ+Cross3D(rmass,q(i,j,k,imom(1:3))/q(i,j,k,1)) 00615 END IF 00616 ELSE 00617 Particle%dQ(1)=Particle%dQ(1)+q(i,j,k,1)*sink_dv*2d0*Pi*xpos(1) 00618 IF (iE /= 0) Particle%dQ(iE)=Particle%dQ(iE)+q(i,j,k,iE)*sink_dv*2d0*Pi*xpos(1) 00619 IF (iAngMom /= 0) Particle%dQ(iAngMom)=Particle%dQ(iAngMom)+q(i,j,k,iAngMom)*sink_dv*2d0*Pi*xpos(1) 00620 END IF 00621 END IF 00622 00623 IF (iE /= 0) THEN 00624 IF (lMHD) q(i,j,k,iE)=q(i,j,k,iE)-half*SUM(q(i,j,k,iBx:iBz)**2) !Subtract off magnetic energy 00625 q(i,j,k,iE)=q(i,j,k,iE)-half*SUM(q(i,j,k,m_low:m_high)**2)/q(i,j,k,1) !Subtract off kinetic energy 00626 END IF 00627 q(i,j,k,1:NrHydroVars)=q(i,j,k,1:NrHydroVars)*(q(i,j,k,1)-drho)/q(i,j,k,1) !Everything else should be accreted proportionally to density... 00628 IF (lMHD) THEN 00629 q(i,j,k,iBx)=.5*(Info%aux(i,j,k,1)+Info%aux(i+1,j,k,1)) 00630 q(i,j,k,iBy)=.5*(Info%aux(i,j,k,2)+Info%aux(i,j+1,k,2)) 00631 IF (nDim == 3) q(i,j,k,iBz)=.5*(Info%aux(i,j,k,3)+Info%aux(i,j,k+1,3)) 00632 END IF 00633 IF (iE /= 0) THEN 00634 q(i,j,k,iE)=q(i,j,k,iE)+half*SUM(q(i,j,k,m_low:m_high)**2)/q(i,j,k,1) !Add in kinetic energy 00635 IF (lMHD) q(i,j,k,iE)=q(i,j,k,iE)+half*SUM(q(i,j,k,iBx:iBz)**2) !Add magnetic energy back in. 00636 END IF 00637 IF (ALL(iipos(1:nDim) >= 1) .AND. ALL(iipos(1:nDim) <= Info%mX(1:nDim))) THEN 00638 IF (iCylindrical == NOCYL) THEN 00639 Particle%dQ(1:NrHydroVars)=Particle%dQ(1:NrHydrovars)-q(i,j,k,1:NrHydroVars)*sink_dv 00640 ELSE 00641 Particle%dQ(1)=Particle%dQ(1)-q(i,j,k,1)*sink_dv*2d0*Pi*xpos(1) 00642 IF (iE /= 0) Particle%dQ(iE)=Particle%dQ(iE)-q(i,j,k,iE)*sink_dv*2d0*Pi*xpos(1) 00643 IF (iAngMom /= 0) Particle%dQ(iAngMom)=Particle%dQ(iAngMom)-q(i,j,k,iAngMom)*sink_dv*2d0*Pi*xpos(1) 00644 END IF 00645 END IF 00646 END IF 00647 END DO 00648 END DO 00649 END DO 00650 DEALLOCATE(count) 00651 END DO 00652 DEALLOCATE(mSs, offsets) 00653 END IF 00654 END SUBROUTINE BondiAccretion 00655 00656 00657 00658 00659 00660 00667 SUBROUTINE GetGasForce(Particle, Info) 00668 INTEGER :: ir(3), il(3), ipos(3) 00669 TYPE(ParticleDef) :: Particle 00670 TYPE(InfoDef) :: Info 00671 REAL(KIND=qPREC), DIMENSION(:,:,:), ALLOCATABLE :: phi 00672 !We can do a third order reconstruction of phi at the cell center and use the coefficients 00673 !to shift the reconstruction to the particle's position and have a 2nd order reconstruction. 00674 REAL(KIND=qPREC) :: phi_x, phi_y, phi_z, phi_xx, phi_xy, phi_xz, phi_yy, phi_yz, phi_zz, 00675 phi_xxx, phi_xxy, phi_xxz, phi_xyy, phi_xyz, phi_xzz, phi_yyy, phi_yyz, phi_yzz, phi_zzz 00676 REAL(KIND=qPREC), PARAMETER, DIMENSION(5,4) :: DC5s = RESHAPE((/ 00677 0.0833333333333333d0, -0.666666666666667d0, 0.00000000000000d0, 0.666666666666667d0, -0.0833333333333333d0, 00678 -0.0416666666666666d0, 0.666666666666666d0, -1.25000000000000d0, 0.666666666666667d0, -0.0416666666666667d0, 00679 -0.0416666666666666d0, 0.666666666666666d0, -1.25000000000000d0, 0.666666666666667d0, -0.0416666666666667d0, 00680 -0.0416666666666666d0, 0.666666666666666d0, -1.25000000000000d0, 0.666666666666667d0, -0.0416666666666667d0/), (/5,4/)) 00681 00682 REAL(KIND=qPREC), PARAMETER, DIMENSION(3,2) :: DC3s = RESHAPE((/ 00683 -.5d0, 0d0, .5d0, 00684 .5d0, -1d0, .5d0/),(/3,2/)) 00685 00686 REAL(KIND=qPREC), DIMENSION(3) :: dpos,f_grav,pos 00687 INTEGER :: i,ip(3,2),j,k,m 00688 REAL(KIND=qPREC) :: weight, total_weight 00689 IF (ParticleInGrid(Particle, Info,ipos)) THEN 00690 00691 !Use surrounding 2**ndim cells to get gas force weighted by proximity 00692 dpos=(Particle%xloc-Info%xBounds(:,1))/sink_dx + half ! 00693 ip=1 00694 ip(1:nDim,1)=floor(dpos) 00695 ip(1:nDim,2)=ip(1:nDim,1)+1 00696 f_grav=0d0 00697 total_weight=0 00698 DO i=ip(1,1),ip(1,2) 00699 DO j=ip(2,1),ip(2,2) 00700 DO k=ip(3,1),ip(3,2) 00701 pos=(/i,j,k/) 00702 ir=pos 00703 il=pos 00704 weight=1d0/(sqrt(sum((real(pos(1:nDim))-dpos(1:nDim))**2))+1e-3) 00705 DO m=1,nDim 00706 ir(m)=pos(m)+1 00707 il(m)=pos(m)-1 00708 f_grav(m)=f_grav(m)+weight*((Info%q(ir(1),ir(2),ir(3),iPhiGas)-Info%q(il(1),il(2),il(3),iPhiGas))) 00709 ir(m)=pos(m) 00710 il(m)=pos(m) 00711 END DO 00712 total_weight=total_weight+weight 00713 END DO 00714 END DO 00715 END DO 00716 Particle%gas_accel=-f_grav/total_weight/sink_dx/2d0 00717 00718 END IF 00719 END SUBROUTINE GetGasForce 00720 00721 00723 00726 00731 SUBROUTINE CalcMomentContributions(Info, Particle, mS, offset) 00732 TYPE(InfoDef) :: Info 00733 TYPE(ParticleDef) :: Particle 00734 INTEGER, DIMENSION(3,2) :: mS 00735 INTEGER :: i,j,k,ipos(3) 00736 REAL(KIND=qpREC) :: r, pos(3), offset(3) 00737 IF (nMoments == 1) THEN 00738 DO i=mS(1,1),mS(1,2) 00739 DO j=mS(2,1),mS(2,2) 00740 DO k=mS(3,1),mS(3,2) 00741 pos=CellPos(Info,i,j,k)+offset-Particle%xloc 00742 r=sqrt(sum(pos(1:nDim)**2)) 00743 IF (r <= r_acc) THEN 00744 ! write(*,*) 'Info%q(',i,',',j,',',k,',1)=',Info%q(i,j,k,1) 00745 Particle%moments(1)=Particle%moments(1)+Info%q(i,j,k,1)*sink_dV 00746 END IF 00747 END DO 00748 END DO 00749 END DO 00750 ELSE 00751 PRINT*, 'error - only 1 moment is currently supported' 00752 STOP 00753 END IF 00754 END SUBROUTINE CalcMomentContributions 00755 00756 00763 SUBROUTINE CalcGridParticleOverlaps(Info, Particle, mTs, nOverlaps, offsets, nGhost) 00764 ! Calculates overlap bounds of cell-centered quantities for 00765 ! overlapping grids. 00766 TYPE(InfoDef) :: Info 00767 INTEGER, DIMENSION(3,2) :: SourcemGlobal 00768 INTEGER, DIMENSION(3,2) :: mO,mGlobal,iOffSet 00769 INTEGER, DIMENSION(3) :: pOffSet 00770 INTEGER, DIMENSION(:,:,:), POINTER :: mTs 00771 REAL(KIND=qPREC), DIMENSION(:,:), POINTER :: offsets 00772 REAL(KIND=qPREC), DIMENSION(27,3) :: MaxOffsets 00773 INTEGER, DIMENSION(27,3,2) :: MaxMTs 00774 INTEGER :: i,j,k,rmbc,nOverlaps,level 00775 INTEGER, OPTIONAL :: nGhost 00776 TYPE(ParticleDef) :: Particle 00777 INTEGER, DIMENSION(3) :: ipos 00778 level=Info%level 00779 ! Calculate particles position in index space 00780 ipos(1:nDim)=floor((Particle%xloc(1:nDim)-GxBounds(1:nDim,1))/levels(level)%dx)+1 00781 ipos(1:nDim)=max(min(ipos(1:nDim),levels(Info%level)%mX(1:nDim)),1) 00782 ! write(*,*) "sink_dx=",sink_dx 00783 ! write(*,*) "ipos=",ipos, info%level 00784 00785 SourcemGlobal(nDim+1:3,:)=1 00786 SourcemGlobal(1:nDim,1)=ipos(1:nDim)-particle%buffer(level) 00787 SourcemGlobal(1:nDim,2)=ipos(1:nDim)+particle%buffer(level) 00788 ! write(*,*) "sourcemglobal=", sourcemglobal 00789 rmbc=nGhost 00790 ! write(*,*) "rmbc=",rmbc 00791 NULLIFY(mTs) 00792 nOverlaps=0 00793 mO(nDim+1:3,:)=1 00794 ioffset=0 00795 WHERE(lHydroPeriodic(1:nDim)) ioffset(1:nDim,2)=1 00796 ioffset(1:nDim,1)=-ioffset(1:nDim,2) 00797 DO i=ioffset(1,1),ioffset(1,2) 00798 DO j=ioffset(2,1),ioffset(2,2) 00799 DO k=ioffset(3,1),ioffset(3,2) 00800 pOffSet=(/i,j,k/)*(GxBounds(:,2)-GxBounds(:,1)) 00801 mGlobal(:,:)=SourcemGlobal(:,:)+SPREAD(pOffSet,2,2) 00802 00803 mO(1:nDim,1)=max(Info%mGlobal(1:nDim,1)-rmbc& 00804 &,mGlobal(1:nDim,1)) 00805 mO(1:nDim,2)=min(Info%mGlobal(1:nDim,2)+rmbc& 00806 &,mGlobal(1:nDim,2)) 00807 00808 IF (ALL(mO(1:nDim,2) >= mO(1:nDim,1))) THEN 00809 nOverlaps=nOverlaps+1 00810 MaxMTs(nOverlaps,:,:)=mO-Spread( Info%mGlobal(:,1)& 00811 &,2,2)+1 00812 MaxOffsets(nOverlaps,:)=pOffSet 00813 END IF 00814 00815 ! write(*,*) "Considering", Info%mGlobal, SourcemGlobal, mO, nOverlaps 00816 END DO 00817 END DO 00818 END DO 00819 IF (nOverlaps > 0) THEN 00820 ALLOCATE(MTs(nOverlaps,3,2), offsets(nOverlaps,3)) 00821 MTs=MaxMTs(1:nOverlaps,:,:) 00822 offsets=MaxOffsets(1:nOverlaps,:) 00823 END IF 00824 00825 END SUBROUTINE CalcGridParticleOverlaps 00826 00827 00832 FUNCTION ParticleInGrid(Particle, Info, ipos) 00833 LOGICAL :: ParticleInGrid 00834 TYPE(InfoDef) :: Info 00835 INTEGER, DIMENSION(:) :: ipos 00836 Type(ParticleDef) :: Particle 00837 ipos=1 00838 ipos=PosCell(Info, Particle%xloc(1), Particle%xloc(2), Particle%xloc(3)) 00839 ParticleInGrid=ALL(ipos(1:nDim) >= 1) .AND. ALL(ipos(1:nDim) <= Info%mX(1:nDim)) 00840 ! write(*,*) 'processor', mpi_id, 'ipos=', ipos, Info%mX 00841 ! IF (ParticleInGrid) write(*,*) MPI_ID, Info%xBounds, Particle%xloc 00842 END FUNCTION ParticleInGrid 00843 00844 00858 00859 00861 00862 END MODULE ParticleInfoOps