/*@@ @file Coord.F77 @date July 2 2000 @author Gabrielle Allen @desc Coordinates for the 2d wave equation solver @enddesc @@*/ #include "cctk.h" #include "cctk_Parameters.h" #include "cctk_Arguments.h" #include "cctk_Functions.h" /*@@ @routine WaveToy2DF77_RegisterCoords @date June 13 2000 @author Gabrielle Allen @desc Routine registers the coordinates for WaveToy2Df77 @enddesc @calls @calledby @history @endhistory @@*/ subroutine WaveToy2DF77_RegisterCoords(CCTK_ARGUMENTS) implicit none DECLARE_CCTK_ARGUMENTS DECLARE_CCTK_PARAMETERS DECLARE_CCTK_FUNCTIONS integer retval,ierr integer one CCTK_INT coord_system_handle retval = 0 ierr= Coord_SystemRegister(cctkGH, 2, "cart2d"); c call CCTK_INFO("Getting System Handle") coord_system_handle = Coord_SystemHandle(cctkGH, "cart2d") c call CCTK_INFO("Setting Coord Handle") ierr = Coord_CoordRegister(cctkGH, coord_system_handle, 1, "x") if (ierr .lt. 0) then call CCTK_WARN(1,"Problem with registering coordinate x") retval = -1 end if ierr = Coord_CoordRegister(cctkGH, coord_system_handle, 2, "y") if (ierr .lt. 0) then call CCTK_WARN(1,"Problem with registering coordinate y") retval = -1 end if c call CCTK_INFO("Setting Coord Type") call Util_TableSetString(ierr, coord_system_handle, "uniform", "TYPE") c call CCTK_INFO("Setting Default System") ierr = Coord_SetDefaultSystem(cctkGH, "cart2d"); return end /*@@ @routine WaveToy2DF77_Coord @date June 13 2000 @author Gabrielle Allen @desc Coordinates for the 2d wave equation @enddesc @calls @calledby @history @endhistory @@*/ subroutine WaveToy2DF77_Coord(CCTK_ARGUMENTS) implicit none DECLARE_CCTK_ARGUMENTS DECLARE_CCTK_FUNCTIONS INTEGER i,j,ierr,iconv CCTK_REAL xmin,xmax,x_origin,lowerx,upperx CCTK_REAL ymin,ymax,y_origin,lowery,uppery CCTK_INT coord_handle INTEGER varindex xmin=0 xmax=1 ymin=0 ymax=1 c Grid spacing on coarsest grid cctk_delta_space(1) = (xmax - xmin) / max(cctk_gsh(1) - 1, 1) cctk_delta_space(2) = (ymax - ymin) / max(cctk_gsh(2) - 1, 1) x_origin = xmin y_origin = ymin do j=1,cctk_lsh(2) do i=1,cctk_lsh(1) x(i,j) = cctk_delta_space(1)*(i-1+cctk_lbnd(1)) + x_origin y(i,j) = cctk_delta_space(2)*(j-1+cctk_lbnd(2)) + y_origin end do end do lowerx = x_origin upperx = x_origin+cctk_delta_space(1)*(cctk_gsh(1)-1) c call CCTK_INFO("Getting Coord Handle") coord_handle = Coord_CoordHandle(cctkGH, "x", "cart2d") if (coord_handle .lt. 0) then call CCTK_WARN(0, "Error retrieving coordinate handle for x of cart2d") endif c call CCTK_INFO("Getting Variable Index") call CCTK_VarIndex(varindex, "wavetoy2df77::x") c call CCTK_INFO("Setting Physical Minimum") call Util_TableSetReal (ierr,coord_handle, lowerx, "PHYSICALMIN") c call CCTK_INFO("Setting Computational Minimum") call Util_TableSetReal (ierr,coord_handle, lowerx, "COMPMIN") c call CCTK_INFO("Setting Physical Maximum") call Util_TableSetReal (ierr,coord_handle, upperx, "PHYSICALMAX") c call CCTK_INFO("Setting Computational Maximum") call Util_TableSetReal (ierr,coord_handle, upperx, "COMPMAX") c call CCTK_INFO("Setting Coordinate Type") call Util_TableSetString (ierr,coord_handle, "uniform", "TYPE") c call CCTK_INFO("Setting Time dependency") call Util_TableSetString (ierr,coord_handle, "no", "TIMEDEPENDENT") c call CCTK_INFO("Setting datatype") call Util_TableSetString (ierr,coord_handle, "CCTK_REAL", "DATATYPE") c call CCTK_INFO("Setting index") call Util_TableSetInt (ierr,coord_handle, varindex, "GAINDEX") c call CCTK_INFO("Setting delta") call Util_TableSetReal (ierr,coord_handle, cctk_delta_space(1), "DELTA") c call CCTK_INFO("Finished Setting up Coordinates") lowery = y_origin uppery = y_origin+cctk_delta_space(2)*(cctk_gsh(2)-1) c call CCTK_INFO("Getting Coord Handle") coord_handle = Coord_CoordHandle(cctkGH, "y", "cart2d") if (coord_handle .lt. 0) then call CCTK_WARN(0, "Error retrieving coordinate handle for y of cart2d") endif c call CCTK_INFO("Getting Variable Index") call CCTK_VarIndex(varindex, "wavetoy2df77::y") c call CCTK_INFO("Setting Physical Minimum") call Util_TableSetReal (ierr,coord_handle, lowery, "PHYSICALMIN") c call CCTK_INFO("Setting Computational Minimum") call Util_TableSetReal (ierr,coord_handle, lowery, "COMPMIN") c call CCTK_INFO("Setting Physical Maximum") call Util_TableSetReal (ierr,coord_handle, uppery, "PHYSICALMAX") c call CCTK_INFO("Setting Computational Maximum") call Util_TableSetReal (ierr,coord_handle, uppery, "COMPMAX") c call CCTK_INFO("Setting Coordinate Type") call Util_TableSetString (ierr,coord_handle, "uniform", "TYPE") c call CCTK_INFO("Setting Time dependency") call Util_TableSetString (ierr,coord_handle, "no", "TIMEDEPENDENT") c call CCTK_INFO("Setting datatype") call Util_TableSetString (ierr,coord_handle, "CCTK_REAL", "DATATYPE") c call CCTK_INFO("Setting index") call Util_TableSetInt (ierr,coord_handle, varindex, "GAINDEX") c call CCTK_INFO("Setting delta") call Util_TableSetReal (ierr,coord_handle, cctk_delta_space(2), "DELTA") c call CCTK_INFO("Finished Setting up Coordinates") return end