v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
CommonData Struct Reference

#include <users_modules/topology_optimization/src/Header.hpp>

Collaboration diagram for CommonData:
[legend]

Public Member Functions

 CommonData (MoFEM::Interface &m_field)
 
MoFEMErrorCode getParameters ()
 
MoFEMErrorCode getSensitivityFromTag (const EntityHandle fe_ent, const int nb_gauss_pts)
 
MoFEMErrorCode setSensitivityToTag (const EntityHandle fe_ent, const double sens)
 
MoFEMErrorCode getDensityFromTag (const EntityHandle fe_ent, const int nb_gauss_pts)
 
MoFEMErrorCode setDensityToTag (const EntityHandle fe_ent, const double rho)
 
MoFEMErrorCode createTags (moab::Interface &moab)
 
MoFEMErrorCode build_main_dm (DM &dm, DM &sub_dm_elastic, DM &sub_dm_rho, MoFEM::Interface &m_field, BitRefLevel &bit_level0, PetscBool is_partitioned)
 
MoFEMErrorCode build_sub_dms (DM &dm, DM &sub_dm_elastic, DM &sub_dm_rho, MoFEM::Interface &m_field, BitRefLevel &bit_level0, PetscBool is_partitioned)
 
MoFEMErrorCode getBoundBox (const double material_coords[3], std::vector< double > &coords, double &boxSize)
 

Public Attributes

MatrixDouble dotEdge
 
MatrixDouble dotEleLeft
 
MatrixDouble dotEleRight
 
MoFEM::InterfacemField
 
Tag thSensitivity
 
Tag thDens
 
Vec massVec
 
Vec changeVec
 
boost::shared_ptr< VectorDouble > dataRhotGaussPtr
 
boost::shared_ptr< MatrixDouble > dispGradAtGaussPts
 
MatrixDouble D
 
VectorDouble sTrain
 
vector< VectorDouble > sTress
 
double volumeFrac
 
MatrixDouble cOmplianceMat
 
VectorDouble sEnsitivityVec
 
VectorDouble rhoVecAtTag
 
MatrixDouble anisLambdaMat
 
double edgeLength
 
double lAmbda
 
double mU
 
double lAmbda0
 
double rHo
 
double yOung
 
double pOisson
 
int penalty
 
double mOve
 
double fIlter_radius
 
double cHange
 
double cUrrent_mass
 
double L_mid
 
bool uPdate_density
 
SmartPetscObj< Mat > M
 Mass matrix. More...
 
SmartPetscObj< KSP > ksp
 Linear solver. More...
 

Detailed Description

Examples
ContactOps.hpp, NavierStokesElement.hpp, PlasticOps.hpp, Remodeling.hpp, SmallStrainPlasticity.hpp, child_and_parent.cpp, continuity_check_on_contact_prism_side_ele.cpp, continuity_check_on_skeleton_3d.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, dg_projection.cpp, dynamic_first_order_con_law.cpp, hanging_node_approx.cpp, higher_derivatives.cpp, mixed_poisson.cpp, photon_diffusion.cpp, and reaction_diffusion.cpp.

Definition at line 22 of file continuity_check_on_skeleton_with_simple_2d_for_h1.cpp.

Constructor & Destructor Documentation

◆ CommonData()

CommonData::CommonData ( MoFEM::Interface m_field)
inline

Definition at line 56 of file Header.hpp.

56 : mField(m_field) {
57
58 dataRhotGaussPtr = boost::shared_ptr<VectorDouble>(new VectorDouble());
59 dispGradAtGaussPts = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
60 volumeFrac = 0.5;
61 yOung = 10;
62 fIlter_radius = 0;
63 penalty = 3;
64 cHange = 0;
65 edgeLength = 1.;
66 lAmbda0 = 0.;
67 cUrrent_mass = 0;
68 L_mid = 0; // langange multiplier
69 uPdate_density = false;
70
71 anisLambdaMat.resize(3,3,false);
72 anisLambdaMat.clear();
73
74 ierr = createTags(m_field.get_moab());
75 }
static PetscErrorCode ierr
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
boost::shared_ptr< VectorDouble > dataRhotGaussPtr
Definition: Header.hpp:28
MoFEM::Interface & mField
Definition: Header.hpp:21
double volumeFrac
Definition: Header.hpp:34
double L_mid
Definition: Header.hpp:53
MatrixDouble anisLambdaMat
Definition: Header.hpp:38
double yOung
Definition: Header.hpp:45
double cHange
Definition: Header.hpp:50
bool uPdate_density
Definition: Header.hpp:54
boost::shared_ptr< MatrixDouble > dispGradAtGaussPts
Definition: Header.hpp:29
double cUrrent_mass
Definition: Header.hpp:52
double lAmbda0
Definition: Header.hpp:43
double fIlter_radius
Definition: Header.hpp:49
int penalty
Definition: Header.hpp:47
double edgeLength
Definition: Header.hpp:40
MoFEMErrorCode createTags(moab::Interface &moab)
Definition: Header.hpp:232
virtual moab::Interface & get_moab()=0

