v0.14.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
MoFEM::OperatorsTester Struct Reference

Calculate directional derivative of the right hand side and compare it with tangent matrix derivative. More...

#include <src/interfaces/OperatorsTester.hpp>

Inheritance diagram for MoFEM::OperatorsTester:
[legend]
Collaboration diagram for MoFEM::OperatorsTester:
[legend]

Public Types

using RandomFieldData = std::pair< std::string, std::array< double, 2 > >
 

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 OperatorsTester (const MoFEM::Core &core)
 
virtual ~OperatorsTester ()=default
 
SmartPetscObj< Vec > setRandomFields (SmartPetscObj< DM > dm, std::vector< RandomFieldData > random_fields, boost::shared_ptr< Range > ents=nullptr)
 Generate random fileds. More...
 
SmartPetscObj< Vec > assembleVec (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr)
 Assemble the right hand side vector. More...
 
SmartPetscObj< Mat > assembleMat (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr)
 Assemble the left hand side vector. More...
 
SmartPetscObj< Vec > directionalCentralFiniteDifference (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, SmartPetscObj< Vec > diff_x, double time, double delta_t, double eps, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
 Calculate directional directive using finite difference. More...
 
SmartPetscObj< Vec > checkCentralFiniteDifference (SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline_rhs, boost::shared_ptr< FEMethod > pipeline_lhs, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, SmartPetscObj< Vec > diff_x, double time, double delta_t, double eps, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
 Check consistency between directional derivative with matrix. More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 

Private Member Functions

std::pair< SmartPetscObj< Vec >, SmartPetscObj< Vec > > setPipelineX (boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double delta_t)
 Set vectors, x, x_t, and x_tt to finite element instance. More...
 

Private Attributes

MoFEM::CorecOre
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 

Detailed Description

Calculate directional derivative of the right hand side and compare it with tangent matrix derivative.

Examples
free_surface.cpp, level_set.cpp, operators_tests.cpp, and tensor_divergence_operator.cpp.

Definition at line 21 of file OperatorsTester.hpp.

Member Typedef Documentation

◆ RandomFieldData

using MoFEM::OperatorsTester::RandomFieldData = std::pair<std::string, std::array<double, 2> >

Definition at line 29 of file OperatorsTester.hpp.

Constructor & Destructor Documentation

◆ OperatorsTester()

MoFEM::OperatorsTester::OperatorsTester ( const MoFEM::Core core)

Definition at line 13 of file OperatorsTester.cpp.

14 : cOre(const_cast<Core &>(core)) {}
CoreTmp< 0 > Core
Definition: Core.hpp:1094

◆ ~OperatorsTester()

virtual MoFEM::OperatorsTester::~OperatorsTester ( )
virtualdefault

Member Function Documentation

◆ assembleMat()

SmartPetscObj< Mat > MoFEM::OperatorsTester::assembleMat ( SmartPetscObj< DM >  dm,
std::string  fe_name,
boost::shared_ptr< FEMethod pipeline,
SmartPetscObj< Vec >  x,
SmartPetscObj< Vec >  delta_x,
SmartPetscObj< Vec >  delta2_x,
double  time,
double  delta_t,
CacheTupleWeakPtr  cache_ptr 
)

Assemble the left hand side vector.

Parameters
dm
fe_name// fe name
pipeline// pipeline, i.e. fe instance
x// problem (dm) vector
delta_x// vector for x rate, can be null, i.e, SmartPetscObj<Vec>()
delta2_x// vector for x second rate, i.e. acceleration
time// time
delta_t// time increment
cache_ptr// finite element data cache, can be null
Returns
SmartPetscObj<Vec>

Definition at line 120 of file OperatorsTester.cpp.

124 {
125 auto m = createDMMatrix(dm);
126
127 pipeline->data_ctx =
129 pipeline->A = m;
130 pipeline->B = m;
131
132 pipeline->ts_t = time;
133
134 pipeline->ts_dt = delta_t;
135 pipeline->ts_a = 1. / delta_t;
136 pipeline->ts_aa = 1. / pow(delta_t, 2);
137
138 auto [v, a] = setPipelineX(pipeline, x, delta_x, delta2_x, delta_t);
139 std::ignore = a;
140 std::ignore = v;
141
143 DMoFEMMeshToLocalVector(dm, x, INSERT_VALUES, SCATTER_REVERSE),
144 "loc to vec");
145 CHK_THROW_MESSAGE(DMoFEMLoopFiniteElements(dm, fe_name, pipeline, cache_ptr),
146 "Run assemble pipeline");
147 CHK_THROW_MESSAGE(MatAssemblyBegin(m, MAT_FINAL_ASSEMBLY), "Asmb mat");
148 CHK_THROW_MESSAGE(MatAssemblyEnd(m, MAT_FINAL_ASSEMBLY), "Asmb mat");
149
150 return m;
151}
constexpr double a
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:595
FTensor::Index< 'm', SPACE_DIM > m
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:509
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:572
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition: DMMoFEM.hpp:988
const double v
phase velocity of light in medium (cm/ns)
std::pair< SmartPetscObj< Vec >, SmartPetscObj< Vec > > setPipelineX(boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double delta_t)
Set vectors, x, x_t, and x_tt to finite element instance.
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38
static constexpr Switches CtxSetTime
Definition: LoopMethods.hpp:42

◆ assembleVec()

SmartPetscObj< Vec > MoFEM::OperatorsTester::assembleVec ( SmartPetscObj< DM >  dm,
std::string  fe_name,
boost::shared_ptr< FEMethod pipeline,
SmartPetscObj< Vec >  x,
SmartPetscObj< Vec >  delta_x,
SmartPetscObj< Vec >  delta2_x,
double  time,
double  delta_t,
CacheTupleWeakPtr  cache_ptr 
)

Assemble the right hand side vector.

Parameters
dm
fe_name// fe name
pipeline// pipeline, i.e. fe instance
x// problem (dm) vector
delta_x// vector for x rate, can be null, i.e, SmartPetscObj<Vec>()
delta2_x// vector for x second rate, i.e. acceleration
time// time
delta_t// time increment
cache_ptr// finite element data cache, can be null
Returns
SmartPetscObj<Vec>

Definition at line 83 of file OperatorsTester.cpp.

87 {
88
89 pipeline->data_ctx = PetscData::CtxSetF | PetscData::CtxSetTime;
90 auto f = createDMVector(dm);
91 pipeline->f = f;
92 pipeline->ts_t = time;
93 pipeline->ts_dt = delta_t;
94
95 auto [v, a] = setPipelineX(pipeline, x, delta_x, delta2_x, delta_t);
96 std::ignore = a;
97 std::ignore = v;
98
100 DMoFEMMeshToLocalVector(dm, x, INSERT_VALUES, SCATTER_REVERSE),
101 "loc to vec");
102 CHK_THROW_MESSAGE(DMoFEMLoopFiniteElements(dm, fe_name, pipeline, cache_ptr),
103 "Run assemble pipeline");
104 CHK_THROW_MESSAGE(VecAssemblyBegin(f), "Asmb vector");
105 CHK_THROW_MESSAGE(VecAssemblyEnd(f), "Asmb vector");
106
107 CHK_THROW_MESSAGE(VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE),
108 "Update vec");
109 CHK_THROW_MESSAGE(VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE),
110 "Update vec");
111 CHK_THROW_MESSAGE(VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD),
112 "Update vec");
113 CHK_THROW_MESSAGE(VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD),
114 "Update vec");
115
116 return f;
117}
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1003
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:36

