v0.14.0
Functions | Variables
testing_jacobian_of_hook_element.cpp File Reference
#include <BasicFiniteElements.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "\n"
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
testing_jacobian_of_hook_element.cpp.

Definition at line 18 of file testing_jacobian_of_hook_element.cpp.

18  {
19 
20  // Initialize MoFEM
21  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
22 
23  // Create mesh database
24  moab::Core mb_instance; // create database
25  moab::Interface &moab = mb_instance; // create interface to database
26 
27  try {
28  // Create MoFEM database and link it to MoAB
29  MoFEM::Core core(moab);
30  MoFEM::Interface &m_field = core;
31 
32  PetscBool ale = PETSC_FALSE;
33  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-ale", &ale, PETSC_NULL);
34  PetscBool test_jacobian = PETSC_FALSE;
35  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-test_jacobian", &test_jacobian,
36  PETSC_NULL);
37 
38  CHKERR DMRegister_MoFEM("DMMOFEM");
39 
40  Simple *si = m_field.getInterface<MoFEM::Simple>();
41 
42  CHKERR si->getOptions();
43  CHKERR si->loadFile();
45  const int order = 2;
46  CHKERR si->setFieldOrder("x", order);
47 
48  if (ale == PETSC_TRUE) {
50  CHKERR si->setFieldOrder("X", 2);
51  }
52 
53  CHKERR si->setUp();
54 
55  // create DM
56  DM dm;
57  CHKERR si->getDM(&dm);
58 
59  // Projection on "x" field
60  {
61  Projection10NodeCoordsOnField ent_method(m_field, "x");
62  CHKERR m_field.loop_dofs("x", ent_method);
63  // CHKERR m_field.getInterface<FieldBlas>()->fieldScale(2, "x");
64  }
65 
66  // Project coordinates on "X" field
67  if (ale == PETSC_TRUE) {
68  Projection10NodeCoordsOnField ent_method(m_field, "X");
69  CHKERR m_field.loop_dofs("X", ent_method);
70  // CHKERR m_field.getInterface<FieldBlas>()->fieldScale(2, "X");
71  }
72 
73  boost::shared_ptr<ForcesAndSourcesCore> fe_lhs_ptr(
75  boost::shared_ptr<ForcesAndSourcesCore> fe_rhs_ptr(
77  struct VolRule {
78  int operator()(int, int, int) const { return 2 * (order - 1); }
79  };
80  fe_lhs_ptr->getRuleHook = VolRule();
81  fe_rhs_ptr->getRuleHook = VolRule();
82 
83  CHKERR DMMoFEMSNESSetJacobian(dm, si->getDomainFEName(), fe_lhs_ptr,
84  nullptr, nullptr);
85  CHKERR DMMoFEMSNESSetFunction(dm, si->getDomainFEName(), fe_rhs_ptr,
86  nullptr, nullptr);
87 
89  boost::shared_ptr<map<int, BlockData>> block_sets_ptr =
90  boost::make_shared<map<int, BlockData>>();
91  (*block_sets_ptr)[0].iD = 0;
92  (*block_sets_ptr)[0].E = 1;
93  (*block_sets_ptr)[0].PoissonRatio = 0.25;
95  si->getDomainFEName(), 3, (*block_sets_ptr)[0].tEts);
96 
97  CHKERR HookeElement::setOperators(fe_lhs_ptr, fe_rhs_ptr, block_sets_ptr,
98  "x", "X", ale, false);
99 
100  Vec x, f;
101  CHKERR DMCreateGlobalVector(dm, &x);
102  CHKERR VecDuplicate(x, &f);
103  CHKERR DMoFEMMeshToLocalVector(dm, x, INSERT_VALUES, SCATTER_FORWARD);
104 
105  // CHKERR VecDuplicate(x, &dx);
106  // PetscRandom rctx;
107  // PetscRandomCreate(PETSC_COMM_WORLD, &rctx);
108  // VecSetRandom(dx, rctx);
109  // PetscRandomDestroy(&rctx);
110  // CHKERR DMoFEMMeshToGlobalVector(dm, x, INSERT_VALUES, SCATTER_REVERSE);
111 
112  Mat A, fdA;
113  CHKERR DMCreateMatrix(dm, &A);
114  CHKERR MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &fdA);
115 
116  if (test_jacobian == PETSC_TRUE) {
117  char testing_options[] =
118  "-snes_test_jacobian -snes_test_jacobian_display "
119  "-snes_no_convergence_test -snes_atol 0 -snes_rtol 0 -snes_max_it 1 "
120  "-pc_type none";
121  CHKERR PetscOptionsInsertString(NULL, testing_options);
122  } else {
123  char testing_options[] = "-snes_no_convergence_test -snes_atol 0 "
124  "-snes_rtol 0 -snes_max_it 1 -pc_type none";
125  CHKERR PetscOptionsInsertString(NULL, testing_options);
126  }
127 
128  SNES snes;
129  CHKERR SNESCreate(PETSC_COMM_WORLD, &snes);
130  MoFEM::SnesCtx *snes_ctx;
131  CHKERR DMMoFEMGetSnesCtx(dm, &snes_ctx);
132  CHKERR SNESSetFunction(snes, f, SnesRhs, snes_ctx);
133  CHKERR SNESSetJacobian(snes, A, A, SnesMat, snes_ctx);
134  CHKERR SNESSetFromOptions(snes);
135 
136  CHKERR SNESSolve(snes, NULL, x);
137 
138  if (test_jacobian == PETSC_FALSE) {
139  double nrm_A0;
140  CHKERR MatNorm(A, NORM_INFINITY, &nrm_A0);
141 
142  char testing_options_fd[] = "-snes_fd";
143  CHKERR PetscOptionsInsertString(NULL, testing_options_fd);
144 
145  CHKERR SNESSetFunction(snes, f, SnesRhs, snes_ctx);
146  CHKERR SNESSetJacobian(snes, fdA, fdA, SnesMat, snes_ctx);
147  CHKERR SNESSetFromOptions(snes);
148 
149  CHKERR SNESSolve(snes, NULL, x);
150  CHKERR MatAXPY(A, -1, fdA, SUBSET_NONZERO_PATTERN);
151 
152  double nrm_A;
153  CHKERR MatNorm(A, NORM_INFINITY, &nrm_A);
154  PetscPrintf(PETSC_COMM_WORLD, "Matrix norms %3.4e %3.4e\n", nrm_A,
155  nrm_A / nrm_A0);
156  nrm_A /= nrm_A0;
157 
158  const double tol = 1e-5;
159  if (nrm_A > tol) {
160  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
161  "Difference between hand-calculated tangent matrix and finite "
162  "difference matrix is too big");
163  }
164  }
165 
166  CHKERR VecDestroy(&x);
167  CHKERR VecDestroy(&f);
168  CHKERR MatDestroy(&A);
169  CHKERR MatDestroy(&fdA);
170  CHKERR SNESDestroy(&snes);
171 
172  // destroy DM
173  CHKERR DMDestroy(&dm);
174  }
175  CATCH_ERRORS;
176 
177  // finish work cleaning memory, getting statistics, etc
179 
180  return 0;
181 }

