v0.10.0
hcurl_curl_operator.cpp

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.

/**
* \file hcurl_curl_operator.cpp
* \brief Testich curl-curl operator by applying Stokes theorem
* \example hcurl_curl_operator.cpp
*
* 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.
*
*/
/* This file is part of MoFEM.
* MoFEM 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.
*
* MoFEM 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 MoFEM. If not, see <http://www.gnu.org/licenses/>. */
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
"HCURL", UserDataOperator::OPROW),
tCurl(t_curl) {}
MoFEMErrorCode doWork(int side, EntityType type,
};
"HCURL", UserDataOperator::OPROW),
tCurl(t_curl) {}
MoFEMErrorCode doWork(int side, EntityType type,
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
enum bases { AINSWORTH, DEMKOWICZ, LASTOP };
const char *list[] = {"ainsworth", "demkowicz"};
PetscBool flg;
PetscInt choise_value = AINSWORTH;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list, LASTOP,
&choise_value, &flg);
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "base not set");
}
PetscBool ho_geometry = PETSC_FALSE;
CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-ho_geometry", &ho_geometry,
PETSC_NULL);
DMType dm_name = "DMMOFEM";
// Create MoAB
moab::Core mb_instance; ///< database
moab::Interface &moab = mb_instance; ///< interface
// Create MoFEM
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
Simple *simple_interface = m_field.getInterface<Simple>();
PipelineManager *pipeline_mng = m_field.getInterface<PipelineManager>();
CHKERR simple_interface->getOptions();
CHKERR simple_interface->loadFile("");
// fields
switch (choise_value) {
case AINSWORTH:
CHKERR simple_interface->addDomainField("HCURL", HCURL,
CHKERR simple_interface->addBoundaryField("HCURL", HCURL,
break;
case DEMKOWICZ:
CHKERR simple_interface->addDomainField("HCURL", HCURL,
CHKERR simple_interface->addBoundaryField("HCURL", HCURL,
break;
}
if (ho_geometry == PETSC_TRUE)
CHKERR simple_interface->addDataField("MESH_NODE_POSITIONS", H1,
constexpr int order = 5;
CHKERR simple_interface->setFieldOrder("HCURL", order);
if (ho_geometry == PETSC_TRUE)
CHKERR simple_interface->setFieldOrder("MESH_NODE_POSITIONS", 2);
CHKERR simple_interface->setUp();
auto integration_rule = [](int, int, int p_data) { return 2 * p_data; };
CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpTetCurl(t_curl_vol));
pipeline_mng->getOpBoundaryRhsPipeline().push_back(
new OpFacesRot(t_curl_skin));
t_curl_vol(i) = 0;
t_curl_skin(i) = 0;
// project geometry form 10 node tets on higher order approx. functions
if (ho_geometry == PETSC_TRUE) {
Projection10NodeCoordsOnField ent_method(m_field, "MESH_NODE_POSITIONS");
CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method);
}
// Run pipelines on mesh
CHKERR pipeline_mng->loopFiniteElements();
std::cout.precision(12);
std::cout << "curl_vol " << t_curl_vol << std::endl;
std::cout << "curl_skin " << t_curl_skin << std::endl;
t_curl_vol(i) -= t_curl_skin(i);
double nrm2 = sqrt(t_curl_vol(i) * t_curl_vol(i));
constexpr double eps = 1e-8;
if (fabs(nrm2) > eps)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Curl operator not passed test\n");
}
}
MoFEMErrorCode OpTetCurl::doWork(int side, EntityType type,
if (data.getFieldData().size() == 0)
const unsigned int nb_gauss_pts = data.getDiffN().size1();
const unsigned int nb_dofs = data.getFieldData().size();
MatrixDouble curl_mat;
unsigned int gg = 0;
for (; gg < nb_gauss_pts; gg++) {
double w = getGaussPts()(3, gg) * getVolume();
if (getHoGaussPtsDetJac().size() == nb_gauss_pts) {
// if ho geometry is given
w *= getHoGaussPtsDetJac()(gg);
}
CHKERR getCurlOfHCurlBaseFunctions(side, type, data, gg, curl_mat);
FTensor::Tensor1<double *, 3> t_curl(&curl_mat(0, 0), &curl_mat(0, 1),
&curl_mat(0, 2), 3);
for (unsigned int dd = 0; dd != nb_dofs; dd++) {
tCurl(i) += w * t_curl(i);
++t_curl;
}
}
}
MoFEMErrorCode OpFacesRot::doWork(int side, EntityType type,
int nb_dofs = data.getFieldData().size();
if (nb_dofs == 0)
int nb_gauss_pts = data.getN().size1();
auto t_curl_base = data.getFTensor1N<3>();
// double area = getArea();
double n0 = getNormal()[0] * 0.5;
double n1 = getNormal()[1] * 0.5;
double n2 = getNormal()[2] * 0.5;
for (int gg = 0; gg < nb_gauss_pts; gg++) {
for (int dd = 0; dd < nb_dofs; dd++) {
double w = getGaussPts()(2, gg);
if (getNormalsAtGaussPts().size1() == (unsigned int)nb_gauss_pts) {
n0 = getNormalsAtGaussPts(gg)[0] * 0.5;
n1 = getNormalsAtGaussPts(gg)[1] * 0.5;
n2 = getNormalsAtGaussPts(gg)[2] * 0.5;
}
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.
Definition: definitions.h:441
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:158
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
@ H1
continuous field
Definition: definitions.h:177
@ HCURL
field with continuous tangents
Definition: definitions.h:178
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
@ MOFEM_IMPOSIBLE_CASE
Definition: definitions.h:127
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:132
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define CHKERR
Inline error check.
Definition: definitions.h:604
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:48
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.
int main(int argc, char *argv[])
static char help[]
FTensor::Index< 'j', 3 > j
FTensor::Index< 'i', 3 > i
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)
Definition: ddTensor0.hpp:33
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:76
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
CoreTmp< 0 > Core
Definition: Core.hpp:1129
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1943
double w(const double g, const double t)
DataForcesAndSourcesCore::EntData EntData
Core (interface) class.
Definition: Core.hpp:77
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:60
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:100
Deprecated interface functions.
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
VolumeElementForcesAndSourcesCoreBase::UserDataOperator UserDataOperator
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Operator for linear form, usually to calculate values on right hand side.