v0.15.0
Loading...
Searching...
No Matches
ObjectiveFunctionDataImpl Struct Reference
Inheritance diagram for ObjectiveFunctionDataImpl:
[legend]
Collaboration diagram for ObjectiveFunctionDataImpl:
[legend]

Public Member Functions

 ObjectiveFunctionDataImpl ()=default
 
virtual ~ObjectiveFunctionDataImpl ()=default
 
MoFEMErrorCode initPython (const std::string py_file)
 
MoFEMErrorCode evalObjectiveFunction (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< VectorDouble > o_ptr)
 
MoFEMErrorCode evalObjectiveGradientStress (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< MatrixDouble > o_ptr)
 
MoFEMErrorCode evalObjectiveGradientStrain (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< MatrixDouble > o_ptr)
 
MoFEMErrorCode evalObjectiveGradientU (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< MatrixDouble > o_ptr)
 
MoFEMErrorCode numberOfModes (int block_id, int &modes)
 
MoFEMErrorCode blockModes (int block_id, MatrixDouble &coords, std::array< double, 3 > &centroid, std::array< double, 6 > &bbodx, MatrixDouble &o_ptr)
 
- Public Member Functions inherited from ObjectiveFunctionData
virtual ~ObjectiveFunctionData ()=default
 

Private Member Functions

MoFEMErrorCode objectiveFunctionImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 
MoFEMErrorCode objectiveGradientStressImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 
MoFEMErrorCode objectiveGradientStrainImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 
MoFEMErrorCode objectiveGradientUImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 
MoFEMErrorCode blockModesImpl (int block_id, np::ndarray coords, np::ndarray centroid, np::ndarray bbodx, np::ndarray &o_ptr)
 
np::ndarray convertToNumPy (std::vector< double > &data, int rows, int nb_gauss_pts)
 Converts a std::vector<double> to a NumPy ndarray.
 
np::ndarray convertToNumPy (double *ptr, int s)
 
MatrixDouble copyToFull (MatrixDouble &s)
 
void copyToSymmetric (double *ptr, MatrixDouble &s)
 

Private Attributes

bp::object mainNamespace
 
bp::object objectiveFunction
 
bp::object objectiveGradientStress
 
bp::object objectiveGradientStrain
 
bp::object objectiveGradientU
 
bp::object objectNumberOfModes
 
bp::object objectBlockModes
 

Detailed Description

Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2087 of file adjoint.cpp.

Constructor & Destructor Documentation

◆ ObjectiveFunctionDataImpl()

ObjectiveFunctionDataImpl::ObjectiveFunctionDataImpl ( )
default

◆ ~ObjectiveFunctionDataImpl()

virtual ObjectiveFunctionDataImpl::~ObjectiveFunctionDataImpl ( )
virtualdefault

Member Function Documentation

◆ blockModes()

MoFEMErrorCode ObjectiveFunctionDataImpl::blockModes ( int block_id,
MatrixDouble & coords,
std::array< double, 3 > & centroid,
std::array< double, 6 > & bbodx,
MatrixDouble & o_ptr )
virtual

Implements ObjectiveFunctionData.

Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2391 of file adjoint.cpp.

2393 {
2395 try {
2396
2397 int nb_modes = bp::extract<int>(objectNumberOfModes(block_id));
2398
2399 auto np_coords =
2400 convertToNumPy(coords.data(), coords.size1(), coords.size2());
2401
2402 auto np_centroid = convertToNumPy(centroid.data(), 3);
2403 auto np_bbodx = convertToNumPy(bbodx.data(), 6);
2404
2405 np::ndarray np_output =
2406 np::empty(bp::make_tuple(coords.size1(), nb_modes, 3),
2407 np::dtype::get_builtin<double>());
2408 CHKERR blockModesImpl(block_id, np_coords, np_centroid, np_bbodx,
2409 np_output);
2410
2411 o_ptr.resize(nb_modes, coords.size1() * coords.size2(), false);
2412 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
2413 std::copy(val_ptr, val_ptr + coords.size1() * coords.size2() * nb_modes,
2414 o_ptr.data().begin());
2415
2416 } catch (bp::error_already_set const &) {
2417 PyErr_Print();
2419 }
2421}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MoFEMErrorCode blockModesImpl(int block_id, np::ndarray coords, np::ndarray centroid, np::ndarray bbodx, np::ndarray &o_ptr)
Definition adjoint.cpp:2524
np::ndarray convertToNumPy(std::vector< double > &data, int rows, int nb_gauss_pts)
Converts a std::vector<double> to a NumPy ndarray.
Definition adjoint.cpp:2563

