KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
scheme.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Riccardo Rossi
11 //
12 
13 #pragma once
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "includes/model_part.h"
21 #include "utilities/openmp_utils.h" //TODO: SOME FILES INCLUDING scheme.h RELY ON THIS. LEAVING AS FUTURE TODO.
25 
26 namespace Kratos
27 {
42 
52 template<class TSparseSpace,
53  class TDenseSpace //= DenseSpace<double>
54  >
55 class Scheme
56 {
57 public:
60 
63 
66 
68  using TDataType = typename TSparseSpace::DataType;
69 
72 
75 
78 
81 
84 
87 
90 
93 
97 
102  explicit Scheme()
103  {
104  mSchemeIsInitialized = false;
105  mElementsAreInitialized = false;
107  }
111  explicit Scheme(Parameters ThisParameters)
112  {
113  // Validate default parameters
114  ThisParameters = this->ValidateAndAssignParameters(ThisParameters, this->GetDefaultParameters());
115  this->AssignSettings(ThisParameters);
116 
117  mSchemeIsInitialized = false;
118  mElementsAreInitialized = false;
120  }
121 
124  explicit Scheme(Scheme& rOther)
128  {
129  }
130 
133  virtual ~Scheme()
134  {
135  }
136 
140 
144 
149  virtual typename ClassType::Pointer Create(Parameters ThisParameters) const
150  {
151  return Kratos::make_shared<ClassType>(ThisParameters);
152  }
153 
158  virtual Pointer Clone()
159  {
160  return Kratos::make_shared<Scheme>(*this) ;
161  }
162 
168  virtual void Initialize(ModelPart& rModelPart)
169  {
170  KRATOS_TRY
171  mSchemeIsInitialized = true;
172  KRATOS_CATCH("")
173  }
174 
180  {
181  return mSchemeIsInitialized;
182  }
183 
188  void SetSchemeIsInitialized(bool SchemeIsInitializedFlag = true)
189  {
190  mSchemeIsInitialized = SchemeIsInitializedFlag;
191  }
192 
198  {
200  }
201 
206  void SetElementsAreInitialized(bool ElementsAreInitializedFlag = true)
207  {
208  mElementsAreInitialized = ElementsAreInitializedFlag;
209  }
210 
216  {
218  }
219 
224  void SetConditionsAreInitialized(bool ConditionsAreInitializedFlag = true)
225  {
226  mConditionsAreInitialized = ConditionsAreInitializedFlag;
227  }
228 
234  virtual void InitializeElements( ModelPart& rModelPart)
235  {
236  KRATOS_TRY
237 
238  EntitiesUtilities::InitializeEntities<Element>(rModelPart);
239 
241 
242  KRATOS_CATCH("")
243  }
244 
250  virtual void InitializeConditions(ModelPart& rModelPart)
251  {
252  KRATOS_TRY
253 
254  KRATOS_ERROR_IF_NOT(mElementsAreInitialized) << "Before initializing Conditions, initialize Elements FIRST" << std::endl;
255 
256  EntitiesUtilities::InitializeEntities<Condition>(rModelPart);
257 
259 
260  KRATOS_CATCH("")
261  }
262 
273  ModelPart& rModelPart,
275  TSystemVectorType& Dx,
277  )
278  {
279  KRATOS_TRY
280 
281  // Initializes solution step for all of the elements, conditions and constraints
283 
284  KRATOS_CATCH("")
285  }
286 
294  virtual void FinalizeSolutionStep(
295  ModelPart& rModelPart,
297  TSystemVectorType& Dx,
299  {
300  KRATOS_TRY
301 
302  // Finalizes solution step for all of the elements, conditions and constraints
304 
305  KRATOS_CATCH("")
306  }
307 
308  /************************ BEGIN FRACTIONAL STEP METHODS ****************************/
309  /********************* TODO: DECIDE IF NECESSARY TO DEFINE *************************/
310  /***********************************************************************************/
311 
312 // /**
313 // * @brief Initializes solution step, to be used when system is not explicitly defined
314 // * @details For example for fractional step strategies
315 // * @warning Must be defined in derived classes
316 // * @param rModelPart The model part of the problem to solve
317 // */
318 // virtual void InitializeSolutionStep(ModelPart& rModelPart)
319 // {
320 // KRATOS_TRY
321 // KRATOS_CATCH("")
322 // }
323 //
324 // /**
325 // * @brief Finalizes solution step, to be used when system is not explicitly defined
326 // * @details For example for fractional step strategies
327 // * @warning Must be defined in derived classes
328 // * @param rModelPart The model part of the problem to solve
329 // */
330 // virtual void FinalizeSolutionStep(ModelPart& rModelPart)
331 // {
332 // KRATOS_TRY
333 // KRATOS_CATCH("")
334 // }
335 //
336 // /**
337 // * @brief Executed before each fractional step
338 // * @warning Must be defined in derived classes
339 // * @param rModelPart The model part of the problem to solve
340 // */
341 // virtual void InitializeFractionalSolutionStep(ModelPart& rModelPart)
342 // {
343 // KRATOS_TRY
344 // KRATOS_CATCH("")
345 // }
346 //
347 // /**
348 // * @brief Executed after each fractional step
349 // * @warning Must be defined in derived classes
350 // * @param rModelPart The model part of the problem to solve
351 // */
352 // virtual void FinalizeFractionalSolutionStep(ModelPart& rModelPart)
353 // {
354 // KRATOS_TRY
355 // KRATOS_CATCH("")
356 // }
357 
358  /************************ END FRACTIONAL STEP METHODS ****************************/
359  /***********************************************************************************/
360 
372  ModelPart& rModelPart,
374  TSystemVectorType& Dx,
376  )
377  {
378  KRATOS_TRY
379 
380  // Finalizes non-linear iteration for all of the elements, conditions and constraints
382 
383  KRATOS_CATCH("")
384  }
385 
394  ModelPart& rModelPart,
396  TSystemVectorType& Dx,
398  )
399  {
400  KRATOS_TRY
401 
402  // Finalizes non-linear iteration for all of the elements, conditions and constraints
404 
405  KRATOS_CATCH("")
406  }
407 
416  virtual void Predict(
417  ModelPart& rModelPart,
418  DofsArrayType& rDofSet,
420  TSystemVectorType& Dx,
422  )
423  {
424  KRATOS_TRY
425  KRATOS_CATCH("")
426  }
427 
437  virtual void Update(
438  ModelPart& rModelPart,
439  DofsArrayType& rDofSet,
441  TSystemVectorType& Dx,
443  )
444  {
445  KRATOS_TRY
446  KRATOS_CATCH("")
447  }
448 
458  virtual void CalculateOutputData(
459  ModelPart& rModelPart,
460  DofsArrayType& rDofSet,
462  TSystemVectorType& Dx,
464  )
465  {
466  KRATOS_TRY
467  KRATOS_CATCH("")
468  }
469 
474  virtual void CleanOutputData()
475  {
476  KRATOS_TRY
477  KRATOS_CATCH("")
478  }
479 
484  virtual void Clean()
485  {
486  KRATOS_TRY
487  KRATOS_CATCH("")
488  }
489 
494  virtual void Clear()
495  {
496  KRATOS_TRY
497  KRATOS_CATCH("")
498  }
499 
508  virtual int Check(const ModelPart& rModelPart) const
509  {
510  KRATOS_TRY
511  return 0;
512  KRATOS_CATCH("");
513  }
514 
515  virtual int Check(ModelPart& rModelPart)
516  {
517  // calling the const version for backward compatibility
518  const Scheme& r_const_this = *this;
519  const ModelPart& r_const_model_part = rModelPart;
520  return r_const_this.Check(r_const_model_part);
521  }
522 
533  Element& rElement,
534  LocalSystemMatrixType& LHS_Contribution,
535  LocalSystemVectorType& RHS_Contribution,
536  Element::EquationIdVectorType& rEquationIdVector,
537  const ProcessInfo& rCurrentProcessInfo
538  )
539  {
540  rElement.CalculateLocalSystem(LHS_Contribution, RHS_Contribution, rCurrentProcessInfo);
541  }
542 
552  Condition& rCondition,
553  LocalSystemMatrixType& LHS_Contribution,
554  LocalSystemVectorType& RHS_Contribution,
555  Element::EquationIdVectorType& rEquationIdVector,
556  const ProcessInfo& rCurrentProcessInfo
557  )
558  {
559  rCondition.CalculateLocalSystem(LHS_Contribution, RHS_Contribution, rCurrentProcessInfo);
560  }
561 
570  Element& rElement,
571  LocalSystemVectorType& RHS_Contribution,
572  Element::EquationIdVectorType& rEquationIdVector,
573  const ProcessInfo& rCurrentProcessInfo
574  )
575  {
576  rElement.CalculateRightHandSide(RHS_Contribution, rCurrentProcessInfo);
577  }
578 
587  Condition& rCondition,
588  LocalSystemVectorType& RHS_Contribution,
589  Element::EquationIdVectorType& rEquationIdVector,
590  const ProcessInfo& rCurrentProcessInfo
591  )
592  {
593  rCondition.CalculateRightHandSide(RHS_Contribution, rCurrentProcessInfo);
594  }
595 
604  Element& rElement,
605  LocalSystemMatrixType& LHS_Contribution,
606  Element::EquationIdVectorType& rEquationIdVector,
607  const ProcessInfo& rCurrentProcessInfo
608  )
609  {
610  rElement.CalculateLeftHandSide(LHS_Contribution, rCurrentProcessInfo);
611  }
612 
621  Condition& rCondition,
622  LocalSystemMatrixType& LHS_Contribution,
623  Element::EquationIdVectorType& rEquationIdVector,
624  const ProcessInfo& rCurrentProcessInfo
625  )
626  {
627  rCondition.CalculateLeftHandSide(LHS_Contribution, rCurrentProcessInfo);
628  }
629 
636  virtual void EquationId(
637  const Element& rElement,
638  Element::EquationIdVectorType& rEquationId,
639  const ProcessInfo& rCurrentProcessInfo
640  )
641  {
642  rElement.EquationIdVector(rEquationId, rCurrentProcessInfo);
643  }
644 
651  virtual void EquationId(
652  const Condition& rCondition,
653  Element::EquationIdVectorType& rEquationId,
654  const ProcessInfo& rCurrentProcessInfo
655  )
656  {
657  rCondition.EquationIdVector(rEquationId, rCurrentProcessInfo);
658  }
659 
666  virtual void GetDofList(
667  const Element& rElement,
668  Element::DofsVectorType& rDofList,
669  const ProcessInfo& rCurrentProcessInfo
670  )
671  {
672  rElement.GetDofList(rDofList, rCurrentProcessInfo);
673  }
674 
681  virtual void GetDofList(
682  const Condition& rCondition,
683  Element::DofsVectorType& rDofList,
684  const ProcessInfo& rCurrentProcessInfo
685  )
686  {
687  rCondition.GetDofList(rDofList, rCurrentProcessInfo);
688  }
689 
694  {
695  const Parameters default_parameters = Parameters(R"(
696  {
697  "name" : "scheme"
698  })" );
699  return default_parameters;
700  }
701 
706  static std::string Name()
707  {
708  return "scheme";
709  }
710 
714 
718 
722 
724  virtual std::string Info() const
725  {
726  return "Scheme";
727  }
728 
730  virtual void PrintInfo(std::ostream& rOStream) const
731  {
732  rOStream << Info();
733  }
734 
736  virtual void PrintData(std::ostream& rOStream) const
737  {
738  rOStream << Info();
739  }
740 
744 
746 
747 protected:
750 
754 
758 
762 
766 
774  Parameters ThisParameters,
775  const Parameters DefaultParameters
776  ) const
777  {
778  ThisParameters.ValidateAndAssignDefaults(DefaultParameters);
779  return ThisParameters;
780  }
781 
786  virtual void AssignSettings(const Parameters ThisParameters)
787  {
788  }
789 
793 
797 
801 
803 
804 private:
807 
811 
815 
819 
823 
827 
831 
833 
834 }; // Class Scheme
835 
836 } // namespace Kratos.
Base class for all Conditions.
Definition: condition.h:59
virtual void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:273
virtual void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:426
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:260
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:440
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:408
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
Base class for all Elements.
Definition: element.h:60
virtual void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:423
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:437
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:258
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
virtual void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:271
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:405
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
PointerVectorSet< DofType > DofsArrayType
Definition: model_part.h:115
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
void ValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1306
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
virtual void Predict(ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Performing the prediction of the solution.
Definition: scheme.h:416
virtual void Clean()
This function is intended to be called at the end of the solution step to clean up memory storage not...
Definition: scheme.h:484
virtual void CalculateLHSContribution(Condition &rCondition, LocalSystemMatrixType &LHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: scheme.h:620
virtual void CalculateLHSContribution(Element &rElement, LocalSystemMatrixType &LHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
This function is designed to calculate just the LHS contribution.
Definition: scheme.h:603
virtual void Update(ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Performing the update of the solution.
Definition: scheme.h:437
typename TSparseSpace::MatrixType TSystemMatrixType
Matrix type definition.
Definition: scheme.h:71
virtual void CalculateSystemContributions(Condition &rCondition, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: scheme.h:551
bool SchemeIsInitialized()
This method returns if the scheme is initialized.
Definition: scheme.h:179
virtual Parameters GetDefaultParameters() const
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: scheme.h:693
virtual void InitializeElements(ModelPart &rModelPart)
This is the place to initialize the elements.
Definition: scheme.h:234
void SetSchemeIsInitialized(bool SchemeIsInitializedFlag=true)
This method sets if the elements have been initialized or not (true by default)
Definition: scheme.h:188
bool mConditionsAreInitialized
Flag taking in account if the elements were initialized correctly or not.
Definition: scheme.h:757
virtual int Check(ModelPart &rModelPart)
Definition: scheme.h:515
static std::string Name()
Returns the name of the class as used in the settings (snake_case format)
Definition: scheme.h:706
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: scheme.h:730
virtual void CleanOutputData()
Functions that cleans the results data.
Definition: scheme.h:474
virtual void EquationId(const Element &rElement, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo)
This method gets the eqaution id corresponding to the current element.
Definition: scheme.h:636
typename TSparseSpace::VectorType TSystemVectorType
Vector type definition.
Definition: scheme.h:74
typename TDenseSpace::VectorType LocalSystemVectorType
Local system vector type definition.
Definition: scheme.h:80
KRATOS_CLASS_POINTER_DEFINITION(Scheme)
Pointer definition of Scheme.
virtual void GetDofList(const Condition &rCondition, Element::DofsVectorType &rDofList, const ProcessInfo &rCurrentProcessInfo)
Function that returns the list of Degrees of freedom to be assembled in the system for a Given condit...
Definition: scheme.h:681
virtual void InitializeConditions(ModelPart &rModelPart)
This is the place to initialize the conditions.
Definition: scheme.h:250
virtual void FinalizeSolutionStep(ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Function called once at the end of a solution step, after convergence is reached if an iterative proc...
Definition: scheme.h:294
virtual Parameters ValidateAndAssignParameters(Parameters ThisParameters, const Parameters DefaultParameters) const
This method validate and assign default parameters.
Definition: scheme.h:773
bool ConditionsAreInitialized()
This method returns if the conditions are initialized.
Definition: scheme.h:215
virtual void InitializeSolutionStep(ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Function called once at the beginning of each solution step.
Definition: scheme.h:272
virtual void CalculateSystemContributions(Element &rElement, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
This function is designed to be called in the builder and solver to introduce the selected time integ...
Definition: scheme.h:532
virtual void Initialize(ModelPart &rModelPart)
This is the place to initialize the Scheme.
Definition: scheme.h:168
typename TSparseSpace::DataType TDataType
Data type definition.
Definition: scheme.h:68
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: scheme.h:736
typename TDenseSpace::MatrixType LocalSystemMatrixType
Local system matrix type definition.
Definition: scheme.h:77
Scheme(Parameters ThisParameters)
Constructor with Parameters.
Definition: scheme.h:111
void SetElementsAreInitialized(bool ElementsAreInitializedFlag=true)
This method sets if the elements have been initialized or not (true by default)
Definition: scheme.h:206
virtual void AssignSettings(const Parameters ThisParameters)
This method assigns settings to member variables.
Definition: scheme.h:786
bool ElementsAreInitialized()
This method returns if the elements are initialized.
Definition: scheme.h:197
bool mSchemeIsInitialized
Definition: scheme.h:755
Scheme()
Default Constructor.
Definition: scheme.h:102
virtual void GetDofList(const Element &rElement, Element::DofsVectorType &rDofList, const ProcessInfo &rCurrentProcessInfo)
Function that returns the list of Degrees of freedom to be assembled in the system for a Given elemen...
Definition: scheme.h:666
void SetConditionsAreInitialized(bool ConditionsAreInitializedFlag=true)
This method sets if the conditions have been initialized or not (true by default)
Definition: scheme.h:224
virtual int Check(const ModelPart &rModelPart) const
This function is designed to be called once to perform all the checks needed on the input provided....
Definition: scheme.h:508
virtual ClassType::Pointer Create(Parameters ThisParameters) const
Create method.
Definition: scheme.h:149
virtual void EquationId(const Condition &rCondition, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo)
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: scheme.h:651
virtual void Clear()
Liberate internal storage.
Definition: scheme.h:494
Scheme(Scheme &rOther)
Definition: scheme.h:124
virtual void InitializeNonLinIteration(ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
unction to be called when it is needed to initialize an iteration. It is designed to be called at the...
Definition: scheme.h:371
virtual std::string Info() const
Turn back information as a string.
Definition: scheme.h:724
virtual void CalculateRHSContribution(Element &rElement, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
This function is designed to calculate just the RHS contribution.
Definition: scheme.h:569
bool mElementsAreInitialized
Flag to be used in controlling if the Scheme has been initialized or not.
Definition: scheme.h:756
virtual void CalculateOutputData(ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Functions to be called to prepare the data needed for the output of results.
Definition: scheme.h:458
virtual void CalculateRHSContribution(Condition &rCondition, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: scheme.h:586
virtual void FinalizeNonLinIteration(ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Function to be called when it is needed to finalize an iteration. It is designed to be called at the ...
Definition: scheme.h:393
virtual ~Scheme()
Definition: scheme.h:133
virtual Pointer Clone()
Clone method.
Definition: scheme.h:158
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR_IF_NOT(conditional)
Definition: exception.h:163
void FinalizeNonLinearIterationAllEntities(ModelPart &rModelPart)
This method calls FinalizeNonLinearIteration for all the entities (conditions, elements,...
Definition: entities_utilities.cpp:243
void InitializeSolutionStepAllEntities(ModelPart &rModelPart)
This method calls InitializeSolution for all the entities (conditions, elements, constraints)
Definition: entities_utilities.cpp:210
void FinalizeSolutionStepAllEntities(ModelPart &rModelPart)
This method calls FinalizeSolutionStep for all the entities (conditions, elements,...
Definition: entities_utilities.cpp:221
void InitializeNonLinearIterationAllEntities(ModelPart &rModelPart)
This method calls InitializeNonLinearIteration for all the entities (conditions, elements,...
Definition: entities_utilities.cpp:232
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
b
Definition: generate_total_lagrangian_mixed_volumetric_strain_element.py:31
A
Definition: sensitivityMatrix.py:70