ParPyDTK2 C++ API

Common Definitions

group common

Defines

handle_moab_error(__ret)

macro to handle moab error

throw_error(__msg)

throw runtime_error exception

throw_error_if(__cond, __msg)

conditionally error throw

throw_noimpl(__what)

throw not implemented feature error

throw_noimpl_if(__cond, __what)

throw not implemented feature error with condition

show_warning(__msg)

log warning message in stderr

show_warning_if(__cond, __msg)

log warning with condition

show_experimental(__msg)

log experimental warning in stderr

show_experimental_if(__cond, __msg)

log experimental warning in stderr with condition

show_info(__msg, __rank)

show information in parallel

show_info_master(__msg, __rank)

show information only on master rank

streamer(__rank)

streaming message with specific rank

streamer_master(__rank)

streaming messages only on the master process

Typedefs

typedef entity_t

MOAB entity handle.

Enums

enum [anonymous]

root set

Values:

root_set = 0

Field Variables

class FieldData

a representation of MOAB tag for field data

Public Functions

FieldData(moab::Core &mdb, const std::string &field_name, int dim = 1)

constructor with moab instance

Note
set is not the same as mesh set in MOAB
Parameters
  • mdb: moab instance
  • field_name: field name
  • dim: field dimension

void assign(const moab::Range &range, const double *values)

assign values

Parameters
  • range: entity ranges, for this work, it should be vertices
  • values: data values, for vector/tensor, C order is expected

void assign_1st(const moab::Range &range, const double *values)

assign to first node

This function is used by Python for handling empty partitions

Parameters
  • range: entity ranges
  • values: values for the first node, at least size of dim

void extract(const moab::Range &range, double *values) const

extract values

Parameters
  • range: entity ranges, for this work, it should be vertices
  • values: data values, for vector/tensor, C order is expected

void extract_1st(const moab::Range &range, double *values) const

extract the value from the first node

Parameters
  • range: entity ranges
  • values: values for the first node, at least size of dim

operator const std::string&() const

brief implicitly cast to string

int dim() const

check the dimension

int set() const

check the set ID

const moab::Tag &tag() const

get MOAB tag

Protected Attributes

moab::Core &mdb_

reference to moab instance

std::string fn_

field name

int set_

set count

int dim_

field dimension

moab::Tag tag_

moab tag

class FieldDataSet

a set of field data

Public Types

typedef base_t::iterator iterator

interator type

typedef base_t::const_iterator const_iterator

constant iterator

Public Functions

virtual ~FieldDataSet()

destructor

bool has_field(const std::string &fn) const

check if a field exist

Parameters
  • fn: field name

void create(moab::Core &mdb, const std::string &field_name, int dim = 1)

create an data field

Note
set is not the same as mesh set in MOAB
Parameters
  • mdb: moab data base
  • field_name: field name
  • dim: field dimension

FieldData &operator[](const std::string &fn)

get a reference to a field data

Note
this overloads the base operator[]
Parameters
  • fn: field name

const FieldData &operator[](const std::string &fn) const

get a const reference to a field data

Parameters
  • fn: field name

iterator begin()

get the first iterator

iterator end()

get the end iterator

const_iterator begin() const

get the constant iterator

const_iterator end() const

get the constant end iterator

const_iterator cbegin() const

get the constant iterator

const_iterator cend() const

get the constant end iterator

Protected Attributes

base_t fs_

fields

Private Types

typedef std::unordered_map<std::string, FieldData *> base_t

data structure

Interface Mesh Database

class IMeshDB

interface mesh database, build on top of MOAB

imesh_py_interface

bool created() const

check created mesh database

void begin_create()

begin to create mesh

void create_vset()

create a new vertex set

void create_vertices(int nv, const double *coords)

create vertices

Parameters
  • nv: number of vertices
  • coords: coords values

void extract_vertices(double *coords) const

extract assigned coordinates

Note
coords must be at least n*3 where n is the size of the mesh
Parameters
  • coords: coordinates

void assign_gids(int nv, const int *gids)

assign global IDs

Note
gids should be one-based indices
Parameters
  • nv: number of local vertices
  • gids: global IDs

void extract_gids(int *gids) const

extract global IDs

Parameters
  • gids: global IDs

void finish_create(bool trivial_gid = true)

finish manupilating the mesh

NOTE that if your input mesh is element-based partition and the vertices you create are mesh nodes, then you have to specify the correct global IDs in parallel. However, if the coordinates are face centres, then the global IDs can be trivially computed by MOAB since there are no shared entities cross different processes.

Parameters
  • trivial_gid: true if we let MOAB to compute the GID

bool empty() const

check if empty partition

bool has_empty() const

check if any of the process has an empty partition

const std::vector<int> &_m2s() const

get a reference to the m2s pattern

Note
This is used in Python level as “private” thus having “_”

int size() const

check mesh size

int gsize() const

