REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORAFillPatcher Class Reference

#include <REMORA_FillPatcher.H>

Collaboration diagram for REMORAFillPatcher:

Public Member Functions

 REMORAFillPatcher (amrex::BoxArray const &fba, amrex::DistributionMapping const &fdm, amrex::Geometry const &fgeom, amrex::BoxArray const &cba, amrex::DistributionMapping const &cdm, amrex::Geometry const &cgeom, int nghost, int nghost_set, int ncomp, amrex::InterpBase *interp)
 Fill valid and ghost data with the "state data" at the given time.
 
void Define (amrex::BoxArray const &fba, amrex::DistributionMapping const &fdm, amrex::Geometry const &fgeom, amrex::BoxArray const &cba, amrex::DistributionMapping const &cdm, amrex::Geometry const &cgeom, int nghost, int nghost_set, int ncomp, amrex::InterpBase *interp)
 Redefine the coarse and fine patch MultiFabs.
 
void BuildMask (amrex::BoxArray const &fba, int nghost, int nghost_set)
 Generate masking array.
 
void RegisterCoarseData (amrex::Vector< amrex::MultiFab const * > const &crse_data, amrex::Vector< amrex::Real > const &crse_time)
 Register the coarse data to be used by the REMORAFillPatcher.
 
void InterpFace (amrex::MultiFab &fine, amrex::MultiFab const &crse, int mask_val)
 Interpolate to cell faces.
 
void InterpCell (amrex::MultiFab &fine, amrex::MultiFab const &crse, amrex::Vector< amrex::BCRec > const &bcr, int mask_val)
 Interpolate to cell centers.
 
int GetSetMaskVal ()
 
int GetRelaxMaskVal ()
 
amrex::iMultiFab * GetMask ()
 
