SNOWPACK  SNOWPACK-3.6.0
SnowStation Class Reference

Detailed Description

Station data including all information on snowpack layers (elements and nodes) and on canopy
This is the PRIMARY data structure of the SNOWPACK program
It is used extensively not only during the finite element solution but also to control the post-processing writes. It is initialized from SN_SNOWSOIL_DATA (at present).

#include <DataClasses.h>

Public Member Functions

 SnowStation (const bool &i_useCanopyModel=true, const bool &i_useSoilLayers=true)
 
 SnowStation (const SnowStation &c)
 
 ~SnowStation ()
 
SnowStationoperator= (const SnowStation &)
 Assignement operator. More...
 
void initialize (const SN_SNOWSOIL_DATA &SSdata, const size_t &i_sector)
 This routine initializes the snow cover structure which contains all information about a station including element, nodal and canopy data
Because you are working with layers, the first node is a special case; the rest can be generated in a loop .... The bottom temperature at the beginning of the computation is given by the temperature at the top of the lowest soil or snow layer
IMPORTANT: it is very important for Alpine3D that Cdata.height is initialized even if CANOPY = 0, otherwise SnowInterface will not recognize the canopy grids (David 2007-06-25). More...
 
void resize (const size_t &number_of_elements)
 Reallocate element and node data
Edata and Ndata as well as nElems and nNodes are reallocated or reset, respectively. More...
 
void reduceNumberOfElements (const size_t &rnE)
 Remove the upper "marked" element of two (snow only)

  1. Merging two elements:
    1. density is undefined
    2. take the uppermost node of both
  2. Removing melted or thin elements
    1. density is undefined AND length negative (*= -1.) as the latter will be used!
    2. keep upper node of lowest element.
More...
 
void combineElements (const size_t &number_top_elements, const bool &reduce_n_elements)
 If more than NUMBER_TOP_ELEMENTS snow elements exist, attempt to reduce their number in the FEM mesh, leaving NUMBER_TOP_ELEMENTS surface elements untouched
Pairs of elements within the snow cover satisfying the conditions of combineCondition() are combined by placing everything in the lower element, setting the density of upper element to Constants::undefined, and getting rid of node in between.
The elements being very similar and thus the microstructure parameters being approximately equal as defined in combineCondition(), simply average the microstructure properties
NOTE that the condense element check is placed at the end of a time step, allowing elements do develop on their own. More...
 
void combineElements (const size_t &number_top_elements, const bool &reduce_n_elements, const size_t &cond)
 
void splitElement (const size_t &e)
 Split the element provided as first argument. More...
 
void splitElements ()
 Split elements when they are near the top of the snowpack, when REDUCE_N_ELEMENTS is used. More...
 
void compSnowpackMasses ()
 Computes the actual total masses of the snowpack (kg m-2) More...
 
void compSnowpackInternalEnergyChange (const double &sn_dt)
 Computes the internal energy change of the snowpack during one computation time step (J m-2) More...
 
void compSoilInternalEnergyChange (const double &sn_dt)
 Computes the internal energy change of the soil during one computation time step (J m-2) More...
 
double getLiquidWaterIndex () const
 Computes the liquid water index defined as the ratio of total liquid water content (in mm w.e.) to calculated snow depth (in mm) divided by 0.03. Unit: (1) More...
 
double getModelledTemperature (const double &z) const
 Returns modelled internal snow or/and soil temperature (instantaneous value; degC), at a given position z perpendicular to slope (m)
z must be less than computed height (Xdata->cH), otherwise modeled temperature is set to Constants::undefined. More...
 
size_t getNumberOfElements () const
 
size_t getNumberOfNodes () const
 
bool isGlacier (const bool &hydro=false) const
 returns if a snow profile can be considered as a glacier. Practically, the hydrological criteria is that if a pixel contains more than 2 m of pure ice anywhere, it is considered to be glaciated. The standard criteria is that if the top 5 layers are made of pure ice, the pixel is glaciated. Therefore, a glacier covered by seasonal snow is glaciated in regard to the hydrological criteria but non-glaciated in regard to the standard criteria. More...
 
bool hasSoilLayers () const
 
size_t find_tag (const size_t &tag) const
 Find element with corresponding tag or return IOUtils::npos if not found. More...
 
const std::string toString () const
 

Static Public Member Functions

static bool combineCondition (const ElementData &Edata0, const ElementData &Edata1, const double &depth, const bool &reduce_n_elements)
 Boolean routine to check whether two snow elements can be combined. More...
 
