aboutsummaryrefslogtreecommitdiff
path: root/src/MoL_Init.F90
blob: 4683c0911ed4cff95e249de9fc5e374919cba519 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
! Registration of variables with MoL
! $Header$

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"

subroutine EHFinder_MoLRegister(CCTK_ARGUMENTS)

  implicit none

  CCTK_INT :: ierr, ierr_cum, varindex, rhsindex
  CCTK_INT :: MoLRegisterEvolved, MoLRegisterEvolvedGroup
  CCTK_INT :: l

  character(len=15) :: vname

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  ierr_cum = 0
  do l = 0, eh_number_level_sets - 1
    write(vname,'(a12,i1,a1)') 'ehfinder::f[', l, ']'
    call CCTK_VarIndex(varindex, vname )
    write(vname,'(a13,i1,a1)') 'ehfinder::sf[', l, ']'
    call CCTK_VarIndex(rhsindex, vname)
    ierr_cum = ierr_cum + MoLRegisterEvolved(varindex, rhsindex)
  end do

  if ( evolve_generators .gt. 0 ) then

    call CCTK_GroupIndex (varindex, 'ehfinder::xg')
    call CCTK_GroupIndex(rhsindex, 'ehfinder::dxg')

    ierr_cum = ierr_cum + MoLRegisterEvolvedGroup(varindex, rhsindex)

    call CCTK_GroupIndex (varindex, 'ehfinder::yg')
    call CCTK_GroupIndex(rhsindex, 'ehfinder::dyg')

    ierr_cum = ierr_cum + MoLRegisterEvolvedGroup(varindex, rhsindex)

    call CCTK_GroupIndex (varindex, 'ehfinder::zg')
    call CCTK_GroupIndex(rhsindex, 'ehfinder::dzg')

    ierr_cum = ierr_cum + MoLRegisterEvolvedGroup(varindex, rhsindex)

  end if
 
!  call CCTK_VarIndex(varindex, 'admbase::gxx')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::gxy')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::gxz')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::gyy')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::gyz')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::gzz')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::alp')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::betax')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::betay')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
!  call CCTK_VarIndex(varindex, 'admbase::betaz')
!  call MoL_RegisterPrimitive(ierr, varindex)
!  ierr_cum = ierr_cum + ierr
  print*,'MoLRegister ', ierr_cum
  if ( ierr_cum .gt. 0 ) then
    call CCTK_WARN(0,'Problems registering variables with MoL')
  end if
end subroutine EHFinder_MoLRegister