22 namespace bio = boost::iostreams;
24 using bio::tee_device;
26 using namespace MoFEM;
28 static char help[] =
"...\n\n";
30 constexpr
double eps = 1e-6;
33 if (std::abs(v) <
eps)
37 int main(
int argc,
char *argv[]) {
46 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
48 PetscBool flg = PETSC_TRUE;
50 #if PETSC_VERSION_GE(3, 6, 4)
57 if (flg != PETSC_TRUE) {
58 SETERRQ(PETSC_COMM_SELF, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
65 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab,
MYPCOMM_INDEX);
67 pcomm =
new ParallelComm(&moab, PETSC_COMM_WORLD);
79 EntityHandle meshset_level0;
80 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
107 EntityHandle root_set = moab.get_root_set();
153 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
163 : ofs(
"forces_and_sources_getting_mult_H1_H1_atom_test.txt"),
164 my_tee(std::cout, ofs), my_split(my_tee){};
166 ~my_mult_H1_H1() { my_split.close(); }
168 ublas::matrix<FieldData> NN;
170 MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
178 int nb_row_dofs = row_data.
getN().size2();
179 int nb_col_dofs = col_data.
getN().size2();
181 my_split << row_side <<
" " << col_side <<
" " << row_type <<
" "
182 << col_type << std::endl;
183 my_split <<
"nb_row_dofs " << nb_row_dofs <<
" nb_col_dofs "
184 << nb_col_dofs << std::endl;
185 NN.resize(nb_row_dofs, nb_col_dofs);
190 my_split << std::setprecision(3);
191 my_split << std::fixed;
192 my_split << row_data.
getN() << std::endl;
193 my_split << col_data.
getN() << std::endl;
195 for (
unsigned int gg = 0; gg < row_data.
getN().size1(); gg++) {
197 bzero(&*NN.data().begin(),
198 nb_row_dofs * nb_col_dofs *
sizeof(
FieldData));
201 &row_data.
getN()(gg, 0), 1, &col_data.
getN()(gg, 0), 1,
202 &*NN.data().begin(), nb_col_dofs);
204 my_split <<
"gg " << gg <<
" : ";
205 my_split << std::setprecision(3);
206 my_split << std::fixed;
209 NN - outer_prod(row_data.
getN(gg), col_data.
getN(gg));
212 my_split << difference << std::endl;
213 if (row_type != MBVERTEX) {
214 my_split << row_data.
getDiffN(gg) << std::endl;
217 if (row_type == MBVERTEX) {
218 my_split << row_data.
getDiffN() << std::endl;
220 typedef ublas::array_adaptor<FieldData> storage_t;
221 storage_t st(nb_row_dofs * 3, &row_data.
getDiffN()(gg, 0));
222 ublas::matrix<FieldData, ublas::row_major, storage_t>
223 digNatGaussPt(nb_row_dofs, 3, st);
224 my_split << std::endl << digNatGaussPt << std::endl;
228 my_split << std::endl;
249 CHKERR getSpacesAndBaseOnEntities(data_row);
250 CHKERR getSpacesAndBaseOnEntities(data_col);
252 CHKERR getEntitySense<MBEDGE>(data_row);
253 CHKERR getEntitySense<MBTRI>(data_row);
254 CHKERR getEntitySense<MBEDGE>(data_col);
255 CHKERR getEntitySense<MBTRI>(data_col);
257 CHKERR getEntityDataOrder<MBEDGE>(data_row,
H1);
258 CHKERR getEntityDataOrder<MBEDGE>(data_col,
H1);
259 CHKERR getEntityDataOrder<MBTRI>(data_row,
H1);
260 CHKERR getEntityDataOrder<MBTRI>(data_col,
H1);
261 CHKERR getEntityDataOrder<MBTET>(data_row,
H1);
262 CHKERR getEntityDataOrder<MBTET>(data_col,
H1);
265 CHKERR getEntityFieldData(data_row,
"FIELD1", MBEDGE);
268 CHKERR getEntityFieldData(data_col,
"FIELD2", MBEDGE);
269 CHKERR getRowNodesIndices(data_row,
"FIELD1");
270 CHKERR getColNodesIndices(data_col,
"FIELD2");
271 CHKERR getEntityRowIndices(data_row,
"FIELD1", MBEDGE);
272 CHKERR getEntityColIndices(data_col,
"FIELD2", MBEDGE);
273 CHKERR getFaceTriNodes(data_row);
274 CHKERR getFaceTriNodes(data_col);
277 for (
int gg = 0; gg < 4; gg++) {
292 CHKERR op.opLhs(data_row, data_col);
304 ForcesAndSourcesCore_TestFE fe1(m_field);
void cblas_dger(const enum CBLAS_ORDER order, const int M, const int N, const double alpha, const double *X, const int incX, const double *Y, const int incY, double *A, const int lda)
#define BARRIER_PCOMM_RANK_END(PCMB)
set barrier start Run code in sequence, starting from process 0, and ends on last process.
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#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 ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
#define BARRIER_PCOMM_RANK_START(PCMB)
set barrier start Run code in sequence, starting from process 0, and ends on last process.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
int main(int argc, char *argv[])
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string &name_row)=0
set field row which finite element use
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string &name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string &name_row)=0
set field col which finite element use
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string &name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
double FieldData
Field data type.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
implementation of Data Operators for Forces and Sources
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
DeprecatedCoreInterface Interface
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
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.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
FieldApproximationBase & getBase()
Get approximation base.
data structure for finite element entity
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
base operator to do operations at Gauss Pt. level
Deprecated interface functions.
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
Class used to pass element data to calculate base functions on tet,triangle,edge.
structure to get information form mofem into DataForcesAndSourcesCore
Problem manager is used to build and partition problems.
Calculate base functions on tetrahedral.
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.