◆ blockModesImpl()

MoFEMErrorCode ObjectiveFunctionDataImpl::blockModesImpl ( int block_id,
np::ndarray coords,
np::ndarray centroid,
np::ndarray bbodx,
np::ndarray & o_ptr )
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2524 of file adjoint.cpp.

2528 {
2530 try {
2531 // call python function
2532 o = bp::extract<np::ndarray>(
2533 objectBlockModes(block_id, coords, centroid, bbodx));
2534 } catch (bp::error_already_set const &) {
2535 // print all other errors to stderr
2536 PyErr_Print();
2538 }
2540}

◆ convertToNumPy() [1/2]

np::ndarray ObjectiveFunctionDataImpl::convertToNumPy ( double * ptr,
int s )
inlineprivate

Definition at line 2571 of file adjoint.cpp.

2572 {
2573 auto dtype = np::dtype::get_builtin<double>();
2574 auto size = bp::make_tuple(s);
2575 auto stride = bp::make_tuple(sizeof(double));
2576 return (np::from_data(ptr, dtype, size, stride, bp::object()));
2577}

◆ convertToNumPy() [2/2]

np::ndarray ObjectiveFunctionDataImpl::convertToNumPy ( std::vector< double > & data,
int rows,
int nb_gauss_pts )
inlineprivate

Converts a std::vector<double> to a NumPy ndarray.

This function wraps the given vector data into a NumPy array with the specified number of rows and Gauss points. The resulting ndarray shares memory with the input vector, so changes to one will affect the other.

Parameters
dataReference to the vector containing double values to be converted.
rowsNumber of rows in the resulting NumPy array.
nb_gauss_ptsNumber of Gauss points (columns) in the resulting NumPy array.
Returns
np::ndarray NumPy array view of the input data.
Note
  • size specifies the shape of the resulting ndarray as a tuple (rows, nb_gauss_pts).
  • stride specifies the step size in bytes to move to the next element in memory. Here, it is set to sizeof(double), indicating contiguous storage for each element.
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2563 of file adjoint.cpp.

2564 {
2565 auto dtype = np::dtype::get_builtin<double>();
2566 auto size = bp::make_tuple(rows, nb_gauss_pts);
2567 auto stride = bp::make_tuple(nb_gauss_pts * sizeof(double), sizeof(double));
2568 return (np::from_data(data.data(), dtype, size, stride, bp::object()));
2569}

◆ copyToFull()

MatrixDouble ObjectiveFunctionDataImpl::copyToFull ( MatrixDouble & s)
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2216 of file adjoint.cpp.

2216 {
2217 MatrixDouble f(9, s.size2());
2218 f.clear();
2221 auto t_f = getFTensor2FromMat<3, 3>(f);
2223 for (int ii = 0; ii != s.size2(); ++ii) {
2224 t_f(i, j) = t_s(i, j);
2225 ++t_f;
2226 ++t_s;
2227 }
2228 return f;
2229};
#define FTENSOR_INDEX(DIM, I)
constexpr int SPACE_DIM
[Define dimension]
Definition adjoint.cpp:22
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, Tensor_Dim1, Tensor_Dim2 > getFTensor2FromMat(MatrixDouble &data)
Get tensor rank 2 (matrix) form data matrix.
static auto getFTensor2SymmetricFromMat(ublas::matrix< T, L, A > &data)
Get symmetric tensor rank 2 (matrix) form data matrix.

