FronTier++
Softerware for front tracking method
 All Files Functions Groups
Functions
fapi.h File Reference

The fapi.h contains the functions used to operate the interface. More...

#include <front/fdecs.h>

Go to the source code of this file.

Functions

IMPORT void FT_Init (int argc, char **argv, F_BASIC_DATA *f_basic)
 Process the command line arguments, set IO handles, and store initialization parameters in f_basic, Including read the number of processors, the dimesension of the problem, partition of the computational domain, restart message, the input file and the output directory. -d dim -p nx [ny] [nz] -i input-name -o output-dir-name -r restart-dir-name -t restart-step. More...
 
IMPORT void FT_ReadSpaceDomain (char *in_name, F_BASIC_DATA *f_basic)
 Read from input file the computational domain information of the problem, including the domain limits, computational grid, and the types of the rectangular boundaries. The information is stored in the structure f_basic. More...
 
IMPORT void FT_ReadComparisonDomain (char *in_name, F_BASIC_DATA *f_basic)
 Read from input file the comparison domain information of the problem, including the domain limits, computational grid, and the types of the rectangular boundaries. The information is stored in the structure f_basic. More...
 
IMPORT void FT_StartUp (Front *front, F_BASIC_DATA *f_basic)
 Initialize front computational grid, interface and function hooks, default function calls, and if restart, read interface from restart directory. Constructor for Front object. More...
 
IMPORT void FT_InitDebug (char *inname)
 Initialize strings for debugging option of the problem, read the input file. More...
 
IMPORT void FT_InitIntfc (Front *front, LEVEL_FUNC_PACK *level_func_pack)
 Initialize the interface interface curves (2D) and surfaces (3D) using level function and parameters, or in some cases, point set. Install boundary curves (2D) or surfaces (3D). For 2D, the boundary is installed by default. For 3D, the boundary is installed on request. More...
 
IMPORT void FT_ClipIntfcToSubdomain (Front *front)
 Clip the Initial interface of the front to a parallel subdomain. More...
 
IMPORT void FT_InitFrontVeloFunc (Front *front, VELO_FUNC_PACK *velo_func_pack)
 Initialize front velocity function for front point propagation. The velocity function use point and other related structures as input, must also supply parameters needed for the velocity function. More...
 
IMPORT void FT_ReadTimeControl (char *in_name, Front *front)
 Read the time domain and control information from input file, including maximum time, maximum step, restart print interval, movie frame output interval, CFL factor, and redistribution step interval. More...
 
IMPORT void FT_ResetTime (Front *front)
 Reset the time to 0.0, time step to 0, print interval index to 0, movie interval index to 0. More...
 
IMPORT POINTER * FT_CreateLevelHyperSurfs (RECT_GRID *rgr, INTERFACE *intfc, int neg_comp, int pos_comp, double(*func)(POINTER, double *), POINTER func_params, int w_type, int *num_hs)
 This function creates a set of hypersurfaces (curves in 2D and surfaces in 3D) using a level function (provided by the caller). The function return the handle for the array of hyper surfaces as (POINTER*). More...
 
IMPORT void FT_Propagate (Front *front)
 Propagate the Front for one time step. The process includes advancing the front forward by one step to new positions, redistributing new interface mesh and resolving physical and topological bifurcations. The interface in the front is replaced by a new and propagated one and the old one is freed. More...
 
IMPORT void FT_RedistMesh (Front *front)
 This is an independent call for redistribution and optimization of the interface mesh. A parallel communication of front should be called following this call to ensure that the redistribution is consistent globally. More...
 
IMPORT void FT_OptimizeMesh (Front *front, SCALED_REDIST_PARAMS params)
 This is an independent call for redistribution and optimization of the interface mesh. A parallel communication of front should be called following this call to ensure that the redistribution is consistent globally. The function will redistribute all curves and surfaces (3D). It will recursively do either 10 times or when nothing to be done. More...
 
IMPORT boolean FT_OptimizeSurfMesh (Front *front, SURFACE *surf, SCALED_REDIST_PARAMS params)
 This is an independent call for redistribution and optimization of a surf mesh. No parallel communication is called after redistribution. Return YES (nothing_done == YES) if no more triangle to be redistributed. More...
 
