PRISMS-PF  v2.1
test_setRigidBodyModeConstraints.h
Go to the documentation of this file.
1 // Unit test(s) for the method "setRigidBodyModeConstraints"
2 template <int dim, int degree>
4 {
5  public:
7  //init the MatrixFreePDE class for testing
8  this->initForTests();
9  };
11  this->matrixFreeObject.clear();
12  };
13 
14  void call_setRigidBodyModeConstraints(std::vector<int> rigidBodyModeComponents, unsigned int & num_constraints){
15 
16  this->setRigidBodyModeConstraints(rigidBodyModeComponents,this->constraintsOtherSet_nonconst[0],this->dofHandlersSet_nonconst[0]);
17 
18  // Calculate the number of constraints that were set
19  num_constraints = this->constraintsOtherSet_nonconst[0]->n_constraints();
20  };
21 
22 
23  void setBCs(){};
24 
25  // Function to set the initial conditions (in ICs_and_BCs.h)
26  void setInitialCondition(const dealii::Point<dim> &p, const unsigned int index, double & scalar_IC, dealii::Vector<double> & vector_IC){};
27 
28  // Function to set the non-uniform Dirichlet boundary conditions (in ICs_and_BCs.h)
29  void setNonUniformDirichletBCs(const dealii::Point<dim> &p, const unsigned int index, const unsigned int direction, const double time, double & scalar_BC, dealii::Vector<double> & vector_BC){};
30 
31  private:
32  //RHS implementation for explicit solve
33  void getRHS(const MatrixFree<dim,double> &data,
34  std::vector<vectorType*> &dst,
35  const std::vector<vectorType*> &src,
36  const std::pair<unsigned int,unsigned int> &cell_range) const{};
37 
38  // Function to set the RHS of the governing equations for explicit time dependent equations (in equations.h)
39  void explicitEquationRHS(variableContainer<dim,degree,dealii::VectorizedArray<double> > & variable_list,
40  dealii::Point<dim, dealii::VectorizedArray<double> > q_point_loc) const {};
41 
42  // Function to set the RHS of the governing equations for all other equations (in equations.h)
43  void nonExplicitEquationRHS(variableContainer<dim,degree,dealii::VectorizedArray<double> > & variable_list,
44  dealii::Point<dim, dealii::VectorizedArray<double> > q_point_loc) const {};
45 
46  // Function to set the LHS of the governing equations (in equations.h)
47  void equationLHS(variableContainer<dim,degree,dealii::VectorizedArray<double> > & variable_list,
48  dealii::Point<dim, dealii::VectorizedArray<double> > q_point_loc) const {};
49 
50 };
51 
52 template <int dim,typename T>
53 bool unitTest<dim,T>::test_setRigidBodyModeConstraints(std::vector<int> rigidBodyModeComponents, userInputParameters<dim> userInputs){
54 
55  bool pass = false;
56  std::cout << "\nTesting 'setRigidBodyModeConstraints' with " << rigidBodyModeComponents.size() << " component(s) needing a constraint...'" << std::endl;
57 
58  //create test problem class object
60  unsigned int num_constraints;
61  test.call_setRigidBodyModeConstraints(rigidBodyModeComponents,num_constraints);
62 
63  // Add up the total number of constraints across all processors
64  unsigned int global_num_constraints = rigidBodyModeComponents.size();
65  Utilities::MPI::sum(global_num_constraints,MPI_COMM_WORLD);
66 
67  // Check if calculated value equals expected value
68  if (global_num_constraints == rigidBodyModeComponents.size()){
69  pass = true;
70  }
71 
72  char buffer[100];
73  sprintf (buffer, "Test result for 'setRigidBodyModeConstraints' with %lu component(s) needing a constraint: %u\n", rigidBodyModeComponents.size(), pass);
74  std::cout << buffer;
75 
76  return pass;
77 }
void equationLHS(variableContainer< dim, degree, dealii::VectorizedArray< double > > &variable_list, dealii::Point< dim, dealii::VectorizedArray< double > > q_point_loc) const
void setRigidBodyModeConstraints(const std::vector< int >, ConstraintMatrix *, const DoFHandler< dim > *) const
void explicitEquationRHS(variableContainer< dim, degree, dealii::VectorizedArray< double > > &variable_list, dealii::Point< dim, dealii::VectorizedArray< double > > q_point_loc) const
setRigidBodyModeConstraintsTest(userInputParameters< dim > _userInputs)
bool test_setRigidBodyModeConstraints(std::vector< int >, userInputParameters< dim > userInputs)
MatrixFree< dim, double > matrixFreeObject
void nonExplicitEquationRHS(variableContainer< dim, degree, dealii::VectorizedArray< double > > &variable_list, dealii::Point< dim, dealii::VectorizedArray< double > > q_point_loc) const
std::vector< DoFHandler< dim > * > dofHandlersSet_nonconst
std::vector< ConstraintMatrix * > constraintsOtherSet_nonconst
void getRHS(const MatrixFree< dim, double > &data, std::vector< vectorType *> &dst, const std::vector< vectorType *> &src, const std::pair< unsigned int, unsigned int > &cell_range) const
void setNonUniformDirichletBCs(const dealii::Point< dim > &p, const unsigned int index, const unsigned int direction, const double time, double &scalar_BC, dealii::Vector< double > &vector_BC)
void initForTests()
Definition: initForTests.cc:7
void setInitialCondition(const dealii::Point< dim > &p, const unsigned int index, double &scalar_IC, dealii::Vector< double > &vector_IC)
void call_setRigidBodyModeConstraints(std::vector< int > rigidBodyModeComponents, unsigned int &num_constraints)