ply_sampled_tracking_module Module

Tracking offers the possibility to extract data from a simulation for given subsections of the mesh and specific points in time.

When the data is given in form of polynomials, it usually is required to perform some kind of sampling to obtain a visualization. The ply_sampled_tracking implements this functionality. See ply_sampling_module for details on the configuration of the sampling strategy and tem_tracking_module for details on the general configuration of trackings (description of mesh subsections to extract data for).

An important special case is the tracking of a single point. In this case no sampling has to be done, instead the polynomial representing the state in the element that contains the requested point needs to be evaluated at this point. This can be configured by setting use_get_point=true in the output subtable of the tracking object.

A complete example looks like this:

tracking = {
  label = 'track_shearhat2D',
  folder = './',
  variable = {'momentum','density','energy'},
  shape = {
    kind = 'canoND',
    object= {
      origin ={0.01*dx, 0., 0.}
    }
  },
  time_control = {
    min = 0,
    max = sim_control.time_control.max,
    interval = {iter = 10}
  },
  output = {
    format = 'ascii',
    use_get_point = true
  }
}

This tracks the state variables (momentum, density and energy) in a single point (0.01*dx, 0, 0) and writes them every ten iteration to an ASCII file. Each point in time gets written on a new line in the same file. If you do not use the use_get_point option (its false by default), and no sampling is active, all degrees of freedom of the field in the element, that contains the point, will be written. You can limit the number of degrees of freedom by setting ndofs to some value. The first mode of the Legendre series is the integral mean, so this is usually the value you want to get. Thus, setting ndofs=1 gives you the averaged value of the element the point is found in. The according output table would then look as follows:

  output = {
    format = 'ascii',
    ndofs = 1
  }

Of course, more complex shapes may be tracked, in that case it usually is not sensible to use ascii output anymore. Instead you are than likely to want to view data later on and accordingly write it in VTK format. For this, a sampling table (see ply_sampling_module) should be considered.


Uses

Used by

  • module~~ply_sampled_tracking_module~~UsedByGraph module~ply_sampled_tracking_module ply_sampled_tracking_module program~sdr_harvesting sdr_harvesting program~sdr_harvesting->module~ply_sampled_tracking_module module~sdr_hvs_config_module sdr_hvs_config_module program~sdr_harvesting->module~sdr_hvs_config_module module~sdr_hvs_config_module->module~ply_sampled_tracking_module

Contents


Derived Types

type, public :: ply_sampled_tracking_type

Components

TypeVisibilityAttributesNameInitial
type(tem_tracking_type), private :: tracking

Contains all tracking headers, control and tracking entities active on local process

type(treelmesh_type), private, allocatable:: mesh(:)

Subsampled mesh for each tracking.

Read more…
type(tem_varSys_type), private, allocatable:: varsys(:)

Variable system description after subsampling.

Read more…
type(ply_sampling_type), private :: sampling

Configuration of the subsampling (applied to all trackings).

integer, private :: ndims

Dimensionality of the data to sample.


Subroutines

public subroutine ply_sampled_tracking_load(me, conf)

Load the configuration of sampled tracking objects.

Arguments

TypeIntentOptionalAttributesName
type(ply_sampled_tracking_type), intent(out) :: me

Sampled tracking data to load from the config

type(flu_State) :: conf

Lua config to load the tracking from

public subroutine ply_sampled_track_init(me, mesh, solver, varSys, bc, stencil, proc, nDofs, nDims)

Initialize the sampled tracking entities.

Read more…

Arguments

TypeIntentOptionalAttributesName
type(ply_sampled_tracking_type), intent(inout) :: me

Sampled tracking variable to initialize. It has to be configured by ply_sampled_tracking_load beforehand.

type(treelmesh_type), intent(in) :: mesh

The global mesh.

type(tem_solveHead_type), intent(in) :: solver

Information about the solver (used to construct file name strings).

type(tem_varSys_type), intent(in) :: varSys

Global variable system with description of the data to get the tracking variables from.

type(tem_BC_prop_type), intent(in) :: bc

Boundary condition properties, used to identify elements close to the boundary.

type(tem_stencilHeader_type), intent(in), optional :: stencil

Description of the stencil in the numerical scheme.

This is needed to describe elements adjacent to specific boundary labels.

type(tem_comm_env_type), intent(in) :: proc

General communication environment

integer, intent(in) :: nDofs

Number of degrees of freedom to use in the output.

integer, intent(in) :: nDims

Number of dimensions in the polynomial representations.

public subroutine ply_sampled_track_output(me, mesh, bc, solver, proc, varSys, var_degree, lvl_degree, var_space, simControl, time)

Output sampled tracking data.

Read more…

Arguments

TypeIntentOptionalAttributesName
type(ply_sampled_tracking_type), intent(inout) :: me

Sampled tracking instances.

type(treelmesh_type), intent(in) :: mesh

Global mesh, required for the sampling.

type(tem_BC_prop_type), intent(in) :: bc

Boundary properties, needed to inherit boundary information to refined meshes and allow the extraction of boundary shape geometries.

type(tem_solveHead_type), intent(in) :: solver

Information about the solver, needed for the output file name.

type(tem_comm_env_type), intent(in) :: proc

General communication environment

type(tem_varSys_type), intent(in) :: varSys

Original variable system

integer, intent(in) :: var_degree(:)

Maximal polynomial degree for each variable

Needs to match the size of the variable system.

integer, intent(in) :: lvl_degree(:)

Maximal polynomial degree for each level

integer, intent(in) :: var_space(:)

Maximal polynomial space for each variable

Needs to match the size of the variable system.

type(tem_simControl_type), intent(in), optional :: simControl

Simulation control to determine, whether trackings should be written

If not provided, all trackings will be written unconditionally.

type(tem_time_type), intent(in), optional :: time

Provide a time for the current data set to write in tracking.

This only is respected if no simControl is provided. If simControl is present the time information from it will be used instead.