IMPORT boolean FT_OptimizeCurveMesh (Front *front, CURVE *curve, SCALED_REDIST_PARAMS params)
 This is an independent call for redistribution and optimization of a curve mesh. No parallel communication is called after redistribution. Return YES (nothing_done == YES) if no more bond to be redistributed. More...
 
IMPORT void FT_SetCurveSpacing (Front *front, double scaled_spacing)
 This function set the optimal spacing for curve redistribution. The scaled spacing is in the unit of regular grid spacing h. More...
 
IMPORT void FT_OptimizeCurveMeshWithEqualBonds (Front *front, CURVE *curve)
 This is an independent call for redistribution and optimization of curve. A parallel communication of front should be called following this call to ensure that the redistribution is consistent globally. More...
 
IMPORT void FT_SetTimeStep (Front *front)
 Calculate front->dt for next time step using the recorded maximum speed from previous time step, it step is reduced by a CFL factor to to ensure numerical stability. More...
 
IMPORT void FT_SetOutputCounter (Front *front)
 This function is used in restart to set the printing index and movie output index to appropiate number according to the restart time. More...
 
IMPORT void FT_TimeControlFilter (Front *front)
 To further reduce time step if the printing or movie output time is smaller than the time after next time step. Increment priting or/and movie output index if either of both of them are met. More...
 
IMPORT boolean FT_IsSaveTime (Front *front)
 Signals that time is reached for printing restart files. returns YES to indicate that restart file should be generated. This function does not write the output. More...
 
IMPORT boolean FT_IsDrawTime (Front *front)
 Signals that time is reached for output of a movie frame. returns YES to indicate that a movie frame should be generated. This function does not write the movie frame. More...
 
IMPORT boolean FT_TimeLimitReached (Front *front)
 Signals that time is reached for termination of the run. Return YES either maximum time has been reached or maximum time step has been reached. More...
 
IMPORT void FT_RecordMaxFrontSpeed (int dir, double speed, POINTER state, double *coords, Front *front)
 This function compare and record maximum speed of the front. The recorded speed will be used to determine the time step which must satisfy the CFL condition. More...
 
IMPORT void FT_AddTimeStepToCounter (Front *front)
 Add front->dt to front->time after propagation. More...
 
IMPORT void FT_Save (Front *front)
 Output front geometric data to the directory of out_name. The data can be used for restart of the run. More...
 
IMPORT void FT_Draw (Front *front)
 Output a movie frame, currently includes GD, hdf, vtk formats. More...
 
IMPORT void FT_XgraphSampleLine (char *dirname, char *varname, boolean data_in_domain, int size, double *x, double *var)
 This function output variable sample along a grid line as xgraph data file. It considers parallelization of subdomains. More...
 
IMPORT void FT_MakeGridIntfc (Front *front)
 Make a duplicate interface whose topological grid is the expanded dual grid, currently with buffer of 4h for PERIODIC and SUBDOMAIN boundary, 1h for other boundaries. Install crossings of the interface and the expanded dual grid and store them in the interface table. These crossings can be used to interact with various PDE solvers. More...
 
IMPORT void FT_FreeGridIntfc (Front *front)
 Delete and free space of grid crossing interface made by the function FT_MakeGridIntfc(). More...
 
IMPORT void FT_MakeCompGridIntfc (Front *front)
 Make a duplicate interface whose topological grid is the expanded comp grid, currently with buffer of 4h for PERIODIC and SUBDOMAIN boundary, 1h for other boundaries. Install crossings of the interface and the expanded dual grid and store them in the interface table. These crossings can be used to interact with various PDE solvers. More...
 
IMPORT void FT_FreeCompGridIntfc (Front *front)
 
IMPORT void FT_FreeOldGridIntfc (Front *front)
 Delete and free space of grid crossing interface of front->old_grid_intfc. More...
 
IMPORT void FT_FreeFront (Front *front)
 Delete and free space occupied by the front including grid_intfc if still there, and interf. More...
 
IMPORT void FT_FreeMainIntfc (Front *front)
 Delete and free space occupied by the primary interface of the front. More...
 
