Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #pragma once 11 : 12 : // MOOSE includes 13 : #include "Constraint.h" 14 : #include "NeighborCoupleableMooseVariableDependencyIntermediateInterface.h" 15 : 16 : namespace libMesh 17 : { 18 : template <typename T> 19 : class SparseMatrix; 20 : } 21 : 22 : /** 23 : * A NodeElemConstraintBase is used when you need to create constraints between 24 : * a secondary node and a primary element. It works by allowing you to modify the 25 : * residual and jacobian entries on the secondary node and the primary element. 26 : */ 27 : class NodeElemConstraintBase 28 : : public Constraint, 29 : public NeighborCoupleableMooseVariableDependencyIntermediateInterface, 30 : public NeighborMooseVariableInterface<Real> 31 : { 32 : public: 33 : static InputParameters validParams(); 34 : 35 : NodeElemConstraintBase(const InputParameters & parameters); 36 : ~NodeElemConstraintBase(); 37 : 38 : /// Compute the value the secondary node should have at the beginning of a timestep. 39 : void computeSecondaryValue(NumericVector<Number> & current_solution); 40 : 41 : /** 42 : * Whether or not this constraint should be applied. 43 : * @return bool true if this constraint is active, false otherwise 44 : */ 45 0 : virtual bool shouldApply() { return true; } 46 : 47 : /** 48 : * Whether or not the secondary's residual should be overwritten. 49 : * @return bool When this returns true the secondary's residual as computed by the constraint will 50 : * _replace_ the residual previously at that node for that variable. 51 : */ 52 : virtual bool overwriteSecondaryResidual() const; 53 : 54 : /** 55 : * Whether or not the secondary's Jacobian row should be overwritten. 56 : * @return bool When this returns true the secondary's Jacobian row as computed by the constraint 57 : * will _replace_ the residual previously at that node for that variable. 58 : */ 59 4970 : virtual bool overwriteSecondaryJacobian() const { return overwriteSecondaryResidual(); }; 60 : 61 : /** 62 : * The variable on the primary elem. 63 : * @return MooseVariable & a reference to the primary variable 64 : */ 65 9940 : virtual MooseVariable & primaryVariable() const { return _primary_var; } 66 : 67 : /** 68 : * The variable number that this object operates on. 69 : */ 70 27650 : const MooseVariable & variable() const override { return _var; } 71 : 72 : protected: 73 : /// Compute the value the secondary node should have at the beginning of a timestep. 74 : virtual Real computeQpSecondaryValue() = 0; 75 : 76 : /// secondary node variable 77 : MooseVariable & _var; 78 : /// Primary side variable 79 : MooseVariable & _primary_var; 80 : 81 : /// secondary block id 82 : unsigned short _secondary; 83 : /// primary block id 84 : unsigned short _primary; 85 : 86 : /// current node being processed 87 : const Node * const & _current_node; 88 : /// current element being processed 89 : const Elem * const & _current_elem; 90 : 91 : /// Shape function on the secondary side. 92 : VariablePhiValue _phi_secondary; 93 : /// Shape function on the secondary side. This will always only have one entry and that entry will always be "1" 94 : VariableTestValue _test_secondary; 95 : 96 : /// Side shape function. 97 : const VariablePhiValue & _phi_primary; 98 : 99 : /// Side test function 100 : const VariableTestValue & _test_primary; 101 : 102 : /// MooseMesh map of current nodes to the connected elements 103 : const std::map<dof_id_type, std::vector<dof_id_type>> & _node_to_elem_map; 104 : 105 : /// maps secondary node ids to primary element ids 106 : std::map<dof_id_type, dof_id_type> _secondary_to_primary_map; 107 : 108 : /** 109 : * Whether or not the secondary's residual should be overwritten. 110 : * 111 : * When this is true the secondary's residual as computed by the constraint will _replace_ 112 : * the residual previously at that node for that variable. 113 : */ 114 : bool _overwrite_secondary_residual; 115 : 116 : public: 117 : /// system Jacobian, provides pre-constraint Jacobian for nonAD kinematic constraints 118 : const SparseMatrix<Number> * _jacobian; 119 : /// dofs connected to the secondary node 120 : std::vector<dof_id_type> _connected_dof_indices; 121 : /// stiffness matrix holding primary-secondary jacobian 122 : DenseMatrix<Number> _Kne; 123 : /// stiffness matrix holding secondary-secondary jacobian 124 : DenseMatrix<Number> _Kee; 125 : };