Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! particle_declarations.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 00030 00033 MODULE ParticleDeclarations 00034 USE GlobalDeclarations 00035 USE PhysicsDeclarations 00036 USE PointGravitySrc 00037 USE Outflows 00038 IMPLICIT NONE 00039 SAVE 00040 PUBLIC 00041 INTEGER, PARAMETER :: nMoments=1 00042 INTEGER, PARAMETER :: MaxHydroVars = 20 00043 00044 INTEGER, PARAMETER :: NOACCRETION=0, FEDERRATH_ACCRETION=1, KRUMHOLZ_ACCRETION=2 00045 00046 TYPE, PUBLIC :: ParticleDef 00047 ! REAL(KIND=qPREC) :: mass !particle mass 00048 REAL(KIND=qPREC), DIMENSION(MaxHydroVars) :: Q !Conserved quantities 00049 REAL(KIND=qPREC), DIMENSION(3) :: xloc !location of particle 00050 REAL(KIND=qPREC), DIMENSION(3) :: J 00051 REAL(KIND=qPREC), DIMENSION(MaxHydroVars) :: dQ !Accreted contributions to quantities 00052 REAL(KIND=qPREC), DIMENSION(3) :: dJ 00053 REAL(KIND=qPREC), DIMENSION(3) :: drmass 00054 REAL(KIND=qPREC), DIMENSION(nMoments) :: moments 00055 REAL(KIND=qPREC), DIMENSION(3) :: gas_accel 00056 REAL(KIND=qPREC) :: rho_inf !Used for Bondi Accretion 00057 INTEGER :: radius=4 00058 INTEGER :: Buffer(0:MaxDepth)=0 !number of extra cells of buffer around sink region to set err flags 00059 LOGICAL :: lFixed=.false. !Particle is free to move 00060 INTEGER :: iAccrete=NOACCRETION !particle accretion routine 00061 REAL(KIND=qPREC) :: AccretionRate=0d0 00062 REAL(KIND=qPREC) :: Bondi_kernel=0d0 00063 TYPE(PointGravityDef), POINTER :: PointGravityObj ! Potential Source object 00064 TYPE(OutflowDef), POINTER :: OutflowObj ! Potential Source object 00065 END TYPE ParticleDef 00066 00068 TYPE, PUBLIC :: ParticleListDef 00069 TYPE(ParticleDef), POINTER :: self=>Null() 00070 TYPE(ParticleListDef), POINTER :: next=>Null() 00071 END TYPE ParticleListDef 00072 00073 00074 PUBLIC :: CreateParticle 00075 PUBLIC :: IO_TRACKED_COMPONENTS, Particle_InitChomboDatasets, Particle_WriteObjectToChomboFile 00076 00077 00078 REAL(KIND=qPREC), PARAMETER :: JEAN_CELLS=4 !Number of cells needed to refine jeans length 00079 INTEGER, PARAMETER :: IR_ACC = 4 !Accretion radius in cells 00080 REAL(KIND=qpREC), PARAMETER :: CFL_LEAPFROG=.1 !CFL criteria for leap frog method 00081 REAL(KIND=qPREC) :: r_acc, r_acc2, r_inner_acc, r_soft!Accretion radius and inner accretion radius in physical units 00082 REAL(KIND=qPREC) :: JeansFact !Factor involving scalegrav, pi, gamma, and the fine grid dx 00083 REAL(KIND=qPREC) :: sink_dv, sink_dx !volume of fine grid cells, size of fine grid cells 00084 REAL(KIND=qPREC) :: ScaleGravdV !Product of ScaleGrav and cell volume dV 00085 INTEGER :: nParticleFields 00086 LOGICAL, DIMENSION(:,:,:), ALLOCATABLE :: lControlVolume 00087 TYPE(ParticleListDef), POINTER :: SinkParticles, LastSinkParticle, NewSinkParticles, LastNewSinkParticle, BackupParticles 00088 LOGICAL :: HaveSinkParticles=.false. 00089 INTEGER :: NrSinkParticles 00090 INTEGER :: nAngularMomentum 00091 INTEGER :: nDerivatives 00092 REAL(KIND=qPREC) :: ParticleMaxSpeed 00093 REAL(KIND=qPREC) :: DefaultParticleRadius = IR_ACC 00094 INTEGER, PARAMETER :: IO_TRACKED_COMPONENTS = 11 00095 00096 CONTAINS 00097 00100 SUBROUTINE CreateParticle(particle) 00101 TYPE(ParticleDef), POINTER :: particle 00102 IF (Associated(particle)) THEN 00103 PRINT*, 'particle already associataed in CreateParticle' 00104 STOP 00105 END IF 00106 ALLOCATE(particle) 00107 IF (NrHydroVars > MaxHydroVars) THEN ! 00108 PRINT*, 'err increase maxhydrovars in particle_declarations.' 00109 STOP 00110 END IF 00111 ! write(*,*) 'made a particle' 00112 Particle%J=0 00113 Particle%dQ=0 00114 Particle%dJ=0 00115 Particle%drmass=0 00116 particle%radius=DefaultParticleRadius 00117 particle%Buffer=0 00118 particle%lfixed=.false. 00119 particle%Q=0 00120 particle%xloc=0 00121 particle%moments=0 00122 particle%gas_accel=0 00123 00124 NULLIFY(particle%PointGravityObj) 00125 NULLIFY(particle%OutflowObj) 00126 particle%iaccrete=0 00127 particle%buffer(0:MaxLevel)=ceiling(real(max(particle%radius, IR_ACC)*sink_dx)/levels(0:MaxLevel)%dx) 00128 00129 ! [BDS][20110114]: On a restart, the point gravity object will be read in from the data files. 00130 ! Consequently, There is no need to create the particle here. The point gravity object will be 00131 ! associated in the ParticleDeclarations::Particle_ReadObjectFromChombo() subroutine. 00132 ! IF (.NOT. lRestart) THEN 00133 ! CALL CreatePointGravityObject(Particle%PointGravityObj) 00134 ! Particle%PointGravityObj%soft_length=Particle%radius*sink_dx 00135 ! Particle%PointGravityObj%soft_function=SPLINESOFT 00136 ! ELSE 00137 ! NULLIFY(Particle%PointGravityObj) 00138 ! END IF 00139 00140 END SUBROUTINE CreateParticle 00141 00146 SUBROUTINE AddParticleToList(particle, particlelist, lastparticle) 00147 TYPE(Particlelistdef), POINTER :: lastparticle 00148 TYPE(Particlelistdef), POINTER, OPTIONAL :: particlelist 00149 TYPE(ParticleDef), POINTER :: particle 00150 INTEGER :: iErr 00151 IF (.NOT. ASSOCIATED(LastParticle)) THEN 00152 ALLOCATE(LastParticle, STAT=iErr) 00153 NULLIFY(LastParticle%next) 00154 NULLIFY(LastParticle%self) 00155 IF (iErr /= 0) THEN 00156 PRINT *, "AddParticle() error: unable to allocate LastParticle list object." 00157 STOP 00158 END IF 00159 IF (present(particlelist)) particlelist=>lastparticle 00160 ELSE 00161 IF (ASSOCIATED(LastParticle%next)) THEN 00162 PRINT *, "Error - last particle next allocated" 00163 STOP 00164 END IF 00165 ALLOCATE(LastParticle%next, STAT=iErr) 00166 IF (iErr /= 0) THEN 00167 PRINT *, "AddParticle() error: unable to allocate LastParticle%next list object." 00168 STOP 00169 END IF 00170 LastParticle=>LastParticle%next 00171 NULLIFY(LastParticle%next) 00172 NULLIFY(LastParticle%self) 00173 00174 END IF 00175 IF (ASSOCIATED(particle)) THEN 00176 LastParticle%self=>particle 00177 ELSE 00178 ALLOCATE(LastParticle%self, STAT=iErr) 00179 PRINT*, 'Expecting particle to be associated in AddParticleToList' 00180 00181 IF (iErr /= 0) THEN 00182 PRINT *, "AddParticle() error: unable to allocate LastParticle%self object." 00183 STOP 00184 END IF 00185 particle=>LastParticle%self 00186 END IF 00187 END SUBROUTINE AddParticleToList 00188 00191 RECURSIVE SUBROUTINE ClearParticleList(particlelist) 00192 TYPE(Particlelistdef), POINTER :: particlelist 00193 IF (.NOT. ASSOCIATED(particlelist)) RETURN 00194 IF (ASSOCIATED(particlelist%next)) CALL ClearParticleList(particlelist%next) 00195 DEALLOCATE(particlelist) 00196 NULLIFY(particlelist) 00197 END SUBROUTINE ClearParticleList 00198 00199 00202 SUBROUTINE DestroyParticle(particle) 00203 TYPE(ParticleDef), POINTER :: particle 00204 IF (.NOT. ASSOCIATED(particle)) THEN 00205 PRINT*, 'asked to destroy unassociated particle in DestroyParticle' 00206 ELSE 00207 ! IF (ASSOCIATED(Particle%PointGravityObj)) CALL DestroyPointGravityObject(Particle%PointGravityObj) 00208 ! IF (ASSOCIATED(Particle%OutflowObj)) CALL DestroyOutflowObject(Particle%OutflowObj) 00209 DEALLOCATE(particle) 00210 NULLIFY(particle) 00211 END IF 00212 END SUBROUTINE DestroyParticle 00213 00214 00217 RECURSIVE SUBROUTINE DestroyParticleList(particlelist) 00218 TYPE(Particlelistdef), POINTER :: particlelist 00219 IF (.NOT. ASSOCIATED(particlelist)) RETURN 00220 IF (ASSOCIATED(particlelist%next)) CALL DestroyParticleList(particlelist%next) 00221 IF (ASSOCIATED(particlelist%self)) CALL DestroyParticle(particlelist%self) 00222 DEALLOCATE(particlelist) 00223 NULLIFY(particlelist) 00224 END SUBROUTINE DestroyParticleList 00225 00228 FUNCTION ParticleCount(particlelist) 00229 INTEGER :: ParticleCount 00230 TYPE(ParticleListDef), POINTER :: particlelist, temp 00231 ParticleCount=0 00232 00233 temp=>particlelist 00234 DO WHILE (ASSOCIATED(temp)) 00235 IF (ASSOCIATED(temp%self)) ParticleCount = ParticleCount + 1 00236 temp => temp%next 00237 END DO 00238 END FUNCTION ParticleCount 00239 00242 SUBROUTINE AddNewSinkParticle(Particle) 00243 TYPE(ParticleDef), POINTER :: Particle 00244 CALL AddParticleToList(Particle, NewSinkParticles, LastNewSinkParticle) 00245 ! CALL PrintParticle(Particle) 00246 END SUBROUTINE AddNewSinkParticle 00247 00250 SUBROUTINE AddSinkParticle(Particle) 00251 TYPE(ParticleDef), POINTER :: Particle 00252 CALL AddParticleToList(Particle, SinkParticles, LastSinkParticle) 00253 NrSinkParticles=NrSinkParticles+1 00254 CALL PrintParticle(Particle) 00255 lParticles=.true. 00256 lSinkParticles=.true. 00257 END SUBROUTINE AddSinkParticle 00258 00261 SUBROUTINE PrintParticleList(ParticleList) 00262 TYPE(ParticleListDef), POINTER :: ParticleList 00263 TYPE(ParticleListDef), POINTER :: TempParticleList 00264 TempParticleList=>ParticleList 00265 00266 DO WHILE (ASSOCIATED(TempParticleList)) 00267 CALL PrintParticle(TempParticleList%self) 00268 TempParticleList=>TempParticleList%next 00269 END DO 00270 00271 END SUBROUTINE PrintParticleList 00272 00275 SUBROUTINE PrintParticle(Particle) 00276 TYPE(ParticleDef) :: Particle 00277 IF (MPI_ID == 0) THEN 00278 write(*,'(A,E18.6)') "Mass = ", Particle%Q(1) 00279 write(*,'(A,3E18.6)') "Position = ", Particle%xloc(1:nDim) 00280 write(*,'(A,3E18.6)') "Velocity = ", Particle%Q(imom(1:nDim)) 00281 ! write(*,'(A,3E18.6)') "Angular Momentum = ", Particle%J(1:nAngularMomentum) 00282 ! write(*,'(A,20E18.6)') "Moments = ", Particle%moments 00283 ! write(*,'(A,3E18.6)') "Gas acceleration = ", Particle%gas_accel(1:nDim) 00284 ! write(*,'(A,20E18.6)') "Particle%Q", Particle%Q(1:NrHydroVars) 00285 END IF 00286 END SUBROUTINE PrintParticle 00287 00288 00294 FUNCTION GetBindingEnergy(Particle, pos, q, r) 00295 REAL(KIND=qPREC) :: GetBindingEnergy 00296 TYPE(ParticleDef) :: Particle 00297 REAL(KIND=qPREC), DIMENSION(:) :: pos 00298 REAL(KIND=qPREC), DIMENSION(:) :: q 00299 REAL(KIND=qPREC) :: r 00300 ! r=sqrt(sum((pos(1:nDim)-Particle%xloc(1:nDim))**2)) 00301 IF (nDim == 3 .OR. iCylindrical /= NOCYL) THEN 00302 GetBindingEnergy=SUM((q(m_low:m_high)/q(1)-Particle%Q(m_low:m_high))**2)-ScaleGrav*(Particle%Q(1)+Particle%moments(1))/r 00303 ELSE 00304 GetBindingEnergy=SUM((q(m_low:m_high)/q(1)-Particle%Q(m_low:m_high))**2)+2d0*ScaleGrav*(Particle%Q(1)+Particle%moments(1))*log(r/R2DEff) 00305 END IF 00306 END FUNCTION GetBindingEnergy 00307 00308 00310 INTEGER FUNCTION Particle_CountSinkParticles() 00311 00312 TYPE(ParticleListDef), POINTER :: part_list 00313 INTEGER :: working_count 00314 00315 00316 part_list => SinkParticles 00317 working_count = 0 00318 00319 DO WHILE (ASSOCIATED(part_list)) 00320 working_count = working_count + 1 00321 part_list => part_list%next 00322 END DO 00323 00324 Particle_CountSinkParticles = working_count 00325 00326 END FUNCTION Particle_CountSinkParticles 00327 00328 00331 SUBROUTINE Particle_InitChomboDatasets(chandle) 00332 00333 USE ChomboDeclarations, ONLY: ChomboHandle 00334 USE HDF5Declarations, ONLY: Initialize_HDF5_Dataset_Int, Initialize_HDF5_Dataset_Double 00335 00336 TYPE(ChomboHandle), POINTER :: chandle 00337 INTEGER :: i 00338 CHARACTER(LEN=30) :: parttracername 00339 IF (.NOT. ASSOCIATED(chandle)) THEN 00340 PRINT *, "Particle_InitChomboDatasets error::chombo handle not associated." 00341 STOP 00342 END IF 00343 00344 ! Position coordinates 00345 CALL Initialize_HDF5_Dataset_Double("position_x", chandle%particle_group_id, NrSinkParticles) 00346 CALL Initialize_HDF5_Dataset_Double("position_y", chandle%particle_group_id, NrSinkParticles) 00347 CALL Initialize_HDF5_Dataset_Double("position_z", chandle%particle_group_id, NrSinkParticles) 00348 00349 ! Particle mass 00350 CALL Initialize_HDF5_Dataset_Double("P_mass", chandle%particle_group_id, NrSinkParticles) 00351 00352 ! velocity 00353 CALL Initialize_HDF5_Dataset_Double("P_vx", chandle%particle_group_id, NrSinkParticles) 00354 CALL Initialize_HDF5_Dataset_Double("P_vy", chandle%particle_group_id, NrSinkParticles) 00355 IF (ivz /= 0) CALL Initialize_HDF5_Dataset_Double("P_vz", chandle%particle_group_id, NrSinkParticles) 00356 00357 IF (iE /= 0) CALL Initialize_HDF5_Dataset_Double("P_E", chandle%particle_group_id, NrSinkParticles) 00358 ! Creates tags for tracer variables. 00359 DO i = 1, NrTracerVars !NrCons+1,NrHydroVars 00360 write(parttracername,'(A2,A)') 'P_', trim(TracerName(i)) 00361 CALL Initialize_HDF5_Dataset_Double(trim(parttracername),chandle%particle_group_id, NrSinkParticles) 00362 END DO 00363 00364 ! Angular momentum 00365 CALL Initialize_HDF5_Dataset_Double("P_Jx", chandle%particle_group_id, NrSinkParticles) 00366 CALL Initialize_HDF5_Dataset_Double("P_Jy", chandle%particle_group_id, NrSinkParticles) 00367 CALL Initialize_HDF5_Dataset_Double("P_Jz", chandle%particle_group_id, NrSinkParticles) 00368 00369 00370 ! Particle radius in cells 00371 CALL Initialize_HDF5_Dataset_Int("P_radius", chandle%particle_group_id, NrSinkParticles) 00372 00373 ! A converted logical flag, indicating that the particle is free to move. 00374 CALL Initialize_HDF5_Dataset_Int("P_lFixed", chandle%particle_group_id, NrSinkParticles) 00375 00376 ! A converted logical flag, indicating that the particle is free to move. 00377 CALL Initialize_HDF5_Dataset_Int("P_iAccrete", chandle%particle_group_id, NrSinkParticles) 00378 00379 ! The cell buffer around each level. 00380 CALL Initialize_HDF5_Dataset_Int("Buffer", chandle%particle_group_id, NrSinkParticles * (MaxDepth + 1)) 00381 00382 ! A tag for the particle acceleration routine. 00383 CALL Initialize_HDF5_Dataset_Int("point_gravity_id", chandle%particle_group_id, NrSinkParticles) 00384 00385 END SUBROUTINE Particle_InitChomboDatasets 00386 00390 SUBROUTINE Particle_WriteObjectToChomboFile(chandle, particle) 00391 00392 USE ChomboDeclarations, ONLY: ChomboHandle 00393 USE PointGravitySrc, ONLY: PointGravityDef 00394 USE HDF5Declarations, ONLY: Write_Slab_To_Dataset_Int, Write_Slab_To_Dataset_Double 00395 00396 TYPE(ChomboHandle), POINTER :: chandle 00397 TYPE(ParticleDef), POINTER :: particle 00398 INTEGER :: i 00399 CHARACTER(LEN=30) :: parttracername 00400 00401 IF (.NOT. ASSOCIATED(chandle)) THEN 00402 PRINT *, "Particle_WriteObjectToChomboFile error::Chombo handle not associated." 00403 STOP 00404 END IF 00405 00406 IF (.NOT. ASSOCIATED(particle)) THEN 00407 PRINT *, "Particle_WriteObjectToChomboFile error::particle not associated." 00408 STOP 00409 END IF 00410 00411 ! Write position components to datasets. Be sure to index the xloc array using the i:i notation, as 00412 ! this passes a 1-element array to the subroutine instead of a scalar. 00413 CALL Write_Slab_To_Dataset_Double("position_x", chandle%particle_group_id, particle%xloc(1:1), chandle%particle_offset) 00414 CALL Write_Slab_To_Dataset_Double("position_y", chandle%particle_group_id, particle%xloc(2:2), chandle%particle_offset) 00415 CALL Write_Slab_To_Dataset_Double("position_z", chandle%particle_group_id, particle%xloc(3:3), chandle%particle_offset) 00416 00417 ! Write mass to dataset. 00418 CALL Write_Slab_To_Dataset_Double("P_mass", chandle%particle_group_id, (/ particle%Q(1) /), chandle%particle_offset) 00419 00420 ! Write velocity components to datasets. 00421 CALL Write_Slab_To_Dataset_Double("P_vx", chandle%particle_group_id, particle%Q(ivx:ivx), chandle%particle_offset) 00422 CALL Write_Slab_To_Dataset_Double("P_vy", chandle%particle_group_id, particle%Q(ivy:ivy), chandle%particle_offset) 00423 IF (ivz /= 0) THEN 00424 CALL Write_Slab_To_Dataset_Double("P_vz", chandle%particle_group_id, particle%Q(ivz:ivz), chandle%particle_offset) 00425 END IF 00426 00427 IF (iE /= 0) CALL Write_Slab_To_Dataset_Double("P_E", chandle%particle_group_id, particle%Q(iE:iE),chandle%particle_offset) 00428 ! Creates tags for tracer variables. 00429 DO i = 1, NrTracerVars !NrCons+1,NrHydroVars 00430 write(parttracername,'(A2,A)') 'P_', trim(TracerName(i)) 00431 CALL Write_Slab_To_Dataset_Double(trim(partTracerName),chandle%particle_group_id, particle%Q(nTracerLo-1+i:nTracerLo-1+i),chandle%particle_offset) 00432 END DO 00433 00434 ! Write angular momentum components to datasets. 00435 CALL Write_Slab_To_Dataset_Double("P_Jx", chandle%particle_group_id, particle%J(1:1), chandle%particle_offset) 00436 CALL Write_Slab_To_Dataset_Double("P_Jy", chandle%particle_group_id, particle%J(2:2), chandle%particle_offset) 00437 CALL Write_Slab_To_Dataset_Double("P_Jz", chandle%particle_group_id, particle%J(3:3), chandle%particle_offset) 00438 00439 CALL Write_Slab_To_Dataset_Int("P_radius", chandle%particle_group_id, (/ particle%radius /), chandle%particle_offset) 00440 00441 CALL Write_Slab_To_Dataset_Int("P_lFixed", & 00442 chandle%particle_group_id, & 00443 (/ BoolToInt(particle%lFixed) /), & 00444 chandle%particle_offset) 00445 00446 CALL Write_Slab_To_Dataset_Int("P_iAccrete", & 00447 chandle%particle_group_id, & 00448 (/ (particle%iAccrete) /), & 00449 chandle%particle_offset) 00450 00451 00452 ! Write the whole buffer array to the appropriate dataset. 00453 CALL Write_Slab_To_Dataset_Int("Buffer", chandle%particle_group_id, particle%Buffer, chandle%particle_offset * (MaxDepth + 1)) 00454 00455 IF (ASSOCIATED(Particle%PointGravityObj)) THEN 00456 cALL Write_Slab_To_Dataset_Int("point_gravity_id", & 00457 chandle%particle_group_id, & 00458 (/ particle%PointGravityObj%id /), & 00459 chandle%particle_offset) 00460 ELSE 00461 cALL Write_Slab_To_Dataset_Int("point_gravity_id", & 00462 chandle%particle_group_id, & 00463 (/ -1 /), & 00464 chandle%particle_offset) 00465 END IF 00466 chandle%particle_offset = chandle%particle_offset + 1 00467 00468 END SUBROUTINE Particle_WriteObjectToChomboFile 00469 00473 SUBROUTINE Particle_ReadObjectFromChombo(chandle, particle) 00474 00475 USE ChomboDeclarations, ONLY: ChomboHandle 00476 USE HDF5Declarations, ONLY: Read_Slab_From_Dataset_Int, Read_Slab_From_Dataset_Double 00477 00478 TYPE(ChomboHandle), POINTER :: chandle 00479 TYPE(ParticleDef), POINTER :: particle 00480 00481 INTEGER, DIMENSION(1), TARGET :: int_buffer_array 00482 REAL(KIND=qPrec), DIMENSION(1), TARGET :: dbl_buffer_array 00483 INTEGER, DIMENSION(:), POINTER :: int_buffer 00484 REAL(KIND=qPrec), DIMENSION(:), POINTER :: dbl_buffer 00485 INTEGER :: i 00486 CHARACTER(LEN=30) :: parttracername 00487 00488 ! Fail if a null pointer is passed in for a ChomboHandle. 00489 IF (.NOT. ASSOCIATED(chandle)) THEN 00490 PRINT "('Particle_ReadObjectFromChombo() error::particle number ', i4, ' object not created.')", & 00491 chandle%particle_offset 00492 STOP 00493 END IF 00494 00495 ! Fail if a null pointer is passed in for a particle object. 00496 IF (.NOT. ASSOCIATED(particle)) THEN 00497 PRINT "('Particle_ReadObjectFromChombo() error::particle number ', i4, ' object not created.')", & 00498 chandle%particle_offset 00499 STOP 00500 END IF 00501 00502 int_buffer => int_buffer_array 00503 dbl_buffer => dbl_buffer_array 00504 00505 int_buffer = 0 00506 dbl_buffer = 0.d0 00507 00508 ! Read in the position coordinates for this particle. This takes three reads from the file, because 00509 ! the coordinates stored in three separate arrays. 00510 CALL Read_Slab_From_Dataset_Double("position_x", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00511 particle%xloc(1) = dbl_buffer(1) 00512 00513 CALL Read_Slab_From_Dataset_Double("position_y", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00514 particle%xloc(2) = dbl_buffer(1) 00515 00516 CALL Read_Slab_From_Dataset_Double("position_z", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00517 particle%xloc(3) = dbl_buffer(1) 00518 00519 ! Read in the particle mass. 00520 CALL Read_Slab_From_Dataset_Double("P_mass", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00521 particle%Q(1) = dbl_buffer(1) 00522 00523 ! Read in the velocity for this particle. This takes three reads from the file, because 00524 ! the angular momentum components are stored in three separate arrays. 00525 CALL Read_Slab_From_Dataset_Double("P_vx", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00526 particle%Q(ivx) = dbl_buffer(1) 00527 00528 CALL Read_Slab_From_Dataset_Double("P_vy", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00529 particle%Q(ivy) = dbl_buffer(1) 00530 00531 IF (ivz /= 0) THEN 00532 CALL Read_Slab_From_Dataset_Double("P_vz", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00533 particle%Q(ivz) = dbl_buffer(1) 00534 END IF 00535 00536 IF (iE /= 0) THEN 00537 CALL Read_Slab_From_Dataset_Double("P_E", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00538 particle%Q(iE) = dbl_buffer(1) 00539 END IF 00540 00541 ! Creates tags for tracer variables. 00542 DO i = 1, NrTracerVars !NrCons+1,NrHydroVars 00543 write(parttracername,'(A2,A)') 'P_', trim(TracerName(i)) 00544 CALL Read_Slab_From_Dataset_Double(trim(partTracerName), chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00545 particle%Q(nTracerlo-1+i)=dbl_buffer(1) 00546 END DO 00547 00548 00549 ! Read in the angular momentum for this particle. This takes three reads from the file, because 00550 ! the angular momentum components are stored in three separate arrays. 00551 CALL Read_Slab_From_Dataset_Double("P_Jx", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00552 particle%J(1) = dbl_buffer(1) 00553 00554 CALL Read_Slab_From_Dataset_Double("P_Jy", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00555 particle%J(2) = dbl_buffer(1) 00556 00557 CALL Read_Slab_From_Dataset_Double("P_Jz", chandle%particle_group_id, dbl_buffer, chandle%particle_offset) 00558 particle%J(3) = dbl_buffer(1) 00559 00560 CALL Read_Slab_From_Dataset_Int("P_radius", chandle%particle_group_id, int_buffer, chandle%particle_offset) 00561 particle%radius = int_buffer(1) 00562 00563 ! The ALL mask converts the integer value from the chombo file into a logical value. 00564 CALL Read_Slab_From_Dataset_Int("P_lFixed", chandle%particle_group_id, int_buffer, chandle%particle_offset) 00565 particle%lFixed = ALL(int_buffer /= 0) 00566 00567 CALL Read_Slab_From_Dataset_Int("P_iAccrete", chandle%particle_group_id, int_buffer, chandle%particle_offset) 00568 particle%iAccrete = int_buffer(1) 00569 00570 00571 ! Most of the objects in particle are scalars, or their components are stored as scalars. Particle%Buffer, 00572 ! on the other hand, is an array whose components are stored as a contiguous array. Consequently, we 00573 ! can pass particle%Buffer directly to the Read_Slab function. The buffer is of size MaxDepth + 1, hence the 00574 ! multiplier for the particle offset. 00575 CALL Read_Slab_From_Dataset_Int("Buffer", chandle%particle_group_id, particle%Buffer, & 00576 chandle%particle_offset * (MaxDepth + 1)) 00577 00578 ! Retrieve the ID of the point gravity object associated with this particle, then find that object 00579 ! and point the particle's PointGravityObj toward it. 00580 CALL Read_Slab_From_Dataset_Int("point_gravity_id", chandle%particle_group_id, int_buffer, chandle%particle_offset) 00581 00582 00583 IF (int_buffer(1) > -1) THEN 00584 CALL FindPointGravityObject(int_buffer(1), particle%PointGravityObj) 00585 END IF 00586 00587 chandle%particle_offset = chandle%particle_offset + 1 00588 00589 END SUBROUTINE Particle_ReadObjectFromChombo 00590 00591 FUNCTION GetParticleBounds(Particle) 00592 TYPE(ParticleDef) :: Particle 00593 REAL(KIND=qPREC), DIMENSION(3,2) :: GetParticleBounds 00594 GetParticleBounds=GxBounds 00595 GetParticleBounds(1:nDim,1)=Particle%xloc(1:nDim)-Particle%radius*sink_dx 00596 GetParticleBounds(1:nDim,2)=Particle%xloc(1:nDim)+Particle%radius*sink_dx 00597 END FUNCTION GetParticleBounds 00598 00599 00600 FUNCTION GetParticlePos(Particle, t) 00601 TYPE(ParticleDef) :: Particle 00602 REAL(KIND=qPREC) :: t 00603 REAL(KIND=qPREC) :: GetParticlePos(3) 00604 GetParticlePos=Particle%PointGravityObj%x0+Particle%PointGravityObj%v0*(t-Particle%PointGravityObj%t0) 00605 00606 END FUNCTION GetParticlePos 00607 00608 END MODULE ParticleDeclarations 00609 00610