aboutsummaryrefslogtreecommitdiff
path: root/src/MoL_Init.F90
diff options
context:
space:
mode:
authordiener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c>2003-06-26 12:09:25 +0000
committerdiener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c>2003-06-26 12:09:25 +0000
commitffb807f5d7b17fe526de381d63098349199f2900 (patch)
tree39a9801a2339d9076af472692bfd90ac6e74ea3e /src/MoL_Init.F90
parent88a32348f469a016ee14a18ea4f52f3817943a6c (diff)
Change from registering the individual generator coordinates separately with
MoL to register the whole group at once. From Ian Hawke. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/EHFinder/trunk@120 2a26948c-0e4f-0410-aee8-f1d3e353619c
Diffstat (limited to 'src/MoL_Init.F90')
-rw-r--r--src/MoL_Init.F9027
1 files changed, 15 insertions, 12 deletions
diff --git a/src/MoL_Init.F90 b/src/MoL_Init.F90
index bccc870..6a43ddd 100644
--- a/src/MoL_Init.F90
+++ b/src/MoL_Init.F90
@@ -10,7 +10,7 @@ subroutine EHFinder_MoLRegister(CCTK_ARGUMENTS)
implicit none
CCTK_INT :: ierr, ierr_cum, varindex, rhsindex
- CCTK_INT :: MoLRegisterEvolved
+ CCTK_INT :: MoLRegisterEvolved, MoLRegisterEvolvedGroup
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_ARGUMENTS
@@ -25,20 +25,23 @@ subroutine EHFinder_MoLRegister(CCTK_ARGUMENTS)
if ( evolve_generators .gt. 0 ) then
-! call CCTK_GroupIndex (varindex, 'ehfinder::generators')
-! call CCTK_GroupIndex(rhsindex, 'ehfinder::generators_rhs')
+ call CCTK_GroupIndex (varindex, 'ehfinder::generators')
+ call CCTK_GroupIndex(rhsindex, 'ehfinder::generators_rhs')
- call CCTK_VarIndex(varindex, 'ehfinder::xg')
- call CCTK_VarIndex(rhsindex, 'ehfinder::dxg')
- ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
+!!$ call CCTK_VarIndex(varindex, 'ehfinder::xg')
+!!$ call CCTK_VarIndex(rhsindex, 'ehfinder::dxg')
+!!$ ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
+!!$
+!!$ call CCTK_VarIndex(varindex, 'ehfinder::yg')
+!!$ call CCTK_VarIndex(rhsindex, 'ehfinder::dyg')
+!!$ ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
+!!$
+!!$ call CCTK_VarIndex(varindex, 'ehfinder::zg')
+!!$ call CCTK_VarIndex(rhsindex, 'ehfinder::dzg')
+!!$ ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
- call CCTK_VarIndex(varindex, 'ehfinder::yg')
- call CCTK_VarIndex(rhsindex, 'ehfinder::dyg')
- ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
+ ierr_cum = ierr_cum + MoLRegisterEvolvedGroup(varindex, rhsindex)
- call CCTK_VarIndex(varindex, 'ehfinder::zg')
- call CCTK_VarIndex(rhsindex, 'ehfinder::dzg')
- ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
end if
! call CCTK_VarIndex(varindex, 'admbase::gxx')