Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! particle_advance.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 ParticleAdvance 00033 USE GlobalDeclarations 00034 USE ParticleDeclarations 00035 IMPLICIT NONE 00036 00037 PUBLIC AdvanceParticles 00038 PRIVATE 00039 00040 CONTAINS 00042 SUBROUTINE AdvanceParticles 00043 !Subroutine that advances the particles positions using the particles positions, velocities, and gas_accel 00044 REAL(KIND=qPREC), DIMENSION(:,:), ALLOCATABLE :: vel, pos, gas_accel, sink_accel 00045 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: gmasses 00046 LOGICAL, DIMENSION(:), ALLOCATABLE :: lFixed 00047 TYPE(ParticleListDef), POINTER :: particlelist 00048 TYPE(ParticleDef), POINTER :: particle 00049 INTEGER :: i,j, iters 00050 REAL(KIND=qPREC) :: dtmax, dt_remaining, dt_next, temp 00051 REAL(KIND=qPREC), DIMENSION(3) :: accumgas = 0 00052 INTEGER, PARAMETER :: maxiters=1e4 00053 ParticleMaxSpeed=0 00054 IF (NrSinkParticles == 0) RETURN 00055 00056 ! temp=SinkParticles%self%q(1)*SinkParticles%self%q(ivx) 00057 ALLOCATE(gmasses(NrSinkParticles), vel(NrSinkParticles,nDim), pos(NrSinkParticles,nDim), gas_accel(NrSinkParticles,nDerivatives), sink_accel(NrSinkParticles, nDim), lFixed(NrSinkParticles)) 00058 particleList=>SinkParticles 00059 i=0 00060 DO WHILE (ASSOCIATED(particlelist)) 00061 i=i+1 00062 particle=>particlelist%self 00063 gmasses(i)=particle%Q(1)*ScaleGrav 00064 vel(i,1:nDim)=particle%Q(imom(1:nDim)) 00065 pos(i,1:nDim)=particle%xloc(1:nDim) 00066 lFixed(i)=particle%lFixed 00067 particlelist=>particlelist%next 00068 00069 00070 END DO 00071 dt_remaining=levels(MaxLevel)%dt 00072 ! IF (MPI_ID == 0) write(*,*) 'preadvance pos=', pos 00073 ! IF (MPI_ID == 0) write(*,*) 'preadvance vel=', vel 00074 00075 ! sink_accel=gas_accel 00076 ! CALL GetGasAcceleration(particlelist, pos, sink_accel)!Start with gas acceleration 00077 CALL GetGasAcceleration(SinkParticles, pos, sink_accel) 00078 ! write(*,*) 'gas_acceleration=', sink_accel 00079 ! write(*,*) 'preadvance accel=', sink_accel 00080 CALL CalcSinkAcc(gmasses, pos, sink_accel, dtmax) !Add in particle-particle acceleration 00081 iters=0 00082 DO WHILE (dt_remaining > 0d0) 00083 iters=iters+1 00084 dt_next=min(dtmax,dt_remaining) 00085 FORALL(i=1:NrSinkParticles, .NOT. lFixed(i)) 00086 pos(i,1:nDim)=pos(i,1:nDim)+vel(i,1:nDim)*dt_next+half*(sink_accel(i,1:nDim))*dt_next**2 00087 vel(i,1:nDim)=vel(i,1:nDim)+half*sink_accel(i,1:nDim)*dt_next 00088 END FORALL 00089 ! sink_accel=gas_accel 00090 ! particleList=>SinkParticles - commented out by Erica 1/16/2012 00091 ! sink_accel=0 00092 CALL GetGasAcceleration(SinkParticles, pos, sink_accel)!Recalculate gas force at new position 00093 CALL CalcSinkAcc(gmasses, pos, sink_accel, dtmax) 00094 FORALL(i=1:NrSinkParticles, .NOT. lfixed(i)) 00095 vel(i,1:nDim)=vel(i,1:nDim)+half*sink_accel(i,1:nDim)*dt_next 00096 END FORALL 00097 dt_remaining=dt_remaining-dt_next 00098 IF (iters > MaxIters) THEN 00099 lRequestRestart=.true. 00100 IF (MPI_ID == 0) THEN 00101 write(*,*) 'too many iterations required ' 00102 write(*,*) 'pos=', pos 00103 END IF 00104 EXIT 00105 END IF 00106 END DO 00107 particleList=>SinkParticles 00108 i=0 00109 DO WHILE (ASSOCIATED(particlelist)) 00110 i=i+1 00111 particle=>particlelist%self 00112 ! gmasses(i)=particle%Q(1)*ScaleGrav 00113 particle%Q(imom(1:nDim))=vel(i,1:nDim) 00114 ParticleMaxSpeed=max(ParticleMaxSpeed, sqrt(sum(vel(i,1:nDim)**2))) 00115 particle%xloc(1:nDim)=pos(i,1:nDim) 00116 DO j=1,nDim 00117 ! write(*,*) particle%xloc(j), floor((particle%xloc(j)-GxBounds(j,1))/(GxBounds(j,2)-GxBounds(j,1))) 00118 IF (lHydroPeriodic(j)) particle%xloc(j)=particle%xloc(j)-(GxBounds(j,2)-GxBounds(j,1))*floor((particle%xloc(j)-GxBounds(j,1))/(GxBounds(j,2)-GxBounds(j,1))) 00119 IF (lHydroPeriodic(j) .AND. (particle%xloc(j) < GxBounds(j,1) .OR. particle%xloc(j) > GxBounds(j,2))) THEN 00120 write(*,*) 'whoops in particleadvance', j,particle%xloc(j), GxBounds(j,:) 00121 ! STOP 00122 END IF 00123 END DO 00124 00125 ! gas_accel(i,1:nDim)=particle%gas_accel(1:nDim) 00126 particlelist=>particlelist%next 00127 00128 ! IF (MPI_ID == 0) THEN 00129 ! write(*,*) 'pos= ', particle%xloc 00130 ! write(*,*) 'vel= ', particle%vel 00131 ! write(*,*) 'fgas= ', particle%Q(1)*particle%gas_accel(1:nDim)*levels(MaxLevel)%dt 00132 ! accumgas=accumgas+particle%Q(1)*particle%gas_accel(1:nDim)*levels(MaxLevel)%dt 00133 ! write(*,*) 'accumgas=', accumgas 00134 ! END IF 00135 END DO 00136 ! write(*,*) 'momenta', SinkParticles%self%q(1)*SinkParticles%self%q(ivx), temp, temp+sink_accel(1,1)*levels(MaxLevel)%dt*SinkParticles%self%q(1) 00137 ! accumgas(1)=accumgas(1)+sink_accel(1,1)*levels(MaxLevel)%dt*SinkParticles%self%q(1) 00138 ! write(*,*) 'accumulated gas back force', accumgas(1) 00139 DEALLOCATE(gmasses, vel, pos, gas_accel, sink_accel, lFixed) 00140 00141 END SUBROUTINE AdvanceParticles 00142 00145 00150 SUBROUTINE GetGasAcceleration(ParticleList, pos, sink_accel) 00151 TYPE(ParticleListDef), POINTER :: ParticleList, TempParticleList 00152 TYPE(ParticleDef), POINTER :: Particle 00153 REAL(KIND=qPREC), DIMENSION(:,:) :: pos, sink_accel 00154 INTEGER :: i, j 00155 i=0 00156 TempParticleList=>ParticleList 00157 DO WHILE (ASSOCIATED(TempParticleList)) 00158 i=i+1 00159 Particle=>TempParticleList%self 00160 IF (ASSOCIATED(Particle%pointgravityobj) .and. Particle%q(1)>0d0& 00161 .AND. .NOT.Particle%lFixed) & !jonathan added 8mar '12 00162 THEN 00163 sink_accel(i,1:ndim)=(/(SUM(Particle%pointgravityobj%dmom(0:maxlevel,j)/(levels(0:maxlevel)%dt)), j=1, ndim)/)/particle%q(1) 00164 ELSE 00165 Sink_accel(i,1:ndim) = 0d0 00166 END IF 00167 00168 00169 ! sink_accel(i,:)=Particle%gas_accel(nDim+1:2*nDim) !First derivatives at particles initial position (phi_x, phi_y, phi_z) 00170 ! IF (nDim == 2) THEN 00171 ! sink_accel=sink_accel+(pos(i,nDim:1:-1)-Particle%xloc(nDim:1:-1))*Particle%gas_accel(2*nDim+1) !Cross Derivatives (phi_xy) 00172 ! ELSE 00173 ! sink_accel=sink_accel+(/ & 00174 ! (pos(i,(/2,3/))-Particle%xloc((/2,3/)))*Particle%gas_accel((/7,9/)), & !Cross Derivatives (phi_xy phi_xz) 00175 ! (pos(i,(/1,3/))-Particle%xloc((/1,3/)))*Particle%gas_accel((/7,8/)), & !Cross Derivatives (phi_yx phi_yz) 00176 ! (pos(i,(/1,2/))-Particle%xloc((/1,2/)))*Particle%gas_accel((/9,8/))/) !Cross Derivatives (phi_zx phi_zx) 00177 ! END IF 00178 TempParticleList=>TempParticleList%next 00179 END DO 00180 END SUBROUTINE GetGasAcceleration 00181 00182 00188 SUBROUTINE CalcSinkAcc(gmasses, pos, sink_accel,dtmax) 00189 USE CommonFunctions 00190 REAL(KIND=qPREC), DIMENSION(:,:) :: pos, sink_accel 00191 REAL(KIND=qPREC), DIMENSION(:) :: gmasses 00192 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: minradii 00193 REAL(KIND=qpREC) :: dtmax, radius, gacc(3),poffset(3) 00194 INTEGER :: i, j, nParticles, ii, jj, kk, ioffset(3,2) 00195 nParticles=size(gmasses) 00196 ALLOCATE(minradii(nParticles)) 00197 minradii=huge(minradii) 00198 dtmax=huge(dtmax) 00199 ioffset=0 00200 WHERE(lHydroPeriodic(1:nDim)) ioffset(1:nDim,2)=1 00201 ioffset(1:nDim,1)=-ioffset(1:nDim,2) 00202 DO i=1,nParticles 00203 DO j=i+1,nParticles 00204 DO ii=ioffset(1,1),ioffset(1,2) 00205 DO jj=ioffset(2,1),ioffset(2,2) 00206 DO kk=ioffset(3,1),ioffset(3,2) 00207 pOffSet=(/ii,jj,kk/)*(GxBounds(:,2)-GxBounds(:,1)) 00208 radius=sqrt(SUM((pos(i,:)+pOffSet-pos(j,:))**2)) 00209 minradii(i)=min(radius,minradii(i)) 00210 minradii(j)=min(radius,minradii(j)) 00211 gacc(:)=SplineSoftening(pos(j,:)-pos(i,:)-pOffset,r_soft) 00212 sink_accel(i,:)=sink_accel(i,:)+gacc(1:nDim)*gmasses(j) 00213 sink_accel(j,:)=sink_accel(j,:)-gacc(1:nDim)*gmasses(i) 00214 END DO 00215 END DO 00216 END DO 00217 END DO 00218 dtmax=min(dtmax,CFL_LEAPFROG*sqrt(sink_dx/sqrt(sum(sink_accel(i,:)**2)))) 00219 END DO 00220 DEALLOCATE(minradii) 00221 00222 END SUBROUTINE CalcSinkAcc 00223 00225 00226 END MODULE ParticleAdvance