Loading [MathJax]/extensions/tex2jax.js
OGS
ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >

Definition at line 26 of file SmallDeformationNonlocalProcess.h.

#include <SmallDeformationNonlocalProcess.h>

Inheritance diagram for ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >:
[legend]

Public Member Functions

 SmallDeformationNonlocalProcess (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, SmallDeformationNonlocalProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables)
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 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)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep. More...
 
void postTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep. More...
 
void postNonLinearSolver (GlobalVector const &x, GlobalVector const &xdot, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
 compute secondary variables for the coupled equations or for output. More...
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize ()
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void setCoupledSolutionsForStaggeredScheme (CoupledSolutionsForStaggeredScheme *const coupled_solutions)
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
bool isMonolithicSchemeUsed () const
 
virtual void setCoupledTermForTheStaggeredSchemeToLocalAssemblers (int const)
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
 
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
 
MeshLib::MeshgetMesh () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< IntegrationPointWriter > > const * getIntegrationPointWriter (MeshLib::Mesh const &mesh) const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 

Private Types

using LocalAssemblerInterface = SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >
 

Private Member Functions

void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize(). More...
 
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
 
void preAssembleConcreteProcess (const double t, double const dt, GlobalVector const &x) override
 
void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, int const process_id) override
 
NumLib::IterationResult postIterationConcreteProcess (GlobalVector const &x) override
 

Private Attributes

SmallDeformationNonlocalProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 

Additional Inherited Members

- Public Types inherited from ProcessLib::Process
using NonlinearSolver = NumLib::NonlinearSolverBase
 
using TimeDiscretization = NumLib::TimeDiscretization
 
- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Protected Member Functions inherited from ProcessLib::Process
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
 
virtual void constructDofTable ()
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_prosess_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
CoupledSolutionsForStaggeredScheme_coupled_solutions
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Member Typedef Documentation

◆ LocalAssemblerInterface

Definition at line 79 of file SmallDeformationNonlocalProcess.h.

Constructor & Destructor Documentation

◆ SmallDeformationNonlocalProcess()

template<int DisplacementDim>
ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::SmallDeformationNonlocalProcess ( 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,
SmallDeformationNonlocalProcessData< DisplacementDim > &&  process_data,
SecondaryVariableCollection &&  secondary_variables 
)

Definition at line 26 of file SmallDeformationNonlocalProcess.cpp.