template<typename BC >
void FillSet (amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
 Fill fine data in the set region.
 
template<typename BC >
void FillRelax (amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
 Fill fine data in the relax region.
 
template<typename BC >
void Fill (amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs, int mask_val)
 Fill fine data in the relax region.
 

Private Attributes

amrex::BoxArray m_fba
 
amrex::BoxArray m_cba
 
amrex::DistributionMapping m_fdm
 
amrex::DistributionMapping m_cdm
 
amrex::Geometry m_fgeom
 
amrex::Geometry m_cgeom
 
int m_nghost
 
int m_nghost_subset
 
int m_ncomp
 
amrex::InterpBase * m_interp
 
amrex::IntVect m_ratio
 
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_old
 
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_new
 
std::unique_ptr< amrex::iMultiFab > m_cf_mask
 
amrex::Vector< amrex::Real > m_crse_times
 
amrex::Real m_dt_crse
 
int m_set_mask {2}
 
int m_relax_mask {1}
 

Detailed Description

Definition at line 8 of file REMORA_FillPatcher.H.

Constructor & Destructor Documentation

◆ REMORAFillPatcher()

REMORAFillPatcher::REMORAFillPatcher ( amrex::BoxArray const &  fba,
amrex::DistributionMapping const &  fdm,
amrex::Geometry const &  fgeom,
amrex::BoxArray const &  cba,
amrex::DistributionMapping const &  cdm,
amrex::Geometry const &  cgeom,
int  nghost,
int  nghost_set,
int  ncomp,
amrex::InterpBase *  interp 
)

Fill valid and ghost data with the "state data" at the given time.

Parameters
[in]fbaBoxArray of data to be filled at fine level
[in]fdmDistributionMapping of data to be filled at fine level
[in]fgeomcontainer of geometry information at fine level
[in]cbaBoxArray of data to be filled at coarse level
[in]cdmDistributionMapping of data to be filled at coarse level
[in]cgeomcontainer of geometry information at coarse level
[in]nghostnumber of ghost cells to be filled
[in]ncompnumber of components to be filled
[in]interpinterpolation operator to be used

Definition at line 19 of file REMORA_FillPatcher.cpp.

Here is the call graph for this function:

Member Function Documentation

◆ BuildMask()

void REMORAFillPatcher::BuildMask ( amrex::BoxArray const &  fba,
int  nghost,
int  nghost_set 
)

Generate masking array.

Definition at line 115 of file REMORA_FillPatcher.cpp.

Referenced by Define().

Here is the caller graph for this function:

◆ Define()

void REMORAFillPatcher::Define ( amrex::BoxArray const &  fba,
amrex::DistributionMapping const &  fdm,
amrex::Geometry const &  fgeom,
amrex::BoxArray const &  cba,
amrex::DistributionMapping const &  cdm,
amrex::Geometry const &  cgeom,
int  nghost,
int  nghost_set,
int  ncomp,
amrex::InterpBase *  interp 
)

Redefine the coarse and fine patch MultiFabs.

Parameters
[in]fbaBoxArray of data to be filled at fine level
[in]fdmDistributionMapping of data to be filled at fine level
[in]fgeomcontainer of geometry information at fine level
[in]cbaBoxArray of data to be filled at coarse level
[in]cdmDistributionMapping of data to be filled at coarse level
[in]cgeomcontainer of geometry information at coarse level
[in]nghostnumber of ghost cells to be filled
[in]ncompnumber of components to be filled
[in]interpinterpolation operator to be used

Definition at line 53 of file REMORA_FillPatcher.cpp.

Referenced by REMORAFillPatcher().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Fill()

template<typename BC >
void REMORAFillPatcher::Fill ( amrex::MultiFab &  mf,
amrex::Real  time,
BC &  cbc,
amrex::Vector< amrex::BCRec > const &  bcs,
int  mask_val 
)

Fill fine data in the relax region.

Parameters
[out]mfMultiFab to be filled
[in]timeTime at which to fill data
[in]cbcCoarse boundary condition
[in]bcsVector of boundary conditions
[in]mask_valValue to assign mask array

Definition at line 124 of file REMORA_FillPatcher.H.

Referenced by FillRelax(), and FillSet().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ FillRelax()

template<typename BC >
void REMORAFillPatcher::FillRelax ( amrex::MultiFab &  mf,
amrex::Real  time,
BC &  cbc,
amrex::Vector< amrex::BCRec > const &  bcs 
)

Fill fine data in the relax region.

Parameters
[out]mfMultiFab to be filled
[in]timeTime at which to fill data
[in]cbcCoarse boundary condition
[in]bcsVector of boundary conditions

Definition at line 109 of file REMORA_FillPatcher.H.

Here is the call graph for this function:

◆ FillSet()

template<typename BC >
void REMORAFillPatcher::FillSet ( amrex::MultiFab &  mf,
amrex::Real  time,
BC &  cbc,
amrex::Vector< amrex::BCRec > const &  bcs 
)

Fill fine data in the set region.

Parameters
[out]mfMultiFab to be filled
[in]timeTime at which to fill data
[in]cbcCoarse boundary condition
[in]bcsVector of boundary conditions

Definition at line 95 of file REMORA_FillPatcher.H.

Here is the call graph for this function:

◆ GetMask()

amrex::iMultiFab * REMORAFillPatcher::GetMask ( )
inline

Definition at line 49 of file REMORA_FillPatcher.H.

◆ GetRelaxMaskVal()

int REMORAFillPatcher::GetRelaxMaskVal ( )
inline

Definition at line 47 of file REMORA_FillPatcher.H.

◆ GetSetMaskVal()

int REMORAFillPatcher::GetSetMaskVal ( )
inline

Definition at line 45 of file REMORA_FillPatcher.H.

◆ InterpCell()

void REMORAFillPatcher::InterpCell ( amrex::MultiFab &  fine,
amrex::MultiFab const &  crse,
amrex::Vector< amrex::BCRec > const &  bcr,
int  mask_val 
)

Interpolate to cell centers.

Parameters
[in,out]finefine level data
[in]crsecoarse level data
[in]bcrboundary condition type
[in]mask_valmasked value

Definition at line 333 of file REMORA_FillPatcher.cpp.

Referenced by Fill().

Here is the caller graph for this function:

◆ InterpFace()

void REMORAFillPatcher::InterpFace ( amrex::MultiFab &  fine,
amrex::MultiFab const &  crse,
int  mask_val 
)

Interpolate to cell faces.

Parameters
[in,out]finefine level data
[in]crsecoarse level data
[in]mask_valmasked value

Definition at line 209 of file REMORA_FillPatcher.cpp.

Referenced by Fill().

Here is the caller graph for this function:

◆ RegisterCoarseData()

void REMORAFillPatcher::RegisterCoarseData ( amrex::Vector< amrex::MultiFab const * > const &  crse_data,
amrex::Vector< amrex::Real > const &  crse_time 
)

Register the coarse data to be used by the REMORAFillPatcher.

Definition at line 180 of file REMORA_FillPatcher.cpp.

Member Data Documentation

◆ m_cba

amrex::BoxArray REMORAFillPatcher::m_cba
private

Definition at line 68 of file REMORA_FillPatcher.H.

Referenced by Define().

◆ m_cdm

amrex::DistributionMapping REMORAFillPatcher::m_cdm
private

Definition at line 70 of file REMORA_FillPatcher.H.

Referenced by Define().

◆ m_cf_crse_data_new

std::unique_ptr<amrex::MultiFab> REMORAFillPatcher::m_cf_crse_data_new
private

Definition at line 79 of file REMORA_FillPatcher.H.

Referenced by Define(), Fill(), and RegisterCoarseData().

◆ m_cf_crse_data_old

std::unique_ptr<amrex::MultiFab> REMORAFillPatcher::m_cf_crse_data_old
private

Definition at line 78 of file REMORA_FillPatcher.H.

Referenced by Define(), Fill(), and RegisterCoarseData().

◆ m_cf_mask

std::unique_ptr<amrex::iMultiFab> REMORAFillPatcher::m_cf_mask
private

Definition at line 80 of file REMORA_FillPatcher.H.

Referenced by BuildMask(), Define(), GetMask(), InterpCell(), and InterpFace().

◆ m_cgeom

amrex::Geometry REMORAFillPatcher::m_cgeom
private

Definition at line 72 of file REMORA_FillPatcher.H.

Referenced by Define(), InterpCell(), InterpFace(), and RegisterCoarseData().

◆ m_crse_times

amrex::Vector<amrex::Real> REMORAFillPatcher::m_crse_times
private

Definition at line 81 of file REMORA_FillPatcher.H.

Referenced by Fill(), RegisterCoarseData(), and REMORAFillPatcher().

◆ m_dt_crse

amrex::Real REMORAFillPatcher::m_dt_crse
private

Definition at line 82 of file REMORA_FillPatcher.H.

Referenced by Fill(), and RegisterCoarseData().

◆ m_fba

amrex::BoxArray REMORAFillPatcher::m_fba
private

Definition at line 67 of file REMORA_FillPatcher.H.

Referenced by Define().

◆ m_fdm

amrex::DistributionMapping REMORAFillPatcher::m_fdm
private

Definition at line 69 of file REMORA_FillPatcher.H.

Referenced by Define().

◆ m_fgeom

amrex::Geometry REMORAFillPatcher::m_fgeom
private

Definition at line 71 of file REMORA_FillPatcher.H.

Referenced by Define(), and InterpCell().

◆ m_interp

amrex::InterpBase* REMORAFillPatcher::m_interp
private

Definition at line 76 of file REMORA_FillPatcher.H.

Referenced by Define(), and InterpCell().

◆ m_ncomp

int REMORAFillPatcher::m_ncomp
private

Definition at line 75 of file REMORA_FillPatcher.H.

Referenced by Define(), Fill(), InterpCell(), InterpFace(), and RegisterCoarseData().

◆ m_nghost

int REMORAFillPatcher::m_nghost
private

Definition at line 73 of file REMORA_FillPatcher.H.

Referenced by Define().

◆ m_nghost_subset

int REMORAFillPatcher::m_nghost_subset
private

Definition at line 74 of file REMORA_FillPatcher.H.

Referenced by Define().

◆ m_ratio

amrex::IntVect REMORAFillPatcher::m_ratio
private

Definition at line 77 of file REMORA_FillPatcher.H.

Referenced by Define(), InterpCell(), and InterpFace().

◆ m_relax_mask

int REMORAFillPatcher::m_relax_mask {1}
private

Definition at line 84 of file REMORA_FillPatcher.H.

Referenced by Define(), FillRelax(), and GetRelaxMaskVal().

◆ m_set_mask

int REMORAFillPatcher::m_set_mask {2}
private

Definition at line 83 of file REMORA_FillPatcher.H.

Referenced by Define(), FillSet(), and GetSetMaskVal().


The documentation for this class was generated from the following files: