Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! source_control.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 00028 00031 00035 00038 MODULE SourceControl 00039 00040 USE DataDeclarations 00041 USE PhysicsDeclarations 00042 USE SourceDeclarations 00043 USE EOS 00044 USE CoolingSrc 00045 USE CylindricalSrc 00046 USE UniformGravitySrc 00047 USE PointGravitySrc 00048 USE SelfGravitySrc 00049 USE RotatingSrc 00050 ! USE OutflowSrc 00051 00052 IMPLICIT NONE 00053 PRIVATE 00054 PUBLIC Src, SrcInit, SrcFinalizeInit, SrcInitTracers, SrcInitElliptics, ReadSourceObjects, SourceCell !CreateSrc, CoolingDef, CreateCoolingObject 00055 SAVE 00056 TYPE(SourcesDef),POINTER :: sources 00057 TYPE(LevelDef) :: currlevel 00058 00059 ! sources verbosity shorthand 00060 INTEGER :: vb=0 00061 00063 INTERFACE CreateSrc 00064 MODULE PROCEDURE CreateCoolingObject, CreatePointGravityObject!, CreateCylindricalObject, CreateOutflowObject & 00065 !, CreateUniformGravityObject 00066 END INTERFACE 00067 00068 00069 CONTAINS 00070 00072 SUBROUTINE SrcInit(isrcsolvetype,iverbosity) 00073 INTEGER,OPTIONAL :: isrcsolvetype,iverbosity 00074 ! allocate sources & zero out components 00075 CALL CreateSourcesObj(sources,isrcsolvetype,iverbosity) 00076 CALL PointGravityInit() 00077 IF (lRestart) CALL ReadSourceObjects(restart_frame) 00078 00079 END SUBROUTINE SrcInit 00080 00086 SUBROUTINE SrcFinalizeInit 00087 00088 sources%lCooling=CoolingCheck() 00089 IF(sources%lCooling) CALL CoolingFinalizeInit 00090 00091 END SUBROUTINE SrcFinalizeInit 00092 00099 SUBROUTINE Src(Info, mb, tnow, hdt) 00100 TYPE(InfoDef) :: Info 00101 INTEGER,INTENT(IN) :: mb(3,2) 00102 REAL(KIND=qPrec),INTENT(IN) :: hdt, tnow 00103 ! 00104 REAL(KIND=xPrec) :: mx(3,2) 00105 LOGICAL :: lsrc 00106 IF(vb>0) PRINT*,'::Src begin' 00107 IF(vb>1) PRINT*,' sources%iSrcSolveType',sources%iSrcSolveType 00108 00109 IF(vb>0) THEN 00110 PRINT*,'hydro dt is',hdt 00111 PRINT*,'extents are ',mb 00112 END IF 00113 00114 currlevel=levels(Info%level) 00115 !mx(:,1)=()/currlevel%dx 00116 !mx(:,2)=()/currlevel%dx 00117 00118 CALL Protectq(Info,mb, 'before source') 00119 00120 SELECT CASE(sources%iSrcSolveType) 00121 CASE(ExplicitSource) 00122 CALL ExplicitSrc(Info, mb, tnow, hdt) 00123 CASE(ImplicitSource) 00124 CALL ImplicitSrc(Info, mb, tnow, hdt) 00125 CASE(ExactSource) 00126 CALL ExactSrc(Info, mb, tnow, hdt) 00127 CASE(NoSource) 00128 CASE DEFAULT 00129 PRINT*,'source_control.f90::Src error -- Unknown source solve type requested. Halting.' 00130 STOP 00131 END SELECT 00132 CALL Protectq(Info,mb, 'after source') 00133 IF(vb>0) PRINT*,'::Src complete' 00134 END SUBROUTINE Src 00135 00136 00137 ! ================================================================== 00138 ! = Explicit Scheme Section = 00139 ! ================================================================== 00140 00141 00150 SUBROUTINE ExplicitSrc(Info, mb, tnow, hdt) 00151 ! Interface declarations 00152 TYPE(InfoDef) :: Info 00153 INTEGER,INTENT(IN) :: mb(3,2) 00154 REAL(KIND=qPrec),INTENT(IN) :: hdt 00155 ! Internal declarations 00156 REAL(KIND=qPrec) :: pos(3) 00157 REAL(KIND=qPrec) :: dq(NrHydroVars),tnow,tnext, dv 00158 INTEGER :: lvl,i,j,k,ip(3) 00159 REAL(KIND=qPrec) :: error, volumefactor 00160 LOGICAL :: success, ghost, ghostz, ghostyz 00161 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: q 00162 dv = Levels(Info%level)%dx**ndim 00163 ALLOCATE(q(NrHydroVars)) 00164 DO k=mb(3,1),mb(3,2) 00165 Ghostz = .NOT. (k >= 1 .AND. k <= Info%mx(3)) 00166 DO j=mb(2,1),mb(2,2) 00167 Ghostyz = Ghostz .OR. .NOT. (j >= 1 .AND. j <= Info%mx(2)) 00168 DO i=mb(1,1),mb(1,2) 00169 Ghost = Ghostyz .OR. .NOT. (i >= 1 .AND. i <= Info%mx(1)) 00170 IF (Ghost) THEN 00171 VolumeFactor = 0d0 00172 ELSE 00173 IF (Info%Level==MaxLevel) THEN 00174 VolumeFactor = 1d0 00175 ELSE 00176 IF (Info%ChildMask(i,j,k) >= 1) THEN 00177 VolumeFactor = 0d0 00178 ELSE 00179 VolumeFactor = 1d0 00180 END IF 00181 END IF 00182 END IF 00183 00184 ! Method: 00185 ! 1) Try 2nd order Runge-Kutta to get estimate of error for timestep. 00186 ! Should speed things up where source terms are weak. 00187 ! 2) If error is above tolerance, reduce timestep and switch to 00188 ! 5th/4th order RK to do the integration. 00189 ! 3) Iterate at 5th/4th until full timestep reached. 00190 00191 q=Info%q(i,j,k,1:NrHydroVars) 00192 CALL Cons_to_prim(q) 00193 ip=(/i,j,k/) 00194 00195 00196 00197 CALL SourceCell(q, Info, ip, tnow, hdt, dv*VolumeFactor, PRIMITIVE)!where it gets complicated 00198 CALL Prim_to_cons(q) 00199 IF(ANY(q(:).ne.q(:) .OR. abs(q(:)).gt.huge(abs(q(:)))) .AND. .NOT. lRequestRestart) THEN 00200 PRINT*, 'Src routine produced a NAN' 00201 CALL PrintQ(Info, q, tnow, i, j, k) 00202 lRequestRestart=.true. 00203 ! STOP 00204 DEALLOCATE(q) 00205 RETURN 00206 END IF 00207 Info%q(i,j,k,1:NrHydroVars)=q 00208 00209 END DO 00210 END DO 00211 END DO 00212 DEALLOCATE(q) 00213 !IF(vb>0) THEN 00214 00215 !END IF 00216 END SUBROUTINE ExplicitSrc 00217 00218 SUBROUTINE SourceCell(q, Info, ip, tnow, hdt, dv, lForm) 00219 ! Interface declarations 00220 TYPE(InfoDef) :: Info 00221 REAL(KIND=qPrec),INTENT(IN) :: hdt 00222 ! Internal declarations 00223 REAL(KIND=qPrec) :: q(:) 00224 REAL(KIND=qPrec) :: pos(3) 00225 REAL(KIND=qPrec) :: dq(NrHydroVars),tnow,tnext, dv 00226 INTEGER :: ip(3) 00227 REAL(KIND=qPrec) :: error 00228 LOGICAL :: success 00229 LOGICAL :: lForm 00230 00231 tnext=tnow+hdt 00232 pos(1:nDim)=Info%xBounds(1:nDim,1)+(REAL(ip(1:nDim))-half)*levels(Info%level)%dx 00233 IF (nDim == 2) pos(3)=Info%xBounds(3,1) 00234 IF (iCylindrical == NOCYL) CALL PointGravity(q,hdt,pos,tnow,dv,info%Level,lform) 00235 ! CALL PointGravity(q,hdt,pos,tnow,dv,info%Level,lform) 00236 ! 2nd order RK, uses hdt to get full timestep error 00237 CALL RKOrder2(q,dq,pos,tnow,hdt,error,Info,ip,lform) 00238 IF(vb>1) PRINT*,'::ExplicitSrc:RKOrder2 -- dq,error',dq,error 00239 00240 IF(Error > sources%SrcTol) THEN 00241 IF(vb>1) PRINT*,'::ExplicitSrc:RKOrder45, qbefore, tnow, tnext',q,tnow,tnext 00242 ! iterate with 5th/4th order RK to get error & adjust timestep 00243 ! until full timestep is completed 00244 ! The subroutine updates q 00245 CALL RKOrder45(q,pos,hdt,tnow,tnext,Info,ip,success,lform) 00246 IF (.NOT. success .AND. .NOT. lRequestRestart) THEN 00247 write(*,*) 'RK Failed' 00248 CALL PrintQ(Info, q, tnow, ip(1), ip(2), ip(3)) 00249 lRequestRestart=.true. 00250 ! CALL OutputDoubleArray(Info%q(i-1:i+1,j-1:j+1,k-1,iPhi)) 00251 ! CALL OutputDoubleArray(Info%q(i-1:i+1,j-1:j+1,k,iPhi)) 00252 ! CALL OutputDoubleArray(Info%q(i-1:i+1,j-1:j+1,k+1,iPhi)) 00253 00254 ! CALL OutputDoubleArray(Info%q(i-1:i+1,j-1:j+1,k-1,iPhiDot)) 00255 ! CALL OutputDoubleArray(Info%q(i-1:i+1,j-1:j+1,k,iPhiDot)) 00256 ! CALL OutputDoubleArray(Info%q(i-1:i+1,j-1:j+1,k+1,iPhiDot)) 00257 ! STOP 00258 RETURN 00259 END IF 00260 IF(vb>1) PRINT*,'::ExplicitSrc:RKOrder45, qafter, tnow, tnext',q,tnow,tnext 00261 00262 ELSE 00263 ! good enough at 2nd order, update q 00264 q=q+dq 00265 END IF 00266 00267 END SUBROUTINE SourceCell 00268 00269 00270 00276 SUBROUTINE RKOrder2(q,dq,pos,t,dt,error,Info,ip,lform) 00277 ! Interface declarations 00278 TYPE(InfoDef) :: Info 00279 INTEGER :: ip(3) 00280 REAL(KIND=qPrec) :: q(:) 00281 REAL(KIND=qPrec) :: dq(NrHydroVars), pos(3), dt, error,t 00282 ! Internal declarations 00283 REAL(KIND=qPrec), ALLOCATABLE :: qstar(:) 00284 REAL(KIND=qPrec) :: k1(NrHydroVars),k2(NrHydroVars), qerr(NrHydroVars), minscales(NrHydroVars) 00285 LOGICAL :: lform 00286 dq=0d0 00287 00288 ALLOCATE(qStar(NrHydroVars)) 00289 ! Runge-Kutta 2nd order: for a function y(t), going from t=t0:t1 with dt=t1-t0, 00290 ! the second order integration is 00291 ! k1 = d(y0)/dt 00292 ! y* = y0 + dt k1 00293 ! k2 = d(y*)/dt 00294 ! y1 = y0 + dt/2 (k1+k2) 00295 ! where in our case y(t) -> q(t) and dq/dt is calculated from the source terms. 00296 ! The error is 00297 ! error = |y*/y1| 00298 00299 k1 = SrcDerivs(q,pos,t,Info,ip, lform) 00300 WHERE(k1/=k1)k1=0d0 00301 qstar = q + dt*k1 00302 k2 = SrcDerivs(qstar,pos,t+dt,Info,ip, lform) 00303 WHERE(k2/=k2)k2=0d0 00304 00305 dq=5d-1*dt*(k1+k2) 00306 00307 IF(ALL(ABS(dq(:))<TINY(dq(:)))) THEN 00308 dq=0d0 00309 error=0d0 00310 ELSE 00311 minscales(1)=minDensity 00312 minscales(m_low:m_high)=max(max(minDensity, q(1))*Iso_Speed,sqrt(SUM(q(m_low:m_high)**2))) 00313 IF (iE .ne. 0) minscales(iE)=max(minDensity, q(1))*Iso_Speed2 00314 IF (lMHD) minscales(iBx:iBz) = max(max(minDensity, q(1))*Iso_Speed2,sqrt(SUM(q(iBx:iBz)**2))) 00315 IF (nTracerLo /= 0) minscales(nTracerLo:nTracerHi) = minDensity 00316 qerr=q+dq-qstar 00317 error=MAXVAL(ABS(qErr(:))/max(abs(q(:)), minscales(:))) 00318 ! error = MAXVAL(ABS(qstar(:))/ABS(q(:)+dq(:))) !5d-1*dt*(k1+k2))) 00319 END IF 00320 00321 00322 00323 00324 00325 DEALLOCATE(qStar) 00326 END SUBROUTINE RKOrder2 00327 00328 00341 SUBROUTINE RKOrder45(q,pos,hdt,htnow,htnext,Info,ip,success,lform) 00342 ! Interface declarations 00343 TYPE(InfoDef) :: Info 00344 INTEGER :: ip(3) 00345 REAL(KIND=qPrec) :: q(:) 00346 REAL(KIND=qPrec) :: pos(3),hdt,htnow,htnext 00347 ! Internal declarations 00348 REAL(KIND=qPrec) :: dt,tnow,tnext 00349 INTEGER,PARAMETER :: MaxIters=1d4 00350 INTEGER :: niter 00351 REAL(KIND=qPrec),ALLOCATABLE :: qStar(:) 00352 ! CKRK stuff, see above reference for table 00353 REAL(KIND=qPrec) :: k1(NrHydroVars),k2(NrHydroVars),k3(NrHydroVars),k4(NrHydroVars),k5(NrHydroVars),k6(NrHydroVars), qTemp(NrHydroVars), qErr(NrHydroVars), maxerror, minscales(NrHydroVars) 00354 REAL(KIND=qPrec),PARAMETER :: A1=0.0, A2=0.2, A3=0.3, A4=0.6, A5=1.0, A6=0.875 00355 , B21=0.2 00356 , B31=3./40., B32=9./40. 00357 , B41=0.3, B42=-0.9, B43=1.2 00358 , B51=-11./54., B52=2.5, B53=-70./27., B54=35./27. 00359 , B61=1631./55296., B62=175./512., B63=575./13824., B64=44275./110592., B65=253./4096. 00360 , C1=37./378., C3=250./621., C4=125./594., C6=512./1771. 00361 , DC1=C1-2825./27648., DC3=C3-18575./48384., DC4=C4-13525./55296., DC5=-277./14336., DC6=C6-0.25 00362 REAL(KIND=qPrec),PARAMETER :: SAFETY=0.9, PSHRNK=-0.25, PGROW=-0.2 00363 LOGICAL :: success 00364 LOGICAL :: lform 00365 ! 00366 ALLOCATE(qStar(NrHydroVars)) 00367 dt=hdt 00368 tnow=htnow 00369 tnext=htnext 00370 00371 k1=0d0;k2=0d0;k3=0d0;k4=0d0;k5=0d0;k6=0d0 00372 00373 DO niter=1,MaxIters 00374 ! 00375 ! Cash-Karp RK functions, see above reference for details 00376 ! 00377 k1 = SrcDerivs(q,pos,tnow,Info,ip,lform) 00378 WHERE(k1/=k1)k1=0 00379 qStar = q + dt*( B21*k1 ) 00380 k2 = SrcDerivs(qStar,pos,tnow+dt*A2,Info,ip,lform) 00381 WHERE(k2/=k2)k2=0 00382 qStar = q + dt*( B31*k1 + B32*k2 ) 00383 k3 = SrcDerivs(qStar,pos,tnow+dt*(A3),Info,ip,lform) 00384 WHERE(k3/=k3)k3=0 00385 qStar = q + dt*( B41*k1 + B42*k2 + B43*k3 ) 00386 k4 = SrcDerivs(qStar,pos,tnow+dt*(A4),Info,ip,lform) 00387 WHERE(k4/=k4)k4=0 00388 qStar = q + dt*( B51*k1 + B52*k2 + B53*k3 + B54*k4 ) 00389 k5 = SrcDerivs(qStar,pos,tnow+dt*(A5),Info,ip,lform) 00390 WHERE(k5/=k5)k5=0 00391 qStar = q + dt*( B61*k1 + B62*k2 + B63*k3 + B64*k4 + B65*k5 ) 00392 k6 = SrcDerivs(qStar,pos,tnow+dt*(A6),Info,ip,lform) 00393 WHERE(k6/=k6)k6=0 00394 qTemp = q + dt*( C1*k1 + C3*k3 + C4*k4 + C6*k6 ) 00395 qErr = dt*( DC1*k1 + DC3*k3 + DC4*k4 + DC5*k5 + DC6*k6 ) 00396 00397 ! 00398 ! Error evaluation & timestep adjustment 00399 ! 00400 ! IF(vb>1) WRITE(*,'(A,10E20.12)')'qtemp ',qtemp(1) 00401 ! IF(vb>1) WRITE(*,'(A,10E20.12)')'qerr ',qerr(1) 00402 IF(vb>1) PRINT*,'qtemp ',qtemp 00403 IF(vb>1) PRINT*,'qerr ',qerr 00404 minscales(1)=minDensity 00405 minscales(m_low:m_high)=max(minDensity, q(1))*Iso_Speed 00406 IF (iE .ne. 0) minscales(iE)=max(minDensity, q(1))*Iso_Speed2 00407 IF (lMHD) minscales(iBx:iBz) = max(minDensity, q(1))*Iso_Speed2 00408 IF (nTracerLo /= 0) minscales(nTracerLo:nTracerHi) = minDensity 00409 maxerror=MAXVAL(ABS(qErr(:))/max(abs(q(:)), minscales(:))) 00410 IF (vb>1) print*,maxerror 00411 IF(maxerror > sources%SrcTol) THEN 00412 ! reduce timestep and retry 00413 dt = MAX( SAFETY*dt*(maxerror/sources%SrcTol)**PSHRNK, dt/MaxIters) 00414 IF(dt <= TINY(dt) .AND. .NOT. lRequestRestart) THEN 00415 PRINT*,'Src Error::RKOrder45 -- timestep insignificant. Halting.' 00416 ! PRINT*,' Iterations/MaxIters: ',niter, MaxIters 00417 ! PRINT*,' Last timestep, full timestep: ',dt, hdt 00418 ! PRINT*,' Current time, hydro begin/end time: ',tnow,htnow,htnext 00419 ! PRINT*,' Time left to go: ',htnow-tnow 00420 CALL PrintQ(Info, q, htnow, ip(1),ip(2),ip(3)) 00421 lRequestRestart=.true. 00422 success=.false. 00423 EXIT 00424 ! STOP 00425 ! RETURN 00426 00427 END IF 00428 ! tnow remains the same; update tnext 00429 tnext=tnow+dt 00430 ELSE 00431 q = qTemp 00432 tnow=tnow+dt 00433 dt = MAX( MIN( SAFETY*dt*(maxerror/sources%SrcTol)**PGROW, 5*dt, htnext-tnow), 0d0) 00434 tnext=tnow+dt 00435 IF(dt <= TINY(dt)) THEN 00436 ! end of full timestep reached 00437 success=.true. 00438 EXIT 00439 END IF 00440 END IF 00441 END DO 00442 00443 IF(niter>MaxIters .AND. .NOT. lRequestRestart) THEN 00444 PRINT*,'Src Error::RKOrder45 -- maximum number of iterations reached.' 00445 ! PRINT*,' Max iterations: ',MaxIters 00446 ! PRINT*,' Last timestep, full timestep: ',dt, hdt 00447 ! PRINT*,' Current time, hydro begin/end time: ',tnow,htnow,htnext 00448 ! PRINT*,' Time left to go: ',htnow-tnow 00449 CALL PrintQ(Info, q, htnow, ip(1), ip(2), ip(3)) 00450 ! STOP 00451 lRequestRestart=.true. 00452 success=.false. 00453 END IF 00454 DEALLOCATE(qStar) 00455 END SUBROUTINE RKOrder45 00456 00457 00458 00459 ! ================================================================== 00460 ! = Implicit Scheme Section = 00461 ! ================================================================== 00462 00470 SUBROUTINE ImplicitSrc(Info, mb, tnow,hdt) 00471 TYPE(InfoDef) :: Info 00472 INTEGER,INTENT(IN) :: mb(3,2) 00473 REAL(KIND=qPrec),INTENT(IN) :: hdt,tnow 00474 IF(vb>0) PRINT*,'::ImplicitSrc begin' 00475 IF(vb>0) PRINT*,'::ImplicitSrc complete' 00476 END SUBROUTINE ImplicitSrc 00477 00478 00479 ! ================================================================== 00480 ! = Exact Scheme Section = 00481 ! ================================================================== 00482 00489 SUBROUTINE ExactSrc(Info, mb, tnow,hdt) 00490 TYPE(InfoDef) :: Info 00491 INTEGER,INTENT(IN) :: mb(3,2) 00492 REAL(KIND=qPrec),INTENT(IN) :: tnow,hdt 00493 IF(vb>0) PRINT*,'::ExactSrc begin' 00494 IF(vb>0) PRINT*,'::ExactSrc complete' 00495 END SUBROUTINE ExactSrc 00496 00497 00498 00499 ! ================================================================== 00500 ! = ~~~ Wrapper Subroutines/Functions ~~~ = 00501 ! ================================================================== 00502 00506 FUNCTION SrcDerivs(q,pos,t,Info,ip,lform) 00507 ! Interface declarations 00508 TYPE(InfoDef) :: Info 00509 INTEGER :: ip(3) 00510 REAL(KIND=qPrec) :: q(:) 00511 REAL(KIND=qPrec) :: pos(3), SrcDerivs(NrHydroVars), t 00512 ! Internal declarations 00513 REAL(KIND=qPrec) :: dqdt(NrHydroVars) 00514 LOGICAL :: lform 00515 dqdt=0d0 00516 00517 !IF(sources%lCooling) 00518 CALL Cooling(q,dqdt,pos,currlevel%dx,lform) 00519 00520 IF (iCylindrical.ge.NoCyl .and. iCylindrical.le.WithAngMom) THEN 00521 IF (iCylindrical.ne.NoCyl) CALL Cylindrical(q,dqdt,pos,lform)!,currlevel%dx) 00522 ELSE 00523 print*,'';print*,'iCylindrical= ',iCylindrical 00524 print*,'but should be either 0, 1 or 2. Aborting.' 00525 stop 00526 END IF 00527 00528 IF(luniformgravity) CALL UniformGravity_src(q,dqdt,t,lform) 00529 IF (OmegaRot /= 0d0) CALL Rotating(q, dqdt, pos, lform) 00530 IF (iCylindrical /= NOCYL) CALL PointGravity_inst(q,dqdt,pos,t,lform) !moved to source cell for improved mom. cons. 00531 ! CALL PointGravity_inst(q,dqdt,pos,t,lform) !moved to source cell for improved mom. cons. 00532 00533 ! IF(sources%lOutflows) CALL Outflows(q,dqdt,pos) 00534 00535 ! IF (lSelfGravity .AND. SourceMethod == OPERATORSPLIT) CALL SelfGravity(q,dqdt,pos,currlevel%dx,t,Info,ip,lform) 00536 !Self Gravity now always implemented in hyperbolic solve directly 00537 SrcDerivs=dqdt 00538 END FUNCTION SrcDerivs 00539 00540 00544 SUBROUTINE SrcInitTracers 00545 ! Interface declarations 00546 !IF(sources%lCooling) 00547 CALL InitCoolingTracers(sources) 00548 !IF(sources%lCylindrical) CALL InitCylindricalTracers(NrVars,NrTracerVars,sources) 00549 ! IF(sources%lUniformGravity) CALL InitUniformGravTracers(NrVars,NrTracerVars,sources) 00550 ! CALL InitPointGravTracers(NrVars,NrTracerVars,sources) 00551 ! IF(sources%lOutflows) CALL InitOutflowTracers(NrVars,NrTracerVars,sources) 00552 END SUBROUTINE SrcInitTracers 00553 00557 SUBROUTINE SrcInitElliptics() 00558 ! Interface declarations 00559 !IF(sources%lCooling) 00560 CALL InitCoolingElliptics(sources) 00561 !IF(sources%lCylindrical) CALL InitCylindricalElliptics(NrVars,NrEllipticVars,sources) 00562 ! IF(sources%lUniformGravity) CALL InitUniformGravElliptics(NrVars,NrEllipticVars,sources) 00563 ! IF(sources%lPointGravity) CALL InitPointGravElliptics(NrVars,NrEllipticVars,sources) 00564 ! IF(sources%lOutflows) CALL InitOutflowElliptics(NrVars,NrEllipticVars,sources) 00565 END SUBROUTINE SrcInitElliptics 00566 00567 00570 SUBROUTINE ReadSourceObjects(nframe) 00571 00572 USE ChomboDeclarations, ONLY: ChomboHandle, CreateChomboHandle, CloseChomboHandle, CHOMBO_HANDLE_READ, & 00573 Chombo_OpenSourceGroup, Chombo_CloseSourceGroup 00574 00575 USE CoolingSrc, ONLY: CoolingDef, CreateCoolingObject, Cooling_ReadObjectFromChombo 00576 00577 USE PointGravitySrc, ONLY: PointGravityDef, CreatePointGravityObject, PointGravity_ReadObjectFromChombo 00578 00579 INTEGER :: nframe 00580 00581 ! Variable declarations 00582 TYPE(ChomboHandle), POINTER :: chandle 00583 TYPE(CoolingDef), POINTER :: coolingobj 00584 TYPE(PointGravityDef), POINTER :: gravity_obj 00585 INTEGER :: i_err 00586 INTEGER :: nr_objects 00587 CHARACTER(LEN=23) :: s_filename 00588 00589 00590 ! Open a reading handle for the specified frame. 00591 WRITE(s_filename, '(A10,I5.5,A4)') 'out/chombo', nframe, '.hdf' 00592 CALL CreateChomboHandle(s_filename, chandle, CHOMBO_HANDLE_READ) 00593 00594 !!! Check For Cooling Terms !!! 00595 00596 ! Open the 'cooling_objects' group and save the handle. 00597 nr_objects = Chombo_OpenSourceGroup(chandle, "cooling_objects") 00598 00599 chandle%source_offset = 0 00600 00601 ! Read source objects from the data file until the expected number of objects is read. 00602 DO WHILE (chandle%source_offset < nr_objects) 00603 00604 ! Create a new cooling object. 00605 NULLIFY(coolingobj) 00606 CALL CreateCoolingObject(coolingobj) 00607 00608 ! Read in the cooling data from the Chombo file. This subroutine also advances the 00609 ! chandle%source_offset variable. 00610 CALL Cooling_ReadObjectFromChombo(chandle, coolingobj) 00611 00612 END DO 00613 00614 ! Close the source term group. This also clears the source term offset. 00615 CALL Chombo_CloseSourceGroup(chandle) 00616 00617 !!! End Check for Cooling Terms !!! 00618 00619 !!! Begin Check For Point Gravity Terms !!! 00620 ! Open the 'point_gravity_objects' group and save the handle. 00621 nr_objects = Chombo_OpenSourceGroup(chandle, "point_gravity_objects") 00622 00623 chandle%source_offset = 0 00624 00625 ! Read source objects from the data file until the expected number of objects is read. 00626 DO WHILE (chandle%source_offset < nr_objects) 00627 00628 ! Create a new gravity object. 00629 NULLIFY(gravity_obj) 00630 CALL CreatePointGravityObject(gravity_obj) 00631 00632 ! Read in the point gravity data from the Chombo file. This subroutine also advances the 00633 ! chandle%source_offset variable. 00634 CALL PointGravity_ReadObjectFromChombo(chandle, gravity_obj) 00635 00636 END DO 00637 00638 ! Close the source term group. This also clears the source term offset. 00639 CALL Chombo_CloseSourceGroup(chandle) 00640 00641 00642 00643 00644 CALL CloseChomboHandle(chandle) 00645 00646 END SUBROUTINE ReadSourceObjects 00647 00648 ! ========================================== 00649 ! = Sources creation/destruction = 00650 ! = and list manipulation section = 00651 ! ========================================== 00652 00657 SUBROUTINE CreateSourcesObj(sources,isrcsolvetype,iverbosity) 00658 ! Interface declarations 00659 TYPE(SourcesDef),POINTER :: sources 00660 INTEGER,OPTIONAL :: isrcsolvetype,iverbosity 00661 00662 IF(ASSOCIATED(sources)) THEN 00663 PRINT*,'source_control.f90::CreateSourcesObj error -- sources already associated. Halting.' 00664 STOP 00665 END IF 00666 00667 IF(PRESENT(iverbosity)) THEN 00668 vb=iverbosity 00669 ELSE 00670 vb=0 00671 END IF 00672 00673 ALLOCATE(sources) 00674 ! Initialize everything 00675 sources%iSrcSolveType=explicitSource 00676 sources%level=0 00677 sources%lPrimitive=.FALSE. 00678 sources%SrcTol=SrcPrecision ! tolerance for source terms (may vary with AMR level) 00679 !sources%lCooling=.FALSE. 00680 !sources%lCylindrical=.FALSE. 00681 !sources%lUniformGravity=.FALSE. 00682 !sources%lPointGravity=.FALSE. 00683 !sources%lOutflows=.FALSE. 00684 ! 00685 IF(PRESENT(isrcsolvetype)) sources%iSrcSolveType=isrcsolvetype 00686 sources%iverbosity=vb 00687 00688 IF(vb>0) PRINT*,'::CreateSourcesObj successful' 00689 END SUBROUTINE CreateSourcesObj 00690 00692 SUBROUTINE DestroySourcesObj 00693 END SUBROUTINE DestroySourcesObj 00694 00695 00696 00697 END MODULE SourceControl 00698