diff options
author | rhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45> | 2012-05-14 15:04:12 +0000 |
---|---|---|
committer | rhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45> | 2012-05-14 15:04:12 +0000 |
commit | 1e7919fe644c5591010a9d29bea68a644b504fd9 (patch) | |
tree | 16e21d22eb22ab199e00bf1df9370fd1804a7177 /src/GRHydro_CalcUpdate.F90 | |
parent | ad867c1344d6939d25dbf5c3b85304601ec20a90 (diff) |
GRHydro: fixes for GRMHD
All by Philipp Moesta.
1) Fix parity of psidc and divb
2) Fix a wrong index in the source terms of scon
3) Fix wrong indices of derivatives of space-time metric in the source of the
divergence cleaning scalar.
4) Calculate divergence of B in MoL PseudoEvolution and set its
Prolongation="Restrict".
5) Correct the source terms and fluxes for the Bfield and the divergence
cleaning field when having a non-flat space-time.
6) Make sure alpha factors match between UpdateCalculation and fluxes
definition.
7) Include 1/sqrt(detg) factor in calculation of \epsilon^{\muijk} in the cross
product to obtain the Bfield form the vector potential.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@330 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_CalcUpdate.F90')
-rw-r--r-- | src/GRHydro_CalcUpdate.F90 | 24 |
1 files changed, 12 insertions, 12 deletions
diff --git a/src/GRHydro_CalcUpdate.F90 b/src/GRHydro_CalcUpdate.F90 index 8c634d4..62653e9 100644 --- a/src/GRHydro_CalcUpdate.F90 +++ b/src/GRHydro_CalcUpdate.F90 @@ -39,7 +39,7 @@ subroutine UpdateCalculation(CCTK_ARGUMENTS) DECLARE_CCTK_FUNCTIONS CCTK_INT :: i,j,k,itracer - CCTK_REAL :: idx, alp_l, alp_r, Bvec_l, Bvec_r, alp_tmp + CCTK_REAL :: idx, alp_l, alp_r, Bcons_l, Bcons_r, alp_tmp idx = 1.d0 / CCTK_DELTA_SPACE(flux_direction) @@ -47,7 +47,7 @@ subroutine UpdateCalculation(CCTK_ARGUMENTS) if (use_weighted_fluxes == 0) then - !$OMP PARALLEL DO PRIVATE(i,j,k,itracer,alp_l,alp_r,alp_tmp,Bvec_l,Bvec_r) + !$OMP PARALLEL DO PRIVATE(i,j,k,itracer,alp_l,alp_r,alp_tmp,Bcons_l,Bcons_r) do k = GRHydro_stencil + 1 - transport_constraints, cctk_lsh(3) - GRHydro_stencil ! we need to compute Evec on all faces/edges where the fluxes are defined do j = GRHydro_stencil + 1 - transport_constraints, cctk_lsh(2) - GRHydro_stencil do i = GRHydro_stencil + 1 - transport_constraints, cctk_lsh(1) - GRHydro_stencil @@ -119,11 +119,11 @@ subroutine UpdateCalculation(CCTK_ARGUMENTS) Bcons(i+xoffset,j+1-xoffset,k+1 ,flux_direction)-Bcons(i ,j+zoffset ,k+1-zoffset,flux_direction)+ & Bcons(i+1 ,j+1 ,k+1 ,flux_direction)-Bcons(i+1-xoffset,j+1-yoffset,k+1-zoffset,flux_direction))*idx else - Bvec_l = 0.5d0 * (Bvec(i,j,k,flux_direction) + & - Bvec(i-xoffset,j-yoffset,k-zoffset,flux_direction)) - Bvec_r = 0.5d0 * (Bvec(i,j,k,flux_direction) + & - Bvec(i+xoffset,j+yoffset,k+zoffset,flux_direction)) - divB(i,j,k) = divB(i,j,k) + ( alp_l * Bvec_l - alp_r * Bvec_r ) * idx + Bcons_l = 0.5d0 * (Bcons(i,j,k,flux_direction) + & + Bcons(i-xoffset,j-yoffset,k-zoffset,flux_direction)) + Bcons_r = 0.5d0 * (Bcons(i,j,k,flux_direction) + & + Bcons(i+xoffset,j+yoffset,k+zoffset,flux_direction)) + divB(i,j,k) = divB(i,j,k) + (Bcons_l - Bcons_r ) * idx endif endif @@ -276,11 +276,11 @@ subroutine UpdateCalculation(CCTK_ARGUMENTS) psidcflux(i,j,k)) * idx endif if(track_divB.ne.0) then - Bvec_l = 0.5d0 * (Bvec(i,j,k,flux_direction) + & - Bvec(i-xoffset,j-yoffset,k-zoffset,flux_direction)) - Bvec_r = 0.5d0 * (Bvec(i,j,k,flux_direction) + & - Bvec(i+xoffset,j+yoffset,k+zoffset,flux_direction)) - divB(i,j,k) = divB(i,j,k) + ( Bvec_l - Bvec_r ) * idx + Bcons_l = 0.5d0 * (Bcons(i,j,k,flux_direction) + & + Bcons(i-xoffset,j-yoffset,k-zoffset,flux_direction)) + Bcons_r = 0.5d0 * (Bcons(i,j,k,flux_direction) + & + Bcons(i+xoffset,j+yoffset,k+zoffset,flux_direction)) + divB(i,j,k) = divB(i,j,k) + ( Bcons_l - Bcons_r ) * idx endif endif |