IMPORT boolean FT_NormalAtGridCrossing (Front *front, int *icoords, GRID_DIRECTION dir, int comp, double *nor, HYPER_SURF **hs, double *crx_coords)
 Standing at grid icoords, looking to the direction dir, this function looks for the nearest interface cross on the grid line segment. The function returns YES if the crossing exists, in such case, the crossing coordinates are copied to crx_coords, the corresponding hyper surface (curce in 2D and surface in 3D) is assigned to hs, and the normal vector to the side of comp. If no crossing exists, the function return NO;. More...
 
IMPORT boolean FT_StateStructAtGridCrossing (Front *front, INTERFACE *grid_intfc, int *icoords, GRID_DIRECTION dir, int comp, POINTER *state, HYPER_SURF **hs, double *crx_coords)
 
IMPORT boolean FT_StateVarAtGridCrossing (Front *front, int *icoords, GRID_DIRECTION dir, int comp, double(*state_func)(Locstate), double *ans, double *crx_coords)
 
IMPORT HYPER_SURF * FT_HyperSurfAtGridCrossing (Front *front, int *icoords, GRID_DIRECTION dir, int wave_type)
 Sitting at icoords and look to the direction dir, this function detects the nearest hyper surface (curve in 2D and surface in 3D) on the grid segment. Return pointer to hyper surface if there is one, return NULL if no crossing hyper surface is found. More...
 
IMPORT boolean FT_IntrpStateVarAtCoords (Front *front, int comp, double *coords, double *var_array, double(*state_func)(POINTER), double *ans, double *default_ans)
 Interpolate a state variable at a space point with coords. If comp == NO_COMP, it interpolates with no regard of interface. Otherwise it will interpolate in the subdomain of comp. The state_func() is needed to tell the function how to retrieve the variable from the interface state. The interpolated variable is assigned in ans. Return YES if the interpolation is successful. More...
 
IMPORT boolean FT_CompGridIntrpStateVarAtCoords (Front *front, int comp, double *coords, double *var_array, double(*state_func)(POINTER), double *ans, double *default_ans)
 Interpolate a state variable at a space point with coords on computational grid. If comp == NO_COMP, it interpolates with no regard of interface. Otherwise it will interpolate in the subdomain of comp. The state_func() is needed to tell the function how to retrieve the variable from the interface state. The interpolated variable is assigned in ans. Return YES if the interpolation is successful. More...
 
IMPORT boolean FT_NearestRectGridVarInRange (Front *front, int comp, double *coords, double *var_array, int range, double *ans)
 Find the state variable on rectangular grid point which has the same component as the input and is nearest to the input coordinate. Return YES if such point is found, and NO if no such point is found. In the latter case, the value of the ans is set to zero. More...
 
IMPORT boolean FT_FindNearestIntfcPointInRange (Front *front, int comp, double *coords, USE_BOUNDARIES bdry, double *intfc_point, double *t, HYPER_SURF_ELEMENT **hse, HYPER_SURF **hs, int range)
 
IMPORT void FT_NormalAtPoint (POINT *p, Front *front, double *normal, int comp)
 Get the normal vector at the point p. More...
 
IMPORT void FT_CurvatureAtPoint (POINT *p, Front *front, double *curvature)
 Get the curvature at the point p. More...
 
IMPORT double FT_GridSizeInDir (double *dir, Front *front)
 The grid size in direction dir is (dir dot h). If h are the same in all directions (square mesh), it returns the h. More...
 
IMPORT Nor_stencil * FT_CreateNormalStencil (Front *front, POINT *p, int comp, int num_pts)
 This function create a normal stencil at the interface point p with size (number of points) num_pts. The normal stencil in in the ambient with component comp. More...
 
IMPORT boolean FT_ReflectPointThroughBdry (Front *front, HYPER_SURF *hs, double *coords, int comp, double *coords_bdry, double *coords_ref, double *normal)
 Given the coordinates coords, this function find the reflected coordinates coordsrefl through the hypersurface hs, it also provide the normal vector nor at the reflection. Return NO if conditions not satisfied. More...
 
IMPORT void FT_SetDirichletBoundary (Front *front, void(*state_func)(double *, HYPER_SURF *, Front *, POINTER, POINTER), const char *state_func_name, POINTER state_func_params, POINTER state, HYPER_SURF *hs)
 This function sets state of a Dirichlet boundary as a hyper surface. It provide two methods to set the boundary: (a) set the boundary with a constant state, in this case, the address of the state as a pointer must be supplied as a pointer while state_func, state_func_name and state_func_params must be set to NULL; (b) set the boundary state through a state function, in this case, state must be set to NULL while state_func, state_func_name and state_func_params must be supplied. Currently only flow-through and time-dependent functions are available. More...
 
