Using PipelineManager interface calculate the curl of base functions, and integral of the vector tangent vector with normal on the boundary. Since the h-curl space is used, volume integral and boundary integral should give the same result, as a result, as we are applying Stokes theorem on h-curl space.
static char help[] =
"...\n\n";
"HCURL", UserDataOperator::OPROW),
};
};
int main(
int argc,
char *argv[]) {
try {
enum bases { AINSWORTH, DEMKOWICZ,
LASTOP };
const char *list[] = {"ainsworth", "demkowicz"};
PetscBool flg;
PetscInt choise_value = AINSWORTH;
&choise_value, &flg);
if (flg != PETSC_TRUE) {
}
PetscBool ho_geometry = PETSC_FALSE;
PETSC_NULL);
DMType dm_name = "DMMOFEM";
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
switch (choise_value) {
case AINSWORTH:
break;
case DEMKOWICZ:
break;
}
if (ho_geometry == PETSC_TRUE)
if (ho_geometry == PETSC_TRUE)
auto material_grad_mat = boost::make_shared<MatrixDouble>();
auto material_det_vec = boost::make_shared<VectorDouble>();
auto material_inv_grad_mat = boost::make_shared<MatrixDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
auto det_ptr = boost::make_shared<VectorDouble>();
boost::dynamic_pointer_cast<VolumeElementForcesAndSourcesCore>(
pipeline_mng->getDomainRhsFE())
->meshPositionsFieldName = "none";
boost::dynamic_pointer_cast<PipelineManager::FaceEle>(
pipeline_mng->getBoundaryRhsFE())
->meshPositionsFieldName = "none";
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
if (ho_geometry) {
pipeline_mng->getOpDomainRhsPipeline().push_back(
material_grad_mat));
material_grad_mat, material_det_vec, material_inv_grad_mat));
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
}
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpVolCurl(t_curl_vol));
pipeline_mng->getOpBoundaryRhsPipeline().push_back(
pipeline_mng->getOpBoundaryRhsPipeline().push_back(
pipeline_mng->getOpBoundaryRhsPipeline().push_back(
if (ho_geometry == PETSC_TRUE) {
}
CHKERR pipeline_mng->loopFiniteElements();
std::cout.precision(12);
t_curl_vol(
i) -= t_curl_skin(
i);
double nrm2 = sqrt(t_curl_vol(
i) * t_curl_vol(
i));
constexpr double eps = 1e-8;
"Curl operator not passed test\n");
}
}
const unsigned int nb_gauss_pts = data.
getDiffN().size1();
unsigned int gg = 0;
for (; gg < nb_gauss_pts; gg++) {
double w = getGaussPts()(3, gg) * getVolume();
for (
unsigned int dd = 0;
dd != nb_dofs;
dd++) {
++t_curl_base;
}
}
}
if (nb_dofs == 0)
int nb_gauss_pts = data.
getN().size1();
for (int gg = 0; gg < nb_gauss_pts; gg++) {
for (
int dd = 0;
dd < nb_dofs;
dd++) {
double w = getGaussPts()(2, gg);
const double n0 = getNormalsAtGaussPts(gg)[0];
const double n1 = getNormalsAtGaussPts(gg)[1];
const double n2 = getNormalsAtGaussPts(gg)[2];
if (getFEType() == MBTRI) {
}
tCurl(0) += (n1 * t_curl_base(2) - n2 * t_curl_base(1)) * w;
tCurl(1) += (n2 * t_curl_base(0) - n0 * t_curl_base(2)) * w;
tCurl(2) += (n0 * t_curl_base(1) - n1 * t_curl_base(0)) * w;
++t_curl_base;
}
}
}
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
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
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
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.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
FTensor::Tensor1< double, 3 > & tCurl
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
FTensor::Tensor1< double, 3 > & tCurl