PRISMS-PF  v2.1
computeRHS.cc
Go to the documentation of this file.
1 //computeRHS() method for MatrixFreePDE class
2 
3 #include "../../include/matrixFreePDE.h"
4 #include "../../include/variableContainer.h"
5 
6 //update RHS of each field
7 template <int dim, int degree>
9  //log time
10  computing_timer.enter_section("matrixFreePDE: computeRHS");
11 
12  //clear residual vectors before update
13  for(unsigned int fieldIndex=0; fieldIndex<fields.size(); fieldIndex++){
14  if (userInputs.var_eq_type[fieldIndex] == EXPLICIT_TIME_DEPENDENT){
15  (*residualSet[fieldIndex])=0.0;
16  }
17  }
18 
19  //call to integrate and assemble
20  matrixFreeObject.cell_loop (&MatrixFreePDE<dim,degree>::getExplicitRHS, this, residualSet, solutionSet);
21 
22  //end log
23  computing_timer.exit_section("matrixFreePDE: computeRHS");
24 }
25 
26 template <int dim, int degree>
27 void MatrixFreePDE<dim,degree>::getExplicitRHS(const MatrixFree<dim,double> &data,
28  std::vector<vectorType*> &dst,
29  const std::vector<vectorType*> &src,
30  const std::pair<unsigned int,unsigned int> &cell_range) const{
31 
32  variableContainer<dim,degree,dealii::VectorizedArray<double> > variable_list(data,userInputs.varInfoListExplicitRHS);
33 
34  //loop over cells
35  for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell){
36 
37  // Initialize, read DOFs, and set evaulation flags for each variable
38  variable_list.reinit_and_eval(src, cell);
39 
40  unsigned int num_q_points = variable_list.get_num_q_points();
41 
42  //loop over quadrature points
43  for (unsigned int q=0; q<num_q_points; ++q){
44  variable_list.q_point = q;
45 
46  dealii::Point<dim, dealii::VectorizedArray<double> > q_point_loc = variable_list.get_q_point_location();
47 
48  // Calculate the residuals
49  explicitEquationRHS(variable_list,q_point_loc);
50  }
51 
52  variable_list.integrate_and_distribute(dst);
53  }
54 }
55 
56 //update RHS of each field
57 template <int dim, int degree>
59  //log time
60  computing_timer.enter_section("matrixFreePDE: computeRHS");
61 
62  //clear residual vectors before update
63  for(unsigned int fieldIndex=0; fieldIndex<fields.size(); fieldIndex++){
64  if (userInputs.var_eq_type[fieldIndex] != EXPLICIT_TIME_DEPENDENT){
65  (*residualSet[fieldIndex])=0.0;
66  }
67  }
68 
69  //call to integrate and assemble
70  matrixFreeObject.cell_loop (&MatrixFreePDE<dim,degree>::getNonexplicitRHS, this, residualSet, solutionSet);
71 
72  //end log
73  computing_timer.exit_section("matrixFreePDE: computeRHS");
74 }
75 
76 template <int dim, int degree>
77 void MatrixFreePDE<dim,degree>::getNonexplicitRHS(const MatrixFree<dim,double> &data,
78  std::vector<vectorType*> &dst,
79  const std::vector<vectorType*> &src,
80  const std::pair<unsigned int,unsigned int> &cell_range) const{
81 
82  variableContainer<dim,degree,dealii::VectorizedArray<double> > variable_list(data,userInputs.varInfoListNonexplicitRHS);
83 
84  //loop over cells
85  for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell){
86 
87  // Initialize, read DOFs, and set evaulation flags for each variable
88  variable_list.reinit_and_eval(src, cell);
89 
90  unsigned int num_q_points = variable_list.get_num_q_points();
91 
92  //loop over quadrature points
93  for (unsigned int q=0; q<num_q_points; ++q){
94  variable_list.q_point = q;
95 
96  dealii::Point<dim, dealii::VectorizedArray<double> > q_point_loc = variable_list.get_q_point_location();
97 
98  // Calculate the residuals
99  nonExplicitEquationRHS(variable_list,q_point_loc);
100  }
101 
102  variable_list.integrate_and_distribute(dst);
103  }
104 }
105 
106  #include "../../include/matrixFreePDE_template_instantiations.h"
void integrate_and_distribute(std::vector< vectorType *> &dst)
unsigned int q_point
void getExplicitRHS(const MatrixFree< dim, double > &data, std::vector< vectorType *> &dst, const std::vector< vectorType *> &src, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition: computeRHS.cc:27
void reinit_and_eval(const std::vector< vectorType *> &src, unsigned int cell)
unsigned int get_num_q_points()
void getNonexplicitRHS(const MatrixFree< dim, double > &data, std::vector< vectorType *> &dst, const std::vector< vectorType *> &src, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition: computeRHS.cc:77
void computeExplicitRHS()
Definition: computeRHS.cc:8
dealii::Point< dim, T > get_q_point_location()
void computeNonexplicitRHS()
Definition: computeRHS.cc:58