inciter::KozCG class

KozCG Charm++ chare array used to advance PDEs in time with KozCG.

Public static functions

static void registerReducers()
Configure Charm++ custom reduction types initiated from this chare array.

Constructors, destructors, conversion operators

KozCG(CkMigrateMessage* m) explicit
Migrate constructor.

Public functions

auto KozCG(const CProxy_Discretization& disc, const std::map<int, std::vector<std::size_t>>& bface, const std::map<int, std::vector<std::size_t>>& bnode, const std::vector<std::size_t>& triinpoel) -> KozCG_SDAG_CODE explicit
Constructor.
void ResumeFromSync() override
Return from migration.
void setup(tk::real v)
Start setup for solution.
void advance(tk::real newdt)
Advance equations to next time step.
void feop()
Start (re-)computing domain and boundary integrals.
void comnorm(const std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>>& inbnd)
Receive contributions to boundary point normals on chare-boundaries.
void comrhs(const std::unordered_map<std::size_t, std::vector<tk::real>>& inrhs)
Receive contributions to right-hand side vector on chare-boundaries.
void comaec(const std::unordered_map<std::size_t, std::vector<tk::real>>& inaec)
Receive antidiffusive and low-order contributions on chare-boundaries.
void comalw(const std::unordered_map<std::size_t, std::vector<tk::real>>& inalw)
Receive allowed limits contributions on chare-boundaries.
void comlim(const std::unordered_map<std::size_t, std::vector<tk::real>>& inlim)
Receive limited antidiffusive contributions on chare-boundaries.
void evalres(const std::vector<tk::real>& l2res)
Evaluate residuals.
void resizePostAMR(const std::vector<std::size_t>& ginpoel, const tk::UnsMesh::Chunk& chunk, const tk::UnsMesh::Coords& coord, const std::unordered_map<std::size_t, tk::UnsMesh::Edge>& addedNodes, const std::unordered_map<std::size_t, std::size_t>& addedTets, const std::set<std::size_t>& removedNodes, const std::unordered_map<int, std::unordered_set<std::size_t>>& nodeCommMap, const std::map<int, std::vector<std::size_t>>& bface, const std::map<int, std::vector<std::size_t>>& bnode, const std::vector<std::size_t>& triinpoel)
Receive new mesh from Refiner.
auto solution() const -> const tk::Fields&
void integrals()
Compute integral quantities for output.
void step()
Evaluate whether to continue with next time step.
void evalLB(int nrestart)
Evaluate whether to do load balancing.

Charm++ pack/unpack serializer member functions

CProxy_Discretization m_disc
Discretization proxy.
std::size_t m_nrhs
Counter for right-hand side vector nodes updated.
std::size_t m_nnorm
Counter for receiving boundary point normals.
std::size_t m_naec
Counter for receiving antidiffusive contributions.
std::size_t m_nalw
Counter for receiving allowed limits.
std::size_t m_nlim
Counter for receiving limited antidiffusive contributions.
std::map<int, std::vector<std::size_t>> m_bnode
Boundary node lists mapped to side set ids used in the input file.
std::map<int, std::vector<std::size_t>> m_bface
Boundary face lists mapped to side set ids used in the input file.
std::vector<std::size_t> m_triinpoel
Boundary triangle face connecitivity where BCs are set by user.
tk::Fields m_u
Unknown/solution vector at mesh nodes.
tk::Fields m_p
Max/min antidiffusive edge contributions at mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_pc
tk::Fields m_q
Max/min allowed limits at mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_qc
tk::Fields m_a
Limited antidiffusive contributions at mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_ac
tk::Fields m_rhs
Right-hand side vector (for the high order system)
std::unordered_map<std::size_t, std::vector<tk::real>> m_rhsc
NodeDiagnostics m_diag
Diagnostics object.
std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>> m_bnorm
decltype(m_bnorm) m_bnormc
std::unordered_map<std::size_t, std::array<tk::real, 3>> m_bndpoinint
std::vector<std::size_t> m_dirbcmasks
Nodes and their Dirichlet BC masks.
std::vector<std::size_t> m_prebcnodes
Nodes at pressure BCs.
std::vector<tk::real> m_prebcvals
Density and pressure values at pressure BCs.
std::set<std::size_t> m_symbcnodeset
Unique set of ordered nodes at which symmetry BCs are set.
std::vector<std::size_t> m_symbcnodes
Streamable nodes at which symmetry BCs are set.
std::vector<tk::real> m_symbcnorms
Streamable normals at nodes at which symmetry BCs are set.
std::set<std::size_t> m_farbcnodeset
Unique set of ordered nodes at which farfield BCs are set.
std::vector<std::size_t> m_farbcnodes
Streamable nodes at which farfield BCs are set.
std::vector<tk::real> m_farbcnorms
Streamable normals at nodes at which farfield BCs are set.
std::map<int, std::pair<std::vector<std::size_t>, std::vector<tk::real>>> m_surfint
Streamable surface integral nodes and normals * dA on surfaces.
std::vector<tk::real> m_dtp
Time step size for each mesh node.
std::vector<tk::real> m_tp
Physical time for each mesh node.
int m_finished
True in the last time step.
tk::real m_freezeflow
dt multiplier after flow no longer updated
void pup(PUP::er& p) override
Pack/Unpack serialize member function.
void operator|(PUP::er& p, KozCG& i)
Pack/Unpack serialize operator|.

Function documentation

static void inciter::KozCG::registerReducers()

Configure Charm++ custom reduction types initiated from this chare array.

Since this is a [initnode] routine, the runtime system executes the routine exactly once on every logical node early on in the Charm++ init sequence. Must be static as it is called without an object. See also: Section "Initializations at Program Startup" at in the Charm++ manual http://charm.cs.illinois.edu/manuals/html/charm++/manual.html.