Variable Documentation

◆ help

char help[] = "\n"
static
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreInterface::loop_dofs
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.
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
VolRule
Set integration rule.
Definition: simple_interface.cpp:88
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:527
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::Simple::addDomainField
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.
Definition: Simple.cpp:264
MoFEM::DMMoFEMGetSnesCtx
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition: DMMoFEM.cpp:1098
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
MoFEM::Simple::getDomainFEName
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition: Simple.hpp:341
MoFEM::DMMoFEMSNESSetJacobian
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition: DMMoFEM.cpp:763
VolRule::operator()
int operator()(int, int, int) const
Definition: simple_interface.cpp:89
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
help
static char help[]
Definition: testing_jacobian_of_hook_element.cpp:16
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEM::CoreInterface::get_finite_element_entities_by_dimension
virtual MoFEMErrorCode get_finite_element_entities_by_dimension(const std::string name, int dim, Range &ents) const =0
get entities in the finite element by dimension
MoFEM::SnesRhs
PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition: SnesCtx.cpp:27
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
NonlinearElasticElement::BlockData::iD
int iD
Definition: HookeElement.hpp:33
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::BlockData
Definition: MeshsetsManager.cpp:752
setOperators
static MoFEMErrorCode setOperators(boost::shared_ptr< ForcesAndSourcesCore > fe_lhs_ptr, boost::shared_ptr< ForcesAndSourcesCore > fe_rhs_ptr, boost::shared_ptr< map< int, BlockData >> block_sets_ptr, const std::string x_field, const std::string X_field, const bool ale, const bool field_disp, const EntityType type=MBTET, boost::shared_ptr< DataAtIntegrationPts > data_at_pts=nullptr)
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::SnesCtx
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:13
MoFEM::SnesMat
PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition: SnesCtx.cpp:139
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
BlockData
NonlinearElasticElement::BlockData BlockData
Definition: elasticity.cpp:74
tol
double tol
Definition: mesh_smoothing.cpp:27
MoFEM::DMMoFEMSNESSetFunction
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition: DMMoFEM.cpp:722
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182