OGS
ThermoHydroMechanicsProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <cassert>
7
12#include "ProcessLib/Process.h"
17
18namespace ProcessLib
19{
21{
22template <int DisplacementDim>
24 std::string name, MeshLib::Mesh& mesh,
25 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
26 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
27 unsigned const integration_order,
28 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
29 process_variables,
31 SecondaryVariableCollection&& secondary_variables,
32 bool const use_monolithic_scheme, bool const is_linear)
33 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
34 integration_order, std::move(process_variables),
35 std::move(secondary_variables), use_monolithic_scheme),
37 *_jacobian_assembler, is_linear, use_monolithic_scheme},
38 _process_data(std::move(process_data))
39
40{
41 // For numerical Jacobian
42 if (this->_jacobian_assembler->isPerturbationEnabled())
43 {
45 "Numerical Jacobian is not supported for the "
46 "ThermoHydroMechanicsProcess yet.");
47 }
48
49 _integration_point_writer.emplace_back(
50 std::make_unique<MeshLib::IntegrationPointWriter>(
51 "sigma_ip",
52 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
53 integration_order, local_assemblers_,
55
56 _integration_point_writer.emplace_back(
57 std::make_unique<MeshLib::IntegrationPointWriter>(
58 "sigma_ice_ip",
59 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
60 integration_order, local_assemblers_,
62
63 _integration_point_writer.emplace_back(
64 std::make_unique<MeshLib::IntegrationPointWriter>(
65 "ice_volume_fraction_ip", 1, integration_order, local_assemblers_,
67 _integration_point_writer.emplace_back(
68 std::make_unique<MeshLib::IntegrationPointWriter>(
69 "epsilon_m_ip",
70 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
71 integration_order, local_assemblers_,
73
74 _integration_point_writer.emplace_back(
75 std::make_unique<MeshLib::IntegrationPointWriter>(
76 "epsilon_ip",
77 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
78 integration_order, local_assemblers_,
80}
81
82template <int DisplacementDim>
88
89template <int DisplacementDim>
92 const int process_id) const
93{
94 // For the monolithic scheme or the M process (deformation) in the staggered
95 // scheme.
96 if (_use_monolithic_scheme || process_id == 2)
97 {
98 auto const& l = *_local_to_global_index_map;
99 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
100 &l.getGhostIndices(), &this->_sparsity_pattern};
101 }
102
103 // For staggered scheme and T or H process (pressure).
105 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
106 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
107}
108
109template <int DisplacementDim>
111{
112 // Create single component dof in every of the mesh's nodes.
114 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
115 // Create single component dof in the mesh's base nodes.
116 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
118 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
119
120 // TODO move the two data members somewhere else.
121 // for extrapolation of secondary variables of stress or strain
122 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
125 std::make_unique<NumLib::LocalToGlobalIndexMap>(
126 std::move(all_mesh_subsets_single_component),
127 // by location order is needed for output
129
131 {
132 // For temperature, which is the first
133 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
135
136 // For pressure, which is the second
137 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
138
139 // For displacement.
140 const int monolithic_process_id = 0;
141 std::generate_n(std::back_inserter(all_mesh_subsets),
142 getProcessVariables(monolithic_process_id)[2]
143 .get()
144 .getNumberOfGlobalComponents(),
145 [&]() { return *_mesh_subset_all_nodes; });
146
147 std::vector<int> const vec_n_components{1, 1, DisplacementDim};
149 std::make_unique<NumLib::LocalToGlobalIndexMap>(
150 std::move(all_mesh_subsets), vec_n_components,
153 }
154 else
155 {
156 // For displacement equation.
157 const int process_id = 2;
158 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
159 std::generate_n(std::back_inserter(all_mesh_subsets),
160 getProcessVariables(process_id)[0]
161 .get()
162 .getNumberOfGlobalComponents(),
163 [&]() { return *_mesh_subset_all_nodes; });
164
165 std::vector<int> const vec_n_components{DisplacementDim};
167 std::make_unique<NumLib::LocalToGlobalIndexMap>(
168 std::move(all_mesh_subsets), vec_n_components,
170
171 // For pressure equation or temperature equation.
172 // Collect the mesh subsets with base nodes in a vector.
173 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
176 std::make_unique<NumLib::LocalToGlobalIndexMap>(
177 std::move(all_mesh_subsets_base_nodes),
178 // by location order is needed for output
180
183
186 }
187}
188
189template <int DisplacementDim>
191 NumLib::LocalToGlobalIndexMap const& dof_table,
192 MeshLib::Mesh const& mesh,
193 unsigned const integration_order)
194{
197 mesh.getElements(), dof_table, local_assemblers_,
198 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
200
201 auto add_secondary_variable = [&](std::string const& name,
202 int const num_components,
203 auto get_ip_values_function)
204 {
205 _secondary_variables.addSecondaryVariable(
206 name,
207 makeExtrapolator(num_components, getExtrapolator(),
209 std::move(get_ip_values_function)));
210 };
211
212 add_secondary_variable(
213 "sigma",
215 DisplacementDim>::RowsAtCompileTime,
217
218 add_secondary_variable(
219 "sigma_ice",
221 DisplacementDim>::RowsAtCompileTime,
223
224 add_secondary_variable(
225 "epsilon_0",
227 DisplacementDim>::RowsAtCompileTime,
229
230 add_secondary_variable(
231 "epsilon_m",
233 DisplacementDim>::RowsAtCompileTime,
235
236 add_secondary_variable(
237 "epsilon",
239 DisplacementDim>::RowsAtCompileTime,
241
242 add_secondary_variable(
243 "ice_volume_fraction", 1,
245
246 add_secondary_variable(
247 "velocity", mesh.getDimension(),
249
250 add_secondary_variable(
251 "fluid_density", 1,
253
254 add_secondary_variable(
255 "viscosity", 1,
257
259 const_cast<MeshLib::Mesh&>(mesh), "ice_volume_fraction_avg",
261
262 _process_data.element_fluid_density =
264 const_cast<MeshLib::Mesh&>(mesh), "fluid_density_avg",
266
268 const_cast<MeshLib::Mesh&>(mesh), "viscosity_avg",
270
272 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
275 DisplacementDim>::RowsAtCompileTime);
276
277 _process_data.element_ice_stresses =
279 const_cast<MeshLib::Mesh&>(mesh), "sigma_ice_avg",
282 DisplacementDim>::RowsAtCompileTime);
283
284 _process_data.pressure_interpolated =
286 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
288
289 _process_data.temperature_interpolated =
291 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
293
294 //
295 // enable output of internal variables defined by material models
296 //
299 add_secondary_variable);
300
303 _process_data.solid_materials, local_assemblers_,
304 _integration_point_writer, integration_order);
305
308
309 // Initialize local assemblers after all variables have been set.
313}
314
315template <int DisplacementDim>
317 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
318{
320 {
321 const int process_id_of_thermohydromechancs = 0;
323 *_local_to_global_index_map, process_id_of_thermohydromechancs,
324 media);
325 return;
326 }
327
328 // Staggered scheme:
329 // for the equations of heat transport
330 const int thermal_process_id = 0;
332 *_local_to_global_index_map_with_base_nodes, thermal_process_id, media);
333
334 // for the equations of mass balance
335 const int hydraulic_process_id = 1;
337 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id,
338 media);
339
340 // for the equations of deformation.
341 const int mechanical_process_id = 2;
343 *_local_to_global_index_map, mechanical_process_id, media);
344}
345
346template <int DisplacementDim>
348 setInitialConditionsConcreteProcess(std::vector<GlobalVector*>& x,
349 double const t,
350 int const process_id)
351{
352 DBUG("SetInitialConditions ThermoHydroMechanicsProcess.");
353
356 local_assemblers_, getDOFTables(x.size()), x, t, process_id);
357}
358
359template <int DisplacementDim>
361 const double t, double const dt, std::vector<GlobalVector*> const& x,
362 std::vector<GlobalVector*> const& x_prev, int const process_id,
364{
365 DBUG("Assemble the equations for ThermoHydroMechanics");
366
368 t, dt, x, x_prev, process_id, M, K, b);
369}
370
371template <int DisplacementDim>
374 const double t, double const dt, std::vector<GlobalVector*> const& x,
375 std::vector<GlobalVector*> const& x_prev, int const process_id,
376 GlobalVector& b, GlobalMatrix& Jac)
377{
378 // For the monolithic scheme
380 {
381 DBUG(
382 "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
383 "scheme.");
384 }
385 else
386 {
387 // For the staggered scheme
388 if (process_id == 0)
389 {
390 DBUG(
391 "Assemble the Jacobian equations of heat transport process in "
392 "ThermoHydroMechanics for the staggered scheme.");
393 }
394 else if (process_id == 1)
395 {
396 DBUG(
397 "Assemble the Jacobian equations of liquid fluid process in "
398 "ThermoHydroMechanics for the staggered scheme.");
399 }
400 else
401 {
402 DBUG(
403 "Assemble the Jacobian equations of mechanical process in "
404 "ThermoHydroMechanics for the staggered scheme.");
405 }
406 }
407
409 assembleWithJacobian(t, dt, x, x_prev, process_id, b, Jac);
410}
411
412template <int DisplacementDim>
414 std::vector<GlobalVector*> const& x, double const t, double const dt,
415 const int process_id)
416{
417 DBUG("PreTimestep ThermoHydroMechanicsProcess.");
418
419 if (hasMechanicalProcess(process_id))
420 {
424 dt);
425
427 updateActiveElements();
428 }
429}
430
431template <int DisplacementDim>
433 std::vector<GlobalVector*> const& x,
434 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
435 const int process_id)
436{
437 if (process_id != 0)
438 {
439 return;
440 }
441
442 DBUG("PostTimestep ThermoHydroMechanicsProcess.");
443
447 x_prev, t, dt, process_id);
448}
449
450template <int DisplacementDim>
451std::vector<std::vector<std::string>>
453 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
454{
455 INFO("ThermoHydroMechanicsProcess process initializeSubmeshOutput().");
456 std::vector<std::vector<std::string>> per_process_residuum_names;
457
458 std::string const flow_rate_name =
459 _process_data.is_volume_balance_equation_type ? "VolumetricFlowRate"
460 : "MassFlowRate";
461
462 if (_process_variables.size() == 1) // monolithic
463 {
464 per_process_residuum_names = {
465 {"HeatFlowRate", flow_rate_name, "NodalForces"}};
466 }
467 else // staggered
468 {
469 per_process_residuum_names = {
470 {"HeatFlowRate"}, {flow_rate_name}, {"NodalForces"}};
471 }
472
474 initializeAssemblyOnSubmeshes(meshes, per_process_residuum_names);
475
476 return per_process_residuum_names;
477}
478
479template <int DisplacementDim>
481 computeSecondaryVariableConcrete(double const t, double const dt,
482 std::vector<GlobalVector*> const& x,
483 GlobalVector const& x_prev,
484 const int process_id)
485{
486 if (process_id != 0)
487 {
488 return;
489 }
490
491 DBUG("Compute the secondary variables for ThermoHydroMechanicsProcess.");
492
496 x, x_prev, process_id);
497}
498
499template <int DisplacementDim>
500std::tuple<NumLib::LocalToGlobalIndexMap*, bool> ThermoHydroMechanicsProcess<
501 DisplacementDim>::getDOFTableForExtrapolatorData() const
502{
503 const bool manage_storage = false;
504 return std::make_tuple(_local_to_global_index_map_single_component.get(),
505 manage_storage);
506}
507
508template <int DisplacementDim>
511 const int process_id) const
512{
513 if (hasMechanicalProcess(process_id))
514 {
516 }
517
518 // For the equation of pressure
520}
521
522template class ThermoHydroMechanicsProcess<2>;
523template class ThermoHydroMechanicsProcess<3>;
524
525} // namespace ThermoHydroMechanics
526} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
bool isAxiallySymmetric() const
Definition Mesh.h:130
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:101
unsigned getDimension() const
Definition Mesh.h:80
Properties & getProperties()
Definition Mesh.h:127
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void setInitialConditions(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
std::string const name
Definition Process.h:361
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:389
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:83
void assemble(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
Definition Process.cpp:258
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:37
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:399
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:149
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:391
const bool _use_monolithic_scheme
Definition Process.h:378
Handles configuration of several secondary variables from the project file.
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
THERMOHYDROMECHANICS_EXPORT NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
THERMOHYDROMECHANICS_EXPORT ThermoHydroMechanicsProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, ThermoHydroMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, bool const is_linear)
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) override
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
THERMOHYDROMECHANICS_EXPORT MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > local_assemblers_
ThermoHydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.
void solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, std::vector< std::unique_ptr< LocalAssemblerInterface > > const &local_assemblers, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writer, int const integration_order)
void solidMaterialInternalToSecondaryVariables(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
void createLocalAssemblersHM(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtFluidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon0(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaIce(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtIceVolume(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtViscosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonM(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0