v0.14.0
Functions | Variables
testing_jacobian_of_navier_stokes_element.cpp File Reference
#include <BasicFiniteElements.hpp>
#include <NavierStokesElement.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[] 
)

Definition at line 20 of file testing_jacobian_of_navier_stokes_element.cpp.

20  {
21 
22  // Initialise MoFEM
23  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
24 
25  try {
26 
27  PetscBool test_jacobian = PETSC_FALSE;
28 
29  // Create mesh database
30  moab::Core mb_instance; // create database
31  moab::Interface &moab = mb_instance; // create interface to database
32 
33  // Create moab communicator
34  // Create separate MOAB communicator, it is duplicate of PETSc communicator.
35  // NOTE That this should eliminate potential communication problems between
36  // MOAB and PETSC functions.
37  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
38  auto moab_comm_wrap =
39  boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
40  if (pcomm == NULL)
41  pcomm =
42  new ParallelComm(&moab, moab_comm_wrap->get_comm());
43 
44  // Get command line options
45  char mesh_file_name[255];
46  PetscBool flg_file;
47  int order_p = 1; // default approximation order_p
48  int order_u = 2; // default approximation order_u
49  PetscBool is_partitioned = PETSC_FALSE;
50  PetscBool flg_test = PETSC_FALSE; // true check if error is numerical error
51 
52  PetscBool only_stokes = PETSC_FALSE;
53  PetscBool discont_pressure = PETSC_FALSE;
54 
55  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "TEST_NAVIER_STOKES problem",
56  "none");
57 
58  CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
59  mesh_file_name, 255, &flg_file);
60  // Set approximation order
61  CHKERR PetscOptionsInt("-my_order_p", "approximation order_p", "", order_p,
62  &order_p, PETSC_NULL);
63  CHKERR PetscOptionsInt("-my_order_u", "approximation order_u", "", order_u,
64  &order_u, PETSC_NULL);
65 
66  CHKERR PetscOptionsBool("-is_partitioned", "is_partitioned?", "",
67  is_partitioned, &is_partitioned, PETSC_NULL);
68 
69  CHKERR PetscOptionsBool("-only_stokes", "only stokes", "", only_stokes,
70  &only_stokes, PETSC_NULL);
71 
72  CHKERR PetscOptionsBool("-discont_pressure", "discontinuous pressure", "",
73  discont_pressure, &discont_pressure, PETSC_NULL);
74 
75  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-test_jacobian", &test_jacobian,
76  PETSC_NULL);
77  ierr = PetscOptionsEnd();
78 
79  if (flg_file != PETSC_TRUE) {
80  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
81  }
82 
83  // Read whole mesh or part of it if partitioned
84  if (is_partitioned == PETSC_TRUE) {
85  // This is a case of distributed mesh and algebra. In that case each
86  // processor keeps only part of the problem.
87  const char *option;
88  option = "PARALLEL=READ_PART;"
89  "PARALLEL_RESOLVE_SHARED_ENTS;"
90  "PARTITION=PARALLEL_PARTITION;";
91  CHKERR moab.load_file(mesh_file_name, 0, option);
92  } else {
93  // In this case, we have distributed algebra, i.e. assembly of vectors and
94  // matrices is in parallel, but whole mesh is stored on all processors.
95  // snes and matrix scales well, however problem set-up of problem is
96  // not fully parallel.
97  const char *option;
98  option = "";
99  CHKERR moab.load_file(mesh_file_name, 0, option);
100  }
101 
102  // Create MoFEM database and link it to MoAB
103  MoFEM::Core core(moab);
104  MoFEM::Interface &m_field = core;
105 
106  // Print boundary conditions and material parameters
107  MeshsetsManager *meshsets_mng_ptr;
108 
109  BitRefLevel bit_level0;
110  bit_level0.set(0);
111  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
112  0, 3, bit_level0);
113 
114  // **** ADD FIELDS **** //
115 
116  CHKERR m_field.add_field("U", H1, AINSWORTH_LEGENDRE_BASE, 3);
117  if (discont_pressure) {
118  CHKERR m_field.add_field("P", L2, AINSWORTH_LEGENDRE_BASE, 1);
119  } else {
120  CHKERR m_field.add_field("P", H1, AINSWORTH_LEGENDRE_BASE, 1);
121  }
122  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
123  3);
124 
125  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "U");
126  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "P");
127  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "MESH_NODE_POSITIONS");
128 
129  CHKERR m_field.set_field_order(0, MBVERTEX, "U", 1);
130  CHKERR m_field.set_field_order(0, MBEDGE, "U", order_u);
131  CHKERR m_field.set_field_order(0, MBTRI, "U", order_u);
132  CHKERR m_field.set_field_order(0, MBTET, "U", order_u);
133 
134  CHKERR m_field.set_field_order(0, MBVERTEX, "P", 1);
135  CHKERR m_field.set_field_order(0, MBEDGE, "P", order_p);
136  CHKERR m_field.set_field_order(0, MBTRI, "P", order_p);
137  CHKERR m_field.set_field_order(0, MBTET, "P", order_p);
138 
139  // Set 2nd order of approximation of geometry
140  auto setting_second_order_geometry = [&m_field]() {
142  Range tets, edges;
143  CHKERR m_field.get_moab().get_entities_by_type(0, MBTET, tets);
144  CHKERR m_field.get_moab().get_adjacencies(tets, 1, false, edges,
145  moab::Interface::UNION);
146 
147  CHKERR m_field.set_field_order(edges, "MESH_NODE_POSITIONS", 2);
149  };
150 
151  CHKERR m_field.set_field_order(0, MBVERTEX, "MESH_NODE_POSITIONS", 1);
152  CHKERR setting_second_order_geometry();
153 
154  CHKERR m_field.build_fields();
155 
156  Projection10NodeCoordsOnField ent_method_material(m_field,
157  "MESH_NODE_POSITIONS");
158  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
159 
160  PetscRandom rctx;
161  PetscRandomCreate(PETSC_COMM_WORLD, &rctx);
162 
163  auto set_velocity = [&](VectorAdaptor &&field_data, double *x, double *y,
164  double *z) {
166  double value;
167  double scale = 2.0;
168  PetscRandomGetValueReal(rctx, &value);
169  field_data[0] = (value - 0.5) * scale;
170  PetscRandomGetValueReal(rctx, &value);
171  field_data[1] = (value - 0.5) * scale;
172  PetscRandomGetValueReal(rctx, &value);
173  field_data[2] = (value - 0.5) * scale;
175  };
176 
177  auto set_pressure = [&](VectorAdaptor &&field_data, double *x, double *y,
178  double *z) {
180  double value;
181  double scale = 1.0;
182  PetscRandomGetValueReal(rctx, &value);
183  field_data[0] = value * scale;
185  };
186 
187  CHKERR m_field.getInterface<FieldBlas>()->setVertexDofs(set_velocity, "U");
188  CHKERR m_field.getInterface<FieldBlas>()->setVertexDofs(set_pressure, "P");
189 
190  PetscRandomDestroy(&rctx);
191 
192  // **** ADD ELEMENTS **** //
193 
194  // Add finite element (this defines element, declaration comes later)
195  CHKERR m_field.add_finite_element("TEST_NAVIER_STOKES");
196 
197  CHKERR m_field.modify_finite_element_add_field_row("TEST_NAVIER_STOKES",
198  "U");
199  CHKERR m_field.modify_finite_element_add_field_col("TEST_NAVIER_STOKES",
200  "U");
201  CHKERR m_field.modify_finite_element_add_field_data("TEST_NAVIER_STOKES",
202  "U");
203 
204  CHKERR m_field.modify_finite_element_add_field_row("TEST_NAVIER_STOKES",
205  "P");
206  CHKERR m_field.modify_finite_element_add_field_col("TEST_NAVIER_STOKES",
207  "P");
208  CHKERR m_field.modify_finite_element_add_field_data("TEST_NAVIER_STOKES",
209  "P");
210 
211  CHKERR m_field.modify_finite_element_add_field_data("TEST_NAVIER_STOKES",
212  "MESH_NODE_POSITIONS");
213 
214  // Add entities to that element
216  "TEST_NAVIER_STOKES");
217  // build finite elements
218  CHKERR m_field.build_finite_elements();
219  // build adjacencies between elements and degrees of freedom
220  CHKERR m_field.build_adjacencies(bit_level0);
221 
222  // **** BUILD DM **** //
223  DM dm;
224  DMType dm_name = "DM_TEST_NAVIER_STOKES";
225  // Register DM problem
226  CHKERR DMRegister_MoFEM(dm_name);
227  CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
228  CHKERR DMSetType(dm, dm_name);
229  // Create DM instance
230  CHKERR DMMoFEMCreateMoFEM(dm, &m_field, dm_name, bit_level0);
231  // Configure DM form line command options (DM itself, sness,
232  // pre-conditioners, ... )
233  CHKERR DMSetFromOptions(dm);
234  // Add elements to dm (only one here)
235  CHKERR DMMoFEMAddElement(dm, "TEST_NAVIER_STOKES");
236 
237  CHKERR DMMoFEMSetIsPartitioned(dm, is_partitioned);
238  // setup the DM
239  CHKERR DMSetUp(dm);
240 
241  boost::shared_ptr<FEMethod> nullFE;
242 
243  boost::shared_ptr<VolumeElementForcesAndSourcesCore> feLhs(
244  new VolumeElementForcesAndSourcesCore(m_field));
245  boost::shared_ptr<VolumeElementForcesAndSourcesCore> feRhs(
246  new VolumeElementForcesAndSourcesCore(m_field));
247 
248  feLhs->getRuleHook = NavierStokesElement::VolRule();
249  feRhs->getRuleHook = NavierStokesElement::VolRule();
250 
251  boost::shared_ptr<NavierStokesElement::CommonData> commonData =
252  boost::make_shared<NavierStokesElement::CommonData>(m_field);
253 
255  if (bit->getName().compare(0, 9, "MAT_FLUID") == 0) {
256  const int id = bit->getMeshsetId();
257  CHKERR m_field.get_moab().get_entities_by_type(
258  bit->getMeshset(), MBTET, commonData->setOfBlocksData[id].eNts,
259  true);
260  std::vector<double> attributes;
261  bit->getAttributes(attributes);
262  if (attributes.size() != 2) {
263  SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
264  "should be 2 attributes but is %d", attributes.size());
265  }
266  commonData->setOfBlocksData[id].iD = id;
267  commonData->setOfBlocksData[id].fluidViscosity = attributes[0];
268  commonData->setOfBlocksData[id].fluidDensity = attributes[1];
269 
270  double reNumber = attributes[1] / attributes[0];
271  commonData->setOfBlocksData[id].inertiaCoef = reNumber;
272  commonData->setOfBlocksData[id].viscousCoef = 1.0;
273  }
274  }
275 
276  CHKERR addHOOpsVol("MESH_NODE_POSITIONS", *feLhs, true, false, false, true);
277  CHKERR addHOOpsVol("MESH_NODE_POSITIONS", *feRhs, true, false, false, true);
278  if (only_stokes) {
279  CHKERR NavierStokesElement::setStokesOperators(feRhs, feLhs, "U", "P",
280  commonData);
281  } else {
283  "P", commonData);
284  }
285 
286  CHKERR DMMoFEMSNESSetJacobian(dm, "TEST_NAVIER_STOKES", feLhs, nullFE,
287  nullFE);
288  CHKERR DMMoFEMSNESSetFunction(dm, "TEST_NAVIER_STOKES", feRhs, nullFE,
289  nullFE);
290 
291  // Vector of DOFs and the RHS
292  auto D = createDMVector(dm);
293  auto F = vectorDuplicate(D);
294 
295  auto D2 = vectorDuplicate(D);
296  auto F2 = vectorDuplicate(D);
297 
298  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
299  CHKERR VecZeroEntries(F);
300 
301  CHKERR DMoFEMMeshToLocalVector(dm, D2, INSERT_VALUES, SCATTER_FORWARD);
302  CHKERR VecZeroEntries(F2);
303 
304  // Stiffness matrix
305  auto A = createDMMatrix(dm);
306 
307  CHKERR MatSetOption(A, MAT_SPD, PETSC_TRUE);
308  CHKERR MatZeroEntries(A);
309 
310  auto fdA = matDuplicate(A, MAT_COPY_VALUES);
311 
312  if (test_jacobian == PETSC_TRUE) {
313  char testing_options[] =
314  "-snes_test_jacobian 1e-7 -snes_test_jacobian_view "
315  "-snes_no_convergence_test -snes_atol 0 -snes_rtol 0 "
316  "-snes_max_it 1 ";
317  CHKERR PetscOptionsInsertString(NULL, testing_options);
318  } else {
319  char testing_options[] = "-snes_no_convergence_test -snes_atol 0 "
320  "-snes_rtol 0 "
321  "-snes_max_it 1 ";
322  CHKERR PetscOptionsInsertString(NULL, testing_options);
323  }
324 
325  auto snes = MoFEM::createSNES(m_field.get_comm());
326  SNESConvergedReason snes_reason;
327  SnesCtx *snes_ctx;
328 
329  // create snes nonlinear solver
330  {
331  CHKERR DMMoFEMGetSnesCtx(dm, &snes_ctx);
332  CHKERR SNESSetFunction(snes, F, SnesRhs, snes_ctx);
333  CHKERR SNESSetJacobian(snes, A, A, SnesMat, snes_ctx);
334  CHKERR SNESSetFromOptions(snes);
335  }
336 
337  CHKERR SNESSolve(snes, PETSC_NULL, D);
338 
339  if (test_jacobian == PETSC_FALSE) {
340  double nrm_A0;
341  CHKERR MatNorm(A, NORM_INFINITY, &nrm_A0);
342 
343  char testing_options_fd[] = "-snes_fd";
344  CHKERR PetscOptionsInsertString(NULL, testing_options_fd);
345 
346  CHKERR SNESSetFunction(snes, F2, SnesRhs, snes_ctx);
347  CHKERR SNESSetJacobian(snes, fdA, fdA, SnesMat, snes_ctx);
348  CHKERR SNESSetFromOptions(snes);
349 
350  CHKERR SNESSolve(snes, NULL, D2);
351 
352  CHKERR MatAXPY(A, -1, fdA, SUBSET_NONZERO_PATTERN);
353 
354  double nrm_A;
355  CHKERR MatNorm(A, NORM_INFINITY, &nrm_A);
356  PetscPrintf(PETSC_COMM_WORLD, "Matrix norms %3.4e %3.4e\n", nrm_A,
357  nrm_A / nrm_A0);
358  nrm_A /= nrm_A0;
359 
360  constexpr double tol = 1e-6;
361  if (nrm_A > tol) {
362  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
363  "Difference between hand-calculated tangent matrix and finite "
364  "difference matrix is too big");
365  }
366  }
367  }
368  CATCH_ERRORS;
369 
370  // finish work cleaning memory, getting statistics, etc
372 
373  return 0;
374 }

