|
| v0.14.0
|
Go to the source code of this file.
|
int | main (int argc, char *argv[]) |
|
|
static char | help [] = "\n" |
|
◆ main()
int main |
( |
int |
argc, |
|
|
char * |
argv[] |
|
) |
| |
Definition at line 20 of file testing_jacobian_of_navier_stokes_element.cpp.
27 PetscBool test_jacobian = PETSC_FALSE;
37 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab,
MYPCOMM_INDEX);
39 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD,
false);
42 new ParallelComm(&moab, moab_comm_wrap->get_comm());
49 PetscBool is_partitioned = PETSC_FALSE;
50 PetscBool flg_test = PETSC_FALSE;
52 PetscBool only_stokes = PETSC_FALSE;
53 PetscBool discont_pressure = PETSC_FALSE;
55 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD,
"",
"TEST_NAVIER_STOKES problem",
58 CHKERR PetscOptionsString(
"-my_file",
"mesh file name",
"",
"mesh.h5m",
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);
66 CHKERR PetscOptionsBool(
"-is_partitioned",
"is_partitioned?",
"",
67 is_partitioned, &is_partitioned, PETSC_NULL);
69 CHKERR PetscOptionsBool(
"-only_stokes",
"only stokes",
"", only_stokes,
70 &only_stokes, PETSC_NULL);
72 CHKERR PetscOptionsBool(
"-discont_pressure",
"discontinuous pressure",
"",
73 discont_pressure, &discont_pressure, PETSC_NULL);
77 ierr = PetscOptionsEnd();
79 if (flg_file != PETSC_TRUE) {
80 SETERRQ(PETSC_COMM_SELF, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
84 if (is_partitioned == PETSC_TRUE) {
88 option =
"PARALLEL=READ_PART;"
89 "PARALLEL_RESOLVE_SHARED_ENTS;"
90 "PARTITION=PARALLEL_PARTITION;";
117 if (discont_pressure) {
140 auto setting_second_order_geometry = [&m_field]() {
145 moab::Interface::UNION);
152 CHKERR setting_second_order_geometry();
157 "MESH_NODE_POSITIONS");
161 PetscRandomCreate(PETSC_COMM_WORLD, &rctx);
163 auto set_velocity = [&](
VectorAdaptor &&field_data,
double *x,
double *y,
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;
177 auto set_pressure = [&](
VectorAdaptor &&field_data,
double *x,
double *y,
182 PetscRandomGetValueReal(rctx, &value);
183 field_data[0] = value *
scale;
190 PetscRandomDestroy(&rctx);
212 "MESH_NODE_POSITIONS");
216 "TEST_NAVIER_STOKES");
224 DMType dm_name =
"DM_TEST_NAVIER_STOKES";
227 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
228 CHKERR DMSetType(dm, dm_name);
233 CHKERR DMSetFromOptions(dm);
241 boost::shared_ptr<FEMethod> nullFE;
243 boost::shared_ptr<VolumeElementForcesAndSourcesCore> feLhs(
245 boost::shared_ptr<VolumeElementForcesAndSourcesCore> feRhs(
251 boost::shared_ptr<NavierStokesElement::CommonData> commonData =
252 boost::make_shared<NavierStokesElement::CommonData>(m_field);
255 if (
bit->getName().compare(0, 9,
"MAT_FLUID") == 0) {
256 const int id =
bit->getMeshsetId();
258 bit->getMeshset(), MBTET, commonData->setOfBlocksData[
id].eNts,
260 std::vector<double> attributes;
261 bit->getAttributes(attributes);
262 if (attributes.size() != 2) {
264 "should be 2 attributes but is %d", attributes.size());
266 commonData->setOfBlocksData[id].iD = id;
267 commonData->setOfBlocksData[id].fluidViscosity = attributes[0];
268 commonData->setOfBlocksData[id].fluidDensity = attributes[1];
270 double reNumber = attributes[1] / attributes[0];
271 commonData->setOfBlocksData[id].inertiaCoef = reNumber;
272 commonData->setOfBlocksData[id].viscousCoef = 1.0;
302 CHKERR VecZeroEntries(F2);
307 CHKERR MatSetOption(
A, MAT_SPD, PETSC_TRUE);
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 "
317 CHKERR PetscOptionsInsertString(NULL, testing_options);
319 char testing_options[] =
"-snes_no_convergence_test -snes_atol 0 "
322 CHKERR PetscOptionsInsertString(NULL, testing_options);
326 SNESConvergedReason snes_reason;
334 CHKERR SNESSetFromOptions(snes);
337 CHKERR SNESSolve(snes, PETSC_NULL,
D);
339 if (test_jacobian == PETSC_FALSE) {
341 CHKERR MatNorm(
A, NORM_INFINITY, &nrm_A0);
343 char testing_options_fd[] =
"-snes_fd";
344 CHKERR PetscOptionsInsertString(NULL, testing_options_fd);
348 CHKERR SNESSetFromOptions(snes);
350 CHKERR SNESSolve(snes, NULL, D2);
352 CHKERR MatAXPY(
A, -1, fdA, SUBSET_NONZERO_PATTERN);
355 CHKERR MatNorm(
A, NORM_INFINITY, &nrm_A);
356 PetscPrintf(PETSC_COMM_WORLD,
"Matrix norms %3.4e %3.4e\n", nrm_A,
360 constexpr
double tol = 1e-6;
363 "Difference between hand-calculated tangent matrix and finite "
364 "difference matrix is too big");
◆ help
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
#define MYPCOMM_INDEX
default communicator number PCOMM
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.
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
virtual MPI_Comm & get_comm() const =0
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
auto createSNES(MPI_Comm comm)
@ L2
field with C-1 continuity
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Projection of edge entities with one mid-node on hierarchical basis.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
auto createDMMatrix(DM dm)
Get smart matrix from DM.
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.
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Deprecated interface functions.
Set integration rule to volume elements.
DeprecatedCoreInterface Interface
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
#define CHKERR
Inline error check.
auto createDMVector(DM dm)
Get smart vector from DM.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual moab::Interface & get_moab()=0
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
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
SmartPetscObj< Mat > matDuplicate(Mat mat, MatDuplicateOption op)
VectorShallowArrayAdaptor< double > VectorAdaptor
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.
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
Volume finite element base.
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.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
#define CATCH_ERRORS
Catch errors.
#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.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Interface for managing meshsets containing materials and boundary conditions.
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
const double D
diffusivity
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
@ MOFEM_ATOM_TEST_INVALID
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
Interface for nonlinear (SNES) solver.
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.
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.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
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.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
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.
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)