check the global mesh size

void get_bbox(double *v) const

get bounding box

Parameters
  • v: values

void get_gbbox(double *v) const

get global bounding box

Parameters
  • v: values

void create_field(const std::string &field_name, int dim = 1)

create a field

Parameters
  • field_name: field name
  • dim: field dimension

bool has_field(const std::string &field_name) const

check if we have a field

Parameters
  • field_name: field name

int field_dim(const std::string &field_name) const

check field dimension

Parameters
  • field_name: field name

void assign_field(const std::string &field_name, const double *values)

assign a value to a field

Parameters
  • field_name: field name
  • values: field data values

void _assign_1st(const std::string &field_name, const double *values)

assign to the first node

This function is used by Python to resolve the issues when empty empty partitions happen. Therefore, this function has an “_” prefix to indicate “private” usage!

Parameters
  • field_name: field name
  • values: field data values, at least size of field dimension

void extract_field(const std::string &field_name, double *values) const

extract value

Parameters
  • field_name: field name
  • values: field data values

void _extract_1st(const std::string &field_name, double *values) const

extract first value

This function is used by Python to resolve the issues when empty empty partitions happen. Therefore, this function has an “_” prefix to indicate “private” usage!

Parameters
  • field_name: field name
  • values: field data values

Public Functions

IMeshDB(MPI_Comm comm = MPI_COMM_WORLD)

constructor with communicator

Parameters
  • comm: communicator

int ranks() const

get total ranks

int rank() const

get my rank

MPI_Comm comm() const

get the communicator

Teuchos::RCP<moab::ParallelComm> pcomm() const

get mesh

std::vector<DataTransferKit::MoabManager> &mangers()

get the manger

void set_dimension(int dim)

set geometry dimension

Parameters
  • dim: dimension

bool ready() const

check if ready

dtk_field_t &dtk_fields()

get the dtk fields

Protected Types

typedef std::unordered_map<std::string, std::pair<int, Teuchos::RCP<Tpetra::MultiVector<double, int, DataTransferKit::SupportId>>>> dtk_field_t

handy typedef

Protected Attributes

moab::Core mdb_

moab instance

Teuchos::RCP<moab::ParallelComm> par_

moab parallel interface

std::vector<entity_t> vsets_

vertex sets

std::vector<moab::Range> locals_

local vertex range

std::vector<std::array<double, 6>> bboxes_

bounding boxes

std::vector<std::array<double, 6>> gbboxes_

global bounding boxes

FieldDataSet fields_

field data set

bool created_

flag to indicate whether users are done with creating mesh

moab::Tag gidtag_

global ID tag

moab::Tag parttag_

partition tag

bool usergid_

flag to indicate if we have user computed global ID

std::vector<DataTransferKit::MoabManager> mngrs_

DTK MOAB manager.

dtk_field_t dtkfields_

DTK multi vector for MOAB tags.

bool empty_

check empty partition

bool has_empty_

check if any process is empty

std::vector<int> m2s_

comm pattern for master2slaves for handling empty partitions

int gsize_

global size

Private Functions

void init_(bool del = false)

helper for clean up mesh

Parameters
  • del: whether or not delete mesh

void reset_vecs_()

handle all vectors

void init_bbox_(int i = -1)

initialize empty bounding boxes

Parameters
  • i: index if < 0 then init all

void cmpt_bboxes_()

compute all bounding box

Solution Mapper

class Mapper

the mapper interface for interface solution transfer

mapper_py_interface

int ranks() const

get the ranks

int rank() const

get my rank

MPI_Comm comm() const

get the communicator

void set_dimension(int dim)

set dimension

Parameters
  • dim: geometry dimension

void use_mmls()

use moving least square, this is the default method

See
use_spline, use_n2n

void use_spline()

use spline interpolation method

See
use_mmls, use_n2n

void use_n2n(bool matching = false)

use node 2 node project

See
use_mmls, use_spline
Parameters
  • matching: are the interfaces mathing?

void use_awls()

use adaptive weighted least square fitting

See
use_mmls

void set_basis(int basis)

set basis function, default is Wendland 4th order

See
BasisFunctions
Parameters
  • basis: basis function and order

void use_knn_b(int)

use knn for blue mesh

void use_knn_g(int)

use knn for green mesh

void use_radius_b(double r)

use radius for blue

See
use_knn
Parameters
  • r: physical domain radius support

void use_radius_g(double r)

use radius for green

See
use_knn
Parameters
  • r: physical domain radius support

int check_method() const

check method

int check_basis() const

check basis

int knn_b() const

check blue knn

Note
if blue does not use knn, then negative value returned

int knn_g() const

check green knn

double radius_b() const

check blue radius

Note
if blue does not use radius, then -1.0 is returned

double radius_g() const

check green radius

int dimension() const

get the dimension

QRCP_ONLY

void set_leaf_b(int size)

set the leaf size