KozCG_SDAG_CODE inciter::KozCG::KozCG(const CProxy_Discretization& disc, const std::map<int, std::vector<std::size_t>>& bface, const std::map<int, std::vector<std::size_t>>& bnode, const std::vector<std::size_t>& triinpoel) explicit

Constructor.

Parameters
disc in Discretization proxy
bface in Boundary-faces mapped to side sets used in the input file
bnode in Boundary-node lists mapped to side sets used in input file
triinpoel in Boundary-face connectivity where BCs set (global ids)

void inciter::KozCG::ResumeFromSync() override

Return from migration.

This is called when load balancing (LB) completes. The presence of this function does not affect whether or not we block on LB.

void inciter::KozCG::setup(tk::real v)

Start setup for solution.

Parameters
in Total volume within user-specified box

void inciter::KozCG::advance(tk::real newdt)

Advance equations to next time step.

Parameters
newdt in The smallest dt across the whole problem

void inciter::KozCG::comnorm(const std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>>& inbnd)

Receive contributions to boundary point normals on chare-boundaries.

Parameters
inbnd in Incoming partial sums of boundary point normals

void inciter::KozCG::comrhs(const std::unordered_map<std::size_t, std::vector<tk::real>>& inrhs)

Receive contributions to right-hand side vector on chare-boundaries.

Parameters
inrhs in Partial contributions of RHS to chare-boundary nodes. Key: global mesh node IDs, value: contributions for all scalar components.

void inciter::KozCG::comaec(const std::unordered_map<std::size_t, std::vector<tk::real>>& inaec)

Receive antidiffusive and low-order contributions on chare-boundaries.

Parameters
inaec in Partial contributions of antidiffusive edge and low-order solution contributions on chare-boundary nodes. Key: global mesh node IDs, value: 0: antidiffusive contributions, 1: low-order solution.

void inciter::KozCG::comalw(const std::unordered_map<std::size_t, std::vector<tk::real>>& inalw)

Receive allowed limits contributions on chare-boundaries.

Parameters
inalw in Partial contributions of allowed limits contributions on chare-boundary nodes. Key: global mesh node IDs, value: allowed limit contributions.

void inciter::KozCG::comlim(const std::unordered_map<std::size_t, std::vector<tk::real>>& inlim)

Receive limited antidiffusive contributions on chare-boundaries.

Parameters
inlim in Partial contributions of limited contributions on chare-boundary nodes. Key: global mesh node IDs, value: limited contributions.

void inciter::KozCG::evalres(const std::vector<tk::real>& l2res)

Evaluate residuals.

Parameters
l2res in L2-norms of the residual for each scalar component computed across the whole problem

void inciter::KozCG::resizePostAMR(const std::vector<std::size_t>& ginpoel, const tk::UnsMesh::Chunk& chunk, const tk::UnsMesh::Coords& coord, const std::unordered_map<std::size_t, tk::UnsMesh::Edge>& addedNodes, const std::unordered_map<std::size_t, std::size_t>& addedTets, const std::set<std::size_t>& removedNodes, const std::unordered_map<int, std::unordered_set<std::size_t>>& nodeCommMap, const std::map<int, std::vector<std::size_t>>& bface, const std::map<int, std::vector<std::size_t>>& bnode, const std::vector<std::size_t>& triinpoel)

Receive new mesh from Refiner.

Parameters
ginpoel in Mesh connectivity with global node ids
chunk in New mesh chunk (connectivity and global<->local id maps)
coord in New mesh node coordinates
addedNodes in Newly added mesh nodes and their parents (local ids)
addedTets in Newly added mesh cells and their parents (local ids)
removedNodes in Newly removed mesh node local ids
nodeCommMap in New node communication map
bface in Boundary-faces mapped to side set ids
bnode in Boundary-node lists mapped to side set ids
triinpoel in Boundary-face connectivity

const tk::Fields& inciter::KozCG::solution() const

Returns Const-ref to current solution

Const-ref access to current solution

void inciter::KozCG::evalLB(int nrestart)

Evaluate whether to do load balancing.

Parameters
nrestart in Number of times restarted

void inciter::KozCG::pup(PUP::er& p) override

Pack/Unpack serialize member function.

Parameters
in/out Charm++'s PUP::er serializer object reference

void inciter::KozCG::operator|(PUP::er& p, KozCG& i)

Pack/Unpack serialize operator|.

Parameters
in/out Charm++'s PUP::er serializer object reference
in/out KozCG object reference

Variable documentation

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::KozCG::m_pc

Receive buffer for max/min antidiffusive edge contributions

Key: global node id, value: max/min antidiff edge contributions in nodes.

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::KozCG::m_qc

Receive buffer for max/min allowed limits

Key: global node id, value: max/min allowed limits in nodes.

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::KozCG::m_ac

Receive buffer for limited antidiffusive contributions

Key: global node id, value: limited antidiffusive contributions in nodes.

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::KozCG::m_rhsc

Receive buffer for communication of the right hand side

Key: global node id, value: rhs for all scalar components per node.

std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>> inciter::KozCG::m_bnorm

Boundary point normals

Outer key: side set id. Inner key: global node id of boundary point, value: weighted normals, inverse distance square, nodal area.

decltype(m_bnorm) inciter::KozCG::m_bnormc

Boundary point normals receive buffer

Outer key: side set id. Inner key: global node id of boundary point, value: weighted normals and inverse distance square.

std::unordered_map<std::size_t, std::array<tk::real, 3>> inciter::KozCG::m_bndpoinint

Boundary point integrals

Key: global node id of boundary point, value: boundary point integral contributions.