Page MenuHomec4science

ref_node_paradis.hh
No OneTemporary

File Metadata

Created
Wed, Jun 12, 19:36

ref_node_paradis.hh

/**
* @file ref_node_paradis.hh
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
* @author Max Hodapp <max.hodapp@epfl.ch>
* @author Jaehyun Cho <jaehyun.cho@epfl.ch>
*
* @date Mon Oct 28 19:23:14 2013
*
* @brief Nodal reference of ParaDiS
*
* @section LICENSE
*
* Copyright (©) 2010-2011 EPFL (Ecole Polytechnique Fédérale de Lausanne)
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* LibMultiScale is free software: you can redistribute it and/or modify it
* under the
* terms of the GNU Lesser General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option) any
* later version.
*
* LibMultiScale is distributed in the hope that it will be useful, but
* WITHOUT ANY
* WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
* A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
* details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with LibMultiScale. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef __LIBMULTISCALE_REF_NODE_PARADIS_HH__
#define __LIBMULTISCALE_REF_NODE_PARADIS_HH__
/* -------------------------------------------------------------------------- */
// ParaDis includes
#define PARALLEL
#include <mpi.h>
extern "C" {
#include "Home.h"
#include "Typedefs.h"
}
#undef PARALLEL
/* -------------------------------------------------------------------------- */
#include "ref_node_dd.hh"
/* -------------------------------------------------------------------------- */
__BEGIN_LIBMULTISCALE__
/* -------------------------------------------------------------------------- */
class DomainParaDiS;
class ComparatorRefNodeParaDiS;
/* -------------------------------------------------------------------------- */
/**
* Class RefNodeParaDiS
*
*/
class RefNodeParaDiS : public RefDDNode<3, RefNodeParaDiS> {
/* ------------------------------------------------------------------------ */
/* Typedefs */
/* ------------------------------------------------------------------------ */
public:
//! generic container of the original model container
using Domain = DomainParaDiS;
//! node comparator class to be usable in maps
using RefComparator = ComparatorRefNodeParaDiS;
static const UInt Dim = 3;
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public:
RefNodeParaDiS() : fakestress(0) {
// Request creation of a new node
Node_t *newNode = GetNewNativeNode(ParaDiSHome::paradis_ptr);
AssignNodeToCell(ParaDiSHome::paradis_ptr, newNode);
this->index_node = newNode->myTag.index;
}
RefNodeParaDiS(int index) : fakestress(0) { index_node = index; }
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
bool operator==(RefNodeParaDiS &ref) {
return (ref.index_node == index_node);
};
VectorProxy<3> position0() { return this->position(); }
VectorProxy<3> position() {
return VectorProxy<3>(
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->x,
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->y,
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->z);
}
VectorProxy<3> velocity() {
return VectorProxy<3>(
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->vX,
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->vY,
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->vZ);
}
VectorProxy<3> force() {
return VectorProxy<3>(
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->fX,
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->fY,
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->fZ);
}
UInt nbArms() const {
return ParaDiSHome::paradis_ptr->nodeKeys[index_node]->numNbrs;
}
// Real &prev_position(unsigned int i) {
// switch (i) {
// case 0:
// return ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldx;
// break;
// case 1:
// return ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldy;
// break;
// case 2:
// return ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldz;
// break;
// }
// Real *p = NULL;
// return *p;
// }
// Real &getNeighborNodeCoord(unsigned int n, unsigned int i) {
// UInt my_domain = ParaDiSHome::paradis_ptr->myDomain;
// Node_t current = *ParaDiSHome::paradis_ptr->nodeKeys[this->index_node];
// Tag_t *tag = current.nbrTag;
// UInt tag_domain = tag[n].domainID;
// int idx_for_neighbor = tag[n].index;
// if (my_domain == tag_domain) {
// switch (i) {
// case 0:
// return ParaDiSHome::paradis_ptr->nodeKeys[idx_for_neighbor]->x;
// break;
// case 1:
// return ParaDiSHome::paradis_ptr->nodeKeys[idx_for_neighbor]->y;
// break;
// case 2:
// return ParaDiSHome::paradis_ptr->nodeKeys[idx_for_neighbor]->z;
// break;
// }
// } else {
// switch (i) {
// case 0:
// return ParaDiSHome::paradis_ptr->remoteDomainKeys[tag_domain]
// ->nodeKeys[idx_for_neighbor]
// ->x;
// break;
// case 1:
// return ParaDiSHome::paradis_ptr->remoteDomainKeys[tag_domain]
// ->nodeKeys[idx_for_neighbor]
// ->y;
// break;
// case 2:
// return ParaDiSHome::paradis_ptr->remoteDomainKeys[tag_domain]
// ->nodeKeys[idx_for_neighbor]
// ->z;
// break;
// }
// }
// Real *p = NULL;
// return *p;
// }
// int &getNeighborNodeIndex0(unsigned int n) {
// Node_t current = *ParaDiSHome::paradis_ptr->nodeKeys[this->index_node];
// Tag_t *tag = current.nbrTag;
// return tag[n].index0;
// }
// int &getNeighborNodeIndex(unsigned int n) {
// Node_t current = *ParaDiSHome::paradis_ptr->nodeKeys[this->index_node];
// Tag_t *tag = current.nbrTag;
// return tag[n].index;
// }
// int &getNeighborNodeSlave(unsigned int n) {
// UInt my_domain = ParaDiSHome::paradis_ptr->myDomain;
// Node_t current = *ParaDiSHome::paradis_ptr->nodeKeys[this->index_node];
// Tag_t *tag = current.nbrTag;
// UInt tag_domain = tag[n].domainID;
// int idx_for_neighbor = tag[n].index;
// if (my_domain == tag_domain)
// return ParaDiSHome::paradis_ptr->nodeKeys[idx_for_neighbor]->slave;
// else
// return ParaDiSHome::paradis_ptr->remoteDomainKeys[tag_domain]
// ->nodeKeys[idx_for_neighbor]
// ->slave;
// }
// int &getNbrNeighborsOfNeighbor(unsigned int n) {
// UInt my_domain = ParaDiSHome::paradis_ptr->myDomain;
// Node_t current = *ParaDiSHome::paradis_ptr->nodeKeys[this->index_node];
// Tag_t *tag = current.nbrTag;
// UInt tag_domain = tag[n].domainID;
// int idx_for_neighbor = tag[n].index;
// if (my_domain == tag_domain)
// return ParaDiSHome::paradis_ptr->nodeKeys[idx_for_neighbor]->numNbrs;
// else
// return ParaDiSHome::paradis_ptr->remoteDomainKeys[tag_domain]
// ->nodeKeys[idx_for_neighbor]
// ->numNbrs;
// }
UInt tag() { return this->index_node; }
int &getIndex0() {
return ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->myTag.index0;
}
int &getIndex() {
return ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->myTag.index;
}
void getIndexes(int &n) { n = index_node; }
void setIndex(UInt index) { index_node = index; }
// // Get the number of neighbors
// UInt getNbrOfNeighs() {
// return
// UInt(ParaDiSHome::ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->numNbrs);
// }
int getSlave() {
return ParaDiSHome::ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]
->slave;
}
//! constraint tag:
int getConstraint() {
return ParaDiSHome::ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]
->constraint;
}
void slaving() {
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->slave = 1;
}
void unslaving() {
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->slave = 0;
}
void constrain() {
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->constraint = 7;
}
void unconstrain() {
ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->constraint = 0;
}
// void setOldVel(unsigned int i, Real oldv) {
// switch (i) {
// case 0:
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldvX = oldv;
// return;
// case 1:
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldvY = oldv;
// return;
// case 2:
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldvZ = oldv;
// return;
// }
// return;
// }
// void release(UInt nbr_index) {
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->constraint = 0;
// UInt index_for_neighbor =
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->nbrTag[0].index;
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldvX =
// 0.5 * (ParaDiSHome::paradis_ptr->nodeKeys[index_for_neighbor]->oldvX
// +
// ParaDiSHome::paradis_ptr->nodeKeys[nbr_index]->oldvX);
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldvY =
// 0.5 * (ParaDiSHome::paradis_ptr->nodeKeys[index_for_neighbor]->oldvY
// +
// ParaDiSHome::paradis_ptr->nodeKeys[nbr_index]->oldvY);
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->oldvZ =
// 0.5 * (ParaDiSHome::paradis_ptr->nodeKeys[index_for_neighbor]->oldvZ
// +
// ParaDiSHome::paradis_ptr->nodeKeys[nbr_index]->oldvZ);
// InsertArm(ParaDiSHome::paradis_ptr,
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node],
// &ParaDiSHome::paradis_ptr->nodeKeys[nbr_index]->myTag,
// -1.0 *
// (ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->burgX[0]),
// -1.0 *
// (ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->burgY[0]),
// -1.0 *
// (ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->burgZ[0]),
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->nx[0],
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->ny[0],
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->nz[0], 0,
// 1);
// ParaDiSHome::paradis_ptr->nodeKeys[this->index_node]->numNbrs = 2;
// }
friend class ContainerNodesParaDiS;
friend class IteratorNodesParaDiS;
friend class IteratorElemsParaDiS;
private:
unsigned int index_node;
unsigned int index0_node;
Real fakestress; // This is a temprary fix to make the dumpers work: FIXMETILL
};
/* -------------------------------------------------------------------------- */
class ComparatorRefNodeParaDiS {
public:
bool operator()(RefNodeParaDiS &r1, RefNodeParaDiS &r2) const {
// return r1.index_atome < r2.index_atome;
LM_TOIMPLEMENT;
}
bool operator()(const RefNodeParaDiS &r1, const RefNodeParaDiS &r2) const {
LM_TOIMPLEMENT;
}
};
/* -------------------------------------------------------------------------- */
inline std::ostream &operator<<(std::ostream &os, RefNodeParaDiS &ref) {
int n;
ref.getIndexes(n);
os << "ParaDiSNode reference : " << n;
return os;
}
/* -------------------------------------------------------------------------- */
__END_LIBMULTISCALE__
#endif /* __LIBMULTISCALE_REF_NODE_ParaDiS_HH__ */

Event Timeline