Variable Documentation

◆ help

char help[] = "\n"
static
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
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::addHOOpsVol
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
Definition: HODataOperators.hpp:674
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
MoFEM::createSNES
auto createSNES(MPI_Comm comm)
Definition: PetscSmartObj.hpp:255
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
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::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
scale
double scale
Definition: plastic.cpp:119
MoFEM::createDMMatrix
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition: DMMoFEM.hpp:1056
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:497
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
NavierStokesElement::VolRule
Set integration rule to volume elements.
Definition: NavierStokesElement.hpp:263
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
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
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
help
static char help[]
Definition: testing_jacobian_of_navier_stokes_element.cpp:18
MoFEM::DMMoFEMGetSnesCtx
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition: DMMoFEM.cpp:1094
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
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:759
MoFEM::matDuplicate
SmartPetscObj< Mat > matDuplicate(Mat mat, MatDuplicateOption op)
Definition: PetscSmartObj.hpp:234
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
MoFEM::DMMoFEMCreateMoFEM
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMoFEM.cpp:114
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
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
Range
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::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::SnesCtx
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:13
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
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
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(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_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
NavierStokesElement::setNavierStokesOperators
static MoFEMErrorCode setNavierStokesOperators(boost::shared_ptr< VolumeElementForcesAndSourcesCore > feRhs, boost::shared_ptr< VolumeElementForcesAndSourcesCore > feLhs, const std::string velocity_field, const std::string pressure_field, boost::shared_ptr< CommonData > common_data, const EntityType type=MBTET)
Setting up operators for solving Navier-Stokes equations.
Definition: NavierStokesElement.cpp:22
NavierStokesElement::setStokesOperators
static MoFEMErrorCode setStokesOperators(boost::shared_ptr< VolumeElementForcesAndSourcesCore > feRhs, boost::shared_ptr< VolumeElementForcesAndSourcesCore > feLhs, const std::string velocity_field, const std::string pressure_field, boost::shared_ptr< CommonData > common_data, const EntityType type=MBTET)
Setting up operators for solving Stokes equations.
Definition: NavierStokesElement.cpp:79
tol
double tol
Definition: mesh_smoothing.cpp:26
F
@ F
Definition: free_surface.cpp:394
MoFEM::DMMoFEMSetIsPartitioned
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1123
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:718
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182