v0.14.0
Loading...
Searching...
No Matches
hertz_surface.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

static char help [] = "...\n\n"
 
static int debug = 1
 

Function Documentation

◆ main()

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

Definition at line 32 of file hertz_surface.cpp.

32 {
33
34 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
35
36 try {
37
38 moab::Core mb_instance;
39 moab::Interface &moab = mb_instance;
40 int rank;
41 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
42
43 // Read parameters from line command
44 PetscBool flg_mesh_file = PETSC_TRUE;
45 PetscBool flg_out_file = PETSC_TRUE;
46 char mesh_file_name[255];
47 char out_file_name[255] = "out.h5m";
48
49 double radius = 1.0;
50 double height = 1.0;
51 int dim = 2;
52
53 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "ADD_PRISMS_LAYER", "none");
54
55 CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
56 mesh_file_name, 255, &flg_mesh_file);
57 CHKERR PetscOptionsString("-my_out_file", "out file name", "", "out.h5m",
58 out_file_name, 255, &flg_out_file);
59
60 CHKERR PetscOptionsInt("-my_dim", "dimension (2 or 3)", "", dim, &dim,
61 PETSC_NULL);
62
63 CHKERR PetscOptionsReal("-my_radius", "roughness wavelength", "", radius,
64 &radius, PETSC_NULL);
65 CHKERR PetscOptionsReal("-my_height", "vertical dimension of the mesh", "",
66 height, &height, PETSC_NULL);
67
68 int ierr = PetscOptionsEnd();
70
71 if (flg_mesh_file != PETSC_TRUE) {
72 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
73 }
74
75 // Read mesh to MOAB
76 const char *option;
77 option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
78 CHKERR moab.load_file(mesh_file_name, 0, option);
79 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
80 if (pcomm == NULL)
81 pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
82
83 // Create MoFEM (Joseph) databas
84 MoFEM::Core core(moab);
85 MoFEM::Interface &m_field = core;
86
87 Range tets, tris, nodes, all_nodes;
89 if (bit->getName().compare(0, 11, "MAT_ELASTIC") == 0) {
90 const int id = bit->getMeshsetId();
91 if (id == 1) {
92 CHKERR m_field.get_moab().get_entities_by_dimension(bit->getMeshset(),
93 3, tets, true);
94 CHKERR m_field.get_moab().get_entities_by_dimension(bit->getMeshset(),
95 2, tris, true);
96 }
97 }
98 }
99
100 CHKERR m_field.get_moab().get_connectivity(tets, nodes);
101 all_nodes.merge(nodes);
102 nodes.clear();
103 CHKERR m_field.get_moab().get_connectivity(tris, nodes);
104 all_nodes.merge(nodes);
105
106 double coords[3];
107 for (Range::iterator nit = all_nodes.begin(); nit != all_nodes.end();
108 nit++) {
109 CHKERR moab.get_coords(&*nit, 1, coords);
110 double x = coords[0];
111 double y = coords[1];
112 double z = coords[2];
113 double r = sqrt(x * x + y * y);
114 double coef = (height + z) / height;
115 switch (dim) {
116 case 2:
117 coords[2] -= coef * 0.5 * x * x / radius;
118 break;
119 case 3:
120 coords[2] -= coef * 0.5 * r * r / radius;
121 break;
122 default:
123 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
124 "Wrong dimension = %d", dim);
125 }
126
127 CHKERR moab.set_coords(&*nit, 1, coords);
128 }
129
130 EntityHandle rootset = moab.get_root_set();
131 CHKERR moab.write_file(out_file_name, "MOAB", "", &rootset, 1);
132 MOFEM_LOG_C("WORLD", Sev::inform, "Write file %s\n",
133 string(out_file_name).c_str());
134 }
136
138 return 0;
139}
#define MOFEM_LOG_C(channel, severity, format,...)
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ BLOCKSET
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define CHKERR
Inline error check.
#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.
auto bit
set bit
static char help[]
char out_file_name[255]
char mesh_file_name[255]
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
int r
Definition sdf.py:8
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.

Variable Documentation

◆ debug

int debug = 1
static

Definition at line 30 of file hertz_surface.cpp.

◆ help

char help[] = "...\n\n"
static

Definition at line 29 of file hertz_surface.cpp.