![]() |
Block-Structured AMR Software Framework
|
#include <AMReX_Config.H>#include <AMReX_MultiFab.H>#include <AMReX_Array.H>#include <AMReX_Vector.H>#include <AMReX_BCRec.H>Go to the source code of this file.
Namespaces | |
| namespace | amrex |
Functions | |
| void | amrex::EB_set_covered (MultiFab &mf, Real val) |
Fill all covered cells with a single value val. | |
| void | amrex::EB_set_covered (MultiFab &mf, int icomp, int ncomp, int ngrow, Real val) |
Fill ncomp components (starting at icomp) of covered cells with val. | |
| void | amrex::EB_set_covered (MultiFab &mf, int icomp, int ncomp, const Vector< Real > &vals) |
Fill covered cells with per-component constants provided in vals. | |
| void | amrex::EB_set_covered (MultiFab &mf, int icomp, int ncomp, int ngrow, const Vector< Real > &vals) |
Fill covered cells (including ngrow ghost cells) with per-component values. | |
| void | amrex::EB_set_covered_faces (const Array< MultiFab *, 3 > &umac, Real val) |
Fill covered faces in all directions with constant val. | |
| void | amrex::EB_set_covered_faces (const Array< MultiFab *, 3 > &umac, int scomp, int ncomp, const Vector< Real > &vals) |
| Fill covered faces with per-component values. | |
| void | amrex::EB_average_down (const MultiFab &S_fine, MultiFab &S_crse, const MultiFab &vol_fine, const MultiFab &vfrac_fine, int scomp, int ncomp, const IntVect &ratio) |
| Volume-weighted average-down from fine to coarse using EB volume fractions. | |
| void | amrex::EB_average_down (const MultiFab &S_fine, MultiFab &S_crse, int scomp, int ncomp, int ratio) |
| Convenience average-down that assumes unit volumes (with EB volume fraction weights). | |
| void | amrex::EB_average_down (const MultiFab &S_fine, MultiFab &S_crse, int scomp, int ncomp, const IntVect &ratio) |
| Same as above but with anisotropic ratio. | |
| void | amrex::EB_average_down_faces (const Array< const MultiFab *, 3 > &fine, const Array< MultiFab *, 3 > &crse, int ratio, int ngcrse) |
| Average face-centered data from fine to coarse grids weighted by area fractions. | |
| void | amrex::EB_average_down_faces (const Array< const MultiFab *, 3 > &fine, const Array< MultiFab *, 3 > &crse, const IntVect &ratio, int ngcrse) |
| Variant that accepts anisotropic refinement ratio. | |
| void | amrex::EB_average_down_faces (const Array< const MultiFab *, 3 > &fine, const Array< MultiFab *, 3 > &crse, const IntVect &ratio, const Geometry &crse_geom) |
Average face-centered data weighted by area fractions while respecting periodicity in crse_geom. | |
| void | amrex::EB_average_down_boundaries (const MultiFab &fine, MultiFab &crse, int ratio, int ngcrse) |
| Average EB boundary data from fine to coarse faces. | |
| void | amrex::EB_average_down_boundaries (const MultiFab &fine, MultiFab &crse, const IntVect &ratio, int ngcrse) |
| Anisotropic refinement version of EB_average_down_boundaries(). | |
| void | amrex::EB_computeDivergence (MultiFab &divu, const Array< MultiFab const *, 3 > &umac, const Geometry &geom, bool already_on_centroids) |
Compute divergence of face-centered velocity umac. | |
| void | amrex::EB_computeDivergence (MultiFab &divu, const Array< MultiFab const *, 3 > &umac, const Geometry &geom, bool already_on_centroids, const MultiFab &vel_eb) |
Compute divergence including EB-normal velocities stored in vel_eb. | |
| void | amrex::EB_average_face_to_cellcenter (MultiFab &ccmf, int dcomp, const Array< MultiFab const *, 3 > &fmf) |
| Average face-centered values to cell centers. | |
| void | amrex::EB_interp_CC_to_Centroid (MultiFab ¢, const MultiFab &cc, int scomp, int dcomp, int ncomp, const Geometry &geom) |
| Interpolate cell-centered data to cell centroids. | |
| void | amrex::EB_interp_CC_to_FaceCentroid (const MultiFab &cc, MultiFab &fc_x, MultiFab &fc_y, MultiFab &fc_z, int scomp, int dcomp, int nc, const Geometry &geom, const amrex::Vector< amrex::BCRec > &a_bcs) |
| Interpolate cell-centered data to face centroids for each direction. | |
| void | amrex::EB_interp_CellCentroid_to_FaceCentroid (const MultiFab &phi_centroid, const Array< MultiFab *, 3 > &phi_faces, int scomp, int dcomp, int nc, const Geometry &geom, const amrex::Vector< amrex::BCRec > &a_bcs) |
| Interpolate cell centroids to face centroids (array overload). | |
| void | amrex::EB_interp_CellCentroid_to_FaceCentroid (const MultiFab &phi_centroid, const Vector< MultiFab * > &phi_faces, int scomp, int dcomp, int nc, const Geometry &geom, const amrex::Vector< amrex::BCRec > &a_bcs) |
| Same as above but using Vector container for faces. | |
| void | amrex::EB_interp_CellCentroid_to_FaceCentroid (const MultiFab &phi_centroid, MultiFab &phi_xface, MultiFab &phi_yface, MultiFab &phi_zface, int scomp, int dcomp, int nc, const Geometry &geom, const amrex::Vector< amrex::BCRec > &a_bcs) |
| Component-wise overload taking explicit face references. | |
Utility routines for manipulating MultiFabs that carry EB data.