Member Function Documentation

◆ build_main_dm()

MoFEMErrorCode CommonData::build_main_dm ( DM &  dm,
DM &  sub_dm_elastic,
DM &  sub_dm_rho,
MoFEM::Interface m_field,
BitRefLevel &  bit_level0,
PetscBool  is_partitioned 
)
inline

Definition at line 244 of file Header.hpp.

246 {
248 DMType dm_name = "DM_TOPO";
249 // Register DM problem
250 CHKERR DMRegister_MoFEM(dm_name);
251 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
252 CHKERR DMSetType(dm, dm_name);
253 CHKERR DMMoFEMCreateMoFEM(dm, &m_field, dm_name, bit_level0);
254 CHKERR DMSetFromOptions(dm);
255 CHKERR DMMoFEMAddElement(dm, "ELASTIC");
256 CHKERR DMMoFEMAddElement(dm, "RHO_APPROX");
257 if (m_field.check_finite_element("PRESSURE_FE"))
258 CHKERR DMMoFEMAddElement(dm, "PRESSURE_FE");
259 if (m_field.check_finite_element("FORCE_FE"))
260 CHKERR DMMoFEMAddElement(dm, "FORCE_FE");
261 CHKERR DMMoFEMSetIsPartitioned(dm, is_partitioned);
262 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
263 CHKERR DMSetUp(dm);
264 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
265
266
268 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1109
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:483
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:118
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ build_sub_dms()

MoFEMErrorCode CommonData::build_sub_dms ( DM &  dm,
DM &  sub_dm_elastic,
DM &  sub_dm_rho,
MoFEM::Interface m_field,
BitRefLevel &  bit_level0,
PetscBool  is_partitioned 
)
inline

Definition at line 270 of file Header.hpp.

273 {
274
276 DMType dm_name_sub_elast = "DM_ELASTIC";
277 CHKERR DMRegister_MoFEM(dm_name_sub_elast);
278 CHKERR DMCreate(PETSC_COMM_WORLD, &sub_dm_elastic);
279 CHKERR DMSetType(sub_dm_elastic, dm_name_sub_elast);
280 // m_field.getInterface<ProblemsManager>()->buildProblemFromFields =
281 // PETSC_TRUE;
282 CHKERR DMMoFEMCreateSubDM(sub_dm_elastic, dm, "ELASTIC_PROB");
283 CHKERR DMMoFEMSetDestroyProblem(sub_dm_elastic, PETSC_TRUE);
284 CHKERR DMMoFEMAddSubFieldRow(sub_dm_elastic, "DISPLACEMENTS");
285 CHKERR DMMoFEMAddSubFieldCol(sub_dm_elastic, "DISPLACEMENTS");
286 CHKERR DMMoFEMAddElement(sub_dm_elastic, "ELASTIC");
287 if (m_field.check_finite_element("PRESSURE_FE"))
288 CHKERR DMMoFEMAddElement(sub_dm_elastic, "PRESSURE_FE");
289 if (m_field.check_finite_element("FORCE_FE"))
290 CHKERR DMMoFEMAddElement(sub_dm_elastic, "FORCE_FE");
291
292 CHKERR DMMoFEMSetSquareProblem(sub_dm_elastic, PETSC_TRUE);
293 CHKERR DMMoFEMSetIsPartitioned(sub_dm_elastic, is_partitioned);
294 CHKERR DMSetUp(sub_dm_elastic);
295
296 DMType dm_name_sub_rho = "DM_DENSITY";
297 CHKERR DMRegister_MoFEM(dm_name_sub_rho);
298 CHKERR DMCreate(PETSC_COMM_WORLD, &sub_dm_rho);
299 CHKERR DMSetType(sub_dm_rho, dm_name_sub_rho);
300 // m_field.getInterface<ProblemsManager>()->buildProblemFromFields =
301 // PETSC_TRUE;
302 CHKERR DMMoFEMCreateSubDM(sub_dm_rho, dm, "DENSITY_PROB");
303 CHKERR DMMoFEMSetDestroyProblem(sub_dm_rho, PETSC_TRUE);
304 CHKERR DMMoFEMAddSubFieldRow(sub_dm_rho, "RHO");
305 CHKERR DMMoFEMAddSubFieldCol(sub_dm_rho, "RHO");
306 CHKERR DMMoFEMAddElement(sub_dm_rho, "RHO_APPROX");
307 CHKERR DMMoFEMSetSquareProblem(sub_dm_rho, PETSC_TRUE);
308 CHKERR DMMoFEMSetIsPartitioned(sub_dm_rho, is_partitioned);
309 CHKERR DMSetUp(sub_dm_rho);
311 }
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition: DMMoFEM.cpp:219
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition: DMMoFEM.cpp:442
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:242
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:275
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
Definition: DMMoFEM.cpp:424

◆ createTags()

MoFEMErrorCode CommonData::createTags ( moab::Interface &  moab)
inline

Definition at line 232 of file Header.hpp.

232 {
234 std::vector<double> def_val(9, 0);
235 CHKERR mField.get_moab().tag_get_handle(
236 "SENSITIVITY", 1, MB_TYPE_DOUBLE, thSensitivity,
237 MB_TAG_CREAT | MB_TAG_VARLEN | MB_TAG_SPARSE, &def_val[0]);
238 CHKERR mField.get_moab().tag_get_handle(
239 "RHO", 1, MB_TYPE_DOUBLE, thDens,
240 MB_TAG_CREAT | MB_TAG_VARLEN | MB_TAG_SPARSE, &def_val[0]);
242 }
Tag thDens
Definition: Header.hpp:23
Tag thSensitivity
Definition: Header.hpp:22

◆ getBoundBox()

MoFEMErrorCode CommonData::getBoundBox ( const double  material_coords[3],
std::vector< double > &  coords,
double boxSize 
)
inline

Definition at line 313 of file Header.hpp.

314 {
316 FTensor::Tensor1<double, 3> t_material_coords(
317 material_coords[0], material_coords[1], material_coords[2]);
318 FTensor::Tensor1<double *, 3> t_coords(&coords[0], &coords[1], &coords[2],
319 3);
321 double max_xyz[] = {0, 0, 0};
322 for (int bb = 0; bb != coords.size() / 3; bb++) {
323 t_coords(i) -= t_material_coords(i);
324 for (int dd = 0; dd != 3; dd++) {
325 max_xyz[dd] = (max_xyz[dd] < fabs(t_coords(dd))) ? fabs(t_coords(dd))
326 : max_xyz[dd];
327 }
328 ++t_coords;
329 }
330 //take the max of the box
331 boxSize = (max_xyz[0] > max_xyz[1]) ? max_xyz[0] : max_xyz[1];
332 boxSize = (boxSize > max_xyz[2]) ? boxSize : max_xyz[2];
334 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
FTensor::Index< 'i', SPACE_DIM > i
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33

◆ getDensityFromTag()

MoFEMErrorCode CommonData::getDensityFromTag ( const EntityHandle  fe_ent,
const int  nb_gauss_pts 
)
inline

Definition at line 198 of file Header.hpp.

198 {
200 double *tag_data;
201 int tag_size;
202 CHKERR mField.get_moab().tag_get_by_ptr(
203 thDens, &fe_ent, 1, (const void **)&tag_data, &tag_size);
204
205 if (tag_size == 1) {
206 rhoVecAtTag.resize(nb_gauss_pts , false);
207 rhoVecAtTag.clear();
208 void const *tag_data[] = {&*rhoVecAtTag.data().begin()};
209 const int tag_size = rhoVecAtTag.data().size();
210 CHKERR mField.get_moab().tag_set_by_ptr(thDens, &fe_ent, 1,
211 tag_data, &tag_size);
212 } else if (tag_size != nb_gauss_pts) {
213 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
214 "Wrong size of the tag, data inconsistency");
215 } else {
217 tag_size, ublas::shallow_array_adaptor<double>(tag_size, tag_data));
218 rhoVecAtTag = tag_vec;
219 }
221 }
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
VectorDouble rhoVecAtTag
Definition: Header.hpp:37

