19static char help[] =
"...\n\n";
45int main(
int argc,
char *argv[]) {
51 enum bases { AINSWORTH, DEMKOWICZ, LASTOP };
53 const char *list[] = {
"ainsworth",
"demkowicz"};
56 PetscInt choise_value = AINSWORTH;
59 if (flg != PETSC_TRUE) {
63 PetscBool ho_geometry = PETSC_FALSE;
67 DMType dm_name =
"DMMOFEM";
71 moab::Core mb_instance;
72 moab::Interface &moab = mb_instance;
84 switch (choise_value) {
99 if (ho_geometry == PETSC_TRUE)
103 constexpr int order = 3;
109 if (ho_geometry == PETSC_TRUE)
120 auto material_grad_mat = boost::make_shared<MatrixDouble>();
121 auto material_det_vec = boost::make_shared<VectorDouble>();
122 auto material_inv_grad_mat = boost::make_shared<MatrixDouble>();
123 auto jac_ptr = boost::make_shared<MatrixDouble>();
124 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
125 auto det_ptr = boost::make_shared<VectorDouble>();
127 boost::dynamic_pointer_cast<VolumeElementForcesAndSourcesCore>(
128 pipeline_mng->getDomainRhsFE())
129 ->meshPositionsFieldName =
"none";
130 boost::dynamic_pointer_cast<PipelineManager::FaceEle>(
131 pipeline_mng->getBoundaryRhsFE())
132 ->meshPositionsFieldName =
"none";
134 pipeline_mng->getOpDomainRhsPipeline().push_back(
136 pipeline_mng->getOpDomainRhsPipeline().push_back(
138 pipeline_mng->getOpDomainRhsPipeline().push_back(
140 pipeline_mng->getOpDomainRhsPipeline().push_back(
142 pipeline_mng->getOpDomainRhsPipeline().push_back(
146 pipeline_mng->getOpDomainRhsPipeline().push_back(
150 material_grad_mat, material_det_vec, material_inv_grad_mat));
151 pipeline_mng->getOpDomainRhsPipeline().push_back(
153 pipeline_mng->getOpDomainRhsPipeline().push_back(
155 pipeline_mng->getOpDomainRhsPipeline().push_back(
158 pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpVolCurl(t_curl_vol));
161 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
163 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
165 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
173 if (ho_geometry == PETSC_TRUE) {
179 CHKERR pipeline_mng->loopFiniteElements();
181 std::cout.precision(12);
183 t_curl_vol(
i) -= t_curl_skin(
i);
184 double nrm2 = sqrt(t_curl_vol(
i) * t_curl_vol(
i));
186 constexpr double eps = 1e-8;
187 if (fabs(nrm2) >
eps)
189 "Curl operator not passed test\n");
203 const unsigned int nb_gauss_pts = data.
getDiffN().size1();
204 const unsigned int nb_dofs = data.
getFieldData().size();
214 for (; gg < nb_gauss_pts; gg++) {
215 double w = getGaussPts()(3, gg) * getVolume();
217 for (
unsigned int dd = 0; dd != nb_dofs; dd++) {
218 t_curl(
i) = levi_civita(
j,
i,
k) * t_curl_base(
j,
k);
234 int nb_gauss_pts = data.
getN().size1();
241 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
242 for (
int dd = 0; dd < nb_dofs; dd++) {
243 double w = getGaussPts()(2, gg);
244 const double n0 = getNormalsAtGaussPts(gg)[0];
245 const double n1 = getNormalsAtGaussPts(gg)[1];
246 const double n2 = getNormalsAtGaussPts(gg)[2];
247 if (getFEType() == MBTRI) {
251 tCurl(0) += (n1 * t_curl_base(2) - n2 * t_curl_base(1)) * w;
252 tCurl(1) += (n2 * t_curl_base(0) - n0 * t_curl_base(2)) * w;
253 tCurl(2) += (n0 * t_curl_base(1) - n1 * t_curl_base(0)) * w;
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#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()
@ HCURL
field with continuous tangents
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
virtual bool check_field(const std::string &name) const =0
check if field is in database
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
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.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
get dofs values
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
default operator for TRI element
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Calculate normals at Gauss points of triangle element.
transform local reference derivatives of shape function to global derivatives if higher order geometr...
Set inverse jacobian to base functions.
PipelineManager interface.
Projection of edge entities with one mid-node on hierarchical basis.
Simple interface for fast problem set-up.
MoFEMErrorCode addDataField(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_ZERO, int verb=-1)
Add field on domain.
MoFEMErrorCode addDomainField(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_ZERO, int verb=-1)
Add field on domain.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name)
Load mesh file.
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
MoFEMErrorCode addBoundaryField(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_ZERO, int verb=-1)
Add field on boundary.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Volume finite element base.
OpFacesRot(FTensor::Tensor1< double, 3 > &t_curl)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
FTensor::Tensor1< double, 3 > & tCurl
OpVolCurl(FTensor::Tensor1< double, 3 > &t_curl)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
FTensor::Tensor1< double, 3 > & tCurl