|
| v0.14.0
|
Go to the source code of this file.
|
int | main (int argc, char *argv[]) |
|
◆ main()
int main |
( |
int |
argc, |
|
|
char * |
argv[] |
|
) |
| |
- Examples
- elasticity_mixed_formulation.cpp.
Definition at line 22 of file elasticity_mixed_formulation.cpp.
24 const string default_options =
"-ksp_type gmres \n"
26 "-pc_factor_mat_solver_type mumps \n"
27 "-mat_mumps_icntl_20 0 \n"
30 string param_file =
"param_file.petsc";
31 if (!
static_cast<bool>(ifstream(param_file))) {
32 std::ofstream file(param_file.c_str(), std::ios::ate);
34 file << default_options;
51 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab,
MYPCOMM_INDEX);
53 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD,
false);
55 pcomm =
new ParallelComm(&moab, moab_comm_wrap->get_comm());
62 PetscBool is_partitioned = PETSC_FALSE;
63 PetscBool flg_test = PETSC_FALSE;
65 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD,
"",
"Mix elastic problem",
68 CHKERR PetscOptionsString(
"-my_file",
"mesh file name",
"",
"mesh.h5m",
71 CHKERR PetscOptionsInt(
"-my_order_p",
"approximation order_p",
"", order_p,
72 &order_p, PETSC_NULL);
73 CHKERR PetscOptionsInt(
"-my_order_u",
"approximation order_u",
"", order_u,
74 &order_u, PETSC_NULL);
76 CHKERR PetscOptionsBool(
"-is_partitioned",
"is_partitioned?",
"",
77 is_partitioned, &is_partitioned, PETSC_NULL);
80 CHKERR PetscOptionsBool(
"-test",
"if true is ctest",
"", flg_test,
81 &flg_test, PETSC_NULL);
82 ierr = PetscOptionsEnd();
84 if (flg_file != PETSC_TRUE) {
85 SETERRQ(PETSC_COMM_SELF, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
89 if (is_partitioned == PETSC_TRUE) {
93 option =
"PARALLEL=READ_PART;"
94 "PARALLEL_RESOLVE_SHARED_ENTS;"
95 "PARTITION=PARALLEL_PARTITION;";
149 "MESH_NODE_POSITIONS");
169 "MESH_NODE_POSITIONS");
180 DMType dm_name =
"DM_ELASTIC_MIX";
183 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
184 CHKERR DMSetType(dm, dm_name);
189 CHKERR DMSetFromOptions(dm);
200 DataAtIntegrationPts commonData(m_field);
201 CHKERR commonData.getParameters();
203 boost::shared_ptr<FEMethod> nullFE;
204 boost::shared_ptr<VolumeElementForcesAndSourcesCore> feLhs(
206 boost::shared_ptr<VolumeElementForcesAndSourcesCore> feRhs(
212 for (
auto &sit : commonData.setOfBlocksData) {
214 feLhs->getOpPtrVector().push_back(
216 feLhs->getOpPtrVector().push_back(
218 feLhs->getOpPtrVector().push_back(
221 auto u_ptr = boost::make_shared<MatrixDouble>();
223 post_proc.getOpPtrVector().push_back(
225 post_proc.getOpPtrVector().push_back(
227 post_proc.getOpPtrVector().push_back(
229 commonData.gradDispPtr));
233 post_proc.getOpPtrVector().push_back(
237 post_proc.getPostProcMesh(), post_proc.getMapGaussPts(),
238 {{
"P", commonData.pPtr}},
247 post_proc.getPostProcMesh(), post_proc.getMapGaussPts(), commonData,
257 CHKERR VecDuplicate(
d, &F_ext);
259 CHKERR VecGhostUpdateBegin(
d, INSERT_VALUES, SCATTER_FORWARD);
260 CHKERR VecGhostUpdateEnd(
d, INSERT_VALUES, SCATTER_FORWARD);
262 CHKERR MatZeroEntries(Aij);
267 feLhs->ksp_f = F_ext;
269 boost::shared_ptr<DirichletDisplacementBc> dirichlet_bc_ptr(
271 dirichlet_bc_ptr->snes_ctx = FEMethod::CTX_SNESNONE;
272 dirichlet_bc_ptr->ts_ctx = FEMethod::CTX_TSNONE;
277 boost::ptr_map<std::string, NeumannForcesSurface> neumann_forces;
282 boost::ptr_map<std::string, NeumannForcesSurface>::iterator mit =
283 neumann_forces.begin();
284 for (; mit != neumann_forces.end(); mit++) {
286 &mit->second->getLoopFe());
290 boost::ptr_map<std::string, NodalForce> nodal_forces;
294 boost::ptr_map<std::string, NodalForce>::iterator fit =
295 nodal_forces.begin();
296 for (; fit != nodal_forces.end(); fit++) {
298 &fit->second->getLoopFe());
302 boost::ptr_map<std::string, EdgeForce> edge_forces;
306 boost::ptr_map<std::string, EdgeForce>::iterator fit =
308 for (; fit != edge_forces.end(); fit++) {
310 &fit->second->getLoopFe());
316 CHKERR VecAssemblyBegin(F_ext);
317 CHKERR VecAssemblyEnd(F_ext);
323 CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
324 CHKERR KSPSetFromOptions(solver);
325 CHKERR KSPSetOperators(solver, Aij, Aij);
327 CHKERR KSPSolve(solver, F_ext,
d);
335 PetscPrintf(PETSC_COMM_WORLD,
"Output file: %s\n",
"out.h5m");
336 CHKERR post_proc.writeFile(
"out.h5m");
340 CHKERR VecDestroy(&F_ext);
◆ 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.
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
Get values at integration pts for tensor filed rank 1, i.e. vector field.
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.
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Set Dirichlet boundary conditions on displacements.
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.
DeprecatedCoreInterface Interface
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
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
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
#define CHKERR
Inline error check.
MoFEMErrorCode printForceSet() const
print meshsets with force boundary conditions data structure
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
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.
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Get value at integration points for scalar field.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
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.
MoFEMErrorCode printMaterialsSet() const
print meshsets with material data structure set on it
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define CATCH_ERRORS
Catch errors.
MoFEMErrorCode printDisplacementSet() const
print meshsets with displacement boundary conditions data structure
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
const FTensor::Tensor2< T, Dim, Dim > Vec
Interface for managing meshsets containing materials and boundary conditions.
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
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 DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
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.
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Post post-proc data at points from hash maps.