|
| v0.14.0
|
Calculate directional derivative of the right hand side and compare it with tangent matrix derivative.
More...
#include <src/interfaces/OperatorsTester.hpp>
|
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
|
| OperatorsTester (const MoFEM::Core &core) |
|
virtual | ~OperatorsTester ()=default |
|
SmartPetscObj< Vec > | setRandomFields (SmartPetscObj< DM > dm, std::vector< RandomFieldData > random_fields, boost::shared_ptr< Range > ents=nullptr) |
| Generate random fields. More...
|
|
SmartPetscObj< Vec > | assembleVec (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr) |
| Assemble the right hand side vector. More...
|
|
SmartPetscObj< Mat > | assembleMat (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr) |
| Assemble the left hand side vector. More...
|
|
SmartPetscObj< Vec > | directionalCentralFiniteDifference (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, SmartPetscObj< Vec > diff_x, double time, double delta_t, double eps, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr()) |
| Calculate directional directive using finite difference. More...
|
|
SmartPetscObj< Vec > | checkCentralFiniteDifference (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline_rhs, boost::shared_ptr< FEMethod > pipeline_lhs, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, SmartPetscObj< Vec > diff_x, double time, double delta_t, double eps, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr()) |
| Check consistency between directional derivative with matrix. More...
|
|
template<class IFACE > |
MoFEMErrorCode | registerInterface (bool error_if_registration_failed=true) |
| Register interface. More...
|
|
template<class IFACE > |
MoFEMErrorCode | getInterface (IFACE *&iface) const |
| Get interface reference to pointer of interface. More...
|
|
template<class IFACE > |
MoFEMErrorCode | getInterface (IFACE **const iface) const |
| Get interface pointer to pointer of interface. More...
|
|
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0> |
IFACE | getInterface () const |
| Get interface pointer to pointer of interface. More...
|
|
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0> |
IFACE | getInterface () const |
| Get reference to interface. More...
|
|
template<class IFACE > |
IFACE * | getInterface () const |
| Function returning pointer to interface. More...
|
|
virtual | ~UnknownInterface ()=default |
|
Calculate directional derivative of the right hand side and compare it with tangent matrix derivative.
- Examples
- free_surface.cpp, level_set.cpp, operators_tests.cpp, plastic.cpp, and tensor_divergence_operator.cpp.
Definition at line 21 of file OperatorsTester.hpp.
◆ RandomFieldData
◆ OperatorsTester()
MoFEM::OperatorsTester::OperatorsTester |
( |
const MoFEM::Core & |
core | ) |
|
◆ ~OperatorsTester()
virtual MoFEM::OperatorsTester::~OperatorsTester |
( |
| ) |
|
|
virtualdefault |
◆ assembleMat()
Assemble the left hand side vector.
- Parameters
-
dm | |
fe_name | // fe name |
pipeline | // pipeline, i.e. fe instance |
x | // problem (dm) vector |
delta_x | // vector for x rate, can be null, i.e, SmartPetscObj<Vec>() |
delta2_x | // vector for x second rate, i.e. acceleration |
time | // time |
delta_t | // time increment |
cache_ptr | // finite element data cache, can be null |
- Returns
- SmartPetscObj<Vec>
Definition at line 140 of file OperatorsTester.cpp.
152 pipeline->ts_t = time;
154 pipeline->ts_dt = delta_t;
155 pipeline->ts_a = 1. / delta_t;
156 pipeline->ts_aa = 1. / pow(delta_t, 2);
158 pipeline->ksp_ctx = KspMethod::KSPContext::CTX_OPERATORS;
159 pipeline->snes_ctx = SnesMethod::SNESContext::CTX_SNESSETJACOBIAN;
160 pipeline->ts_ctx = TSMethod::TSContext::CTX_TSSETIJACOBIAN;
162 auto [
v,
a] =
setPipelineX(pipeline, x, delta_x, delta2_x, delta_t);
170 "Run assemble pipeline");
174 auto get_norm = [](
auto m) {
176 MatNorm(
m, NORM_1, &fnorm);
◆ assembleVec()
Assemble the right hand side vector.
- Parameters
-
dm | |
fe_name | // fe name |
pipeline | // pipeline, i.e. fe instance |
x | // problem (dm) vector |
delta_x | // vector for x rate, can be null, i.e, SmartPetscObj<Vec>() |
delta2_x | // vector for x second rate, i.e. acceleration |
time | // time |
delta_t | // time increment |
cache_ptr | // finite element data cache, can be null |
- Returns
- SmartPetscObj<Vec>
Definition at line 90 of file OperatorsTester.cpp.
99 pipeline->ts_t = time;
100 pipeline->ts_dt = delta_t;
102 pipeline->ksp_ctx = KspMethod::KSPContext::CTX_SETFUNCTION;
103 pipeline->snes_ctx = SnesMethod::SNESContext::CTX_SNESSETFUNCTION;
104 pipeline->ts_ctx = TSMethod::TSContext::CTX_TSSETIFUNCTION;
106 auto [
v,
a] =
setPipelineX(pipeline, x, delta_x, delta2_x, delta_t);
114 "Run assemble pipeline");
127 auto get_norm = [](
auto f) {
129 VecNorm(
f, NORM_2, &fnorm);
◆ checkCentralFiniteDifference()
SmartPetscObj< Vec > MoFEM::OperatorsTester::checkCentralFiniteDifference |
( |
SmartPetscObj< DM > |
dm, |
|
|
std::string |
fe_name, |
|
|
boost::shared_ptr< FEMethod > |
pipeline_rhs, |
|
|
boost::shared_ptr< FEMethod > |
pipeline_lhs, |
|
|
SmartPetscObj< Vec > |
x, |
|
|
SmartPetscObj< Vec > |
delta_x, |
|
|
SmartPetscObj< Vec > |
delta2_x, |
|
|
SmartPetscObj< Vec > |
diff_x, |
|
|
double |
time, |
|
|
double |
delta_t, |
|
|
double |
eps, |
|
|
CacheTupleWeakPtr |
cache_ptr = CacheTupleSharedPtr() |
|
) |
| |
Check consistency between directional derivative with matrix.
- Parameters
-
dm | |
fe_name | |
pipeline_rhs | |
pipeline_lhs | |
x | |
delta_x | |
delta2_x | |
diff_x | |
time | |
delta_t | |
eps | |
cache_ptr | |
- Returns
- SmartPetscObj<Vec>
Definition at line 243 of file OperatorsTester.cpp.
252 dm, fe_name, pipeline_rhs, x, delta_x, delta2_x, diff_x, time, delta_t,
255 auto m =
assembleMat(dm, fe_name, pipeline_lhs, x, delta_x, delta2_x, time,
◆ directionalCentralFiniteDifference()
SmartPetscObj< Vec > MoFEM::OperatorsTester::directionalCentralFiniteDifference |
( |
SmartPetscObj< DM > |
dm, |
|
|
std::string |
fe_name, |
|
|
boost::shared_ptr< FEMethod > |
pipeline, |
|
|
SmartPetscObj< Vec > |
x, |
|
|
SmartPetscObj< Vec > |
delta_x, |
|
|
SmartPetscObj< Vec > |
delta2_x, |
|
|
SmartPetscObj< Vec > |
diff_x, |
|
|
double |
time, |
|
|
double |
delta_t, |
|
|
double |
eps, |
|
|
CacheTupleWeakPtr |
cache_ptr = CacheTupleSharedPtr() |
|
) |
| |
Calculate directional directive using finite difference.
- Parameters
-
dm | |
fe_name | |
pipeline | |
x | |
delta_x | |
delta2_x | |
diff_x | // direction of derivative |
time | |
delta_t | |
eps | |
cache_ptr | |
- Returns
- SmartPetscObj<Vec>
Definition at line 186 of file OperatorsTester.cpp.
193 auto axpy = [&](
auto v0,
auto diff_v,
auto eps) {
194 SmartPetscObj<Vec>
v;
195 if (v0.use_count()) {
209 axpy(x, diff_x, -
eps),
211 axpy(delta_x, diff_x, -
eps),
213 axpy(delta2_x, diff_x, -
eps),
215 time, delta_t, cache_ptr);
219 axpy(x, diff_x,
eps),
221 axpy(delta_x, diff_x,
eps),
223 axpy(delta2_x, diff_x,
eps),
225 time, delta_t, cache_ptr);
227 auto calculate_derivative = [
eps](
auto fm1,
auto fp1) {
229 VecGetLocalSize(fp1, &size);
230 VecAXPY(fp1, -1, fm1);
232 VecGetArray(fp1, &array);
233 for (
int j = 0;
j != size; ++
j) {
236 VecRestoreArray(fp1, &array);
240 return calculate_derivative(fm1, fp1);
◆ query_interface()
◆ setPipelineX()
Set vectors, x, x_t, and x_tt to finite element instance.
Finite element instance is a pipeline. x_t and x_tt are evaluated for given delta_x, delta2_x and delta_t.
- Parameters
-
pipeline | |
x | |
delta_x | |
delta2_x | |
delta_t | |
- Returns
- std::pair<SmartPetscObj<Vec>, SmartPetscObj<Vec>>
Definition at line 266 of file OperatorsTester.cpp.
274 SmartPetscObj<Vec> x_t, x_tt;
280 VecCopy(delta_x, x_t);
281 VecScale(x_t, 1. / delta_t);
289 pipeline->x_t = PETSC_NULL;
296 VecCopy(delta2_x, x_tt);
297 VecScale(x_tt, 1. / pow(delta_t, 2));
302 pipeline->x_tt = x_tt;
305 pipeline->x_tt = PETSC_NULL;
308 return std::make_pair(x_t, x_tt);
◆ setRandomFields()
Generate random fields.
Example: generate random vector for DM (problem) from simple interface, where FIELD random values of DOFs are in range from -1 to 1, and FIELD2 random values are in range from 0 to 1.
auto x = opt->setRandomFields(
simple->getDM(),
{{"FIELD1", {-1, 1}}, {"FIELD2", {0,1}}});
TODO: Set random field to specific entities, and potentially order for testing proposes to dissect error in tangent matrix.
- Parameters
-
- Returns
- SmartPetscObj<Vec> smart vector
- Examples
- tensor_divergence_operator.cpp.
Definition at line 24 of file OperatorsTester.cpp.
32 auto get_random_number = [](
auto &range) {
33 const auto r =
static_cast<double>(std::rand()) / RAND_MAX;
34 return range[0] +
r * (range[1] - range[0]);
42 VecGetLocalSize(
v, &size);
46 for (
auto &rf : random_fields) {
49 <<
"Set random field " << rf.first <<
" " <<
static_cast<int>(field_id);
52 auto lo = prb_ptr->getNumeredRowDofsPtr()->lower_bound(
54 const auto hi = prb_ptr->getNumeredRowDofsPtr()->upper_bound(
56 for (; lo != hi; ++lo) {
57 const auto idx = (*lo)->getPetscLocalDofIdx();
58 if (idx >= 0 && idx < size) {
59 array[idx] = get_random_number(rf.second);
63 for (
auto pit = ents_ptr->const_pair_begin();
64 pit != ents_ptr->const_pair_end(); ++pit) {
66 prb_ptr->getNumeredRowDofsPtr()->get<Unique_mi_tag>().lower_bound(
69 prb_ptr->getNumeredRowDofsPtr()->get<Unique_mi_tag>().upper_bound(
71 for (; lo != hi; ++lo) {
72 const auto idx = (*lo)->getPetscLocalDofIdx();
73 if (idx >= 0 && idx < size) {
74 array[idx] = get_random_number(rf.second);
81 CHKERR VecRestoreArray(
v, &array);
83 CHKERR VecGhostUpdateBegin(
v, INSERT_VALUES, SCATTER_FORWARD);
84 CHKERR VecGhostUpdateEnd(
v, INSERT_VALUES, SCATTER_FORWARD);
◆ cOre
The documentation for this struct was generated from the following files:
OperatorsTester(const MoFEM::Core &core)
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
static constexpr Switches CtxSetB
static constexpr Switches CtxSetA
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
std::pair< SmartPetscObj< Vec >, SmartPetscObj< Vec > > setPipelineX(boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double delta_t)
Set vectors, x, x_t, and x_tt to finite element instance.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Deprecated interface functions.
SmartPetscObj< Mat > assembleMat(SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr)
Assemble the left hand side vector.
#define CHKERR
Inline error check.
auto createDMVector(DM dm)
Get smart vector from DM.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
static constexpr Switches CtxSetTime
SmartPetscObj< Vec > directionalCentralFiniteDifference(SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, SmartPetscObj< Vec > diff_x, double time, double delta_t, double eps, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Calculate directional directive using finite difference.
static UId getLoFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
const double v
phase velocity of light in medium (cm/ns)
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
SmartPetscObj< Vec > assembleVec(SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr)
Assemble the right hand side vector.
static UId getHiFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
FTensor::Index< 'j', 3 > j
auto getProblemPtr(DM dm)
get problem pointer from DM
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
FTensor::Index< 'm', 3 > m
static constexpr Switches CtxSetF
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.