Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! EOS.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 EOS 00033 USE GlobalDeclarations 00034 USE PhysicsDeclarations 00035 USE TreeDeclarations 00036 USE DataDeclarations 00037 00038 IMPLICIT NONE 00039 PUBLIC Press, PrimPress, Prim_To_Cons, Cons_To_Prim, ConvertTotalToInternalEnergy, ConvertInternalToTotalEnergy, Protectq, InternalEnergy, PrintQ, soundspeed, primsoundspeed, SetPress, GetInternalEnergy, GetMaxSpeed, calc_eigens 00040 00041 INTERFACE ConvertTotalToInternalEnergy 00042 MODULE PROCEDURE ConvertTotalToInternalEnergy0D, ConvertTotalToInternalEnergy1D, ConvertTotalToInternalEnergy2D, ConvertTotalToInternalEnergy3D 00043 END INTERFACE ConvertTotalToInternalEnergy 00044 00045 INTERFACE ConvertInternalToTotalEnergy 00046 MODULE PROCEDURE ConvertInternalToTotalEnergy0D, ConvertInternalToTotalEnergy1D, ConvertInternalToTotalEnergy2D, ConvertInternalToTotalEnergy3D 00047 END INTERFACE ConvertInternalToTotalEnergy 00048 00049 INTERFACE Cons_to_prim 00050 MODULE PROCEDURE Cons_to_prim_info, cons_to_prim_q 00051 END INTERFACE Cons_to_prim 00052 00053 INTERFACE prim_to_cons 00054 MODULE PROCEDURE prim_to_cons_info, prim_to_cons_q 00055 END INTERFACE prim_to_cons 00056 00057 00058 CONTAINS 00059 00060 00063 SUBROUTINE Protectq(Info,ip,caller, lStopProtect) 00064 TYPE(InfoDef) :: Info 00065 INTEGER :: i,j,k,ii,jj,kk 00066 INTEGER, DIMENSION(3,2) :: ip, mB 00067 INTEGER :: level, rmbc,l, n, m 00068 REAL(KIND=qPREC), PARAMETER :: Precision=1e-14 00069 REAL(KIND=qPREC) :: vel(3),temp 00070 CHARACTER(LEN=*) :: caller 00071 REAL(KIND=qPREC), DIMENSION(:), POINTER :: q 00072 REAL(KIND=qPREC), DIMENSION(:,:,:), POINTER :: pressure, rho 00073 LOGICAL, DIMENSIOn(:,:,:), ALLOCATABLE :: mask 00074 LOGICAL, OPTIONAL :: lStopProtect 00075 00076 level=Info%level 00077 DO i=ip(1,1),ip(1,2) 00078 DO j=ip(2,1),ip(2,2) 00079 DO k=ip(3,1),ip(3,2) 00080 q=>Info%q(i,j,k,:) 00081 IF(ANY(q(1:NrVars).ne.q(1:NrVars) .OR. abs(q(1:NrVars)).gt.huge(abs(q(1:NrVars)))) .AND. .NOT. lRequestRestart) THEN 00082 write(*,*) 'protection found Nans in ', caller 00083 CALL Printq(Info, q, levels(Info%level)%tnow, i,j,k) 00084 PRINT*, 'Processor', MPI_ID, 'requesting restart' 00085 lRequestRestart=.true. 00086 IF (PRESENT(lStopProtect)) THEN 00087 IF (lStopProtect) THEN 00088 WRITE(*,*) "*** Simulation stopped due to protections in the initialization! ***" 00089 STOP 00090 END IF 00091 END IF 00092 RETURN 00093 END IF 00094 END DO 00095 END DO 00096 END DO 00097 00098 ALLOCATE(rho(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2))) 00099 rho=Info%q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),1) 00100 IF (iE /= 0) THEN 00101 ALLOCATE(pressure(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2))) 00102 DO i=ip(1,1),ip(1,2) 00103 DO j=ip(2,1),ip(2,2) 00104 DO k=ip(3,1),ip(3,2) 00105 pressure(i,j,k)=Press(Info%q(i,j,k,:)) 00106 END DO 00107 END DO 00108 END DO 00109 END IF 00110 00111 DO i=ip(1,1),ip(1,2) 00112 DO j=ip(2,1),ip(2,2) 00113 DO k=ip(3,1),ip(3,2) 00114 q=>Info%q(i,j,k,:) 00115 IF (q(1) < MinDensity) THEN 00116 IF (lRestartOnDensityProtections) THEN 00117 IF (.NOT. lRequestRestart) THEN 00118 WRITE(*,'(A,I5,3A)') 'Processor', MPI_ID, ' encountered a density protection from ', caller, ' and lRestartOnDensityProtections is true so requesting restart. Consider setting lRestartOnDensityProtections=.false. in physics.data if you wish to not restart' 00119 CALL Printq(Info, q, levels(Info%level)%tnow, i,j,k) 00120 lRequestRestart=.true. 00121 END IF 00122 ELSE 00123 !First protect momentum 00124 SELECT CASE(iMomentumProtect) 00125 CASE(CONSERVE_MOMENTUM) 00126 CASE(CLEAR_MOMENTUM) 00127 q(m_low:m_high)=0 00128 CASE(AVG_NEARBY_VEL) 00129 DO l=1, 5 00130 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00131 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00132 ALLOCATE(mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00133 mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity 00134 n=count(mask) 00135 IF (n > 0) THEN 00136 DO m=m_low,m_high 00137 q(m)=sum(Info%q(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2),m)/rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)), mask)/REAL(n,8) 00138 END DO 00139 END IF 00140 DEALLOCATE(mask) 00141 IF (n > 0) EXIT 00142 END DO 00143 IF (n == 0) q(m_low:m_high)=0d0 00144 END SELECT 00145 !Then protect density 00146 SELECT CASE(iDensityProtect) 00147 CASE(MIN_DENSITY) 00148 q(1)=MinDensity 00149 CASE(MIN_NEARBY_DENSITY) 00150 DO l=1, 5 00151 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00152 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00153 ALLOCATE(mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00154 mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity 00155 n=count(mask) 00156 IF (n > 0) q(1)=minval(rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)), mask) 00157 DEALLOCATE(mask) 00158 IF (n > 0) EXIT 00159 END DO 00160 IF (n == 0) q(1)=MinDensity 00161 CASE(AVG_NEARBY_DENSITY) 00162 DO l=1, 5 00163 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00164 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00165 ALLOCATE(mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00166 mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity 00167 n=count(mask) 00168 IF (n > 0) q(1)=sum(rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)), mask)/real(n,8) 00169 DEALLOCATE(mask) 00170 IF (n > 0) EXIT 00171 END DO 00172 IF (n == 0) q(1) = MinDensity 00173 END SELECT 00174 IF (iMomentumProtect==AVG_NEARBY_VEL) q(m_low:m_high)=q(m_low:m_high)*q(1) 00175 END IF 00176 END IF 00177 DO l=nTracerLo, nTracerHI 00178 q(l)=max(q(l),0d0) 00179 END DO 00180 IF (iE /= 0) THEN 00181 temp=pressure(i,j,k)/rho(i,j,k) 00182 IF (temp < MinTemp) THEN 00183 IF (lRestartOnPressureProtections) THEN 00184 IF (.NOT. lRequestRestart) THEN 00185 PRINT*, 'Processor', MPI_ID, ' encountered a pressure protection from ', caller, ' and lRestartOnPressureProtections is true so requesting restart. Consider setting lRestartOnPressureProtections=.false. in physics.data if you wish to not restart', temp, mintemp 00186 CALL Printq(Info, q, levels(Info%level)%tnow, i,j,k) 00187 lRequestRestart=.true. 00188 END IF 00189 ELSE 00190 SELECT CASE (iPressureProtect) 00191 CASE(MIN_TEMP) 00192 CALL SetTemp(Info%q(i,j,k,:), MinTemp) 00193 CASE(MIN_NEARBY_PRESS) 00194 DO l=1, 5 00195 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00196 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00197 ALLOCATE(Mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00198 Mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity .AND. pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))*MinTemp 00199 n=count(Mask) 00200 IF (n > 0) CALL SetPress(q(:), minval(pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)),Mask)) 00201 DEALLOCATE(Mask) 00202 IF (n > 0) EXIT 00203 END DO 00204 IF (n == 0) CALL SetPress(q(:), q(1)*MinTemp) 00205 CASE(AVG_NEARBY_PRESS) 00206 DO l=1, 5 00207 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00208 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00209 ALLOCATE(Mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00210 Mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity .AND. pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))*MinTemp 00211 n=count(Mask) 00212 IF (n > 0) CALL SetPress(Info%q(i,j,k,:), sum(pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)),Mask)/REAL(n,8)) 00213 DEALLOCATE(Mask) 00214 IF (n > 0) EXIT 00215 END DO 00216 IF (n == 0) CALL SetPress(q(:), q(1)*MinTemp) 00217 CASE(MIN_NEARBY_TEMP) 00218 DO l=1, 5 00219 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00220 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00221 ALLOCATE(Mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00222 Mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity .AND. pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))*MinTemp 00223 n=count(Mask) 00224 IF (n > 0) CALL SetTemp(Info%q(i,j,k,:), minval(pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))/rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)),Mask)) 00225 DEALLOCATE(Mask) 00226 IF (n > 0) EXIT 00227 END DO 00228 IF (n == 0) CALL SetPress(q(:), q(1)*MinTemp) 00229 CASE(AVG_NEARBY_TEMP) 00230 DO l=1, 5 00231 mB(:,1)=max((/i,j,k/)-l, ip(:,1)) 00232 mB(:,2)=min((/i,j,k/)+l, ip(:,2)) 00233 ALLOCATE(Mask(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))) 00234 Mask=rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > MinDensity .AND. pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)) > rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))*MinTemp 00235 n=count(Mask) 00236 IF (n > 0) CALL SetTemp(Info%q(i,j,k,:), sum(pressure(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2))/rho(mB(1,1):mB(1,2),mB(2,1):mB(2,2),mB(3,1):mB(3,2)),Mask)/REAL(n,8)) 00237 DEALLOCATE(Mask) 00238 IF (n > 0) EXIT 00239 END DO 00240 IF (n == 0) CALL SetPress(q(:), q(1)*MinTemp) 00241 END SELECT 00242 END IF 00243 END IF 00244 END IF 00245 END DO 00246 END DO 00247 END DO 00248 DEALLOCATE(rho) 00249 IF (iE /= 0) DEALLOCATE(pressure) 00250 IF (lRequestRestart .AND. PRESENT(lStopProtect)) THEN 00251 IF (lStopProtect) THEN 00252 WRITE(*,*) "*** Simulation stopped due to protections in the initialization! ***" 00253 STOP 00254 END IF 00255 END IF 00256 END SUBROUTINE Protectq 00257 00260 00263 SUBROUTINE ConvertTotalToInternalEnergy0D(q) 00264 REAL(KIND=qPREC), DIMENSION(:) :: q 00265 IF (iE .ne. 0) THEN 00266 q(iE)=q(iE)-half*sum(q(m_low:m_high)**2)/q(1) 00267 IF (lMHD) q(iE)=q(iE)-half*sum(q(iBx:iBz)**2) 00268 END IF 00269 END SUBROUTINE ConvertTotalToInternalEnergy0D 00270 00273 SUBROUTINE ConvertTotalToInternalEnergy1D(q) 00274 REAL(KIND=qPREC), DIMENSION(:,:) :: q 00275 IF (iE .ne. 0) THEN 00276 q(:,iE)=q(:,iE)-half*sum(q(:,m_low:m_high)**2, 2)/q(:,1) 00277 IF (lMHD) q(:,iE)=q(:,iE)-half*sum(q(:,iBx:iBz)**2,2) 00278 END IF 00279 END SUBROUTINE ConvertTotalToInternalEnergy1D 00280 00283 SUBROUTINE ConvertTotalToInternalEnergy2D(q) 00284 REAL(KIND=qPREC), DIMENSION(:,:,:) :: q 00285 IF (iE .ne. 0) THEN 00286 q(:,:,iE)=q(:,:,iE)-half*sum(q(:,:,m_low:m_high)**2, 3)/q(:,:,1) 00287 IF (lMHD) q(:,:,iE)=q(:,:,iE)-half*sum(q(:,:,iBx:iBz)**2,3) 00288 END IF 00289 END SUBROUTINE ConvertTotalToInternalEnergy2D 00290 00293 SUBROUTINE ConvertTotalToInternalEnergy3D(q) 00294 REAL(KIND=qPREC), DIMENSION(:,:,:,:) :: q 00295 IF (iE .ne. 0) THEN 00296 q(:,:,:,iE)=q(:,:,:,iE)-half*sum(q(:,:,:,m_low:m_high)**2, 4)/q(:,:,:,1) 00297 IF (lMHD) q(:,:,:,iE)=q(:,:,:,iE)-half*sum(q(:,:,:,iBx:iBz)**2,4) 00298 END IF 00299 END SUBROUTINE ConvertTotalToInternalEnergy3D 00300 00302 00303 00306 00309 SUBROUTINE ConvertInternalToTotalEnergy0D(q) 00310 REAL(KIND=qPREC), DIMENSION(:) :: q 00311 IF (iE .ne. 0) THEN 00312 q(iE)=q(iE)+half*sum(q(m_low:m_high)**2)/q(1) 00313 IF (lMHD) q(iE)=q(iE)+half*sum(q(iBx:iBz)**2) 00314 END IF 00315 END SUBROUTINE ConvertInternalToTotalEnergy0D 00316 00319 SUBROUTINE ConvertInternalToTotalEnergy1D(q) 00320 REAL(KIND=qPREC), DIMENSION(:,:) :: q 00321 IF (iE .ne. 0) THEN 00322 q(:,iE)=q(:,iE)+half*sum(q(:,m_low:m_high)**2, 2)/q(:,1) 00323 IF (lMHD) q(:,iE)=q(:,iE)+half*sum(q(:,iBx:iBz)**2,2) 00324 END IF 00325 END SUBROUTINE ConvertInternalToTotalEnergy1D 00326 00329 SUBROUTINE ConvertInternalToTotalEnergy2D(q) 00330 REAL(KIND=qPREC), DIMENSION(:,:,:) :: q 00331 IF (iE .ne. 0) THEN 00332 q(:,:,iE)=q(:,:,iE)+half*sum(q(:,:,m_low:m_high)**2, 3)/q(:,:,1) 00333 IF (lMHD) q(:,:,iE)=q(:,:,iE)+half*sum(q(:,:,iBx:iBz)**2,3) 00334 END IF 00335 END SUBROUTINE ConvertInternalToTotalEnergy2D 00336 00339 SUBROUTINE ConvertInternalToTotalEnergy3D(q) 00340 REAL(KIND=qPREC), DIMENSION(:,:,:,:) :: q 00341 IF (iE .ne. 0) THEN 00342 q(:,:,:,iE)=q(:,:,:,iE)+half*sum(q(:,:,:,m_low:m_high)**2, 4)/q(:,:,:,1) 00343 IF (lMHD) q(:,:,:,iE)=q(:,:,:,iE)+half*sum(q(:,:,:,iBx:iBz)**2,4) 00344 END IF 00345 END SUBROUTINE ConvertInternalToTotalEnergy3D 00346 00348 00349 00355 FUNCTION Press(q,kappa,chi,DPDQ) 00356 ! USE COOL 00357 ! USE TF 00358 ! USE VANDERWAALS 00359 REAL (KIND=qPrec), DIMENSION(:), INTENT(IN) :: q 00360 REAL (KIND=qPrec), OPTIONAL, INTENT(OUT) :: kappa,chi 00361 REAL (KIND=qPrec), DIMENSION(:), OPTIONAL, INTENT(OUT) :: DPDQ 00362 REAL (KIND=qPrec) :: Press 00363 ! Internal variables 00364 REAL (KIND=qPrec) :: kappa1,chi1,ke 00365 ! Caching certain frequently-used derived values. 00366 SELECT CASE(iEOS) 00367 CASE(EOS_ISOTHERMAL) 00368 kappa1 = 0.0d0 00369 chi1 = Iso_Speed2 00370 Press = Iso_Speed2*q(1) 00371 IF (present(dPdQ)) THEN !This is dPdQ keeping everything other conservative quantity constant 00372 !P=rho*IsoSpeed2 = chi1*rho 00373 dPdQ=0 00374 dPdQ(1)=chi1 00375 END IF 00376 CASE(EOS_IDEAL_GAS) 00377 ke=half*SUM(q(m_low:m_high)**2)/q(1) 00378 Press=q(iE)-ke 00379 IF (lMHD) Press=Press-half*SUM(q(iBx:iBz)**2) 00380 Press=Press*gamma1 00381 kappa1=gamma1 !d(P)/d(rho) 00382 chi1=0.d0 00383 IF(PRESENT(DPDQ)) THEN !This is dPdQ keeping everything other conservative quantity constant 00384 !P=kappa1*(E-half*q(m_low:m_high)**2/q(1)-half*q(iBx:iBz)**2)+chi1*rho 00385 DPDQ = zero 00386 DPDQ(1)=chi1+kappa1*ke/q(1) 00387 DPDQ(m_low:m_high)=-kappa1*q(m_low:m_high)/q(1) 00388 DPDQ(iE) = kappa1 00389 IF (lMHD) THEN 00390 DPDQ(iBx:iBz)= -kappa1*q(iBx:iBz) 00391 END IF 00392 END IF 00393 CASE DEFAULT 00394 PRINT*,'unimplemented EOS, iEOS=',iEOS 00395 STOP 00396 END SELECT 00397 ! if present, return chi and kappa to caller 00398 IF(Present(chi)) chi=chi1 00399 IF(Present(kappa)) kappa=kappa1 00400 END FUNCTION Press 00401 00402 00403 FUNCTION Temperature(q) 00404 REAL(KIND=qPREC) :: Temperature 00405 REAL(KIND=qPREC), DIMENSION(:) :: q 00406 Temperature=Press(q)/q(1) 00407 END FUNCTION Temperature 00408 00410 FUNCTION GetInternalEnergy(q,p) 00411 REAL(KIND=qPREC) :: GetInternalEnergy 00412 REAL(KIND=qPREC) :: p 00413 REAL(KIND=qPREC), DIMENSION(:) :: q 00414 SELECT CASE(iEOS) 00415 CASE(EOS_ISOTHERMAL) 00416 GetInternalEnergy=p 00417 CASE(EOS_IDEAL_GAS) 00418 GetInternalEnergy=p*gamma7 00419 CASE DEFAULT 00420 PRINT*,'unimplemented EOS, iEOS=',iEOS 00421 STOP 00422 END SELECT 00423 END FUNCTION GetInternalEnergy 00424 00426 SUBROUTINE SetPress(q, p) 00427 REAL(KIND=qPREC), DIMENSION(:) :: q 00428 REAL(KIND=qPREC) :: p, E 00429 E=half*SUM(q(m_low:m_high)**2)/q(1) 00430 IF (lMHD) E=E+half*SUM(q(iBx:iBz)**2) 00431 q(iE)=GetInternalEnergy(q,p)+E 00432 END SUBROUTINE SetPress 00433 00434 00436 SUBROUTINE SetTemp(q, t) 00437 REAL(KIND=qPREC), DIMENSION(:) :: q 00438 REAL(KIND=qPREC) :: t, E 00439 E=half*SUM(q(m_low:m_high)**2)/q(1) 00440 IF (lMHD) E=E+half*SUM(q(iBx:iBz)**2) 00441 q(iE)=GetInternalEnergy(q,q(1)*t)+E 00442 END SUBROUTINE SetTemp 00443 00446 FUNCTION SoundSpeed(q) 00447 REAL (KIND=qPrec), DIMENSION(:), INTENT(IN) :: q 00448 REAL (KIND=qPrec) :: SoundSpeed, Eth 00449 IF (iEOS == EOS_ISOTHERMAL) THEN 00450 SoundSpeed = sqrt(Iso_Speed2) 00451 ELSE 00452 Eth=q(iE)- half*DOT_PRODUCT(q(m_low:m_high),q(m_low:m_high))/q(1) 00453 IF(lMHD) Eth=Eth - half*DOT_PRODUCT(q(iBx:iBz),q(iBx:iBz)) 00454 SELECT CASE(iEOS) 00455 CASE(EOS_IDEAL_GAS) ! Ideal 00456 SoundSpeed = sqrt(gamma*gamma1*Eth/q(1)) 00457 CASE DEFAULT 00458 PRINT*,'unimplemented EOS, iEOS=',iEOS 00459 STOP 00460 END SELECT 00461 END IF 00462 END FUNCTION SoundSpeed 00463 00464 00467 FUNCTION PrimSoundSpeed(q) 00468 REAL (KIND=qPrec), DIMENSION(:), INTENT(IN) :: q 00469 REAL (KIND=qPrec) :: PrimSoundSpeed 00470 IF (iEOS == EOS_ISOTHERMAL) THEN 00471 PrimSoundSpeed = sqrt(Iso_Speed2) 00472 ELSE 00473 SELECT CASE(iEOS) 00474 CASE(EOS_IDEAL_GAS) ! Ideal 00475 PrimSoundSpeed = sqrt(gamma*q(iE)/q(1)) 00476 CASE DEFAULT 00477 PRINT*,'unimplemented EOS, iEOS=',iEOS 00478 STOP 00479 END SELECT 00480 END IF 00481 END FUNCTION PrimSoundSpeed 00482 00483 00486 PURE FUNCTION InternalEnergy(q) 00487 REAL (KIND=qPrec), DIMENSION(:), INTENT(IN) :: q 00488 REAL (KIND=qPrec) :: InternalEnergy 00489 ! Internal variables 00490 REAL (KIND=qPrec) :: gam,kappa1,invkappa1,chi1,eth,ke,rho,invq1,Bsq,rho_cgs,eth_cgs 00491 ! Caching certain frequently-used derived values. 00492 IF (iEOS == EOS_ISOTHERMAL) THEN 00493 InternalEnergy=1.5d0*q(1)*Iso_Speed2 00494 ELSE 00495 InternalEnergy=q(iE)-half*DOT_PRODUCT(q(m_low:m_high),q(m_low:m_high))/q(1) 00496 IF(lMHD) InternalEnergy=InternalEnergy-half*SUM(q(iBx:iBz)**2) 00497 END IF 00498 END FUNCTION InternalEnergy 00499 00500 00506 FUNCTION PrimPress(q,kappa,chi,DPDQ) 00507 ! USE COOL 00508 ! USE TF 00509 ! USE VANDERWAALS 00510 REAL (KIND=qPrec), DIMENSION(:), INTENT(IN) :: q 00511 REAL (KIND=qPrec), OPTIONAL, INTENT(OUT) :: kappa,chi 00512 REAL (KIND=qPrec), DIMENSION(:), OPTIONAL, INTENT(OUT) :: DPDQ 00513 REAL (KIND=qPrec) :: PrimPress 00514 ! Internal variables 00515 REAL (KIND=qPrec) :: gam,kappa1,invkappa1,chi1,eth,ke,rho,invq1,Bsq,rho_cgs,eth_cgs 00516 00517 00518 ! Caching certain frequently-used derived values. 00519 00520 ! IF(iCooling/=2) THEN 00521 SELECT CASE(iEOS) 00522 CASE(EOS_IDEAL_GAS) ! Ideal 00523 kappa1=gamma1 00524 chi1=0.d0 00525 PrimPress = q(iE) 00526 ! CASE(1) ! analytic corrected-thomas-fermi 00527 ! print *, "Need to rewrite TFEOS for Primitive Variables" 00528 ! rho_cgs=rho*rScale 00529 ! eth_cgs=eth*pScale 00530 ! CALL TFEOS(rho_cgs,eth_cgs,Press,kappa1,chi1) 00531 ! Press=Press/pScale 00532 ! chi1=chi1*rScale/pScale 00533 ! CASE(2) !Sesame 00534 ! set pressure, kappa1, chi1 here. 00535 ! scaling factors (cgs): nScale (# density), rScale (density), pScale, velScale, lScale, runTimeSc 00536 !rho_cgs=rho*rScale 00537 !eth_cgs=eth*pScale 00538 !CALL SesameEOS(rho*rScale,Eth*pScale,Press,kappa1,chi1) 00539 !Press=Press/pScale 00540 !chi=chi1*rScale/pScale 00541 !kappa=kappa1*rScale/pScale 00542 ! CASE(3) ! van der Waals 00543 ! Call VDW_EOS(q(1),q(iE),PrimPress,chi1,kappa1) 00544 CASE(EOS_ISOTHERMAL) 00545 kappa1=0.0 00546 chi1=Iso_Speed2 00547 PrimPress = Iso_Speed2*q(1) 00548 CASE DEFAULT 00549 PRINT*,'unimplemented EOS, iEOS=',iEOS 00550 STOP 00551 END SELECT 00552 ! ELSE 00553 ! CALL EOS_vars(q,gamma=gam) 00554 ! kappa1=gam-1 00555 ! chi1=0.d0 00556 ! PrimPress = kappa1*q(iE) 00557 ! END IF 00558 00559 ! Return derivatives of p w.r.t. the field variables in q. 00560 IF(PRESENT(DPDQ)) THEN 00561 DPDQ = zero 00562 DPDQ(1) = chi1 00563 IF (iE .ne. 0) DPDQ(iE) = 1d0 00564 END IF 00565 ! if present, return chi and kappa to caller 00566 IF(Present(chi)) chi=chi1 00567 IF(Present(kappa)) kappa=kappa1 00568 END FUNCTION PrimPress 00569 00573 SUBROUTINE Prim_To_Cons_info(Info,mB) 00574 TYPE(InfoDef) :: Info 00575 INTEGER, DIMENSION(:,:) :: mB 00576 INTEGER :: i,j,k 00577 REAL(KIND=qpREC), DIMENSION(:,:,:,:), POINTER :: q 00578 q=>Info%q 00579 DO i=mb(1,1),mb(1,2) 00580 DO j=mb(2,1),mb(2,2) 00581 DO k=mb(3,1),mb(3,2) 00582 q(i,j,k,m_low:m_high)=q(i,j,k,m_low:m_high)*q(i,j,k,1) 00583 IF (iE /= 0) THEN 00584 q(i,j,k,iE)=gamma7*q(i,j,k,iE)+half*(SUM(q(i,j,k,m_low:m_high)**2)/q(i,j,k,1)) 00585 IF (lMHD) q(i,j,k,iE)=q(i,j,k,iE)+half*SUM(q(i,j,k,iBx:iBz)**2) 00586 END IF 00587 END DO 00588 END DO 00589 END DO 00590 END SUBROUTINE Prim_To_Cons_info 00591 00592 00593 00597 SUBROUTINE Cons_To_Prim_info(Info,mB) 00598 TYPE(InfoDef) :: Info 00599 INTEGER, DIMENSION(:,:) :: mB 00600 INTEGER :: i,j,k 00601 REAL(KIND=qpREC), DIMENSION(:,:,:,:), POINTER :: q 00602 q=>Info%q 00603 DO i=mb(1,1),mb(1,2) 00604 DO j=mb(2,1),mb(2,2) 00605 DO k=mb(3,1),mb(3,2) 00606 q(i,j,k,m_low:m_high)=q(i,j,k,m_low:m_high)/q(i,j,k,1) 00607 IF (iE /= 0) THEN 00608 q(i,j,k,iE)=gamma1*(q(i,j,k,iE)-half*(SUM(q(i,j,k,m_low:m_high)**2)*q(i,j,k,1))) 00609 IF (lMHD) q(i,j,k,iE)=q(i,j,k,iE)-gamma1*half*SUM(q(i,j,k,iBx:iBz)**2) 00610 END IF 00611 END DO 00612 END DO 00613 END DO 00614 END SUBROUTINE Cons_To_Prim_info 00615 00619 SUBROUTINE Cons_To_Prim_q(q) 00620 TYPE(InfoDef) :: Info 00621 REAL(KIND=qpREC), DIMENSION(:) :: q 00622 q(m_low:m_high)=q(m_low:m_high)/q(1) 00623 if (iE == 0 ) return 00624 q(iE)=gamma1*(q(iE)-half*(SUM(q(m_low:m_high)**2)*q(1))) 00625 IF (lMHD) q(iE)=q(iE)-gamma1*half*SUM(q(iBx:iBz)**2) 00626 END SUBROUTINE Cons_To_Prim_q 00627 00628 00632 SUBROUTINE prim_To_cons_q(q) 00633 REAL(KIND=qpREC), DIMENSION(:) :: q 00634 q(m_low:m_high)=q(1)*q(m_low:m_high) 00635 IF (iE == 0) RETURN 00636 q(iE)=gamma7*q(iE)+half*(SUM(q(m_low:m_high)**2)/q(1)) 00637 IF (lMHD) q(iE)=q(iE)+half*SUM(q(iBx:iBz)**2) 00638 00639 END SUBROUTINE Prim_To_cons_q 00640 00641 00642 00646 SUBROUTINE rprim_to_cons_Iso(w,u) 00647 REAL(KIND=qPrec), DIMENSION(:) :: w,u 00648 u=w 00649 u(m_low:m_high)=w(m_low:m_high)*w(1) 00650 END SUBROUTINE rprim_to_cons_Iso 00651 00655 SUBROUTINE f_prim_Iso(w,f) 00656 REAL(KIND=qPrec), DIMENSION(:) :: w,f 00657 REAL(KIND=qPrec) :: p, pT, B2,vx 00658 INTEGER :: n 00659 n=size(f) 00660 pT=w(1)*Iso_Speed2 00661 f(1)=w(1)*w(2) 00662 f(2)=f(1)*w(2)+pT 00663 f(3:m_high)=f(1)*w(3:m_high) 00664 f(m_high+1:n)=w(m_high+1:n)*w(2) 00665 END SUBROUTINE f_prim_Iso 00666 00669 FUNCTION fast_speed_Iso(w) 00670 REAL(KIND=qPrec), DIMENSION(:) :: w 00671 REAL(KIND=qPrec) :: fast_speed_Iso, a2b2, b2, bx2 00672 bx2=w(5)**2/w(1); 00673 b2=DOT_PRODUCT(w(6:7),w(6:7))/w(1)+bx2 00674 a2b2 = (Iso_Speed2+b2)/2d0 00675 fast_speed_Iso = sqrt((a2b2+sqrt(a2b2**2-Iso_Speed2*bx2))) 00676 END FUNCTION fast_speed_Iso 00677 00679 FUNCTION miniscule(x) 00680 logical :: miniscule 00681 REAL(8) :: x 00682 miniscule = (abs(x)<=1d-100) 00683 END FUNCTION miniscule 00684 00688 SUBROUTINE rprim_to_cons_MHD_Iso(w,u) 00689 REAL(KIND=qPrec), DIMENSION(:) :: w,u 00690 u=w 00691 u(2:4)=w(2:4)*w(1) 00692 END SUBROUTINE rprim_to_cons_MHD_Iso 00693 00697 SUBROUTINE f_prim_MHD_Iso(w,f) 00698 REAL(KIND=qPrec), DIMENSION(:) :: w,f 00699 REAL(KIND=qPrec) :: p, pT, B2,vx 00700 INTEGER :: n 00701 n=size(f) 00702 B2=DOT_PRODUCT(w(5:7),w(5:7)) 00703 pT=w(1)*Iso_Speed2+half*B2 00704 f(1)=w(1)*w(2) 00705 f(2)=f(1)*w(2)+pT-w(5)**2 00706 f(3)=f(1)*w(3)-w(5)*w(6) 00707 f(4)=f(1)*w(4)-w(5)*w(7) 00708 f(5)=0d0 00709 f(6)=w(6)*w(2)-w(5)*w(3) 00710 f(7)=w(7)*w(2)-w(5)*w(4) 00711 f(8:n)=w(8:n)*w(2) 00712 END SUBROUTINE f_prim_MHD_Iso 00713 00714 00717 FUNCTION sound_speed(w) 00718 REAL(KIND=qPrec), DIMENSION(:) :: w 00719 REAL(KIND=qPrec) :: sound_speed 00720 sound_speed=sqrt(gamma*w(2)/w(1)) 00721 END FUNCTION sound_speed 00722 00725 FUNCTION sound_speed_Roe(w) 00726 REAL(KIND=qPrec), DIMENSION(:) :: w 00727 REAL(KIND=qPrec) :: sound_speed_Roe 00728 sound_speed_Roe=sqrt(gamma1*(w(2)-half*SUM(w(3:m_high+1)**2))) 00729 END FUNCTION sound_speed_Roe 00730 00734 SUBROUTINE prim_to_roe(w,r) 00735 REAL(KIND=qPrec), DIMENSION(:) :: w, r 00736 r=w 00737 r(2)=half*(gamma10*w(2)/w(1)+SUM(w(3:m_high+1)**2)) !gamma10 = 2*gamma/(gamma-1) 00738 END SUBROUTINE prim_to_roe 00739 00743 SUBROUTINE prim_to_roe_MHD(w,r) 00744 REAL(KIND=qPrec), DIMENSION(:) :: w, r 00745 r=w 00746 r(2)=half*(gamma10*w(2)/w(1)+SUM(w(3:m_high+1)**2))+SUM(w(6:8)**2)/w(1) !gamma10 = 2*gamma/(gamma-1) 00747 END SUBROUTINE prim_to_roe_MHD 00748 00751 FUNCTION fast_speed(w) 00752 REAL(KIND=qPrec), DIMENSION(:) :: w 00753 REAL(KIND=qPrec) :: fast_speed, B, B2 00754 B2=DOT_PRODUCT(w(6:8),w(6:8)) 00755 B = gamma*w(2)+B2 00756 fast_speed = sqrt(half* (B + sqrt(B**2-4d0*gamma*w(2)*w(6)**2) ) /w(1) ) 00757 END FUNCTION fast_speed 00758 00761 FUNCTION fast_speed_Roe(w) 00762 REAL(KIND=qPrec), DIMENSION(:) :: w 00763 REAL(KIND=qPrec) :: fast_speed_Roe, B, B2, P 00764 B2=DOT_PRODUCT(w(6:8),w(6:8)) 00765 P=gamma14*(w(1)*(w(2)-half*SUM(w(3:m_high+1)**2))-B2) !gamma14 = (gamma-1)/gamma 00766 B = gamma*P+B2 00767 fast_speed_Roe = sqrt(half*(B+sqrt(B**2-4d0*gamma*P*w(6)**2))/w(1)) 00768 END FUNCTION fast_speed_Roe 00769 00773 SUBROUTINE rprim_to_cons(w,u) 00774 REAL(KIND=qPrec), DIMENSION(:) :: w,u 00775 u=w 00776 u(3:m_high+1)=w(3:m_high+1)*w(1) 00777 u(2)=gamma7*w(2)+half*DOT_PRODUCT(u(3:m_high+1),w(3:m_high+1)) 00778 END SUBROUTINE rprim_to_cons 00779 00783 SUBROUTINE rprim_to_cons_MHD(w,u) 00784 REAL(KIND=qPrec), DIMENSION(:) :: w,u 00785 u=w 00786 u(3:5)=w(3:5)*w(1) 00787 u(2)=gamma7*w(2)+half*(DOT_PRODUCT(u(3:5),w(3:5))+DOT_PRODUCT(w(6:8),w(6:8))) 00788 END SUBROUTINE rprim_to_cons_MHD 00789 00793 SUBROUTINE F_prim(w,f) 00794 REAL(KIND=qPrec), DIMENSION(:) :: w,f 00795 REAL(KIND=qPrec) :: p, Energy 00796 INTEGER :: n 00797 n=size(f) 00798 Energy=gamma7*w(2)+half*DOT_PRODUCT(w(3:m_high+1),w(3:m_high+1))*w(1) 00799 f(1)=w(1)*w(3) 00800 f(2)=(Energy+w(2))*w(3) 00801 f(3)=f(1)*w(3)+w(2) 00802 f(4:m_high+1)=f(1)*w(4:m_high+1) 00803 f(m_high+2:n)=w(m_high+2:n)*w(3) 00804 END SUBROUTINE f_prim 00805 00809 SUBROUTINE F_prim_MHD(w,f) 00810 REAL(KIND=qPrec), DIMENSION(:) :: w,f 00811 REAL(KIND=qPrec) :: p, pT, B2,Energy 00812 INTEGER :: n 00813 n=size(f) 00814 B2=DOT_PRODUCT(w(6:8),w(6:8)) 00815 pT=w(2)+half*B2 00816 Energy=gamma7*w(2)+half*(DOT_PRODUCT(w(3:5),w(3:5))*w(1)+DOT_PRODUCT(w(6:8),w(6:8))) 00817 f(1)=w(1)*w(3) 00818 f(2)=(Energy+pT)*w(3)-w(6)*DOT_PRODUCT(w(3:5),w(6:8)) 00819 f(3)=f(1)*w(3)+pT-w(6)**2 00820 f(4)=f(1)*w(4)-w(6)*w(7) 00821 f(5)=f(1)*w(5)-w(6)*w(8) 00822 f(6)=0d0 00823 f(7)=w(7)*w(3)-w(6)*w(4) 00824 f(8)=w(8)*w(3)-w(6)*w(5) 00825 f(9:n)=w(9:n)*w(3) 00826 END SUBROUTINE F_prim_MHD 00827 00828 00829 00832 FUNCTION GetMaxSpeed(q) 00833 REAL (KIND=qPrec), DIMENSION(:,:,:,:), INTENT(IN) :: q 00834 REAL (KIND=qPREC) :: GetMaxSpeed, fast_speed(3), A2, B2, B 00835 INTEGER :: i,j,k 00836 GetMaxSpeed=0d0 00837 DO i=1, size(q,1) 00838 DO j=1, size(q,2) 00839 DO k=1, size(q,3) 00840 IF (lMHD) THEN 00841 B2 = SUM(q(i,j,k,iBx:iBz)**2) 00842 A2=SoundSpeed(q(i,j,k,:))*q(i,j,k,1) 00843 B = A2+B2 00844 fast_speed(1:nDim) = sqrt(half*(B+sqrt(B**2-4d0*A2*q(i,j,k,iB(1:nDim))**2))/q(i,j,k,1)) 00845 GetMaxSpeed=max(GetMaxSpeed, maxval(q(i,j,k,imom(1:nDim))/q(i,j,k,1)+fast_speed(1:nDim))) 00846 ELSE 00847 GetMaxSpeed=max(GetMaxSpeed, sqrt(sum(q(i,j,k,imom(1:nDim))**2))/q(i,j,k,1)+SoundSpeed(q(i,j,k,:))) 00848 END IF 00849 END DO 00850 END DO 00851 END DO 00852 END FUNCTION GetMaxSpeed 00853 00854 00855 SUBROUTINE PrintQ(Info,q, tnow, i, j, k) 00856 TYPE(InfoDef) :: Info 00857 REAL(KIND=qPREC), DIMENSION(:) :: q 00858 REAL(KIND=qPREC) :: tnow 00859 INTEGER, OPTIONAL :: i,j,k 00860 write(*,'(A,20E13.2)') 'q= ',q(:) 00861 write(*,'(A,20E13.2)') 'Info%q= ',Info%q(i,j,k,:) 00862 IF (PRESENT(i) .AND. PRESENT(j) .AND. (Present(K))) THEN 00863 write(*,'(A,20E20.8)') 'pos= ', Info%xBounds(:,1)+(/i,j,k/)*levels(Info%level)%dx-half*levels(Info%level)%dx 00864 write(*,'(A,6I6)') 'i,j,k,mx,my,mz',i,j,k,Info%mX 00865 ! write(*,'(A,6I6)') 'Info%mGlobal= ',Info%mGlobal 00866 END IF 00867 write(*,'(A,6I6)') 'Info%level',Info%level 00868 write(*,'(A,20E20.8)') 'tnow= ', tnow 00869 ! write(*,'(A,6I6)') 'Info%Step',levels(Info%level)%step 00870 00871 ! write(*,'(A,6I6)') 'mb=', mB 00872 END SUBROUTINE PrintQ 00873 00885 SUBROUTINE calc_eigens(request_eigens, prim, req_dim, lambda, n, l, r,i,j,k,level) 00886 REAL(KIND=qPrec), DIMENSION(:) :: prim 00887 REAL(KIND=qPrec), DIMENSION(:,:) :: lambda 00888 REAL(KIND=qPrec), DIMENSION(:,:,:) :: r,l 00889 INTEGER, DIMENSION(:) :: n 00890 INTEGER :: i,j,k,m,dim,ii,jj,kk,level 00891 REAL(KIND=qPrec), DIMENSION(:,:,:), ALLOCATABLE :: A_matrix 00892 REAL(8) :: rho,u,Bx,By,Bz,rho_i,sqrt_rho,sqrt_rho_i,ca,cx,cy,cz, 00893 C_,Cy_,Cz_,a2,a,rhoa2,cf2,cs2,beta_f,beta_s, srho, alphaf, alphas, betay, betaz, cff, css, qf, qs, af, as, nf, ns, cfmcs, bt2, 00894 alpha,cf,cs,bT,ds,df,s,temp,temp1,temp2,A_,B_,D_,E_,F_,G_,hy,hz,hyr,hzr,b2, a_i, a2_i 00895 REAL(8), PARAMETER :: sqrt2_i=.7071067811865476d0, sqrt2=1.4142135623730950 00896 LOGICAL :: req_dim(:) 00897 LOGICAL :: request_eigens 00898 rho=prim(1) 00899 rho_i=1d0/prim(1) 00900 IF (iEOS==EOS_ISOTHERMAL) THEN 00901 a2=Iso_Speed2 00902 a=Iso_Speed 00903 ELSE 00904 a2=gamma*prim(iE)*rho_i 00905 a=sqrt(a2) 00906 END IF 00907 a_i=1d0/a 00908 a2_i=a_i**2 00909 lambda=0d0 00910 r=0d0 00911 l=0d0 00912 n=0 00913 rhoa2=rho*a2 00914 00915 IF (.NOT. lMHD) THEN 00916 n(1:nDim)=NrWaves 00917 lambda(1:nDim,1)=prim(ivx:ivx+nDim-1)-a 00918 lambda(1:nDim,2:NrWaves-1) = SPREAD(prim(ivx:ivx+nDim-1),2,NrWaves-2) 00919 lambda(1:nDim,NrWaves)=prim(ivx:ivx+nDim-1)+a 00920 00921 IF (iEOS==EOS_ISOTHERMAL) THEN 00922 r(1:nDim,1,1:2)=SPREAD((/1d0,-a*rho_i/),1,nDim) 00923 r(1:nDim,NrWaves,1)=r(1:nDim,1,1) 00924 r(1:nDim,NrWaves,2)=-r(1:nDim,1,2) 00925 00926 l(1:nDim,1,1:2)=SPREAD((/half,-half*rho*a_i/),1,nDim) 00927 l(1:nDim,NrWaves,1)=l(1:nDim,1,1) 00928 l(1:nDim,NrWaves,2)=-l(1:nDim,1,2) 00929 00930 DO dim=2,NrWaves-1 00931 r(1:nDim,dim,dim+1)=1d0 00932 l(1:nDim,dim,dim+1)=1d0 00933 END DO 00934 ELSE 00935 r(1:nDim,1,1:3)=SPREAD((/1d0,a2,-a*rho_i/),1,nDim) 00936 r(1:nDim,NrWaves,1:2)=r(1:nDim,1,1:2) 00937 r(1:nDim,NrWaves,3)=-r(1:nDim,1,3) 00938 r(1:nDim,2,1)=1d0 00939 00940 l(1:nDim,1,2:3)=SPREAD((/half*a2_i,-half*rho*a_i/),1,nDim) 00941 l(1:nDim,NrWaves,2)=l(1:nDim,1,2) 00942 l(1:nDim,NrWaves,3)=-l(1:nDim,1,3) 00943 l(1:nDim,2,1:2)=SPREAD((/1d0,-a2_i/),1,nDim) 00944 00945 DO dim=3,NrWaves-1 00946 r(1:nDim,dim,dim+1)=1d0 00947 l(1:nDim,dim,dim+1)=1d0 00948 END DO 00949 END IF 00950 ELSE 00951 sqrt_rho=sqrt(rho) 00952 sqrt_rho_i=1d0/sqrt_rho 00953 Bx=prim(iBx) 00954 By=prim(iBy) 00955 Bz=prim(iBz) 00956 cx=Bx*sqrt_rho_i 00957 cy=By*sqrt_rho_i 00958 cz=Bz*sqrt_rho_i 00959 C_=half*cx**2 00960 Cy_=half*cy**2 00961 Cz_=half*cz**2 00962 00963 A_=half*a2 00964 D_=A_+Cy_+Cz_+C_ 00965 00966 DO dim=1, nDim 00967 IF (.NOT. req_dim(dim)) THEN 00968 temp=cx 00969 cx=cy 00970 cy=cz 00971 cz=temp 00972 temp=C_ 00973 C_=Cy_ 00974 Cy_=Cz_ 00975 Cz_=temp 00976 temp=Bx 00977 Bx=By 00978 By=Bz 00979 Bz=temp 00980 CYCLE 00981 END IF 00982 00983 IF (Bx == 0) THEN 00984 s=0d0 00985 ELSE 00986 s = sign(1d0,Bx) 00987 END IF 00988 00989 B_=Cy_+Cz_ 00990 E_=2d0*sqrt(A_*C_) 00991 F_ = A_+B_-C_ 00992 G_ = 2d0*sqrt(B_*C_) 00993 bt=sqrt(2d0*B_) 00994 IF (my_small(E_,D_)) THEN 00995 cs2 = E_**2/(2d0*D_) 00996 cf2 = 2d0*D_ 00997 ELSE 00998 cs2 = D_ - sqrt(D_**2 - E_**2) 00999 cf2 = D_ + sqrt(D_**2 - E_**2) 01000 END IF 01001 cs = sqrt(cs2) 01002 cf = sqrt(cf2) 01003 ! write(*,*) 'cs = ', cs 01004 ! write(*,*) 'cf = ', cf 01005 ! write(*,*) A_, B_, C_, D_, E_, F_, G_ 01006 ! STOP 01007 u=prim(ivx+dim-1) 01008 ca=abs(cx) 01009 cfmcs=(cf2-cs2) 01010 IF (cfmcs < tiny(cfmcs)) THEN 01011 alphaf=1 01012 alphas=0 01013 ELSE 01014 alphaf=sqrt((a**2-cs**2)/cfmcs) 01015 alphas=sqrt((cf**2-a**2)/cfmcs) 01016 IF (ISNAN(alphas) .OR. ISNAN(alphaf)) THEN 01017 alphas=0 01018 alphaf=1 01019 END IF 01020 END IF 01021 bt2=sqrt(by**2+bz**2) 01022 IF (bt2 < tiny(bt2)) THEN 01023 betay=0d0 01024 betaz=0d0 01025 ELSE 01026 betay=by/bt2 01027 betaz=bz/bt2 01028 END IF 01029 ! write(*,*) bt2, betay, betaz, alphaf, alphas, cfmcs, tiny(cfmcs), tiny(cfmcs) 01030 cff=cf*alphaf 01031 css=cs*alphas 01032 qf=cf*alphaf*s 01033 qs=cs*alphas*s 01034 af=a*alphaf*sqrt(rho) 01035 as=a*alphas*sqrt(rho) 01036 rhoa2=rho*a2 01037 nf=half*a2_i 01038 ns=nf 01039 01040 ! srho=sqrt_rho 01041 01042 ! write(*,*) 'rho=', rho 01043 ! write(*,*) 'Bx=', Bx 01044 ! write(*,*) 'By=', By 01045 ! write(*,*) 'Bz=', Bz 01046 ! write(*,*) 'a=', a 01047 ! write(*,*) 'cf=', cf 01048 ! write(*,*) 'cs=', cs 01049 ! write(*,*) 'cff=',cff 01050 ! write(*,*) 'css=',css 01051 ! write(*,*) 'qf=', qf 01052 ! write(*,*) 'qs=', qs 01053 ! write(*,*) 'af=', af 01054 ! write(*,*) 'as=', as 01055 ! write(*,*) 'rhoa2=', rhoa2 01056 ! write(*,*) 'sqrt_rho=', sqrt_rho 01057 ! write(*,*) 'nf=', nf 01058 ! write(*,*) 'ns=', ns 01059 01060 IF (iEOS==EOS_ISOTHERMAL) THEN 01061 01062 r(dim,1,1:6)=(/rho*alphaf, -cff, qs*betay, qs*betaz, as*betay, as*betaz/) 01063 r(dim,2,1:6)=(/zero, zero, -betaz, betay, -betaz*s*sqrt_rho, betay*s*sqrt_rho/) 01064 r(dim,3,1:6)=(/rho*alphas, -css, -qf*betay, -qf*betaz, -af*betay, -af*betaz/) 01065 r(dim,NrWaves-2,1:6)=(/r(dim,3,1),-r(dim,3,2:4),r(dim,3,5:6)/) 01066 r(dim,NrWaves-1,1:6)=(/r(dim,2,1),-r(dim,2,2:4),r(dim,2,5:6)/) 01067 r(dim,NrWaves,1:6)=(/r(dim,1,1),-r(dim,1,2:4),r(dim,1,5:6)/) 01068 01069 l(dim,1,1:6)=(/ nf*alphaf*a2*rho_i, -nf*cff, nf*qs*betay, nf*qs*betaz, nf*as*betay*rho_i, nf*as*betaz*rho_i/) 01070 l(dim,2,1:6)=(/ zero, zero , -half*betaz, half*betay, -half*betaz*s*sqrt_rho_i, half*betay*s*sqrt_rho_i/) 01071 l(dim,3,1:6)=(/ ns*alphas*a2*rho_i, -ns*css, -ns*qf*betay, -ns*qf*betaz, -ns*af*betay*rho_i, -ns*af*betaz*rho_i/) 01072 l(dim,NrWaves-2,1:6)=(/l(dim,3,1),-l(dim,3,2:4),l(dim,3,5:6)/) 01073 l(dim,NrWaves-1,1:6)=(/l(dim,2,1),-l(dim,2,2:4),l(dim,2,5:6)/) 01074 l(dim,NrWaves,1:6)=(/l(dim,1,1),-l(dim,1,2:4),l(dim,1,5:6)/) 01075 01076 01077 n(dim)=NrWaves 01078 lambda(dim,1) = u - cf 01079 lambda(dim,2) = u - ca 01080 lambda(dim,3) = u - cs 01081 lambda(dim,NrWaves-2) = u + cs 01082 lambda(dim,NrWaves-1) = u + ca 01083 lambda(dim,NrWaves) = u + cf 01084 01085 lambda(dim,5:NrWaves-3) = u 01086 01087 DO ii=5, NrWaves-3 01088 l(dim,ii,NrCons-5+ii)=1d0 01089 r(dim,ii,NrCons-5+ii)=1d0 01090 END DO 01091 01092 ELSE 01093 01094 r(dim,1,1:7)=(/rho*alphaf, rhoa2*alphaf, -cff, qs*betay, qs*betaz, as*betay, as*betaz/) 01095 01096 r(dim,2,1:7)=(/zero, zero, zero, -betaz, betay, -betaz*s*sqrt_rho, betay*s*sqrt_rho/) 01097 r(dim,3,1:7)=(/rho*alphas, rhoa2*alphas, -css, -qf*betay, -qf*betaz, -af*betay, -af*betaz/) 01098 r(dim,4,1:7)=(/one, zero, zero, zero, zero, zero, zero/) 01099 r(dim,NrWaves-2,1:7)=(/r(dim,3,1:2),-r(dim,3,3:5),r(dim,3,6:7)/) 01100 r(dim,NrWaves-1,1:7)=(/r(dim,2,1:2),-r(dim,2,3:5),r(dim,2,6:7)/) 01101 r(dim,NrWaves,1:7)=(/r(dim,1,1:2),-r(dim,1,3:5),r(dim,1,6:7)/) 01102 01103 l(dim,1,1:7)=(/ zero, nf*alphaf*rho_i, -nf*cff, nf*qs*betay, nf*qs*betaz, nf*as*betay*rho_i, nf*as*betaz*rho_i/) 01104 l(dim,2,1:7)=(/ zero, zero, zero , -half*betaz, half*betay, -half*betaz*s*sqrt_rho_i, half*betay*s*sqrt_rho_i/) 01105 l(dim,3,1:7)=(/ zero, ns*alphas*rho_i, -ns*css, -ns*qf*betay, -ns*qf*betaz, -ns*af*betay*rho_i, -ns*af*betaz*rho_i/) 01106 l(dim,4,1:7)=(/ one, -one/a**2, zero, zero, zero, zero, zero/) 01107 l(dim,NrWaves-2,1:7)=(/l(dim,3,1:2),-l(dim,3,3:5),l(dim,3,6:7)/) 01108 l(dim,NrWaves-1,1:7)=(/l(dim,2,1:2),-l(dim,2,3:5),l(dim,2,6:7)/) 01109 l(dim,NrWaves,1:7)=(/l(dim,1,1:2),-l(dim,1,3:5),l(dim,1,6:7)/) 01110 01111 01112 n(dim)=NrWaves 01113 lambda(dim,1) = u - cf 01114 lambda(dim,2) = u - ca 01115 lambda(dim,3) = u - cs 01116 lambda(dim,4) = u 01117 lambda(dim,NrWaves-2) = u + cs 01118 lambda(dim,NrWaves-1) = u + ca 01119 lambda(dim,NrWaves) = u + cf 01120 01121 lambda(dim,5:NrWaves-3) = u 01122 01123 DO ii=5, NrWaves-3 01124 l(dim,ii,NrCons-5+ii)=1d0 01125 r(dim,ii,NrCons-5+ii)=1d0 01126 END DO 01127 01128 END IF 01129 temp=cx 01130 cx=cy 01131 cy=cz 01132 cz=temp 01133 temp=C_ 01134 C_=Cy_ 01135 Cy_=Cz_ 01136 Cz_=temp 01137 temp=Bx 01138 Bx=By 01139 By=Bz 01140 Bz=temp 01141 01142 END DO 01143 END IF 01144 END SUBROUTINE calc_eigens 01145 01146 01150 FUNCTION my_small(x,y) 01151 REAL(8), PARAMETER:: cutoff=1d-4 01152 logical :: my_small 01153 REAL(8) :: x,y 01154 my_small = (y .ne. 0d0 .AND. abs(x/y) <=cutoff) 01155 END FUNCTION my_small 01156 01160 FUNCTION really_small(x,y) 01161 logical :: really_small 01162 REAL(8) :: x,y 01163 really_small = (y .ne. 0d0 .AND. .1d0*x/y <= 1d-16) 01164 END FUNCTION really_small 01165 01166 END MODULE EOS 01167 01168