Page MenuHomec4science

domain.hh
No OneTemporary

File Metadata

Created
Sun, Oct 20, 20:58

domain.hh

/**
* @file domain.hh
* @author Till Junge <junge@lsmspc42.epfl.ch>
* @date Thu Apr 19 17:34:16 2012
*
* @brief defines the domain class.
*
* @section LICENSE
*
* <insert lisence here>
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __CADD_MESHER_DOMAIN_HH__
#define __CADD_MESHER_DOMAIN_HH__
#include "common.hh"
#include "cadd_mesh.hh"
#include "lattice.hh"
#include "point_container.hh"
#include <string>
#include <vector>
#include <forward_list>
template <Uint DIM> class Tree;
template <Uint DIM>
class Domain {
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public:
/*! \param overall_ computational domain
* \param refined_ subset of \a overall_ where atomistic resolution is required
* in the case of CADD, this includes pad and MD
* \param atomistic_ subset of \a refined_. Contains all points which will be
* real atoms, i.e., interface and free
* \param atomistic_skinned_ subset of \a atomistic_. contains all free non-
* free interface atoms
* \param lattice_
* \param centre The origin from where the domain is to be filled
* \param boundaries_special_treatment whether the boundaries should be meshed
* seperately. This is necessary if the overall sim box size cannot be a
* power-of-2-multiple of a lattice or if the lattice is not aligned with
* the overall geometry
*/
Domain(const Geometry<DIM> & overall_, const Geometry<DIM> & refined_,
const Geometry<DIM> & atomistic_,
const Geometry<DIM> & atomistic_skinned_,
const Lattice<DIM> & lattice_,
bool boundaries_special_treatment = false,
const PointRef<DIM> & centre=*reinterpret_cast<PointRef<DIM> *>(NULL));
virtual ~Domain();
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
public:
/// function to print the contain of the class
virtual void printself(std::ostream & stream, int indent = 0) const;
/*! adds a point-refinement data point. Make sure that in the end, the convex
* hull of the points encompasses the entire Geometry \a overall_.
* \param coords
* \param refinement
*/
void setPointRefinement(const Real * coords, Real & refinement);
void setPointRefinement(const PointRef<DIM> & coords, Real & refinement){
this->setPointRefinement(coords.getArray(), refinement);
}
/// fills the domain
void fill_points() throw(std::string);
/* returns whether a point is to be considered for the fem mesh
* \param point*/
bool inFemDomain(PointRef<DIM> & point);
bool inDomain(PointRef<DIM> & point);
/// simple exposure of Geometry::from_border
inline Real from_border(const Real * point) const {
return this->overall->from_border(point);}
/*! returns the interpolated refinement at point \a coords
* \param coords */
inline Real interpolate(const Real * coords) const {
return this->auxiliary_mesh->interpolate(coords, this->refinement);}
inline Real interpolate(PointRef<DIM> point) const {
return this->interpolate(&(point.getComponent(0)));}
/*! dumps the atoms in lammps format
* \param filename */
void dump_atoms_lammps(std::string filename);
/*! dumps the atoms in paraview format
* \param filename */
void dump_atoms_paraview(std::string filename);
void compute_interface_atoms(PointContainer<DIM> & interface,
PointContainer<DIM> & pad,
const Geometry<DIM> & atomistic = *reinterpret_cast<Geometry<DIM>* >(NULL));
/* ------------------------------------------------------------------------ */
/* Accessors */
/* ------------------------------------------------------------------------ */
public:
/// returns a pointer to the lattice
inline Lattice<DIM> * getLattice() {return this->lattice;}
/// returns a reference to the points
inline const PointContainer<DIM> & getPoints() const {return this->points;}
/*! returns a reference to the main mesh. If the mesh does not yet exist, it
* is created and then returned */
Mesh<DIM> & getMesh() throw(std::string);
/*! returns a reference to the auxiliary mesh. Throws an error \a std::string
* if the mesh does not yet exist */
Mesh<DIM> & getAuxiliaryMesh() throw(std::string);
Tree<DIM> & getTree();
const PointContainer<DIM> & getAtoms();
void splitDegenerates(const Real & radius_threshold);
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
private:
/// Creates the auxiliary Mesh using the points defined in setPointRefinement
void initialiseAuxiliaryMesh();
void fill_atoms();
// geometric definitions
/// computational domain
Geometry<DIM> * overall;
/// subdomain of \a overall, atomistic refinement here
Geometry<DIM> * refined;
/// subdomain of \a refined. all real atoms, no pad
Geometry<DIM> * atomistic;
/// subdomain of \a atomistic. real atoms without interface atoms
Geometry<DIM> * atomistic_skinned;
// lattice
Lattice<DIM> * lattice;
// refinement definition
/// points used as nodes for the auxiliary mesh
PointContainer<DIM> refinement_points;
/// nodal refinement values to be interpolated
std::vector<Real> refinement;
/// auxiliary mesh, used to compute refinement in any point of the computational domain
Mesh<DIM> * auxiliary_mesh;
/// quadtree (2D) or octtree (3D) used to make a continuously refining mesh
Tree<DIM> * tree;
// points definition
/// all points, including pure atoms and pure fem nodes
PointContainer<DIM> points; //this is the main container with nodes and atoms
/// atoms, including pad
PointContainer<DIM> atoms;
/// pad and interface atom indices. this is filles when the atoms are
/// computed, so that we can separate interface and pad
std::forward_list<Uint> pad_atoms_indices;
/// nodes, including interface atoms
PointContainer<DIM> nodes;
/// main mesh object
Mesh<DIM> * main_mesh;
/// whether \a points has been filled
bool filled;
/// whether \a main_mesh has been meshed
bool main_meshed;
/*! whether the boundaries should be meshed seperately (important if the
* overall size cannot be a multiple of 2 lattices in every direction */
bool boundaries_special_treatment;
};
/* -------------------------------------------------------------------------- */
/* inline functions */
/* -------------------------------------------------------------------------- */
/// standard output stream operator
template <Uint DIM>
inline std::ostream & operator <<(std::ostream & stream, const Domain<DIM> & _this)
{
_this.printself(stream);
return stream;
}
#endif /* __CADD_MESHER_DOMAIN_HH__ */

Event Timeline