◆ getParameters()

MoFEMErrorCode CommonData::getParameters ( )
inline

Definition at line 76 of file Header.hpp.

76 {
78 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Topo Config", "none");
79 CHKERR(ierr);
80
81 CHKERR PetscOptionsInt("-penalty", "default penalty ", "", penalty, &penalty,
82 PETSC_NULL);
83
84 CHKERR PetscOptionsScalar("-young", "default yOung modulus", "", yOung, &yOung,
85 PETSC_NULL);
86
87 CHKERR PetscOptionsScalar("-poisson", "poisson ratio", "", 0.1, &pOisson,
88 PETSC_NULL);
89
90 CHKERR PetscOptionsScalar("-move", "max move", "", 0.1, &mOve, PETSC_NULL);
91 // CHKERR PetscOptionsScalar("-lambda", "lambda smooth", "", lAmbda0, &lAmbda0, PETSC_NULL);
92
93 double lambda_coef[3] = {0.,0.,0.};
94 int nmax = 3;
95 PetscBool is_lambda_defined = PETSC_FALSE;
96 PetscOptionsGetRealArray(NULL, "", "-lambda", lambda_coef,
97 &nmax, &is_lambda_defined);
98
99 switch (nmax) {
100
101 case 0:
102 case 1:
103 anisLambdaMat(0, 0) = anisLambdaMat(1, 1) = anisLambdaMat(2, 2) =
104 lambda_coef[0];
105 break;
106 case 3:
107 anisLambdaMat(0, 0) = lambda_coef[0];
108 anisLambdaMat(1, 1) = lambda_coef[1];
109 anisLambdaMat(2, 2) = lambda_coef[2];
110 break;
111 default:
112
113 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
114 "Wrong number of lambda parameters, should be 1 or 3, data "
115 "inconsistency");
116 }
117
118 CHKERR PetscOptionsScalar("-volume_frac", "volume frac", "", volumeFrac,
119 &volumeFrac, PETSC_NULL);
120
121 // CHKERR PetscOptionsScalar("-filter_r", "default filter radius", "", 0,
122 // &fIlter_radius, PETSC_NULL);
123 ierr = PetscOptionsEnd();
124 CHKERRG(ierr);
125
126 lAmbda = (yOung * pOisson) / ((1. + pOisson) * (1. - 2. * pOisson));
127 mU = yOung / (2. * (1. + pOisson));
128
129 const double coefficient = yOung / ((1 + pOisson) * (1 - 2 * pOisson));
130 D.resize(36, 1, false);
131 D.clear();
132
134
139
140 t_D(i, j, k, l) = 0.;
141
142 t_D(0, 0, 0, 0) = 1 - pOisson;
143 t_D(1, 1, 1, 1) = 1 - pOisson;
144 t_D(2, 2, 2, 2) = 1 - pOisson;
145
146 t_D(0, 1, 0, 1) = 0.5 * (1 - 2 * pOisson);
147 t_D(0, 2, 0, 2) = 0.5 * (1 - 2 * pOisson);
148 t_D(1, 2, 1, 2) = 0.5 * (1 - 2 * pOisson);
149
150 t_D(0, 0, 1, 1) = pOisson;
151 t_D(1, 1, 0, 0) = pOisson;
152 t_D(0, 0, 2, 2) = pOisson;
153 t_D(2, 2, 0, 0) = pOisson;
154 t_D(1, 1, 2, 2) = pOisson;
155 t_D(2, 2, 1, 1) = pOisson;
156
157 t_D(i, j, k, l) *= coefficient;
158
160 }
#define MAT_TO_DDG(SM)
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
double lAmbda
Definition: Header.hpp:41
double mOve
Definition: Header.hpp:48
double mU
Definition: Header.hpp:42
double pOisson
Definition: Header.hpp:46
MatrixDouble D
Definition: Header.hpp:31

