v0.14.0
Functions | Variables
gel_jacobian_test.cpp File Reference

Atom test testing calculation of element residual vectors and tangent matrices. More...

#include <MoFEM.hpp>
#include <PostProcOnRefMesh.hpp>
#include <Projection10NodeCoordsOnField.hpp>
#include <boost/numeric/ublas/vector_proxy.hpp>
#include <boost/numeric/ublas/matrix.hpp>
#include <boost/numeric/ublas/matrix_proxy.hpp>
#include <boost/numeric/ublas/vector.hpp>
#include <adolc/adolc.h>
#include <Gels.hpp>
#include <boost/iostreams/tee.hpp>
#include <boost/iostreams/stream.hpp>
#include <fstream>
#include <iostream>

Go to the source code of this file.

Functions

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

Variables

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

Detailed Description

Atom test testing calculation of element residual vectors and tangent matrices.

Definition in file gel_jacobian_test.cpp.

Function Documentation

◆ main()

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

Definition at line 43 of file gel_jacobian_test.cpp.

43  {
44 
45  MoFEM::Core::Initialize(&argc,&argv,(char *)0,help);
46 
47  try {
48 
49  moab::Core mb_instance;
50  moab::Interface& moab = mb_instance;
51 
52  {
53  PetscBool flg = PETSC_TRUE;
54  char mesh_file_name[255];
55  ierr = PetscOptionsGetString(PETSC_NULL,PETSC_NULL,"-my_file",mesh_file_name,255,&flg); CHKERRQ(ierr);
56  if(flg != PETSC_TRUE) {
57  SETERRQ(PETSC_COMM_SELF,1,"*** ERROR -my_file (MESH FILE NEEDED)");
58  }
59  const char *option;
60  option = "";//"PARALLEL=BCAST;";//;DEBUG_IO";
61  rval = moab.load_file(mesh_file_name, 0, option); CHKERRQ_MOAB(rval);
62  ParallelComm* pcomm = ParallelComm::get_pcomm(&moab,MYPCOMM_INDEX);
63  if(pcomm == NULL) pcomm = new ParallelComm(&moab,PETSC_COMM_WORLD);
64  }
65 
66  MoFEM::Core core(moab);
67  MoFEM::Interface& m_field = core;
68  BitRefLevel bit_level0;
69  bit_level0.set(0);
70  ierr = m_field.seed_ref_level_3D(0,bit_level0); CHKERRQ(ierr);
71 
72  // Define fields and finite elements
73  {
74 
75  // Set approximation fields
76  {
77  // Seed all mesh entities to MoFEM database, those entities can be potentially used as finite elements
78  // or as entities which carry some approximation field.
79 
80  bool check_if_spatial_field_exist = m_field.check_field("SPATIAL_POSITION");
81  ierr = m_field.add_field("SPATIAL_POSITION",H1,AINSWORTH_LEGENDRE_BASE,3,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
82  ierr = m_field.add_field("CHEMICAL_LOAD",H1,AINSWORTH_LEGENDRE_BASE,1,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
83  ierr = m_field.add_field("HAT_EPS",L2,AINSWORTH_LEGENDRE_BASE,6,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
84  ierr = m_field.add_field("SPATIAL_POSITION_DOT",H1,AINSWORTH_LEGENDRE_BASE,3,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
85  ierr = m_field.add_field("CHEMICAL_LOAD_DOT",H1,AINSWORTH_LEGENDRE_BASE,1,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
86  ierr = m_field.add_field("HAT_EPS_DOT",L2,AINSWORTH_LEGENDRE_BASE,6,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
87 
88  //Add field H1 space rank 3 to approximate geometry using hierarchical basis
89  //For 10 node tetrahedral, before use, geometry is projected on that field (see below)
90  ierr = m_field.add_field("MESH_NODE_POSITIONS",H1,AINSWORTH_LEGENDRE_BASE,3,MB_TAG_SPARSE,MF_ZERO); CHKERRQ(ierr);
91 
92  //meshset consisting all entities in mesh
93  EntityHandle root_set = moab.get_root_set();
94  //add entities to field (root_mesh, i.e. on all mesh entities fields are approx.)
95  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"SPATIAL_POSITION"); CHKERRQ(ierr);
96  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"CHEMICAL_LOAD"); CHKERRQ(ierr);
97  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"HAT_EPS"); CHKERRQ(ierr);
98  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"SPATIAL_POSITION_DOT"); CHKERRQ(ierr);
99  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"CHEMICAL_LOAD_DOT"); CHKERRQ(ierr);
100  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"HAT_EPS_DOT"); CHKERRQ(ierr);
101 
102  PetscBool flg;
103  PetscInt order;
104  ierr = PetscOptionsGetInt(PETSC_NULL,PETSC_NULL,"-my_order",&order,&flg); CHKERRQ(ierr);
105  if(flg != PETSC_TRUE) {
106  order = 2;
107  }
108  if(order < 2) {
109  //SETERRQ()
110  }
111 
112  ierr = m_field.set_field_order(root_set,MBTET,"SPATIAL_POSITION",order); CHKERRQ(ierr);
113  ierr = m_field.set_field_order(root_set,MBTRI,"SPATIAL_POSITION",order); CHKERRQ(ierr);
114  ierr = m_field.set_field_order(root_set,MBEDGE,"SPATIAL_POSITION",order); CHKERRQ(ierr);
115  ierr = m_field.set_field_order(root_set,MBVERTEX,"SPATIAL_POSITION",1); CHKERRQ(ierr);
116 
117  ierr = m_field.set_field_order(root_set,MBTET,"SPATIAL_POSITION_DOT",order); CHKERRQ(ierr);
118  ierr = m_field.set_field_order(root_set,MBTRI,"SPATIAL_POSITION_DOT",order); CHKERRQ(ierr);
119  ierr = m_field.set_field_order(root_set,MBEDGE,"SPATIAL_POSITION_DOT",order); CHKERRQ(ierr);
120  ierr = m_field.set_field_order(root_set,MBVERTEX,"SPATIAL_POSITION_DOT",1); CHKERRQ(ierr);
121 
122  ierr = m_field.set_field_order(root_set,MBTET,"CHEMICAL_LOAD",order-1); CHKERRQ(ierr);
123  ierr = m_field.set_field_order(root_set,MBTRI,"CHEMICAL_LOAD",order-1); CHKERRQ(ierr);
124  ierr = m_field.set_field_order(root_set,MBEDGE,"CHEMICAL_LOAD",order-1); CHKERRQ(ierr);
125  ierr = m_field.set_field_order(root_set,MBVERTEX,"CHEMICAL_LOAD",1); CHKERRQ(ierr);
126 
127  ierr = m_field.set_field_order(root_set,MBTET,"CHEMICAL_LOAD_DOT",order-1); CHKERRQ(ierr);
128  ierr = m_field.set_field_order(root_set,MBTRI,"CHEMICAL_LOAD_DOT",order-1); CHKERRQ(ierr);
129  ierr = m_field.set_field_order(root_set,MBEDGE,"CHEMICAL_LOAD_DOT",order-1); CHKERRQ(ierr);
130  ierr = m_field.set_field_order(root_set,MBVERTEX,"CHEMICAL_LOAD_DOT",1); CHKERRQ(ierr);
131 
132  ierr = m_field.set_field_order(root_set,MBTET,"HAT_EPS",order-1); CHKERRQ(ierr);
133  ierr = m_field.set_field_order(root_set,MBTET,"HAT_EPS_DOT",order-1); CHKERRQ(ierr);
134 
135  //gemetry approximation is set to 2nd oreder
136  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"MESH_NODE_POSITIONS"); CHKERRQ(ierr);
137  ierr = m_field.set_field_order(0,MBTET,"MESH_NODE_POSITIONS",2); CHKERRQ(ierr);
138  ierr = m_field.set_field_order(0,MBTRI,"MESH_NODE_POSITIONS",2); CHKERRQ(ierr);
139  ierr = m_field.set_field_order(0,MBEDGE,"MESH_NODE_POSITIONS",2); CHKERRQ(ierr);
140  ierr = m_field.set_field_order(0,MBVERTEX,"MESH_NODE_POSITIONS",1); CHKERRQ(ierr);
141 
142  ierr = m_field.build_fields(); CHKERRQ(ierr);
143 
144  // Sett geometry approximation and initial spatial positions
145  // 10 node tets
146  if (!check_if_spatial_field_exist) {
147  Projection10NodeCoordsOnField ent_method_material(m_field, "MESH_NODE_POSITIONS");
148  ierr = m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material); CHKERRQ(ierr);
149  Projection10NodeCoordsOnField ent_method_spatial(m_field, "SPATIAL_POSITION");
150  ierr = m_field.loop_dofs("SPATIAL_POSITION", ent_method_spatial); CHKERRQ(ierr);
151  }
152 
153  }
154 
155  //Set finite elements
156  {
157  ierr = m_field.add_finite_element("GEL_FE",MF_ZERO); CHKERRQ(ierr);
158  ierr = m_field.modify_finite_element_add_field_row("GEL_FE","SPATIAL_POSITION"); CHKERRQ(ierr);
159  ierr = m_field.modify_finite_element_add_field_col("GEL_FE","SPATIAL_POSITION"); CHKERRQ(ierr);
160  ierr = m_field.modify_finite_element_add_field_row("GEL_FE","CHEMICAL_LOAD"); CHKERRQ(ierr);
161  ierr = m_field.modify_finite_element_add_field_col("GEL_FE","CHEMICAL_LOAD"); CHKERRQ(ierr);
162  ierr = m_field.modify_finite_element_add_field_row("GEL_FE","HAT_EPS"); CHKERRQ(ierr);
163  ierr = m_field.modify_finite_element_add_field_col("GEL_FE","HAT_EPS"); CHKERRQ(ierr);
164 
165  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","SPATIAL_POSITION"); CHKERRQ(ierr);
166  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","SPATIAL_POSITION_DOT"); CHKERRQ(ierr);
167  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","CHEMICAL_LOAD"); CHKERRQ(ierr);
168  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","CHEMICAL_LOAD_DOT"); CHKERRQ(ierr);
169  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","HAT_EPS"); CHKERRQ(ierr);
170  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","HAT_EPS_DOT"); CHKERRQ(ierr);
171  ierr = m_field.modify_finite_element_add_field_data("GEL_FE","MESH_NODE_POSITIONS"); CHKERRQ(ierr);
172  EntityHandle root_set = moab.get_root_set();
173  ierr = m_field.add_ents_to_finite_element_by_type(root_set,MBTET,"GEL_FE"); CHKERRQ(ierr);
174 
175  //build finite elemnts
176  ierr = m_field.build_finite_elements(); CHKERRQ(ierr);
177  //build adjacencies
178  ierr = m_field.build_adjacencies(bit_level0); CHKERRQ(ierr);
179  }
180 
181  }
182 
183  // Create gel instance
184  Gel gel(m_field);
185  {
186 
187  map<int,Gel::BlockMaterialData> &material_data = gel.blockMaterialData;
188  gel.constitutiveEquationPtr = boost::shared_ptr<Gel::ConstitutiveEquation<adouble> >(
189  new Gel::ConstitutiveEquation<adouble>(gel.blockMaterialData)
190  );
191 
192  // Set material parameters
193  material_data[0].gAlpha = 1;
194  material_data[0].vAlpha = 0.3;
195  material_data[0].gBeta = 1;
196  material_data[0].vBeta = 0.3;
197  material_data[0].gBetaHat = 1;
198  material_data[0].vBetaHat = 0.3;
199  material_data[0].oMega = 1;
200  material_data[0].vIscosity = 1;
201  material_data[0].pErmeability = 2.;
202 
203  int def_block = 0;
204  Tag th_block_id;
205  rval = moab.tag_get_handle(
206  "BLOCK_ID",1,MB_TYPE_INTEGER,th_block_id,MB_TAG_CREAT|MB_TAG_SPARSE,&def_block
207  ); CHKERRQ_MOAB(rval);
208 
209  Gel::CommonData &common_data = gel.commonData;
210  common_data.spatialPositionName = "SPATIAL_POSITION";
211  common_data.spatialPositionNameDot = "SPATIAL_POSITION_DOT";
212  common_data.muName = "CHEMICAL_LOAD";
213  common_data.muNameDot = "CHEMICAL_LOAD_DOT";
214  common_data.strainHatName = "HAT_EPS";
215  common_data.strainHatNameDot = "HAT_EPS_DOT";
216 
217  Gel::GelFE *fe_ptr[] = { &gel.feRhs, &gel.feLhs };
218  for(int ss = 0;ss<2;ss++) {
219  fe_ptr[ss]->getOpPtrVector().push_back(new Gel::OpGetDataAtGaussPts("SPATIAL_POSITION",common_data,false,true));
220  fe_ptr[ss]->getOpPtrVector().push_back(new Gel::OpGetDataAtGaussPts("SPATIAL_POSITION_DOT",common_data,false,true));
221  fe_ptr[ss]->getOpPtrVector().push_back(new Gel::OpGetDataAtGaussPts("CHEMICAL_LOAD",common_data,true,true));
222  fe_ptr[ss]->getOpPtrVector().push_back(new Gel::OpGetDataAtGaussPts("HAT_EPS",common_data,true,false,MBTET));
223  fe_ptr[ss]->getOpPtrVector().push_back(new Gel::OpGetDataAtGaussPts("HAT_EPS_DOT",common_data,true,false,MBTET));
224  }
225 
226  // attach tags for each recorder
227  vector<int> tags;
228  tags.push_back(1);
229  tags.push_back(2);
230  tags.push_back(3);
231  tags.push_back(4);
232 
233  // Right hand side operators
234  gel.feRhs.getOpPtrVector().push_back(
235  new Gel::OpJacobian("SPATIAL_POSITION", tags,gel.constitutiveEquationPtr,gel.commonData,true,false)
236  );
237  gel.feRhs.getOpPtrVector().push_back(
238  new Gel::OpRhsStressTotal(gel.commonData)
239  );
240  gel.feRhs.getOpPtrVector().push_back(
241  new Gel::OpRhsSolventFlux(gel.commonData)
242  );
243  gel.feRhs.getOpPtrVector().push_back(
244  new Gel::OpRhsSolventConcetrationDot(gel.commonData)
245  );
246  gel.feRhs.getOpPtrVector().push_back(
247  new Gel::OpRhsStrainHat(gel.commonData)
248  );
249  // Left hand side operators
250  gel.feLhs.getOpPtrVector().push_back(
251  new Gel::OpJacobian("SPATIAL_POSITION",tags,gel.constitutiveEquationPtr,gel.commonData,false,true)
252  );
253  gel.feLhs.getOpPtrVector().push_back(
254  new Gel::OpLhsdxdx(gel.commonData)
255  );
256  gel.feLhs.getOpPtrVector().push_back(
257  new Gel::OpLhsdxdMu(gel.commonData)
258  );
259  gel.feLhs.getOpPtrVector().push_back(
260  new Gel::OpLhsdxdStrainHat(gel.commonData)
261  );
262  gel.feLhs.getOpPtrVector().push_back(
263  new Gel::OpLhsdStrainHatdStrainHat(gel.commonData)
264  );
265  gel.feLhs.getOpPtrVector().push_back(
266  new Gel::OpLhsdStrainHatdx(gel.commonData)
267  );
268  gel.feLhs.getOpPtrVector().push_back(
269  new Gel::OpLhsdMudMu(gel.commonData)
270  );
271  gel.feLhs.getOpPtrVector().push_back(
272  new Gel::OpLhsdMudx(gel.commonData)
273  );
274 
275  }
276 
277  // Create dm instance
278  DM dm;
279  DMType dm_name = "DMGEL";
280  {
281  ierr = DMRegister_MoFEM(dm_name); CHKERRQ(ierr);
282  ierr = DMCreate(PETSC_COMM_WORLD,&dm);CHKERRQ(ierr);
283  ierr = DMSetType(dm,dm_name);CHKERRQ(ierr);
284  ierr = DMMoFEMCreateMoFEM(dm,&m_field,dm_name,bit_level0); CHKERRQ(ierr);
285  ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
286  //add elements to dm
287  ierr = DMMoFEMAddElement(dm,"GEL_FE"); CHKERRQ(ierr);
288  ierr = DMSetUp(dm); CHKERRQ(ierr);
289  }
290 
291  // Make calculations
292  Mat M;
293  Vec F,U_t;
294  {
295  ierr = DMCreateGlobalVector_MoFEM(dm,&F); CHKERRQ(ierr);
296  ierr = DMCreateGlobalVector_MoFEM(dm,&U_t); CHKERRQ(ierr);
297  ierr = DMCreateMatrix_MoFEM(dm,&M); CHKERRQ(ierr);
298  ierr = VecZeroEntries(F); CHKERRQ(ierr);
299  ierr = VecGhostUpdateBegin(F,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
300  ierr = VecGhostUpdateEnd(F,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
301  ierr = VecZeroEntries(U_t); CHKERRQ(ierr);
302  ierr = VecGhostUpdateBegin(U_t,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
303  ierr = VecGhostUpdateEnd(U_t,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
304  ierr = MatZeroEntries(M); CHKERRQ(ierr);
305  gel.feRhs.ts_F = F; // Set right hand side vector manually
306  gel.feRhs.ts_u_t = U_t;
307  gel.feRhs.ts_ctx = TSMethod::CTX_TSSETIFUNCTION;
308  ierr = DMoFEMLoopFiniteElements(dm,"GEL_FE",&gel.feRhs); CHKERRQ(ierr);
309  gel.feLhs.ts_B = M; // Set matrix M
310  gel.feLhs.ts_a = 1.0; // Set time step parameter
311  gel.feLhs.ts_ctx = TSMethod::CTX_TSSETIJACOBIAN;
312  ierr = DMoFEMLoopFiniteElements(dm,"GEL_FE",&gel.feLhs); CHKERRQ(ierr);
313  ierr = VecAssemblyBegin(F); CHKERRQ(ierr);
314  ierr = VecAssemblyEnd(F); CHKERRQ(ierr);
315  ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
316  ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
317  }
318 
319  // See results
320  {
321 
322  PetscViewer viewer;
323  ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,"gel_jacobian_test.txt",&viewer); CHKERRQ(ierr);
324  ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
325  //MatView(M,PETSC_VIEWER_DRAW_WORLD);
326  //std::string wait;
327  //std::cin >> wait;
328 
329  }
330 
331  double mnorm;
332  ierr = MatNorm(M,NORM_1,&mnorm); CHKERRQ(ierr);
333  ierr = PetscPrintf(PETSC_COMM_WORLD,"mnorm = %9.8e\n",mnorm); CHKERRQ(ierr);
334 
335  if(fabs(mnorm-6.54443978e+01)>1e-6) {
336  SETERRQ(PETSC_COMM_WORLD,MOFEM_ATOM_TEST_INVALID,"Failed to pass test");
337  }
338 
339  // Clean and destroy
340  {
341  ierr = VecDestroy(&F); CHKERRQ(ierr);
342  ierr = VecDestroy(&U_t); CHKERRQ(ierr);
343  ierr = MatDestroy(&M); CHKERRQ(ierr);
344  ierr = DMDestroy(&dm); CHKERRQ(ierr);
345  }
346 
347  }
348  CATCH_ERRORS;
349 
351 
352  return 0;
353 }

Variable Documentation

◆ help

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

Definition at line 41 of file gel_jacobian_test.cpp.

GelModule::Gel::CommonData::strainHatName
string strainHatName
Definition: Gels.hpp:349
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.
GelModule::Gel::OpLhsdStrainHatdx
Assemble matrix .
Definition: Gels.hpp:1753
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
GelModule::Gel::CommonData::strainHatNameDot
string strainHatNameDot
Definition: Gels.hpp:350
H1
@ H1
continuous field
Definition: definitions.h:85
GelModule::Gel::OpJacobian
Definition: Gels.hpp:541
EntityHandle
CHKERRQ_MOAB
#define CHKERRQ_MOAB(a)
check error code of MoAB function
Definition: definitions.h:467
GelModule::Gel::CommonData::spatialPositionName
string spatialPositionName
Definition: Gels.hpp:347
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
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
help
static char help[]
Definition: gel_jacobian_test.cpp:41
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
GelModule::Gel::CommonData::spatialPositionNameDot
string spatialPositionNameDot
Definition: Gels.hpp:348
GelModule::Gel::OpRhsStrainHat
Residual strain hat.
Definition: Gels.hpp:1311
GelModule::Gel::OpLhsdxdMu
Assemble matrix .
Definition: Gels.hpp:1497
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
order
constexpr int order
Definition: dg_projection.cpp:18
GelModule::Gel
Implementation of Gel constitutive model.
Definition: Gels.hpp:74
GelModule::Gel::CommonData::muName
string muName
Definition: Gels.hpp:351
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
GelModule::Gel::OpLhsdMudx
Assemble matrix .
Definition: Gels.hpp:1925
MoFEM::DMCreateGlobalVector_MoFEM
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1167
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
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
MoFEM::DMCreateMatrix_MoFEM
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1197
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::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.
GelModule::Gel::CommonData
Common data for gel model.
Definition: Gels.hpp:345
GelModule::Gel::CommonData::muNameDot
string muNameDot
Definition: Gels.hpp:352
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
GelModule::Gel::OpGetDataAtGaussPts
Definition: Gels.hpp:444
MoFEM::CoreInterface::check_field
virtual bool check_field(const std::string &name) const =0
check if field is in database
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
GelModule::Gel::OpLhsdMudMu
Assemble matrix .
Definition: Gels.hpp:1842
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
GelModule::Gel::OpRhsSolventConcetrationDot
Calculating right hand side.
Definition: Gels.hpp:1260
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
MF_ZERO
@ MF_ZERO
Definition: definitions.h:111
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
GelModule::Gel::OpLhsdxdStrainHat
Assemble matrix .
Definition: Gels.hpp:1580
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::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::DeprecatedCoreInterface::seed_ref_level_3D
DEPRECATED MoFEMErrorCode seed_ref_level_3D(const EntityHandle meshset, const BitRefLevel &bit, int verb=-1)
seed 2D entities in the meshset and their adjacencies (only TETs adjacencies) in a particular BitRefL...
Definition: DeprecatedCoreInterface.cpp:18
GelModule::Gel::OpRhsStressTotal
Assemble internal force vector.
Definition: Gels.hpp:1154
GelModule::Gel::OpLhsdStrainHatdStrainHat
Assemble matrix .
Definition: Gels.hpp:1671
GelModule::Gel::ConstitutiveEquation
Constitutive model functions.
Definition: Gels.hpp:114
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::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::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
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
GelModule::Gel::OpRhsSolventFlux
Calculate internal forces for solvent flux.
Definition: Gels.hpp:1208
GelModule::Gel::OpLhsdxdx
Assemble matrix .
Definition: Gels.hpp:1402
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.
GelModule::Gel::GelFE
definition of volume element
Definition: Gels.hpp:379
F
@ F
Definition: free_surface.cpp:394