static void mergeElements (ElementData &Edata0, const ElementData &Edata1, const bool &merge, const bool &topElement)
 Merging two elements. More...
 

Public Attributes

mio::StationData meta
 Station meta data. More...
 
double cos_sl
 Cosinus of slope angle, initialized once! More...
 
size_t sector
 current slope sector of width 360./max(1, nSlopes-1) More...
 
CanopyData Cdata
 Pointer to canopy data. More...
 
double pAlbedo
 Parameterized snow albedo. More...
 
double Albedo
 Snow albedo used by the model. More...
 
double SoilAlb
 Soil albedo. More...
 
double BareSoil_z0
 Bare soil roughness in m. More...
 
size_t SoilNode
 The top soil node, 0 in case of SNP_SOIL == 0. More...
 
double Ground
 The ground height – meaning the height of the top soil node. More...
 
double cH
 The CALCULATED height, including soil depth if SNP_SOIL == 1. More...
 
double mH
 The MEASURED height, including soil depth if SNP_SOIL == 1. More...
 
double mass_sum
 Total mass summing mass of snow elements. More...
 
double swe
 Total mass summing snow water equivalent of elements. More...
 
double lwc_sum
 Total liquid water in snowpack. More...
 
double hn
 Depth of new snow to be used on slopes. More...
 
double rho_hn
 Density of new snow to be used on slopes. More...
 
size_t ErosionLevel
 Element where snow erosion stopped previously for the drift index. More...
 
double ErosionMass
 Eroded mass either real or virtually (storage if less than one element) More...
 
char S_class1
 Stability class based on hand hardness, grain class ... More...
 
char S_class2
 Stability class based on hand hardness, grain class ... More...
 
double S_d
 Minimum Direct Action Stability Index ... More...
 
double z_S_d
 Depth of Minimum Direct Action Stability. More...
 
double S_n
 Minimum Natural Stability Index. More...
 
double z_S_n
 Depth of Minimum Natural Stability. More...
 
double S_s
 Minimum Skier Stability Index (SSI) More...
 
double z_S_s
 Depth of Minimum SSI. More...
 
double S_4
 stab_index4 More...
 
double z_S_4
 Depth of stab_index4. More...
 
double S_5
 stab_index5 More...
 
double z_S_5
 Depth of stab_index5. More...
 
std::vector< NodeDataNdata
 pointer to nodal data array (e.g. T, z, u, etc..) More...
 
std::vector< ElementDataEdata
 pointer to element data array (e.g. Te, L, Rho, etc..) More...
 
void * Kt
 Pointer to pseudo-conductivity and stiffnes matrix. More...
 
double ColdContent
 Cold content of snowpack (J m-2) More...
 
double ColdContentSoil
 Cold content of soil (J m-2) More...
 
double dIntEnergy
 Internal energy change of snowpack (J m-2) More...
 
double dIntEnergySoil
 Internal energy change of soil (J m-2) More...
 
double meltFreezeEnergy
 Melt freeze part of internal energy change of snowpack (J m-2) More...
 
double meltFreezeEnergySoil
 Melt freeze part of internal energy change of soil (J m-2) More...
 
double ReSolver_dt
 Last used RE time step in the previous SNOWPACK time step. More...
 
bool windward
 True for windward (luv) slope. More...
 
double WindScalingFactor
 Local scaling factor for wind at drift station. More...
 
double TimeCountDeltaHS
 Time counter tracking erroneous settlement in operational mode. More...
 

Static Public Attributes

static const void * Seaice = NULL
 dummy pointer for compatibility with the dev version More...
 
static const double comb_thresh_l = 0.015
 Both elements must be smaller than COMB_THRESH_L (m) for an action to be taken. More...
 
static const double comb_thresh_ice = 0.05
 Volumetric ice content (1), i.e., about 46 kg m-3. More...
 
static const double comb_thresh_water = 0.01
 Water content (1) More...
 
static const double comb_thresh_dd = 0.2
 Dendricity (1) More...
 
static const double comb_thresh_sp = 0.05
 Sphericity (1) More...
 
static const double comb_thresh_rg = 0.125
 Grain radius (mm) More...
 
static const double thresh_moist_snow = 0.003
 Snow elements with a LWC above this threshold are considered at least to be moist. More...
 
static const double thresh_moist_soil = 0.0001
 
static const size_t number_top_elements = 5
 Number of top elements left untouched by the join functions. More...
 
static unsigned short number_of_solutes = 0
 The model treats that number of solutes. More...
 

Friends

std::ostream & operator<< (std::ostream &os, const SnowStation &data)
 
