aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_MP5Reconstruct_drv.F90
diff options
context:
space:
mode:
authorrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-07-06 18:11:18 +0000
committerrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-07-06 18:11:18 +0000
commit2b6919940238fa6fa03f75237e8d30342f42e799 (patch)
tree442c9256424f48fd2393e60c8ff1ae2d37f9ab92 /src/GRHydro_MP5Reconstruct_drv.F90
parent669d5b37c19d16380ba6c63a7195cbcea8d08f35 (diff)
GRHydro: re-implement prim2con routine in C
* new experimental prim2con routine that is about twice as fast as the old one. * Right now tested only for simple EOS. * Right now handles only prim2con call after reconstruction and must be enabled by changing comments in GRHydro_Reconstruct.F90 * remove duplicate code in all reconstruction *_drv.F90 routines and just have the initialization of plus and minus variables in the main GRHydro_Reconstruct.F90 routine From: Christian Ott <cott@bethe.tapir.caltech.edu> git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@556 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_MP5Reconstruct_drv.F90')
-rw-r--r--src/GRHydro_MP5Reconstruct_drv.F9059
1 files changed, 6 insertions, 53 deletions
diff --git a/src/GRHydro_MP5Reconstruct_drv.F90 b/src/GRHydro_MP5Reconstruct_drv.F90
index fc839c0..bb222dd 100644
--- a/src/GRHydro_MP5Reconstruct_drv.F90
+++ b/src/GRHydro_MP5Reconstruct_drv.F90
@@ -133,65 +133,18 @@ subroutine GRHydro_MP5Reconstruct_drv(CCTK_ARGUMENTS)
ny = cctk_lsh(2)
nz = cctk_lsh(3)
-!!$ Initialize variables that store reconstructed quantities:
-
- !$OMP PARALLEL DO PRIVATE(i,j,k)
- do k=1,cctk_lsh(3)
- do j=1,cctk_lsh(2)
- do i=1,cctk_lsh(1)
+ !initialize trivial_rp to false
+ !$OMP PARALLEL DO
+ do k=1,nz
+ do j=1,ny
+ do i=1,nx
trivial_rp(i,j,k) = .false.
- rhoplus(i,j,k) = 0.0d0
- rhominus(i,j,k)= 0.0d0
- epsplus(i,j,k) = 0.0d0
- epsminus(i,j,k) = 0.0d0
- velxplus(i,j,k) = 0.0d0
- velxminus(i,j,k) = 0.0d0
- velyplus(i,j,k) = 0.0d0
- velyminus(i,j,k) = 0.0d0
- velzplus(i,j,k) = 0.0d0
- velzminus(i,j,k) = 0.0d0
-
- if(evolve_mhd.ne.0) then
- Bvecxplus(i,j,k) = 0.0d0
- Bvecxminus(i,j,k) = 0.0d0
- Bvecyplus(i,j,k) = 0.0d0
- Bvecyminus(i,j,k) = 0.0d0
- Bveczplus(i,j,k) = 0.0d0
- Bveczminus(i,j,k) = 0.0d0
- if(clean_divergence.ne.0) then
- psidcplus(i,j,k) = 0.0d0
- psidcminus(i,j,k) = 0.0d0
- endif
- endif
-
- if (evolve_entropy .ne. 0) then
- entropyplus(i,j,k) = 0.0d0
- entropyminus(i,j,k) = 0.0d0
- endif
-
- if (evolve_tracer .ne. 0) then
- tracerplus(i,j,k,:) = 0.0d0
- tracerminus(i,j,k,:) = 0.0d0
- endif
-
- if (evolve_Y_e .ne. 0) then
- Y_e_plus(i,j,k) = 0.0d0
- Y_e_minus(i,j,k) = 0.0d0
- endif
-
- if (evolve_temper .ne. 0) then
- ! set this to cell center value to have
- ! good initial guesses at interfaces
- ! in case we don't reconstruct temp
- tempplus(i,j,k) = temperature(i,j,k)
- tempminus(i,j,k) = temperature(i,j,k)
- endif
-
enddo
enddo
enddo
!$OMP END PARALLEL DO
+
!!$ MP5 starts:
if (evolve_tracer .ne. 0) then
do itracer=1,number_of_tracers