Warning
This method only works with unifem or chiao45 forked backend
Parameters
  • size: the leaf size in kd-tree

void set_leaf_g(int size)

set the leaf size for green database

Warning
This method only works with unifem or chiao45 forked backend
Parameters
  • size: the leaf size in kd-tree

int leaf_b() const

get the blue leaf size

int leaf_g() const

get the green leaf size

void _set_ind_file(const std::string &fn)

set the indicator tuning filename

Note
internal use
Note
This only works with QRCP implementation, or AWLS method
Parameters
  • fn: filename

void _wipe_ind_file()

wipe indicator file

void set_rho(double rho)

set local scaling rho

Note
This only works with QRCP implementation, or AWLS method
Parameters
  • rho: scaling factor

double rho() const

get local scaling rho

Public Functions

Mapper(std::shared_ptr<IMeshDB> B, std::shared_ptr<IMeshDB> G, const std::string &version = "", const std::string &date = "", bool profiling = true, const std::string &stat_file = "", bool verbose = true)

constructor

Parameters
  • B: input blue mesh
  • G: input green mesh
  • version: passed from Python inteface
  • date: passed from Python interface
  • profiling: whether do simple profiling, i.e. wtime
  • stat_file: statistics filename
  • verbose: doing verbose information printing

std::shared_ptr<IMeshDB> blue_mesh() const

get blue mesh

std::shared_ptr<IMeshDB> green_mesh() const

get green mesh

void begin_initialization()

begin initialization

void register_coupling_fields(const std::string &bf, const std::string &gf, bool direct)

register coupling fields

Parameters
  • bf: blue meshdb field data
  • gf: green meshdb field data
  • direct: true for b->g, false for g->b

bool has_coupling_fields(const std::string &bf, const std::string &gf, bool direct)

check if a coupling data fields exists

Parameters
  • bf: blue meshdb field data
  • gf: green meshdb field data
  • direct: true for b->g, false for g->b

void end_initialization()

end initialization

void begin_transfer()

begin to transfer data

void transfer_data(const std::string &bf, const std::string &gf, bool direct, bool resolve_disc = false, double sigma = -1.0)

transfer data

Notice that internally, the unifem backend occupies the mode parameter to indicate whether or not do post processing to resolve non-smooth solutions. Set mode == Teuchos::TRANS

Parameters
  • bf: blue meshdb field data
  • gf: green meshdb field data
  • direct: true for b->g, false for g->b
  • resolve_disc: (optional) if true, then try to resolve disc
  • sigma: (optional) indicator threshold

void end_transfer()

end transfer

Public Static Functions

static bool is_dtk2_backend()

check the backend

In chiao45 version of DTK2, we modified the exception class to add a prefix of “chiao45”, so it’s feasible to query this information w/o adding a new API

Return
if the underlying DTK2 is using in chiao45 forkedversion, then return true

Protected Attributes

std::shared_ptr<IMeshDB> B_

blue mesh

std::shared_ptr<IMeshDB> G_

green mesh

int dim_

dimension

bool ready_

flag to indicate the mapper is ready for transfering

bool profiling_

whether do simple profiling

bool verbose_

verbose output

std::unique_ptr<Teuchos::ParameterList> opts_[2]

parameter list

std::map<std::pair<std::string, std::string>, Teuchos::RCP<DataTransferKit::MapOperator>> operators_[2]

transfer operators references

Teuchos::RCP<DataTransferKit::MapOperator> optrs_[2]

actual transfer operators

std::map<std::pair<std::string, std::string>, std::vector<StatInfo>> info_[2]

information box for each of the registered pair of fields

std::map<std::pair<std::string, std::string>, unsigned long> counts_[2]

total transfer counter

std::unique_ptr<std::ofstream> stat_

the statistics file handle

int stat_write_freq_counter_

statistics file dump freq counter

Protected Static Attributes

DataTransferKit::MapOperatorFactory factory_

map factory

const std::string title_ = std::string(LEN1, '-')

title of verbose printing

const std::string indentation_ = std::string(LEN2, ' ')

indentation

Private Functions

void init_parlist_()

initialize parameter list

template <bool _Dir, typename _V>
void set_search(const std::string &type, const std::string &value, const _V &v)

helper for set local search

Template Parameters
  • _Dir: direction
  • _V: value type
Parameters
  • type: search type
  • value: tag for value
  • v: actual value

template <bool _Dir, typename _V>
_V get_search(const std::string &type, const std::string &value, const _V &dft) const

helper for get search info

Template Parameters
  • _Dir: direction
  • _V: value type
Parameters
  • type: search type
  • value: tag for value
  • dft: default value

void _dump_stats()

helper function for writing outputs

Private Static Functions

static std::string parse_list_(Teuchos::ParameterList &list)

parse and formatting a parameter list

Parameters
  • list: parameter list

struct StatInfo

a simple structure to store discontinuous points and timing results