std::istream & operator>> (std::istream &is, SnowStation &data)
 

Constructor & Destructor Documentation

◆ SnowStation() [1/2]

SnowStation::SnowStation ( const bool &  i_useCanopyModel = true,
const bool &  i_useSoilLayers = true 
)
explicit

◆ SnowStation() [2/2]

SnowStation::SnowStation ( const SnowStation c)

◆ ~SnowStation()

SnowStation::~SnowStation ( )

Member Function Documentation

◆ combineCondition()

bool SnowStation::combineCondition ( const ElementData Edata0,
const ElementData Edata1,
const double &  depth,
const bool &  reduce_n_elements 
)
static

Boolean routine to check whether two snow elements can be combined.

  • no action will be taken if one of the two elements is
    • a soil element
    • larger than comb_thresh_l
    • tagged
    • dry surface hoar (mk=3)
    • dendritic but not both
  • otherwise we use criteria for dendricity, sphericity, volumetric ice or water content, grain size and marker
  • Whatever type of thin elements are treated in WaterTransport::mergingElements()
Parameters
Edata0Lower element
Edata1Upper element
depthDistance of the element from the snow surface
reduce_n_elementsEnable more "aggressive" combining for layers deeper in the snowpack, to reduce the number of elements and thus the computational load.
Returns
true if the two elements should be combined, false otherwise

◆ combineElements() [1/2]

void SnowStation::combineElements ( const size_t &  i_number_top_elements,
const bool &  reduce_n_elements 
)

If more than NUMBER_TOP_ELEMENTS snow elements exist, attempt to reduce their number in the FEM mesh, leaving NUMBER_TOP_ELEMENTS surface elements untouched
Pairs of elements within the snow cover satisfying the conditions of combineCondition() are combined by placing everything in the lower element, setting the density of upper element to Constants::undefined, and getting rid of node in between.
The elements being very similar and thus the microstructure parameters being approximately equal as defined in combineCondition(), simply average the microstructure properties
NOTE that the condense element check is placed at the end of a time step, allowing elements do develop on their own.

Parameters
i_number_top_elementsThe number of surface elements to be left untouched
reduce_n_elementsEnable more "aggressive" combining for layers deeper in the snowpack

◆ combineElements() [2/2]

void SnowStation::combineElements ( const size_t &  number_top_elements,
const bool &  reduce_n_elements,
const size_t &  cond 
)

◆ compSnowpackInternalEnergyChange()

void SnowStation::compSnowpackInternalEnergyChange ( const double &  sn_dt)

Computes the internal energy change of the snowpack during one computation time step (J m-2)

◆ compSnowpackMasses()

void SnowStation::compSnowpackMasses ( )

Computes the actual total masses of the snowpack (kg m-2)

◆ compSoilInternalEnergyChange()

void SnowStation::compSoilInternalEnergyChange ( const double &  sn_dt)

Computes the internal energy change of the soil during one computation time step (J m-2)

◆ find_tag()

size_t SnowStation::find_tag ( const size_t &  tag) const

Find element with corresponding tag or return IOUtils::npos if not found.

Parameters
tagTag to look for
Returns
Index of tagged element, IOUtils::npos if not found

◆ getLiquidWaterIndex()

double SnowStation::getLiquidWaterIndex ( ) const

Computes the liquid water index defined as the ratio of total liquid water content (in mm w.e.) to calculated snow depth (in mm) divided by 0.03. Unit: (1)

◆ getModelledTemperature()

double SnowStation::getModelledTemperature ( const double &  z) const

Returns modelled internal snow or/and soil temperature (instantaneous value; degC), at a given position z perpendicular to slope (m)
z must be less than computed height (Xdata->cH), otherwise modeled temperature is set to Constants::undefined.

Version
11.03
Parameters
zSensor position perpendicular to slope (m)

◆ getNumberOfElements()

size_t SnowStation::getNumberOfElements ( ) const

◆ getNumberOfNodes()

size_t SnowStation::getNumberOfNodes ( ) const

◆ hasSoilLayers()

bool SnowStation::hasSoilLayers ( ) const

◆ initialize()

void SnowStation::initialize ( const SN_SNOWSOIL_DATA SSdata,
const size_t &  i_sector 
)

This routine initializes the snow cover structure which contains all information about a station including element, nodal and canopy data
Because you are working with layers, the first node is a special case; the rest can be generated in a loop .... The bottom temperature at the beginning of the computation is given by the temperature at the top of the lowest soil or snow layer
IMPORTANT: it is very important for Alpine3D that Cdata.height is initialized even if CANOPY = 0, otherwise SnowInterface will not recognize the canopy grids (David 2007-06-25).

Version
10.02
Parameters
SSdata
i_sectordefines the exposition sector of the slope (width 360./number_slopes)

◆ isGlacier()

bool SnowStation::isGlacier ( const bool &  hydro = false) const

returns if a snow profile can be considered as a glacier. Practically, the hydrological criteria is that if a pixel contains more than 2 m of pure ice anywhere, it is considered to be glaciated. The standard criteria is that if the top 5 layers are made of pure ice, the pixel is glaciated. Therefore, a glacier covered by seasonal snow is glaciated in regard to the hydrological criteria but non-glaciated in regard to the standard criteria.

Parameters
hydroif true, use an hydrologist criteria (default: false)
Returns
true if the profile belongs to a glacier

◆ mergeElements()

void SnowStation::mergeElements ( ElementData EdataLower,
const ElementData EdataUpper,
const bool &  merge,
const bool &  topElement 
)
static

Merging two elements.

  • Joining:
    • Keep the lower element, that is, the lowest snow element always survives!
    • The properties of the upper and lower elements are (depth) averaged in an appropriate way.
    • The new length is the sum of both
    • Keep the birthday of the upper element
    • Keep the tag of the upper element only if the lower is untagged (mk >= 100)
    • Note
      Joining two elements may cause the tag (marker >= 100) to "jump" abruptly
  • Removing:
    • Remaining ice, liquid water, solutes, etc. are added to the lower element
    • The length of the lower element is kept
    • Keep the birthday of the lower element
      Parameters
      EdataLowerProperties of lower element
      EdataUpperProperties of upper element
      mergeTrue if upper element is to be joined with lower one, false if upper element is to be removed
      topElementset to true if the upper element is at the very top of the snow pack

◆ operator=()

SnowStation & SnowStation::operator= ( const SnowStation source)

Assignement operator.

◆ reduceNumberOfElements()

void SnowStation::reduceNumberOfElements ( const size_t &  rnE)

Remove the upper "marked" element of two (snow only)

  1. Merging two elements:
    1. density is undefined
    2. take the uppermost node of both
  2. Removing melted or thin elements
    1. density is undefined AND length negative (*= -1.) as the latter will be used!
    2. keep upper node of lowest element.

Parameters
rnEReduced number of elements

◆ resize()

void SnowStation::resize ( const size_t &  number_of_elements)

Reallocate element and node data
Edata and Ndata as well as nElems and nNodes are reallocated or reset, respectively.

Parameters
number_of_elementsThe new number of elements

◆ splitElement()

void SnowStation::splitElement ( const size_t &  e)

Split the element provided as first argument.

◆ splitElements()

void SnowStation::splitElements ( )

Split elements when they are near the top of the snowpack, when REDUCE_N_ELEMENTS is used.

  • This function split elements when they are getting closer to the top of the snowpack. This is required when using the "aggressive" merging option (REDUCE_N_ELEMENTS). When snow melt brings elements back to the snow surface, smaller layer spacing is required to accurately describe temperature and moisture gradients.

◆ toString()

const std::string SnowStation::toString ( ) const

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const SnowStation data 
)
friend

◆ operator>>

std::istream& operator>> ( std::istream &  is,
SnowStation data 
)
friend

Member Data Documentation

◆ Albedo

double SnowStation::Albedo

Snow albedo used by the model.

◆ BareSoil_z0

double SnowStation::BareSoil_z0

Bare soil roughness in m.

◆ Cdata

CanopyData SnowStation::Cdata

Pointer to canopy data.

◆ cH

double SnowStation::cH

The CALCULATED height, including soil depth if SNP_SOIL == 1.

◆ ColdContent

double SnowStation::ColdContent

Cold content of snowpack (J m-2)

◆ ColdContentSoil

double SnowStation::ColdContentSoil

Cold content of soil (J m-2)

◆ comb_thresh_dd

const double SnowStation::comb_thresh_dd = 0.2
static

Dendricity (1)

◆ comb_thresh_ice

const double SnowStation::comb_thresh_ice = 0.05
static

Volumetric ice content (1), i.e., about 46 kg m-3.

◆ comb_thresh_l

const double SnowStation::comb_thresh_l = 0.015
static

Both elements must be smaller than COMB_THRESH_L (m) for an action to be taken.

◆ comb_thresh_rg

const double SnowStation::comb_thresh_rg = 0.125
static

Grain radius (mm)

◆ comb_thresh_sp

const double SnowStation::comb_thresh_sp = 0.05
static

Sphericity (1)

◆ comb_thresh_water

const double SnowStation::comb_thresh_water = 0.01
static

Water content (1)

◆ cos_sl

double SnowStation::cos_sl

Cosinus of slope angle, initialized once!

◆ dIntEnergy

double SnowStation::dIntEnergy

Internal energy change of snowpack (J m-2)

◆ dIntEnergySoil

double SnowStation::dIntEnergySoil

Internal energy change of soil (J m-2)

◆ Edata

std::vector<ElementData> SnowStation::Edata

pointer to element data array (e.g. Te, L, Rho, etc..)

◆ ErosionLevel

size_t SnowStation::ErosionLevel

Element where snow erosion stopped previously for the drift index.

◆ ErosionMass

double SnowStation::ErosionMass

Eroded mass either real or virtually (storage if less than one element)

◆ Ground

double SnowStation::Ground

The ground height – meaning the height of the top soil node.

◆ hn

double SnowStation::hn

Depth of new snow to be used on slopes.

◆ Kt

void* SnowStation::Kt

Pointer to pseudo-conductivity and stiffnes matrix.

◆ lwc_sum

double SnowStation::lwc_sum

Total liquid water in snowpack.

◆ mass_sum

double SnowStation::mass_sum

Total mass summing mass of snow elements.

◆ meltFreezeEnergy

double SnowStation::meltFreezeEnergy

Melt freeze part of internal energy change of snowpack (J m-2)

◆ meltFreezeEnergySoil

double SnowStation::meltFreezeEnergySoil

Melt freeze part of internal energy change of soil (J m-2)

◆ meta

mio::StationData SnowStation::meta

Station meta data.

◆ mH

double SnowStation::mH

The MEASURED height, including soil depth if SNP_SOIL == 1.

◆ Ndata

std::vector<NodeData> SnowStation::Ndata

pointer to nodal data array (e.g. T, z, u, etc..)

◆ number_of_solutes

unsigned short SnowStation::number_of_solutes = 0
static

The model treats that number of solutes.

◆ number_top_elements

const size_t SnowStation::number_top_elements = 5
static

Number of top elements left untouched by the join functions.

◆ pAlbedo

double SnowStation::pAlbedo

Parameterized snow albedo.

◆ ReSolver_dt

double SnowStation::ReSolver_dt

Last used RE time step in the previous SNOWPACK time step.

◆ rho_hn

double SnowStation::rho_hn

Density of new snow to be used on slopes.

◆ S_4

double SnowStation::S_4

stab_index4

◆ S_5

double SnowStation::S_5

stab_index5

◆ S_class1

char SnowStation::S_class1

Stability class based on hand hardness, grain class ...

◆ S_class2

char SnowStation::S_class2

Stability class based on hand hardness, grain class ...

◆ S_d

double SnowStation::S_d

Minimum Direct Action Stability Index ...

◆ S_n

double SnowStation::S_n

Minimum Natural Stability Index.

◆ S_s

double SnowStation::S_s

Minimum Skier Stability Index (SSI)

◆ Seaice

const void * SnowStation::Seaice = NULL
static

dummy pointer for compatibility with the dev version

◆ sector

size_t SnowStation::sector

current slope sector of width 360./max(1, nSlopes-1)

◆ SoilAlb

double SnowStation::SoilAlb

Soil albedo.

◆ SoilNode

size_t SnowStation::SoilNode

The top soil node, 0 in case of SNP_SOIL == 0.

◆ swe

double SnowStation::swe

Total mass summing snow water equivalent of elements.

◆ thresh_moist_snow

const double SnowStation::thresh_moist_snow = 0.003
static

Snow elements with a LWC above this threshold are considered at least to be moist.

◆ thresh_moist_soil

const double SnowStation::thresh_moist_soil = 0.0001
static

◆ TimeCountDeltaHS

double SnowStation::TimeCountDeltaHS

Time counter tracking erroneous settlement in operational mode.

◆ WindScalingFactor

double SnowStation::WindScalingFactor

Local scaling factor for wind at drift station.

◆ windward

bool SnowStation::windward

True for windward (luv) slope.

◆ z_S_4

double SnowStation::z_S_4

Depth of stab_index4.

◆ z_S_5

double SnowStation::z_S_5

Depth of stab_index5.

◆ z_S_d

double SnowStation::z_S_d

Depth of Minimum Direct Action Stability.

◆ z_S_n

double SnowStation::z_S_n

Depth of Minimum Natural Stability.

◆ z_S_s

double SnowStation::z_S_s

Depth of Minimum SSI.


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