◆ copyToSymmetric()

void ObjectiveFunctionDataImpl::copyToSymmetric ( double * ptr,
MatrixDouble & s )
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2231 of file adjoint.cpp.

2231 {
2235
2236 ptr + 0 * s.size2(), ptr + 1 * s.size2(),
2237 ptr + 2 * s.size2(), ptr + 3 * s.size2(),
2238 ptr + 4 * s.size2(), ptr + 5 * s.size2(),
2239 ptr + 6 * s.size2(), ptr + 7 * s.size2(),
2240 ptr + 8 * s.size2()
2241
2242 };
2244 for (int ii = 0; ii != s.size2(); ++ii) {
2245 t_s(i, j) = (t_f(i, j) || t_f(j, i)) / 2.0;
2246 ++t_f;
2247 ++t_s;
2248 }
2249}

◆ evalObjectiveFunction()

MoFEMErrorCode ObjectiveFunctionDataImpl::evalObjectiveFunction ( MatrixDouble & coords,
boost::shared_ptr< MatrixDouble > u_ptr,
boost::shared_ptr< MatrixDouble > stress_ptr,
boost::shared_ptr< MatrixDouble > strain_ptr,
boost::shared_ptr< VectorDouble > o_ptr )
virtual

Implements ObjectiveFunctionData.

Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2251 of file adjoint.cpp.

2255 {
2257 try {
2258
2259 auto np_coords =
2260 convertToNumPy(coords.data(), coords.size1(), coords.size2());
2261 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
2262
2263 auto full_stress = copyToFull(*(stress_ptr));
2264 auto full_strain = copyToFull(*(strain_ptr));
2265
2266 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
2267 full_stress.size2());
2268 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
2269 full_strain.size2());
2270 np::ndarray np_output = np::empty(bp::make_tuple(strain_ptr->size2()),
2271 np::dtype::get_builtin<double>());
2272 CHKERR objectiveFunctionImpl(np_coords, np_u, np_stress, np_strain,
2273 np_output);
2274 o_ptr->resize(stress_ptr->size2(), false);
2275 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
2276 std::copy(val_ptr, val_ptr + strain_ptr->size2(), o_ptr->data().begin());
2277
2278 } catch (bp::error_already_set const &) {
2279 PyErr_Print();
2281 }
2283}
MatrixDouble copyToFull(MatrixDouble &s)
Definition adjoint.cpp:2216
MoFEMErrorCode objectiveFunctionImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Definition adjoint.cpp:2423

◆ evalObjectiveGradientStrain()

MoFEMErrorCode ObjectiveFunctionDataImpl::evalObjectiveGradientStrain ( MatrixDouble & coords,
boost::shared_ptr< MatrixDouble > u_ptr,
boost::shared_ptr< MatrixDouble > stress_ptr,
boost::shared_ptr< MatrixDouble > strain_ptr,
boost::shared_ptr< MatrixDouble > o_ptr )
virtual

Implements ObjectiveFunctionData.

Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2320 of file adjoint.cpp.

2324 {
2326 try {
2327
2328 auto np_coords =
2329 convertToNumPy(coords.data(), coords.size1(), coords.size2());
2330 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
2331
2332 auto full_stress = copyToFull(*(stress_ptr));
2333 auto full_strain = copyToFull(*(strain_ptr));
2334 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
2335 full_stress.size2());
2336 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
2337 full_strain.size2());
2338
2339 np::ndarray np_output =
2340 np::empty(bp::make_tuple(full_strain.size1(), full_strain.size2()),
2341 np::dtype::get_builtin<double>());
2342 CHKERR objectiveGradientStrainImpl(np_coords, np_u, np_stress, np_strain,
2343 np_output);
2344 o_ptr->resize(strain_ptr->size1(), strain_ptr->size2(), false);
2345 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
2346 copyToSymmetric(val_ptr, *(o_ptr));
2347
2348 } catch (bp::error_already_set const &) {
2349 PyErr_Print();
2351 }
2353}
MoFEMErrorCode objectiveGradientStrainImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Definition adjoint.cpp:2466
void copyToSymmetric(double *ptr, MatrixDouble &s)
Definition adjoint.cpp:2231

◆ evalObjectiveGradientStress()

MoFEMErrorCode ObjectiveFunctionDataImpl::evalObjectiveGradientStress ( MatrixDouble & coords,
boost::shared_ptr< MatrixDouble > u_ptr,
boost::shared_ptr< MatrixDouble > stress_ptr,
boost::shared_ptr< MatrixDouble > strain_ptr,
boost::shared_ptr< MatrixDouble > o_ptr )
virtual

Implements ObjectiveFunctionData.

Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2285 of file adjoint.cpp.

2289 {
2291 try {
2292
2293 auto np_coords =
2294 convertToNumPy(coords.data(), coords.size1(), coords.size2());
2295 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
2296
2297 auto full_stress = copyToFull(*(stress_ptr));
2298 auto full_strain = copyToFull(*(strain_ptr));
2299
2300 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
2301 full_stress.size2());
2302 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
2303 full_strain.size2());
2304 np::ndarray np_output =
2305 np::empty(bp::make_tuple(full_strain.size1(), full_strain.size2()),
2306 np::dtype::get_builtin<double>());
2307 CHKERR objectiveGradientStressImpl(np_coords, np_u, np_stress, np_strain,
2308 np_output);
2309 o_ptr->resize(stress_ptr->size1(), stress_ptr->size2(), false);
2310 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
2311 copyToSymmetric(val_ptr, *(o_ptr));
2312
2313 } catch (bp::error_already_set const &) {
2314 PyErr_Print();
2316 }
2318}
MoFEMErrorCode objectiveGradientStressImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Definition adjoint.cpp:2444

◆ evalObjectiveGradientU()

MoFEMErrorCode ObjectiveFunctionDataImpl::evalObjectiveGradientU ( MatrixDouble & coords,
boost::shared_ptr< MatrixDouble > u_ptr,
boost::shared_ptr< MatrixDouble > stress_ptr,
boost::shared_ptr< MatrixDouble > strain_ptr,
boost::shared_ptr< MatrixDouble > o_ptr )
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2355 of file adjoint.cpp.

2359 {
2361 try {
2362
2363 auto np_coords =
2364 convertToNumPy(coords.data(), coords.size1(), coords.size2());
2365 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
2366
2367 auto full_stress = copyToFull(*(stress_ptr));
2368 auto full_strain = copyToFull(*(strain_ptr));
2369 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
2370 full_stress.size2());
2371 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
2372 full_strain.size2());
2373
2374 np::ndarray np_output =
2375 np::empty(bp::make_tuple(u_ptr->size1(), u_ptr->size2()),
2376 np::dtype::get_builtin<double>());
2377 CHKERR objectiveGradientStrainImpl(np_coords, np_u, np_stress, np_strain,
2378 np_output);
2379 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
2380 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
2381 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
2382 o_ptr->data().begin());
2383
2384 } catch (bp::error_already_set const &) {
2385 PyErr_Print();
2387 }
2389}

◆ initPython()

MoFEMErrorCode ObjectiveFunctionDataImpl::initPython ( const std::string py_file)
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2192 of file adjoint.cpp.

2192 {
2194 try {
2195
2196 // create main module
2197 auto main_module = bp::import("__main__");
2198 mainNamespace = main_module.attr("__dict__");
2199 bp::exec_file(py_file.c_str(), mainNamespace, mainNamespace);
2200 // create a reference to python function
2205 objectNumberOfModes = mainNamespace["number_of_modes"];
2206 objectBlockModes = mainNamespace["block_modes"];
2207
2208 } catch (bp::error_already_set const &) {
2209 // print all other errors to stderr
2210 PyErr_Print();
2212 }
2214}
bp::object objectiveGradientStress
Definition adjoint.cpp:2130
bp::object objectiveGradientStrain
Definition adjoint.cpp:2131

◆ numberOfModes()

MoFEMErrorCode ObjectiveFunctionDataImpl::numberOfModes ( int block_id,
int & modes )
virtual

Implements ObjectiveFunctionData.

Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2509 of file adjoint.cpp.

2510 {
2512 try {
2513
2514 modes = bp::extract<int>(objectNumberOfModes(block_id));
2515
2516 } catch (bp::error_already_set const &) {
2517 // print all other errors to stderr
2518 PyErr_Print();
2520 }
2522}

◆ objectiveFunctionImpl()

MoFEMErrorCode ObjectiveFunctionDataImpl::objectiveFunctionImpl ( np::ndarray coords,
np::ndarray u,
np::ndarray stress,
np::ndarray strain,
np::ndarray & o )
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2423 of file adjoint.cpp.

2429 {
2431 try {
2432
2433 // call python function
2434 o = bp::extract<np::ndarray>(objectiveFunction(coords, u, stress, strain));
2435
2436 } catch (bp::error_already_set const &) {
2437 // print all other errors to stderr
2438 PyErr_Print();
2440 }
2442}

◆ objectiveGradientStrainImpl()

MoFEMErrorCode ObjectiveFunctionDataImpl::objectiveGradientStrainImpl ( np::ndarray coords,
np::ndarray u,
np::ndarray stress,
np::ndarray strain,
np::ndarray & o )
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2466 of file adjoint.cpp.

2472 {
2474 try {
2475
2476 // call python function
2477 o = bp::extract<np::ndarray>(
2478 objectiveGradientStrain(coords, u, stress, strain));
2479
2480 } catch (bp::error_already_set const &) {
2481 // print all other errors to stderr
2482 PyErr_Print();
2484 }
2486}

◆ objectiveGradientStressImpl()

MoFEMErrorCode ObjectiveFunctionDataImpl::objectiveGradientStressImpl ( np::ndarray coords,
np::ndarray u,
np::ndarray stress,
np::ndarray strain,
np::ndarray & o )
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2444 of file adjoint.cpp.

2450 {
2452 try {
2453
2454 // call python function
2455 o = bp::extract<np::ndarray>(
2456 objectiveGradientStress(coords, u, stress, strain));
2457
2458 } catch (bp::error_already_set const &) {
2459 // print all other errors to stderr
2460 PyErr_Print();
2462 }
2464}

◆ objectiveGradientUImpl()

MoFEMErrorCode ObjectiveFunctionDataImpl::objectiveGradientUImpl ( np::ndarray coords,
np::ndarray u,
np::ndarray stress,
np::ndarray strain,
np::ndarray & o )
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2488 of file adjoint.cpp.

2494 {
2496 try {
2497
2498 // call python function
2499 o = bp::extract<np::ndarray>(objectiveGradientU(coords, u, stress, strain));
2500
2501 } catch (bp::error_already_set const &) {
2502 // print all other errors to stderr
2503 PyErr_Print();
2505 }
2507}

Member Data Documentation

◆ mainNamespace

bp::object ObjectiveFunctionDataImpl::mainNamespace
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2128 of file adjoint.cpp.

◆ objectBlockModes

bp::object ObjectiveFunctionDataImpl::objectBlockModes
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2134 of file adjoint.cpp.

◆ objectiveFunction

bp::object ObjectiveFunctionDataImpl::objectiveFunction
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2129 of file adjoint.cpp.

◆ objectiveGradientStrain

bp::object ObjectiveFunctionDataImpl::objectiveGradientStrain
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2131 of file adjoint.cpp.

◆ objectiveGradientStress

bp::object ObjectiveFunctionDataImpl::objectiveGradientStress
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2130 of file adjoint.cpp.

◆ objectiveGradientU

bp::object ObjectiveFunctionDataImpl::objectiveGradientU
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2132 of file adjoint.cpp.

◆ objectNumberOfModes

bp::object ObjectiveFunctionDataImpl::objectNumberOfModes
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 2133 of file adjoint.cpp.


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