/*@@ @file Boundary.F @date 1 Apr 1999 @author Gabrielle Allen @desc Apply the relevant boundary condition to a group of GFs @enddesc @@*/ #include "cctk.h" #include "cctk_Parameters.h" #include "cctk_DefineThorn.h" #ifdef DEVELOPMENT_CARTOON2D #include "Development/Cartoon2D/src/include/cartoon2D_tensortypes.h" #endif /*@@ @routine ADM_Boundary @date 1 April 1999 @author Gabrielle Allen @desc Apply appropriate boundary condition for ADM to one GF @enddesc @calls @calledby @history @endhistory @@*/ subroutine ADM_Boundary(GH,name) implicit none DECLARE_CCTK_PARAMETERS CCTK_POINTER :: GH CHARACTER*(*) :: name integer, dimension(3) :: sw = 1 integer :: CCTK_Equals,ierr CCTK_REAL, PARAMETER :: zero = 0.0 CCTK_REAL, PARAMETER :: one = 1.0 c Apply symmetry boundary conditions. call CartSymGN(ierr,GH,name) c Apply Cartoon boundary conditions: This is a hack for now. #ifdef DEVELOPMENT_CARTOON2D if (cartoon==1) then if (name == "einstein::lapse") then call Cartoon2DBCGroup(ierr,GH,TENSORTYPE_SCALAR,name) else call Cartoon2DBCGroup(ierr,GH,TENSORTYPE_DDSYM,name) end if end if #endif c Apply flat boundaries. if (CCTK_EQUALS(bound,"flat")) then call BndFlatGN(ierr,GH,sw,name) c Apply radiation boundaries. else if (CCTK_EQUALS(bound,"radiation").or. . CCTK_EQUALS(bound,"radiative")) then if (name == "einstein::lapse") then call BndRadiativeVN(ierr,GH,sw,one,one, . "einstein::alp","adm::adm_alp_p") else if (name == "einstein::metric") then call BndRadiativeVN(ierr,GH,sw,one,one, . "einstein::gxx","adm::adm_gxx_p") call BndRadiativeVN(ierr,GH,sw,one,one, . "einstein::gyy","adm::adm_gyy_p") call BndRadiativeVN(ierr,GH,sw,one,one, . "einstein::gzz","adm::adm_gzz_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::gxy","adm::adm_gxy_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::gxz","adm::adm_gxz_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::gyz","adm::adm_gyz_p") else if (name == "einstein::curv") then call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::kxx","adm::adm_kxx_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::kyy","adm::adm_kyy_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::kzz","adm::adm_kzz_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::kxy","adm::adm_kxy_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::kxz","adm::adm_kxz_p") call BndRadiativeVN(ierr,GH,sw,zero,one, . "einstein::kyz","adm::adm_kyz_p") end if end if c Die if boundary condition was not applied. if (ierr < 0) then call CCTK_WARN(0,"Boundary could not be applied - giving up!") end if end subroutine ADM_Boundary