◆ getSensitivityFromTag()

MoFEMErrorCode CommonData::getSensitivityFromTag ( const EntityHandle  fe_ent,
const int  nb_gauss_pts 
)
inline

Definition at line 163 of file Header.hpp.

163 {
165 double *tag_data;
166 int tag_size;
167 CHKERR mField.get_moab().tag_get_by_ptr(
168 thSensitivity, &fe_ent, 1, (const void **)&tag_data, &tag_size);
169
170 if (tag_size == 1) {
171 sEnsitivityVec.resize(nb_gauss_pts , false);
172 sEnsitivityVec.clear();
173 void const *tag_data[] = {&*sEnsitivityVec.data().begin()};
174 const int tag_size = sEnsitivityVec.data().size();
175 CHKERR mField.get_moab().tag_set_by_ptr(thSensitivity, &fe_ent, 1,
176 tag_data, &tag_size);
177 } else if (tag_size != nb_gauss_pts) {
178 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
179 "Wrong size of the tag, data inconsistency");
180 } else {
182 tag_size, ublas::shallow_array_adaptor<double>(tag_size, tag_data));
183 sEnsitivityVec = tag_vec;
184 }
186 }
VectorDouble sEnsitivityVec
Definition: Header.hpp:36

◆ setDensityToTag()

MoFEMErrorCode CommonData::setDensityToTag ( const EntityHandle  fe_ent,
const double  rho 
)
inline

Definition at line 223 of file Header.hpp.