IMPORT void FT_InsertDirichletBoundary (Front *front, void(*state_func)(double *, HYPER_SURF *, Front *, POINTER, POINTER), const char *state_func_name, POINTER state_func_params, POINTER state, HYPER_SURF *hs, int index)
 This function insert state to a Dirichlet boundary. The six (four for 2D) sides of the boundary memory has been allocated, this function fills in the content of the boundary state. It provide two methods to set the boundary: (a) set the boundary with a constant state, in this case, the address of the state as a pointer must be supplied as a pointer while state_func, state_func_name and state_func_params must be set to NULL; (b) set the boundary state through a state function, in this case, state must be set to NULL while state_func, state_func_name and state_func_params must be supplied. Currently only flow-through and time-dependent functions are available. More...
 
IMPORT HYPER_SURF ** FT_MixedBoundaryHypSurfs (INTERFACE *intfc, int idir, int nb, int w_type, int *num_hs)
 This function returns a set of hyper surfaces (curves in 2D and surfaces in 3D) whose wave type matches the input w_type. If the boundary on idir and side nb is not MIXED_TYPE_BOUNDARY, the function returns NULL. Otherwise it will return an array of hyper surfaces. The function also assign num_hs, the total number of hyper surfaces in the returned set. More...
 
IMPORT void FT_PromptSetMixedTypeBoundary2d (char *in_name, Front *front)
 This function sets the node types and curve wave types of nodes and curves at the MIXED_TYPE_BOUNDARY side(s). This function is only for 2D. It will prompt and set nodes in ascending order and then curves between the node pairs. More...
 
IMPORT int FT_Dimension ()
 This function returns the spatial dimension of current run.
 
IMPORT int FT_RectBoundaryType (Front *front, int dir, int side)
 This function returns rectangular boundary type in direction dir and side side. More...
 
IMPORT boolean FT_FrontContainWaveType (Front *front, int w_type)
 This function tells if front contains hyper-surface of certain type. It returns YES if it does, NO if it does not contain. More...
 
IMPORT HYPER_SURF * FT_RectBoundaryHypSurf (INTERFACE *intfc, int wave_type, int dir, int side)
 This function looks for a boundary hyper surface (curve in 2D and surface in 3D) at the rectangular domain border with the given wave type of the hyper surface in direction dir and side side. Returns NULL if no match is found and return pointer to the hyper surface if found. More...
 
IMPORT HYPER_SURF ** FT_InteriorHypSurfs (INTERFACE *intfc, int wave_type, int *num_hs)
 This function looks for interior hyper surfaces (curve in 2D and surface in 3D) with the given wave type. Returns NULL if no match is found and return pointer to an array of hyper surfaces if found. This function allocates the memory for the array of pointers to hyper surfaces. More...
 
IMPORT void FT_ParallelExchIntfcBuffer (Front *front)
 This is a root level parallel communication function for front interface geometry and states. It will cut the old buffer parts of the interface and patch it with new parts received from other subdomains or periodically shifted sides. This is a synchronous function and must be called synchronously by every processor. More...
 
IMPORT void FT_ParallelExchGridArrayBuffer (double *grid_array, Front *front, int *symmetry)
 This is a parallel communication function for a double array on the expanded dual grid of the grid_intfc in front. It will cut the old buffer parts of the array and patch it with new buffer parts received from other subdomains or periodically shifted sides. This is a synchronous function and must be called synchronously by every processor. More...
 
IMPORT void FT_ParallelExchGridIntArrayBuffer (int *iarray, Front *front)
 This is a parallel communication function for a integer array on the expanded dual grid of the grid_intfc in front. It will cut the old buffer parts of the array and patch it with new buffer parts received from other subdomains or periodically shifted sides. This is a synchronous function and must be called synchronously by every processor. More...
 
IMPORT void FT_ParallelExchGridVectorArrayBuffer (double **vec_grid_array, Front *front)
 This is a parallel communication function for a vector double array on the expanded dual grid of the grid_intfc in front. It will cut the old buffer parts of the array and patch it with new buffer parts received from other subdomains or periodically shifted sides. This is a synchronous function and must be called synchronously by every processor. Reflection of vector at the RLECTION_BOUNDARY is considered. More...
 
