33static char help[] =
"...\n\n";
35int main(
int argc,
char *argv[]) {
41 moab::Core mb_instance;
42 moab::Interface &moab = mb_instance;
44 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab,
MYPCOMM_INDEX);
46 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD,
false);
49 new ParallelComm(&moab, moab_comm_wrap->get_comm());
54 PetscOptionsBegin(PETSC_COMM_WORLD,
"",
"Lorenz force configure",
56 CHKERR PetscOptionsString(
"-my_file",
"mesh file name",
"",
"solution.h5m",
75 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
76 CHKERR DMSetType(dm,
"DMMOFEM");
79 CHKERR DMSetFromOptions(dm);
85 using VolOp = VolumeElementForcesAndSourcesCore::UserDataOperator;
92 struct MyOpDebug :
public VolOp {
94 boost::shared_ptr<MatrixDouble>
B;
95 MyOpDebug(
decltype(
B) b) :
VolOp(
"MAGNETIC_POTENTIAL", OPROW),
B(b) {}
100 if (type == MBEDGE && side == 0) {
101 std::cout <<
"found " << (*B) << endl;
109 const double dist = 1e-12;
110 const int nb_random_points = 5000;
111 const int nb_steps = 10000;
112 const int mod_step = 10;
113 const double dt = 1e-5;
114 const double velocity_scale = 0.1;
115 const double magnetic_field_scale = 1;
116 const double scale_box = 0.5;
123 auto vol_ele = data_at_pts->feMethodPtr.lock();
126 "Pointer to element does not exists");
127 auto get_rule = [&](
int order_row,
int order_col,
int order_data) {
142 CHKERR data_at_pts->treePtr->get_bounding_box(box);
144 const double bMin = box.bMin[0];
145 const double bMax = box.bMax[0];
147 auto create_vertices = [nb_random_points, bMin, bMax,
148 scale_box](moab::Interface &moab,
auto &verts,
149 auto &arrays_coord) {
151 ReadUtilIface *iface;
153 CHKERR moab.query_interface(iface);
154 CHKERR iface->get_node_coords(3, nb_random_points, 0, startv,
157 arrays_coord[0], arrays_coord[1], arrays_coord[2]};
159 verts =
Range(startv, startv + nb_random_points - 1);
160 for (
int n = 0;
n != nb_random_points; ++
n) {
162 for (
auto ii : {1, 2}) {
163 t_coords(ii) = scale_box * (bMax - bMin) *
164 (std::rand() /
static_cast<double>(RAND_MAX)) -
172 auto set_positions = [nb_random_points,
dt, velocity_scale](
173 moab::Interface &moab,
auto &arrays_coord) {
176 arrays_coord[0], arrays_coord[1], arrays_coord[2]};
178 &init_pos(0, 0), &init_pos(1, 0), &init_pos(2, 0)};
180 for (
int n = 0;
n != nb_random_points; ++
n) {
183 for (
auto ii : {0, 1, 2})
184 t_velocity(ii) = (rand() /
static_cast<double>(RAND_MAX) - 0.5);
185 t_velocity(
i) /= std::sqrt(t_velocity(
i) * t_velocity(
i));
186 t_init_coords(
i) = t_coords(
i) +
dt * velocity_scale * t_velocity(
i);
194 moab::Core mb_charged_partices;
195 moab::Interface &moab_charged_partices = mb_charged_partices;
196 vector<double *> arrays_coord;
198 CHKERR create_vertices(moab_charged_partices, verts, arrays_coord);
199 auto init_positions = set_positions(moab, arrays_coord);
201 auto get_t_coords = [&]() {
203 arrays_coord[0], arrays_coord[1], arrays_coord[2]);
206 auto get_t_init_coords = [&]() {
208 &init_positions(0, 0), &init_positions(1, 0), &init_positions(2, 0));
211 auto calc_rhs = [&](
auto &t_p,
auto &t_init_p,
auto &t_B) {
216 t_rhs(
k) = (2 * t_p(
k) - t_init_p(
k)) /
dt -
217 levi_civita(
j,
i,
k) * t_init_p(
i) * t_B(
j);
221 auto calc_lhs = [&](
auto &t_B) {
226 t_lhs(
i,
k) = levi_civita(
j,
i,
k) * (-t_B(
j));
227 for (
auto ii : {0, 1, 2})
228 t_lhs(ii, ii) += 1 /
dt;
243 auto is_out = [&](
auto &t_p) {
244 for (
int i : {0, 1, 2})
247 }
else if (t_p(
i) < bMin) {
253 auto cache_ptr = boost::make_shared<CacheTuple>();
256 auto calc_position = [&]() {
257 auto t_p = get_t_coords();
258 auto t_init_p = get_t_init_coords();
262 for (
int n = 0;
n != nb_random_points; ++
n) {
270 std::array<double, 3> point = {t_p(0), t_p(1), t_p(2)};
271 data_at_pts->setEvalPoints(point.data(), 1);
274 point.data(), dist, prb_ptr->
getName(),
"MAGNETIC", data_at_pts,
281 for (
int ii : {0, 1, 2})
282 t_B(ii) = (*B)(ii, 0);
286 t_B(
i) *= magnetic_field_scale * 0.5;
288 auto t_rhs = calc_rhs(t_p, t_init_p, t_B);
289 auto t_lhs = calc_lhs(t_B);
296 t_init_p(
i) = t_p(
i);
297 t_p(
i) = t_inv_lhs(
i,
j) * t_rhs(
j);
306 for (
int t = 0;
t != nb_steps; ++
t) {
308 std::string print_step =
309 "Step : " + boost::lexical_cast<std::string>(
t) +
"\r";
310 std::cout << print_step << std::flush;
314 if ((
t % mod_step) == 0) {
316 CHKERR moab_charged_partices.write_file(
317 (
"step_" + boost::lexical_cast<std::string>(
t / mod_step) +
".vtk")
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
FTensor::Index< 'i', SPACE_DIM > i
const double n
refractive index of diffusive medium
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
constexpr double t
plate stiffness
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual int get_comm_rank() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorDouble & getFieldData() const
get dofs values
Field evaluator interface.
boost::shared_ptr< SPD > getData(const double *ptr=nullptr, const int nb_eval_points=0, const double eps=1e-12, VERBOSITY_LEVELS verb=QUIET)
Get the Data object.
MoFEMErrorCode buildTree(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode evalFEAtThePoint(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetPtsData > data_ptr, int lower_rank, int upper_rank, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
Evaluate field at artbitray position.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
RuleHookFun getRuleHook
Hook to get rule.
Calculate curl of vector field.
keeps basic data about problem
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Volume finite element base.
VolEle::UserDataOperator VolOp