223 {
225 void const *tag_data[] = {&*rhoVecAtTag.data().begin()};
226 const int tag_size = rhoVecAtTag.data().size();
227 CHKERR mField.get_moab().tag_set_by_ptr(thDens, &fe_ent, 1,
228 tag_data, &tag_size);
230 }

◆ setSensitivityToTag()

MoFEMErrorCode CommonData::setSensitivityToTag ( const EntityHandle  fe_ent,
const double  sens 
)
inline

Definition at line 189 of file Header.hpp.

189 {
191 void const *tag_data[] = {&*sEnsitivityVec.data().begin()};
192 const int tag_size = sEnsitivityVec.data().size();
193 CHKERR mField.get_moab().tag_set_by_ptr(thSensitivity, &fe_ent, 1,
194 tag_data, &tag_size);
196 }

Member Data Documentation

◆ anisLambdaMat

MatrixDouble CommonData::anisLambdaMat

Definition at line 38 of file Header.hpp.

◆ cHange

double CommonData::cHange

Definition at line 50 of file Header.hpp.

◆ changeVec

Vec CommonData::changeVec

Definition at line 26 of file Header.hpp.

◆ cOmplianceMat

MatrixDouble CommonData::cOmplianceMat

Definition at line 35 of file Header.hpp.

◆ cUrrent_mass

double CommonData::cUrrent_mass

Definition at line 52 of file Header.hpp.

◆ D

MatrixDouble CommonData::D

Definition at line 31 of file Header.hpp.

◆ dataRhotGaussPtr

boost::shared_ptr<VectorDouble> CommonData::dataRhotGaussPtr

Definition at line 28 of file Header.hpp.

◆ dispGradAtGaussPts

boost::shared_ptr<MatrixDouble> CommonData::dispGradAtGaussPts

Definition at line 29 of file Header.hpp.

◆ dotEdge

MatrixDouble CommonData::dotEdge

◆ dotEleLeft

MatrixDouble CommonData::dotEleLeft

◆ dotEleRight

MatrixDouble CommonData::dotEleRight

◆ edgeLength

double CommonData::edgeLength

Definition at line 40 of file Header.hpp.

◆ fIlter_radius

double CommonData::fIlter_radius

Definition at line 49 of file Header.hpp.

◆ ksp

SmartPetscObj<KSP> CommonData::ksp

Linear solver.

Examples
dynamic_first_order_con_law.cpp.

Definition at line 416 of file dynamic_first_order_con_law.cpp.

◆ L_mid

double CommonData::L_mid

Definition at line 53 of file Header.hpp.

◆ lAmbda

double CommonData::lAmbda

Definition at line 41 of file Header.hpp.

◆ lAmbda0

double CommonData::lAmbda0

Definition at line 43 of file Header.hpp.

◆ M

SmartPetscObj<Mat> CommonData::M

Mass matrix.

Examples
dynamic_first_order_con_law.cpp.

Definition at line 415 of file dynamic_first_order_con_law.cpp.

◆ massVec

Vec CommonData::massVec

Definition at line 25 of file Header.hpp.

◆ mField

MoFEM::Interface& CommonData::mField

Definition at line 21 of file Header.hpp.

◆ mOve

double CommonData::mOve

Definition at line 48 of file Header.hpp.

◆ mU

double CommonData::mU

Definition at line 42 of file Header.hpp.

◆ penalty

int CommonData::penalty

Definition at line 47 of file Header.hpp.

◆ pOisson

double CommonData::pOisson

Definition at line 46 of file Header.hpp.

◆ rHo

double CommonData::rHo

Definition at line 44 of file Header.hpp.

◆ rhoVecAtTag

VectorDouble CommonData::rhoVecAtTag

Definition at line 37 of file Header.hpp.

◆ sEnsitivityVec

VectorDouble CommonData::sEnsitivityVec

Definition at line 36 of file Header.hpp.

◆ sTrain

VectorDouble CommonData::sTrain

Definition at line 32 of file Header.hpp.

◆ sTress

vector<VectorDouble> CommonData::sTress

Definition at line 33 of file Header.hpp.

◆ thDens

Tag CommonData::thDens

Definition at line 23 of file Header.hpp.

◆ thSensitivity

Tag CommonData::thSensitivity

Definition at line 22 of file Header.hpp.

◆ uPdate_density

bool CommonData::uPdate_density

Definition at line 54 of file Header.hpp.

◆ volumeFrac

double CommonData::volumeFrac

Definition at line 34 of file Header.hpp.

◆ yOung

double CommonData::yOung

Definition at line 45 of file Header.hpp.


The documentation for this struct was generated from the following files: