blob: c6d503b0fc62111f6043a2ec232388352ea38357 [file] [log] [blame]
/*BHEADER**********************************************************************
* (c) 1997 The Regents of the University of California
*
* See the file COPYRIGHT_and_DISCLAIMER for a complete copyright
* notice, contact person, and disclaimer.
*
* $Revision$
*********************************************************************EHEADER*/
/******************************************************************************
*
* HYPRE_StructGrid interface
*
*****************************************************************************/
#include "headers.h"
/*--------------------------------------------------------------------------
* HYPRE_StructGridCreate
*--------------------------------------------------------------------------*/
int
HYPRE_StructGridCreate( MPI_Comm comm,
int dim,
HYPRE_StructGrid *grid )
{
int ierr;
ierr = hypre_StructGridCreate(comm, dim, grid);
return ierr;
}
/*--------------------------------------------------------------------------
* HYPRE_StructGridDestroy
*--------------------------------------------------------------------------*/
int
HYPRE_StructGridDestroy( HYPRE_StructGrid grid )
{
return ( hypre_StructGridDestroy(grid) );
}
/*--------------------------------------------------------------------------
* HYPRE_StructGridSetExtents
*--------------------------------------------------------------------------*/
int
HYPRE_StructGridSetExtents( HYPRE_StructGrid grid,
int *ilower,
int *iupper )
{
hypre_Index new_ilower;
hypre_Index new_iupper;
int d;
hypre_ClearIndex(new_ilower);
hypre_ClearIndex(new_iupper);
for (d = 0; d < hypre_StructGridDim((hypre_StructGrid *) grid); d++)
{
hypre_IndexD(new_ilower, d) = ilower[d];
hypre_IndexD(new_iupper, d) = iupper[d];
}
return ( hypre_StructGridSetExtents(grid, new_ilower, new_iupper) );
}
/*--------------------------------------------------------------------------
* HYPRE_SetStructGridPeriodicity
*--------------------------------------------------------------------------*/
int
HYPRE_StructGridSetPeriodic( HYPRE_StructGrid grid,
int *periodic )
{
hypre_Index new_periodic;
int d;
hypre_ClearIndex(new_periodic);
for (d = 0; d < hypre_StructGridDim(grid); d++)
{
hypre_IndexD(new_periodic, d) = periodic[d];
}
return ( hypre_StructGridSetPeriodic(grid, new_periodic) );
}
/*--------------------------------------------------------------------------
* HYPRE_StructGridAssemble
*--------------------------------------------------------------------------*/
int
HYPRE_StructGridAssemble( HYPRE_StructGrid grid )
{
return ( hypre_StructGridAssemble(grid) );
}