v0.14.0
wave_equation.cpp
Go to the documentation of this file.
1 /**
2  * \file wave_equation.cpp
3  * \example wave_equation.cpp
4  *
5  * \brief Solve the time-dependent Wave Equation
6  \f[
7  \begin{aligned}
8 \frac{\partial u(\mathbf{x}, t)}{\partial t}-\Delta u(\mathbf{x}, t)
9 &=f(\mathbf{x}, t) & & \forall \mathbf{x} \in \Omega, t \in(0, T), \\
10 u(\mathbf{x}, 0) &=u_{0}(\mathbf{x}) & & \forall \mathbf{x} \in \Omega, \\
11 u(\mathbf{x}, t) &=g(\mathbf{x}, t) & & \forall \mathbf{x} \in \partial \Omega,
12 t \in(0, T). \end{aligned}
13  \f]
14  **/
15 
16 
17 
18 #include <stdlib.h>
19 #include <cmath>
20 #include <MoFEM.hpp>
21 
22 using namespace MoFEM;
23 
24 static char help[] = "...\n\n";
25 
26 template <int DIM> struct ElementsAndOps {};
27 
28 //! [Define dimension]
29 constexpr int SPACE_DIM = 2; //< Space dimension of problem, mesh
30 //! [Define dimension]
31 
38 
44  PETSC>::LinearForm<GAUSS>::OpBaseTimesScalar<1>;
47 
51  PETSC>::LinearForm<GAUSS>::OpBaseTimesScalar<1>;
53  PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
54 
55 constexpr double wave_speed2 = 1;
56 
57 /**
58  * @brief Monitor solution
59  *
60  * This functions is called by TS solver at the end of each step. It is used
61  * to output results to the hard drive.
62  */
63 struct Monitor : public FEMethod {
64 
65  Monitor(SmartPetscObj<DM> dm, boost::shared_ptr<PostProcEle> post_proc)
66  : dM(dm), postProc(post_proc){};
67 
68  MoFEMErrorCode preProcess() { return 0; }
69  MoFEMErrorCode operator()() { return 0; }
70 
71  static constexpr int saveEveryNthStep = 1;
72 
75  if (ts_step % saveEveryNthStep == 0) {
76  CHKERR DMoFEMLoopFiniteElements(dM, "dFE", postProc);
77  CHKERR postProc->writeFile(
78  "out_level_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
79  }
81  }
82 
83 private:
85  boost::shared_ptr<PostProcEle> postProc;
86 };
87 
88 struct WaveEquation {
89 public:
91 
92  // Declaration of the main function to run analysis
93  MoFEMErrorCode runProgram();
94 
95 private:
96  // Declaration of other main functions called in runProgram()
97  MoFEMErrorCode readMesh();
98  MoFEMErrorCode setupProblem();
99  MoFEMErrorCode setIntegrationRules();
100  MoFEMErrorCode initialCondition();
101  MoFEMErrorCode boundaryCondition();
102  MoFEMErrorCode assembleSystem();
103  MoFEMErrorCode solveSystem();
104  MoFEMErrorCode outputResults();
105 
106  // Main interfaces
108 
109  // Object to mark boundary entities for the assembling of domain elements
110  boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
111 };
112 
113 WaveEquation::WaveEquation(MoFEM::Interface &m_field) : mField(m_field) {}
114 
117 
118  auto *simple = mField.getInterface<Simple>();
120  CHKERR simple->getOptions();
121  CHKERR simple->loadFile();
122 
124 }
125 
128 
129  auto *simple = mField.getInterface<Simple>();
130  CHKERR simple->addDomainField("U", H1, AINSWORTH_LEGENDRE_BASE, 1);
131  CHKERR simple->addBoundaryField("U", H1, AINSWORTH_LEGENDRE_BASE, 1);
132 
133  int order = 3;
134  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
135  CHKERR simple->setFieldOrder("U", order);
136 
137  CHKERR simple->setUp();
138 
140 }
141 
144 
145  auto integration_rule = [](int o_row, int o_col, int approx_order) {
146  return 2 * approx_order;
147  };
148 
149  auto *pipeline_mng = mField.getInterface<PipelineManager>();
150  CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
151  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
152  CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
153  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
154 
156 }
157 
161 }
162 
165 
166  auto bc_mng = mField.getInterface<BcManager>();
167  auto *simple = mField.getInterface<Simple>();
168  CHKERR bc_mng->pushMarkDOFsOnEntities(simple->getProblemName(), "ESSENTIAL",
169  "U", 0, 0);
170 
171  auto &bc_map = bc_mng->getBcMapByBlockName();
172  boundaryMarker = boost::make_shared<std::vector<char unsigned>>();
173  for (auto b : bc_map) {
174  if (std::regex_match(b.first, std::regex("(.*)ESSENTIAL(.*)"))) {
175  boundaryMarker->resize(b.second->bcMarkers.size(), 0);
176  for (int i = 0; i != b.second->bcMarkers.size(); ++i) {
177  (*boundaryMarker)[i] |= b.second->bcMarkers[i];
178  }
179  }
180  }
181 
183 }
184 
187 
188  auto add_domain_base_ops = [&](auto &pipeline) {
189  auto det_ptr = boost::make_shared<VectorDouble>();
190  auto jac_ptr = boost::make_shared<MatrixDouble>();
191  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
192  pipeline.push_back(new OpCalculateHOJac<2>(jac_ptr));
193  pipeline.push_back(new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
194  pipeline.push_back(new OpSetHOInvJacToScalarBases<2>(H1, inv_jac_ptr));
195  pipeline.push_back(new OpSetHOWeightsOnFace());
196  };
197 
198  auto add_domain_lhs_ops = [&](auto &pipeline) {
199  pipeline.push_back(new OpSetBc("U", true, boundaryMarker));
200  pipeline.push_back(new OpDomainGradGrad(
201  "U", "U", [](double, double, double) -> double { return 1; }));
202  auto get_c = [this](const double, const double, const double) {
203  auto pipeline_mng = mField.getInterface<PipelineManager>();
204  auto &fe_domain_lhs = pipeline_mng->getDomainLhsFE();
205  return wave_speed2 * fe_domain_lhs->ts_aa;
206  };
207  pipeline.push_back(new OpDomainMass("U", "U", get_c));
208  pipeline.push_back(new OpUnSetBc("U"));
209  };
210 
211  auto add_domain_rhs_ops = [&](auto &pipeline) {
212  pipeline.push_back(new OpSetBc("U", true, boundaryMarker));
213  auto grad_u_at_gauss_pts = boost::make_shared<MatrixDouble>();
214  auto dot2_u_at_gauss_pts = boost::make_shared<VectorDouble>();
215  pipeline.push_back(new OpCalculateScalarFieldGradient<SPACE_DIM>(
216  "U", grad_u_at_gauss_pts));
217  pipeline.push_back(
218  new OpCalculateScalarFieldValuesDotDot("U", dot2_u_at_gauss_pts));
219  pipeline.push_back(new OpDomainGradTimesVec(
220  "U", grad_u_at_gauss_pts,
221  [](double, double, double) -> double { return 1; }));
222  pipeline.push_back(new OpDomainTimesScalarField(
223  "U", dot2_u_at_gauss_pts,
224  [](const double, const double, const double) { return wave_speed2; }));
225  pipeline.push_back(new OpUnSetBc("U"));
226  };
227 
228  auto add_boundary_base_ops = [&](auto &pipeline) {};
229 
230  auto add_lhs_base_ops = [&](auto &pipeline) {
231  pipeline.push_back(new OpSetBc("U", false, boundaryMarker));
232  pipeline.push_back(new OpBoundaryMass(
233  "U", "U", [](const double, const double, const double) { return 1; }));
234  pipeline.push_back(new OpUnSetBc("U"));
235  };
236 
237  auto add_rhs_base_ops = [&](auto &pipeline) {
238  pipeline.push_back(new OpSetBc("U", false, boundaryMarker));
239  auto u_at_gauss_pts = boost::make_shared<VectorDouble>();
240  auto boundary_function = [&](const double x, const double y,
241  const double z) {
242  auto pipeline_mng = mField.getInterface<PipelineManager>();
243  auto &fe_rhs = pipeline_mng->getBoundaryRhsFE();
244  const double t = fe_rhs->ts_t;
245  if ((t <= 0.5) && (x < 0.) && (y > -1. / 3) && (y < 1. / 3))
246  return sin(4 * M_PI * t);
247  else
248  return 0.;
249  };
250  pipeline.push_back(new OpCalculateScalarFieldValues("U", u_at_gauss_pts));
251  pipeline.push_back(new OpBoundaryTimeScalarField(
252  "U", u_at_gauss_pts,
253  [](const double, const double, const double) { return 1; }));
254  pipeline.push_back(new OpBoundarySource("U", boundary_function));
255  pipeline.push_back(new OpUnSetBc("U"));
256  };
257 
258  auto pipeline_mng = mField.getInterface<PipelineManager>();
259  add_domain_base_ops(pipeline_mng->getOpDomainLhsPipeline());
260  add_domain_base_ops(pipeline_mng->getOpDomainRhsPipeline());
261  add_domain_lhs_ops(pipeline_mng->getOpDomainLhsPipeline());
262  add_domain_rhs_ops(pipeline_mng->getOpDomainRhsPipeline());
263 
264  add_boundary_base_ops(pipeline_mng->getOpBoundaryLhsPipeline());
265  add_boundary_base_ops(pipeline_mng->getOpBoundaryRhsPipeline());
266  add_lhs_base_ops(pipeline_mng->getOpBoundaryLhsPipeline());
267  add_rhs_base_ops(pipeline_mng->getOpBoundaryRhsPipeline());
268 
270 }
271 
274 
275  auto *simple = mField.getInterface<Simple>();
276  auto *pipeline_mng = mField.getInterface<PipelineManager>();
277 
278  auto create_post_process_element = [&]() {
279  auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
280 
281  auto det_ptr = boost::make_shared<VectorDouble>();
282  auto jac_ptr = boost::make_shared<MatrixDouble>();
283  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
284  post_proc_fe->getOpPtrVector().push_back(new OpCalculateHOJac<2>(jac_ptr));
285  post_proc_fe->getOpPtrVector().push_back(
286  new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
287  post_proc_fe->getOpPtrVector().push_back(
288  new OpSetHOInvJacToScalarBases<2>(H1, inv_jac_ptr));
289 
290  auto u_ptr = boost::make_shared<VectorDouble>();
291  post_proc_fe->getOpPtrVector().push_back(
292  new OpCalculateScalarFieldValues("U", u_ptr));
293 
295 
296  post_proc_fe->getOpPtrVector().push_back(
297 
298  new OpPPMap(post_proc_fe->getPostProcMesh(),
299  post_proc_fe->getMapGaussPts(),
300 
301  {{"U", u_ptr}},
302 
303  {}, {}, {}
304 
305  )
306 
307  );
308 
309  return post_proc_fe;
310  };
311 
312  auto set_time_monitor = [&](auto dm, auto solver) {
314  boost::shared_ptr<Monitor> monitor_ptr(
315  new Monitor(dm, create_post_process_element()));
316  boost::shared_ptr<ForcesAndSourcesCore> null;
317  CHKERR DMMoFEMTSSetMonitor(dm, solver, simple->getDomainFEName(),
318  monitor_ptr, null, null);
320  };
321 
322  auto set_fieldsplit_preconditioner = [&](auto solver) {
324  SNES snes;
325  CHKERR TSGetSNES(solver, &snes);
326  KSP ksp;
327  CHKERR SNESGetKSP(snes, &ksp);
328  PC pc;
329  CHKERR KSPGetPC(ksp, &pc);
330  PetscBool is_pcfs = PETSC_FALSE;
331  PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
332 
333  if (is_pcfs == PETSC_TRUE) {
334  auto bc_mng = mField.getInterface<BcManager>();
335  auto name_prb = simple->getProblemName();
336  auto is_all_bc = bc_mng->getBlockIS(name_prb, "ESSENTIAL", "U", 0, 0);
337  int is_all_bc_size;
338  CHKERR ISGetSize(is_all_bc, &is_all_bc_size);
339  MOFEM_LOG("EXAMPLE", Sev::inform)
340  << "Field split block size " << is_all_bc_size;
341  CHKERR PCFieldSplitSetIS(pc, PETSC_NULL,
342  is_all_bc); // boundary block
343  }
345  };
346 
347  auto dm = simple->getDM();
348  auto D = createDMVector(dm);
349  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
350 
351  auto ts = pipeline_mng->createTSIM2();
352 
353  CHKERR TSSetSolution(ts, D);
354  auto DD = vectorDuplicate(D);
355  CHKERR TS2SetSolution(ts, D, DD);
356  CHKERR set_time_monitor(dm, ts);
357  CHKERR TSSetFromOptions(ts);
358  CHKERR set_fieldsplit_preconditioner(ts);
359  CHKERR TSSetUp(ts);
360  CHKERR TSSolve(ts, NULL);
361 
362  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
363  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
364  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
365 
367 }
368 
371 
372  // Processes to set output results are integrated in solveSystem()
373 
375 }
376 
379 
380  CHKERR readMesh();
388 
390 }
391 
392 int main(int argc, char *argv[]) {
393 
394  // Initialisation of MoFEM/PETSc and MOAB data structures
395  const char param_file[] = "param_file.petsc";
396  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
397 
398  // Add logging channel for example
399  auto core_log = logging::core::get();
400  core_log->add_sink(
401  LogManager::createSink(LogManager::getStrmWorld(), "EXAMPLE"));
402  LogManager::setLog("EXAMPLE");
403  MOFEM_LOG_TAG("EXAMPLE", "example")
404 
405  // Error handling
406  try {
407  // Register MoFEM discrete manager in PETSc
408  DMType dm_name = "DMMOFEM";
409  CHKERR DMRegister_MoFEM(dm_name);
410 
411  // Create MOAB instance
412  moab::Core mb_instance; // mesh database
413  moab::Interface &moab = mb_instance; // mesh database interface
414 
415  // Create MoFEM instance
416  MoFEM::Core core(moab); // finite element database
417  MoFEM::Interface &m_field = core; // finite element interface
418 
419  // Run the main analysis
420  WaveEquation heat_problem(m_field);
421  CHKERR heat_problem.runProgram();
422  }
423  CATCH_ERRORS;
424 
425  // Finish work: cleaning memory, getting statistics, etc.
427 
428  return 0;
429 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
WaveEquation::initialCondition
MoFEMErrorCode initialCondition()
Definition: wave_equation.cpp:158
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
Monitor::Monitor
Monitor(SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEle > post_proc)
Definition: wave_equation.cpp:65
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, SPACE_DIM > OpDomainGradTimesVec
Definition: initial_diffusion.cpp:34
Monitor::preProcess
MoFEMErrorCode preProcess()
function is run at the beginning of loop
Definition: wave_equation.cpp:68
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
WaveEquation::setupProblem
MoFEMErrorCode setupProblem()
Definition: wave_equation.cpp:126
OpDomainTimesScalarField
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
Definition: wave_equation.cpp:44
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, SPACE_DIM > OpDomainGradTimesVec
Definition: wave_equation.cpp:46
OpBoundaryTimeScalarField
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpBoundaryTimeScalarField
Definition: photon_diffusion.cpp:47
MoFEM::OpSetHOInvJacToScalarBases< 2 >
Definition: HODataOperators.hpp:78
OpBoundarySource
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpBoundarySource
Definition: wave_equation.cpp:53
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
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:523
MoFEM::PipelineManager::getBoundaryRhsFE
boost::shared_ptr< FEMethod > & getBoundaryRhsFE()
Definition: PipelineManager.hpp:413
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
WaveEquation::mField
MoFEM::Interface & mField
Definition: wave_equation.cpp:107
Monitor::postProcess
MoFEMErrorCode postProcess()
function is run at the end of loop
Definition: wave_equation.cpp:73
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
OpBoundaryMass
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpBoundaryMass
Definition: wave_equation.cpp:49
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Definition: child_and_parent.cpp:53
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpDomainMass
Definition: wave_equation.cpp:40
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::OpCalculateHOJac< 2 >
Definition: HODataOperators.hpp:273
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
OpBoundaryMass
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpBoundaryMass
[Only used with Hencky/nonlinear material]
Definition: seepage.cpp:71
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
WaveEquation::boundaryCondition
MoFEMErrorCode boundaryCondition()
Definition: wave_equation.cpp:163
WaveEquation::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
Definition: wave_equation.cpp:142
MoFEM::PipelineManager::EdgeEle
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
Definition: PipelineManager.hpp:36
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: wave_equation.cpp:42
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
WaveEquation::WaveEquation
WaveEquation(MoFEM::Interface &m_field)
Definition: wave_equation.cpp:113
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
OpDomainTimesScalarField
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
Definition: initial_diffusion.cpp:32
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
WaveEquation::runProgram
MoFEMErrorCode runProgram()
Definition: wave_equation.cpp:377
BoundaryEleOp
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
double
MoFEM::PipelineManager::FaceEle
MoFEM::FaceElementForcesAndSourcesCore FaceEle
Definition: PipelineManager.hpp:35
OpBoundarySource
FormsIntegrators< EdgeEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpBoundarySource
Definition: helmholtz.cpp:31
wave_speed2
constexpr double wave_speed2
Definition: wave_equation.cpp:55
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
WaveEquation::solveSystem
MoFEMErrorCode solveSystem()
Definition: wave_equation.cpp:272
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
WaveEquation::readMesh
MoFEMErrorCode readMesh()
Definition: wave_equation.cpp:115
WaveEquation::outputResults
MoFEMErrorCode outputResults()
Definition: wave_equation.cpp:369
MoFEM::OpSetHOWeightsOnFace
Modify integration weights on face to take in account higher-order geometry.
Definition: HODataOperators.hpp:122
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
BiLinearForm
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
MoFEM::OpCalculateScalarFieldValuesDotDot
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_TT > OpCalculateScalarFieldValuesDotDot
Definition: UserDataOperators.hpp:275
Monitor::operator()
MoFEMErrorCode operator()()
function is run for every finite element
Definition: wave_equation.cpp:69
WaveEquation
Definition: wave_equation.cpp:88
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
MoFEM::OpUnSetBc
Definition: FormsIntegrators.hpp:49
ElementsAndOps
Definition: child_and_parent.cpp:18
main
int main(int argc, char *argv[])
Definition: wave_equation.cpp:392
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
DomainEleOp
help
static char help[]
Definition: wave_equation.cpp:24
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: helmholtz.cpp:25
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Monitor
[Push operators to pipeline]
Definition: dynamic_first_order_con_law.cpp:783
approx_order
int approx_order
Definition: test_broken_space.cpp:50
MoFEM::DMMoFEMTSSetMonitor
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition: DMMoFEM.cpp:1056
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::OpSetBc
Set indices on entities on finite element.
Definition: FormsIntegrators.hpp:38
WaveEquation::boundaryMarker
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
Definition: wave_equation.cpp:110
MoFEM::OpInvertMatrix
Definition: UserDataOperators.hpp:3249
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::SmartPetscObj< DM >
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:586
OpBoundaryTimeScalarField
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpBoundaryTimeScalarField
Definition: wave_equation.cpp:51
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
SPACE_DIM
constexpr int SPACE_DIM
[Define dimension]
Definition: wave_equation.cpp:29
EntData
EntitiesFieldData::EntData EntData
[Define dimension]
Definition: wave_equation.cpp:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
WaveEquation::assembleSystem
MoFEMErrorCode assembleSystem()
Definition: wave_equation.cpp:185
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698