39  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
40  integration_order, std::move(process_variables),
41  std::move(secondary_variables)),
42  _process_data(std::move(process_data))
43 {
44  _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
45  mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
46 
47  _integration_point_writer.emplace_back(
48  std::make_unique<IntegrationPointWriter>(
49  "sigma_ip",
50  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
51  integration_order, _local_assemblers,
53 
54  _integration_point_writer.emplace_back(
55  std::make_unique<IntegrationPointWriter>(
56  "kappa_d_ip", 1 /*n_components*/, integration_order,
58 }
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::string const name
Definition: Process.h:323
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:350
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:22
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
SmallDeformationNonlocalProcessData< DisplacementDim > _process_data

References ProcessLib::Process::_integration_point_writer, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::_local_assemblers, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::_nodal_forces, MeshLib::Mesh::getDimension(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >::getKappaD(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >::getSigma(), and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::assembleConcreteProcess ( const double  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
std::vector< GlobalVector * > const &  xdot,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 207 of file SmallDeformationNonlocalProcess.cpp.

211 {
212  DBUG("Assemble SmallDeformationNonlocalProcess.");
213 
214  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
215  dof_table = {std::ref(*_local_to_global_index_map)};
216 
217  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
218 
219  // Call global assembler for each local assembly item.
222  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
223  b);
224 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( const double  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
std::vector< GlobalVector * > const &  xdot,
const double  dxdot_dx,
const double  dx_dx,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 245 of file SmallDeformationNonlocalProcess.cpp.

251 {
252  DBUG("AssembleWithJacobian SmallDeformationNonlocalProcess.");
253 
254  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
255  dof_table = {std::ref(*_local_to_global_index_map)};
256 
257  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
258 
259  // Call global assembler for each local assembly item.
262  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
263  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
264 
266  std::transform(_nodal_forces->begin(), _nodal_forces->end(),
267  _nodal_forces->begin(), [](double val) { return -val; });
268 }
void copyValues(std::vector< double > &u) const
Copy vector values.
Definition: EigenVector.h:94
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), MathLib::EigenVector::copyValues(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const &  dof_table,
MeshLib::Mesh const &  mesh,
unsigned const  integration_order 
)
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 67 of file SmallDeformationNonlocalProcess.cpp.

71 {
72  // Reusing local assembler creation code.
74  DisplacementDim, SmallDeformationNonlocalLocalAssembler>(
75  mesh.getElements(), dof_table, _local_assemblers,
76  mesh.isAxiallySymmetric(), integration_order, _process_data);
77 
78  // TODO move the two data members somewhere else.
79  // for extrapolation of secondary variables
80  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
83  std::make_unique<NumLib::LocalToGlobalIndexMap>(
84  std::move(all_mesh_subsets_single_component),
85  // by location order is needed for output
87 
89  "sigma",
91  DisplacementDim>::RowsAtCompileTime,
94 
96  "epsilon",
98  DisplacementDim>::RowsAtCompileTime,
101 
103  "eps_p_V",
107  "eps_p_D_xx",
110 
112  "damage",
115 
119 
120  // Set initial conditions for integration point data.
121  for (auto const& ip_writer : _integration_point_writer)
122  {
123  auto const& name = ip_writer->name();
124  // First check the field data, which is used for restart.
125  if (mesh.getProperties().existsPropertyVector<double>(name))
126  {
127  auto const& mesh_property =
128  *mesh.getProperties().template getPropertyVector<double>(name);
129 
130  // The mesh property must be defined on integration points.
131  if (mesh_property.getMeshItemType() !=
133  {
134  continue;
135  }
136 
137  auto const ip_meta_data = getIntegrationPointMetaData(mesh, name);
138 
139  // Check the number of components.
140  if (ip_meta_data.n_components !=
141  mesh_property.getNumberOfGlobalComponents())
142  {
143  OGS_FATAL(
144  "Different number of components in meta data ({:d}) than "
145  "in the integration point field data for '{:s}': {:d}.",
146  ip_meta_data.n_components, name,
147  mesh_property.getNumberOfGlobalComponents());
148  }
149 
150  // Now we have a properly named vtk's field data array and the
151  // corresponding meta data.
152  std::size_t position = 0;
153  for (auto& local_asm : _local_assemblers)
154  {
155  std::size_t const integration_points_read =
156  local_asm->setIPDataInitialConditions(
157  name, &mesh_property[position],
158  ip_meta_data.integration_order);
159  if (integration_points_read == 0)
160  {
161  OGS_FATAL(
162  "No integration points read in the integration point "
163  "initial conditions set function.");
164  }
165  position += integration_points_read * ip_meta_data.n_components;
166  }
167  }
168  else if (mesh.getProperties().existsPropertyVector<double>(name +
169  "_ic"))
170  { // Try to find cell data with '_ic' suffix
171  auto const& mesh_property =
172  *mesh.getProperties().template getPropertyVector<double>(name +
173  "_ic");
174  if (mesh_property.getMeshItemType() != MeshLib::MeshItemType::Cell)
175  {
176  continue;
177  }
178 
179  // Now we have a vtk's cell data array containing the initial
180  // conditions for the corresponding integration point writer.
181 
182  // For each assembler use the single cell value for all
183  // integration points.
184  for (std::size_t i = 0; i < _local_assemblers.size(); ++i)
185  {
186  auto& local_asm = _local_assemblers[i];
187 
188  std::vector<double> value(
189  &mesh_property[i],
190  &mesh_property[i] +
191  mesh_property.getNumberOfGlobalComponents());
192  // TODO (naumov) Check sizes / read size / etc.
193  // OR reconstruct dimensions from size / component =
194  // ip_points
195  local_asm->setIPDataInitialConditionsFromCellData(name, value);
196  }
197  }
198  }
199 
200  // Initialize local assemblers after all variables have been set.
204 }
#define OGS_FATAL(...)
Definition: Error.h:26
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:327
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
@ BY_LOCATION
Ordering data by spatial location.
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface >> &local_assemblers, ExtraCtorArgs &&... extra_ctor_args)
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Mesh const &mesh, std::string const &name)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
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
virtual std::vector< double > const & getIntPtDamage(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual void nonlocal(std::size_t const mesh_item_id, std::vector< std::unique_ptr< SmallDeformationNonlocalLocalAssemblerInterface >> const &local_assemblers)=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 & getIntPtEpsPDXX(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 & getIntPtEpsPV(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0

References ProcessLib::Process::_secondary_variables, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), NumLib::BY_LOCATION, MeshLib::Cell, ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Properties::existsPropertyVector(), MeshLib::Mesh::getElements(), ProcessLib::getIntegrationPointMetaData(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::IntegrationPoint, MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), MaterialPropertyLib::name, and OGS_FATAL.

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::isLinear
override

Definition at line 61 of file SmallDeformationNonlocalProcess.cpp.

62 {
63  return false;
64 }

◆ postIterationConcreteProcess()

template<int DisplacementDim>
NumLib::IterationResult ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::postIterationConcreteProcess ( GlobalVector const &  x)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 291 of file SmallDeformationNonlocalProcess.cpp.

293 {
294  _process_data.crack_volume_old = _process_data.crack_volume;
295  _process_data.crack_volume = 0.0;
296 
297  DBUG("PostNonLinearSolver crack volume computation.");
298 
299  const int process_id = 0;
300  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
301 
305  _process_data.crack_volume);
306 
307  INFO("Integral of crack: {:g}", _process_data.crack_volume);
308 
310 }
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
virtual void computeCrackIntegral(std::size_t mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double &crack_volume)=0

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), INFO(), and NumLib::SUCCESS.

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const &  x,
double const  t,
double const  dt,
int const  process_id 
)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 271 of file SmallDeformationNonlocalProcess.cpp.

276 {
277  DBUG("PostTimestep SmallDeformationNonlocalProcess.");
278  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
279  dof_tables.reserve(x.size());
280  std::generate_n(std::back_inserter(dof_tables), x.size(),
281  [&]() { return _local_to_global_index_map.get(); });
282 
283  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
286  pv.getActiveElementIDs(), dof_tables, x, t, dt);
287 }
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postTimestep().

◆ preAssembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::preAssembleConcreteProcess ( const double  t,
double const  dt,
GlobalVector const &  x 
)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 228 of file SmallDeformationNonlocalProcess.cpp.

231 {
232  DBUG("preAssemble SmallDeformationNonlocalProcess.");
233 
234  const int process_id = 0;
235  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
236 
237  // Call global assembler for each local assembly item.
241  *_local_to_global_index_map, t, dt, x);
242 }
void preAssemble(const std::size_t mesh_item_id, LocalAssemblerInterface &local_assembler, const NumLib::LocalToGlobalIndexMap &dof_table, const double t, double const dt, const GlobalVector &x)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::VectorMatrixAssembler::preAssemble().

Member Data Documentation

◆ _local_assemblers

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 84 of file SmallDeformationNonlocalProcess.h.

◆ _nodal_forces

◆ _process_data

template<int DisplacementDim>
SmallDeformationNonlocalProcessData<DisplacementDim> ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::_process_data
private

Definition at line 77 of file SmallDeformationNonlocalProcess.h.


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