◆ checkCentralFiniteDifference()

SmartPetscObj< Vec > MoFEM::OperatorsTester::checkCentralFiniteDifference ( SmartPetscObj< DM >  dm,
std::string  fe_name,
boost::shared_ptr< FEMethod pipeline_rhs,
boost::shared_ptr< FEMethod pipeline_lhs,
SmartPetscObj< Vec >  x,
SmartPetscObj< Vec >  delta_x,
SmartPetscObj< Vec >  delta2_x,
SmartPetscObj< Vec >  diff_x,
double  time,
double  delta_t,
double  eps,
CacheTupleWeakPtr  cache_ptr = CacheTupleSharedPtr() 
)

Check consistency between directional derivative with matrix.

Parameters
dm
fe_name
pipeline_rhs
pipeline_lhs
x
delta_x
delta2_x
diff_x
time
delta_t
eps
cache_ptr
Returns
SmartPetscObj<Vec>

Definition at line 210 of file OperatorsTester.cpp.

216 {
217
219 dm, fe_name, pipeline_rhs, x, delta_x, delta2_x, diff_x, time, delta_t,
220 eps, cache_ptr);
221
222 auto m = assembleMat(dm, fe_name, pipeline_lhs, x, delta_x, delta2_x, time,
223 delta_t, cache_ptr);
224
225 auto fm = vectorDuplicate(fd_diff);
226 CHK_THROW_MESSAGE(MatMult(m, diff_x, fm), "Mat mult");
227 CHK_THROW_MESSAGE(VecAXPY(fd_diff, -1., fm), "Add");
228
229 return fd_diff;
230}
static const double eps
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
SmartPetscObj< Mat > assembleMat(SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr)
Assemble the left hand side vector.
SmartPetscObj< Vec > directionalCentralFiniteDifference(SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, SmartPetscObj< Vec > diff_x, double time, double delta_t, double eps, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Calculate directional directive using finite difference.

◆ directionalCentralFiniteDifference()

SmartPetscObj< Vec > MoFEM::OperatorsTester::directionalCentralFiniteDifference ( SmartPetscObj< DM >  dm,
std::string  fe_name,
boost::shared_ptr< FEMethod pipeline,
SmartPetscObj< Vec >  x,
SmartPetscObj< Vec >  delta_x,
SmartPetscObj< Vec >  delta2_x,
SmartPetscObj< Vec >  diff_x,
double  time,
double  delta_t,
double  eps,
CacheTupleWeakPtr  cache_ptr = CacheTupleSharedPtr() 
)

Calculate directional directive using finite difference.

Parameters
dm
fe_name
pipeline
x
delta_x
delta2_x
diff_x// direction of derivative
time
delta_t
eps
cache_ptr
Returns
SmartPetscObj<Vec>

Definition at line 153 of file OperatorsTester.cpp.

158 {
159
160 auto axpy = [&](auto v0, auto diff_v, auto eps) {
161 SmartPetscObj<Vec> v;
162 if (v0.use_count()) {
163 v = vectorDuplicate(v0);
164 CHK_THROW_MESSAGE(VecCopy(v0, v), "Cpy");
165 CHK_THROW_MESSAGE(VecAXPY(v, eps, diff_v), "Add");
166 CHK_THROW_MESSAGE(VecGhostUpdateBegin(v, INSERT_VALUES, SCATTER_FORWARD),
167 "ghost");
168 CHK_THROW_MESSAGE(VecGhostUpdateEnd(v, INSERT_VALUES, SCATTER_FORWARD),
169 "ghost");
170 }
171 return v;
172 };
173
174 auto fm1 = assembleVec(dm, fe_name, pipeline,
175 //
176 axpy(x, diff_x, -eps),
177 //
178 axpy(delta_x, diff_x, -eps),
179 //
180 axpy(delta2_x, diff_x, -eps),
181 //
182 time, delta_t, cache_ptr);
183
184 auto fp1 = assembleVec(dm, fe_name, pipeline,
185 //
186 axpy(x, diff_x, eps),
187 //
188 axpy(delta_x, diff_x, eps),
189 //
190 axpy(delta2_x, diff_x, eps),
191 //
192 time, delta_t, cache_ptr);
193
194 auto calculate_derivative = [eps](auto fm1, auto fp1) {
195 int size;
196 VecGetLocalSize(fp1, &size);
197 VecAXPY(fp1, -1, fm1);
198 double *array;
199 VecGetArray(fp1, &array);
200 for (int j = 0; j != size; ++j) {
201 array[j] /= 2 * eps;
202 }
203 VecRestoreArray(fp1, &array);
204 return fp1;
205 };
206
207 return calculate_derivative(fm1, fp1);
208}
FTensor::Index< 'j', 3 > j
SmartPetscObj< Vec > assembleVec(SmartPetscObj< DM > dm, std::string fe_name, boost::shared_ptr< FEMethod > pipeline, SmartPetscObj< Vec > x, SmartPetscObj< Vec > delta_x, SmartPetscObj< Vec > delta2_x, double time, double delta_t, CacheTupleWeakPtr cache_ptr)
Assemble the right hand side vector.

◆ query_interface()

MoFEMErrorCode MoFEM::OperatorsTester::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Implements MoFEM::UnknownInterface.

Definition at line 17 of file OperatorsTester.cpp.

18 {
19 *iface = const_cast<OperatorsTester *>(this);
20 return 0;
21}
OperatorsTester(const MoFEM::Core &core)

◆ setPipelineX()

std::pair< SmartPetscObj< Vec >, SmartPetscObj< Vec > > MoFEM::OperatorsTester::setPipelineX ( boost::shared_ptr< FEMethod pipeline,
SmartPetscObj< Vec >  x,
SmartPetscObj< Vec >  delta_x,
SmartPetscObj< Vec >  delta2_x,
double  delta_t 
)
private

Set vectors, x, x_t, and x_tt to finite element instance.

Finite element instance is a pipeline. x_t and x_tt are evaluated for given delta_x, delta2_x and delta_t.

Parameters
pipeline
x
delta_x
delta2_x
delta_t
Returns
std::pair<SmartPetscObj<Vec>, SmartPetscObj<Vec>>

Definition at line 233 of file OperatorsTester.cpp.

235 {
236
237 // Set dofs vector to finite element instance
238
239 pipeline->x = x;
240 pipeline->data_ctx |= PetscData::CTX_SET_X;
241
242 SmartPetscObj<Vec> x_t, x_tt;
243
244 // Set velocity dofs vector to finite element instance
245
246 if (delta_x.use_count()) {
247 x_t = vectorDuplicate(x);
248 VecCopy(delta_x, x_t);
249 VecScale(x_t, 1. / delta_t);
250 CHK_THROW_MESSAGE(VecGhostUpdateBegin(x_t, INSERT_VALUES, SCATTER_FORWARD),
251 "ghost");
252 CHK_THROW_MESSAGE(VecGhostUpdateEnd(x_t, INSERT_VALUES, SCATTER_FORWARD),
253 "ghost");
254 pipeline->x_t = x_t;
255 pipeline->data_ctx |= PetscData::CTX_SET_X_T;
256 } else {
257 pipeline->x_t = PETSC_NULL;
258 }
259
260 // Set acceleration dofs vector to finite element instance
261
262 if (delta2_x.use_count()) {
263 x_tt = vectorDuplicate(x);
264 VecCopy(delta2_x, x_tt);
265 VecScale(x_tt, 1. / pow(delta_t, 2));
266 CHK_THROW_MESSAGE(VecGhostUpdateBegin(x_tt, INSERT_VALUES, SCATTER_FORWARD),
267 "ghost");
268 CHK_THROW_MESSAGE(VecGhostUpdateEnd(x_tt, INSERT_VALUES, SCATTER_FORWARD),
269 "ghost");
270 pipeline->x_tt = x_tt;
271 pipeline->data_ctx |= PetscData::CTX_SET_X_TT;
272 } else {
273 pipeline->x_tt = PETSC_NULL;
274 }
275
276 return std::make_pair(x_t, x_tt);
277}

◆ setRandomFields()

SmartPetscObj< Vec > MoFEM::OperatorsTester::setRandomFields ( SmartPetscObj< DM >  dm,
std::vector< RandomFieldData random_fields,
boost::shared_ptr< Range ents = nullptr 
)

Generate random fileds.

Example: generate random vector for DM (problem) from simple interface, where FIELD random values of DOFs are in range from -1 to 1, and FIELD2 random values are in range from 0 to 1.

auto x = opt->setRandomFields(simple->getDM(),
{{"FIELD1", {-1, 1}}, {"FIELD2", {0,1}}});
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69

TODO: Set random field to specific entities, and potentially order for testing proposes to dissect error in tangent matrix.

Parameters
dm
random_fieldslook at definition RandomFieldData
Returns
SmartPetscObj<Vec> smart vector
Examples
tensor_divergence_operator.cpp.

Definition at line 24 of file OperatorsTester.cpp.

26 {
27
28 MoFEM::Interface &m_field = cOre;
29
30 auto prb_ptr = getProblemPtr(dm);
31
32 auto get_random_number = [](auto &range) {
33 const auto r = static_cast<double>(std::rand()) / RAND_MAX;
34 return range[0] + r * (range[1] - range[0]);
35 };
36
37 auto v = createDMVector(dm);
38 double *array;
39 CHKERR VecGetArray(v, &array);
40
41 MOFEM_LOG_CHANNEL("WORLD");
42
43 for (auto &rf : random_fields) {
44 const auto field_id = m_field.get_field_bit_number(rf.first);
45 MOFEM_TAG_AND_LOG("WORLD", Sev::noisy, "OperatorsTester")
46 << "Set random field " << rf.first << " " << static_cast<int>(field_id);
47
48 if (!ents_ptr) {
49 auto lo = prb_ptr->getNumeredRowDofsPtr()->lower_bound(
51 const auto hi = prb_ptr->getNumeredRowDofsPtr()->upper_bound(
53 for (; lo != hi; ++lo) {
54 const auto idx = (*lo)->getPetscLocalDofIdx();
55 array[idx] = get_random_number(rf.second);
56 }
57 } else {
58 for (auto pit = ents_ptr->const_pair_begin();
59 pit != ents_ptr->const_pair_end(); ++pit) {
60 auto lo =
61 prb_ptr->getNumeredRowDofsPtr()->get<Unique_mi_tag>().lower_bound(
62 DofEntity::getLoFieldEntityUId(field_id, pit->first));
63 auto hi =
64 prb_ptr->getNumeredRowDofsPtr()->get<Unique_mi_tag>().upper_bound(
65 DofEntity::getHiFieldEntityUId(field_id, pit->second));
66 for (; lo != hi; ++lo) {
67 const auto idx = (*lo)->getPetscLocalDofIdx();
68 array[idx] = get_random_number(rf.second);
69 }
70 }
71 }
72 }
73
74 CHKERR VecRestoreArray(v, &array);
75
76 CHKERR VecGhostUpdateBegin(v, INSERT_VALUES, SCATTER_FORWARD);
77 CHKERR VecGhostUpdateEnd(v, INSERT_VALUES, SCATTER_FORWARD);
78
79 return v;
80};
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
auto getProblemPtr(DM dm)
get problem pointer from DM
Definition: DMMoFEM.hpp:976
int r
Definition: sdf.py:8
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
Deprecated interface functions.
static UId getLoFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
static UId getHiFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)

Member Data Documentation

◆ cOre

MoFEM::Core& MoFEM::OperatorsTester::cOre
private

Definition at line 145 of file OperatorsTester.hpp.


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