IMPORT void FT_ParallelExchCompGridArrayBuffer (double *grid_array, Front *front, int *symmetry)
 This is a parallel communication function for a double array on the expanded comp grid of the grid_intfc in front. It will cut the old buffer parts of the array and patch it with new buffer parts received from other subdomains or periodically shifted sides. This is a synchronous function and must be called synchronously by every processor. More...
 
IMPORT void FT_ParallelExchCellIndex (Front *front, int *lbuf, int *ubuf, POINTER ijk_to_I)
 This is a parallel communication function for the cell index on the expanded dual grid of the grid_intfc in front. The cell index translate the nD (n=2,3) icoordinates to a one dimensional index sequence. The indices are parallely globalized. More...
 
IMPORT void FT_ParallelExchCompGridCellIndex (Front *front, int *lbuf, int *ubuf, POINTER ijk_to_I)
 This is a parallel communication function for the dual cell index on the expanded comp grid of the comp_grid_intfc in front. The cell index translate the nD (n=2,3) icoordinates to a one dimensional index sequence. The indices are parallely globalized. More...
 
IMPORT INTERFACE * FT_CollectHypersurfFromSubdomains (Front *front, int *owner, int w_type)
 This function collect pieces of surface given wave type and patch them to the surface in the owner subdomain. The owner subdomain is given by its rectangular subdomain index owner[] = (ix,iy,iz). The collected surfaces are packed in a interface structure. The function return the interface if the collection is succsessful. The function return NULL if the collection is not succsessful. More...
 
IMPORT void FT_GetStatesAtPoint (POINT *p, HYPER_SURF_ELEMENT *hse, HYPER_SURF *hs, POINTER *sl, POINTER *sr)
 This function retrieves the left and right states at a point. Since a point can be shared by different entities, the associated hyper surface element (bond in 2D and tri in 3D) and hyper surface (curve in 2D and surface in 3D) must be provided as input. More...
 
IMPORT void FT_ScalarMemoryAlloc (POINTER *a, int size)
 This function allocate the memory for a scalar. More...
 
IMPORT void FT_VectorMemoryAlloc (POINTER *a, int n1, int size)
 This function allocate the memory for a vector. More...
 
IMPORT void FT_MatrixMemoryAlloc (POINTER *a, int n1, int n2, int size)
 This function allocate the memory for a matrix. More...
 
IMPORT void FT_TriArrayMemoryAlloc (POINTER *a, int n1, int n2, int n3, int size)
 This function allocate the memory for a tri-array. More...
 
IMPORT void FT_QuadArrayMemoryAlloc (POINTER *a, int n1, int n2, int n3, int n4, int size)
 This function allocate the memory for a quad-array. More...
 
IMPORT void FT_QuinArrayMemoryAlloc (POINTER *a, int n1, int n2, int n3, int n4, int n5, int size)
 This function allocate the memory for a quin-array. More...
 
IMPORT void FT_SexArrayMemoryAlloc (POINTER *a, int n1, int n2, int n3, int n4, int n5, int n6, int size)
 This function allocate the memory for a sex-array. More...
 
IMPORT void FT_SetGlobalIndex (Front *front)
 This function set global index for all point in the interface. More...
 
IMPORT void FT_FreeThese (int n,...)
 This function free memory of items allocated by FT_...MemoryAlloc() functions. The number of arguments is flexible, but needs to equal to the input integer n. More...
 
IMPORT boolean FT_StateStructAtGridCrossing2 (Front *front, int *icoords, GRID_DIRECTION dir, int comp, POINTER *state, HYPER_SURF **hs, HYPER_SURF_ELEMENT **hse, double *crx_coords)
 
IMPORT double FT_ComputeTotalVolumeFraction (Front *front, COMPONENT comp_of_vol)
 This function compute the total volume fraction associated the comp_of_vol based on crossing in front->grid_intfc. It returns the total volume fraction as a double precision value. More...
 
IMPORT void FT_ComputeGridVolumeFraction (Front *front, COMPONENT comp_of_vol, POINTER *grid_vol_frac)
 This function compute the volume fraction on the expanded dual grid associated the comp_of_vol based on crossing in front->grid_intfc. It passes the address of the grid volume fraction to the pointer grid_vol_frac. In 2D, it is a double**, in 3D it is double***. More...
 
IMPORT void FT_CurveSegLengthConstr (CURVE *c, BOND *bs, BOND *be, int nb, double seg_length, REDISTRIBUTION_DIRECTION dir)
 This function set curve to constrained length (seg_length), starting from the bond bs and end at bond be, with total of nb bonds. The cutting direction is the input dir. More...
 
IMPORT CURVE * FT_MakeNodeArrayCurve (Front *front, int num_nodes, double **node_array, COMPONENT neg_comp, COMPONENT pos_comp, boolean is_closed_curve, double scale_factor, int w_type)
 
IMPORT CURVE * FT_MakePointArrayCurve (Front *front, int num_points, double **point_array, COMPONENT neg_comp, COMPONENT pos_comp, boolean is_closed_curve, int w_type)
 
IMPORT void FT_MakeEllipticSurf (Front *front, double *center, double *radius, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, int refinement_level, SURFACE **surf)
 This function inserts an elliptic surface into the front with given information of center, radii, components, and wave type. More...
 
IMPORT void FT_MakeDumbBellSurf (Front *front, double x0, double x1, double y0, double z0, double R, double r, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a dumbbell surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakeProjectileSurf (Front *front, double *center, double R, double r, double h, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a projectile surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT CURVE * FT_MakeParametricCurve (Front *front, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, boolean(*func)(POINTER, double, double *), POINTER func_params, int refinement_level, boolean is_closed)
 This function inserts a parametric curve into the front with given information of its function, parameters, components, and wave type. More...
 
IMPORT void FT_RotateSurface (SURFACE *surf, double *center, double phi, double theta)
 This function rorate surface with azimuthal angle theta and polar angle phi about the given center. More...
 
IMPORT void FT_MakeCuboidSurf (Front *front, double *center, double *edge, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a cuboid surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakeCylinderSurf (Front *front, double *center, double radius, double height, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a cylinder surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakeConeSurf (Front *front, double *center, double slope, double height, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a cone surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakeTetrahedronSurf (Front *front, double *center, double radius, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a tetrahedron surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakePlaneSurf (Front *front, double *plane_nor, double *plane_pt, boolean reset_bdry_comp, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a plane surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_InstallSurfEdge (SURFACE *surf, int hsbdry_type)
 This function install a curve at the surface boundary. The boundary of the surface should contain triangles with NULL side. More...
 
IMPORT void FT_CutSurfBdry (SURFACE *surf, boolean constr_func(POINTER, double *), POINTER func_params, double **insert_coords, int num_pts, int insert_idir)
 This function cut surface at the boundary defined by the constrain function. It first insert significant corner points to make sure it cuts the surface accurately. More...
 
IMPORT void FT_MakeEllipticCurve (Front *front, double *center, double *radius, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, int refinement_level, CURVE **curve)
 This function inserts an elliptic curve into the front with given information of center, radii, components, and wave type. More...
 
IMPORT void FT_PrintWaveType (int w_type)
 This function print wave type as a string. More...
 
IMPORT void FT_PrintBoundaryType (int dir, int side)
 This function print boundary type as a string. More...
 
IMPORT int FT_BoundaryType (int dir, int side)
 This function return boundary type as an enumerated. More...
 
IMPORT double * FT_GridIntfcTopL (Front *)
 This function return lower bounds of grid domain. More...
 
IMPORT double * FT_GridIntfcTopU (Front *)
 This function return upper bounds of grid domain. More...
 
IMPORT double * FT_GridIntfcToph (Front *)
 This function return grid spacing of grid domain. More...
 
IMPORT COMPONENT * FT_GridIntfcTopComp (Front *)
 This function return components of grid domain. More...
 
IMPORT COMPONENT * FT_GridIntfcTopGmax (Front *)
 This function return mesh sizes of grid domain. More...
 
IMPORT RECT_GRID * FT_GridIntfcTopGrid (Front *)
 This function return grid structure of grid domain. More...
 
IMPORT void FT_AddHdfMovieVariable (Front *front, boolean preset_bound, boolean untracked, COMPONENT obst_comp, const char *var_name, int idir, double *var_field, double(*getStateFunc)(POINTER), double max_var, double min_var)
 
IMPORT void FT_AddVtkVectorMovieVariable (Front *front, const char *var_name, double **var_field)
 Initialize a variable and information for vtk vector movie output. More...
 
IMPORT void FT_AddVtkScalarMovieVariable (Front *front, const char *var_name, double *var_field)
 Initialize a variable and information for vtk scalar movie output. More...
 
IMPORT void FT_AddVtkIntfcMovieVariable (Front *front, const char *var_name)
 Initialize the variable name for vtk interface movie output. More...
 
IMPORT void FT_ResetDomainAndGrid (Front *front, double *L, double *U, int *gmax)
 
IMPORT boolean FT_CoordsInSubdomain (Front *front, double *coords)
 This function test if the coordinate is inside my subdomain. The function returns YES if it is in, otherwise it returns NO. More...
 
IMPORT void FT_PrintTimeStamp (Front *front)
 Output time information including current time, step and predicted next time step size. More...
 
IMPORT void FT_MakeCrossCylinderSurf (Front *front, double *center1, double *center2, double radius1, double radius2, double height1, double height2, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a surface in which a cylinder go through another cylinder into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakeBowlSurf (Front *front, double *center, double radius1, double radius2, double radius3, double height1, double height2, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 This function inserts a bowl surface into the front with given information of its parameters, components, and wave type. More...
 
IMPORT void FT_MakePlatformSurf (Front *front, double *center, double radius, double height, double slope, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 
IMPORT void FT_MakeStellatedOctahedronSurf (Front *front, double *center, double edge, COMPONENT neg_comp, COMPONENT pos_comp, int w_type, SURFACE **surf)
 
IMPORT void FT_InitSurfVeloFunc (SURFACE *surf, POINTER vparams, int(*vfunc)(POINTER, Front *, POINT *, HYPER_SURF_ELEMENT *, HYPER_SURF *, double *))
 
IMPORT void FT_InitCurveVeloFunc (CURVE *curve, POINTER vparams, int(*vfunc)(POINTER, Front *, POINT *, HYPER_SURF_ELEMENT *, HYPER_SURF *, double *))
 
IMPORT void FT_InitNodeVeloFunc (NODE *node, POINTER vparams, int(*vfunc)(POINTER, Front *, POINT *, HYPER_SURF_ELEMENT *, HYPER_SURF *, double *))
 
IMPORT boolean FT_CheckSurfCompConsistency (Front *front, SURFACE *surf)
 This function check if the outer (positive) component of a surface is consistent with the ambient component. More...
 

Detailed Description

The fapi.h contains the functions used to operate the interface.

Function Documentation

void FT_CutSurfBdry ( SURFACE *  surf,
boolean   constr_funcPOINTER, double *,
POINTER  func_params,
double **  insert_coords,
int  num_pts,
int  insert_idir 
)

This function cut surface at the boundary defined by the constrain function. It first insert significant corner points to make sure it cuts the surface accurately.

Parameters
surfinout Pointer to the surface to be cut.
constr_funcin constrain function defining the boundary.
func_paramsin anonymous pointer of function parameters.
insert_coordsin coordinates of points to be inserted before cut.
num_ptsin number of points to be inserted before cut.
insert_idirin insert point along this direction.
void FT_InstallSurfEdge ( SURFACE *  surf,
int  hsbdry_type 
)

This function install a curve at the surface boundary. The boundary of the surface should contain triangles with NULL side.

Parameters
surfinout Pointer to the surface where edge to be installed.
hsbdry_typein Boundary curve type.
void FT_MakeBowlSurf ( Front *  front,
double *  center,
double  radius1,
double  radius2,
double  radius3,
double  height1,
double  height2,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a bowl surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center of the sphere.
radius1in radius of the outer sphere.
radius2in radius of the cylinder.
radius3in radius of the inner sphere.
height1in height of body of the bowl.
height2in height of bottom of the bowl.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_MakeConeSurf ( Front *  front,
double *  center,
double  slope,
double  height,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a cone surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin vertex of the cone.
slopein slope of the cone.
heightin height of the cone.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_MakeCrossCylinderSurf ( Front *  front,
double *  center1,
double *  center2,
double  radius1,
double  radius2,
double  height1,
double  height2,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a surface in which a cylinder go through another cylinder into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
center1in center of one cylinder.
center2in center of another cylinder.
radius1in radius of one cylinder.
radius2in radius of another cylinder.
height1in height of one cylinder.
height2in height of another cylinder.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_MakeCuboidSurf ( Front *  front,
double *  center,
double *  edge,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a cuboid surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center of the cuboid.
edgein edge of the cuboid.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_MakeCylinderSurf ( Front *  front,
double *  center,
double  radius,
double  height,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a cylinder surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center of the cylinder.
edgein radius of the cylinder.
heightin height of the cylinder.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_MakeDumbBellSurf ( Front *  front,
double  x0,
double  x1,
double  y0,
double  z0,
double  R,
double  r,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a dumbbell surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
x0in x-coordinate of the left center.
x0in x-coordinate of the right center.
y0in y-coordinate of the axis.
z0in z-coordinate of the axis.
Rin radius of the two end spheres.
rin radius cylinder connecting the two spheres.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function.
void FT_MakeEllipticCurve ( Front *  front,
double *  center,
double *  radius,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
int  refinement_level,
CURVE **  curve 
)

This function inserts an elliptic curve into the front with given information of center, radii, components, and wave type.

Parameters
frontinout Pointer to the front in which curve is inserted.
centerin center of the ellipse.
radiusin radii of the ellipse.
neg_compin index for negative side of the curve (inner side).
pos_compin index for positive side of the curve (outer side).
w_typeint wave type of the curve.
refinement_levelint level of refinement.
curveout curve made by this function.
void FT_MakeEllipticSurf ( Front *  front,
double *  center,
double *  radius,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
int  refinement_level,
SURFACE **  surf 
)

This function inserts an elliptic surface into the front with given information of center, radii, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center of the ellipsoid.
radiusin radii of the ellipsoid.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typeint wave type of the surface.
refinement_levelint refinement level of the surface.
surfout surface made by this function.
CURVE * FT_MakeParametricCurve ( Front *  front,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
boolean(*)(POINTER, double, double *)  func,
POINTER  func_params,
int  refinement_level,
boolean  is_closed 
)

This function inserts a parametric curve into the front with given information of its function, parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which curve is inserted.
neg_compin index for negative side of the curve (left side).
pos_compin index for positive side of the curve (right side).
w_typein wave type of the curve.
funcin parametric function for the curve.
func_paramsin anonymous pointer to parameters for the function.
refinement_levelin refinement level on grid spacing.
is_closedin whether the curve is closed.
void FT_MakePlaneSurf ( Front *  front,
double *  plane_nor,
double *  plane_pt,
boolean  reset_bdry_comp,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a plane surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
plane_norin normal vector of the plane.
plane_ptin coordinates of a point on the plane.
reset_bdry_compin if YES, reset boundary component.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function.
void FT_MakePlatformSurf ( Front *  front,
double *  center,
double  radius,
double  height,
double  slope,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a platform surface derived by cutting the head of the cone curface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center of the bottom of the platform.
radiusin radius of the bottom of the platform.
heightin height of the platform.
slopein slope of the platform.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_MakeProjectileSurf ( Front *  front,
double *  center,
double  R,
double  r,
double  h,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a projectile surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center-coordinate of the projectile.
Rin cylindrical radius of the projectile.
rin head height of the projectile.
hin butt height of the projectile.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function.
void FT_MakeTetrahedronSurf ( Front *  front,
double *  center,
double  radius,
COMPONENT  neg_comp,
COMPONENT  pos_comp,
int  w_type,
SURFACE **  surf 
)

This function inserts a tetrahedron surface into the front with given information of its parameters, components, and wave type.

Parameters
frontinout Pointer to the front in which surface is inserted.
centerin center of the tetrahedron.
radiusin circumcircle of the tetrahedron.
neg_compin index for negative side of the surface (inner side).
pos_compin index for positive side of the surface (outer side).
w_typein wave type of the surface.
surfout surface made by this function
void FT_RotateSurface ( SURFACE *  surf,
double *  center,
double  phi,
double  theta 
)

This function rorate surface with azimuthal angle theta and polar angle phi about the given center.

Parameters
surfinout Pointer to the surface to be rotated.
centerin center-coordinate for the rotation.
phiin polar angle of the rotation.
thetain azimuthal angle of the rotation.