libMesh
Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Static Protected Attributes | List of all members
F0 Struct Reference

#include <assembly.h>

Inheritance diagram for F0:
[legend]

Public Member Functions

virtual void interior_assembly (FEMContext &c)
 Perform the element interior assembly. More...
 
virtual void interior_assembly (FEMContext &c)
 Perform the element interior assembly. More...
 
virtual void interior_assembly (FEMContext &c)
 Perform the element interior assembly. More...
 
virtual void boundary_assembly (FEMContext &c)
 Perform the element boundary assembly. More...
 
virtual void get_nodal_values (std::vector< dof_id_type > &, DenseMatrix< Number > &, DenseVector< Number > &, const System &, const Node &)
 Get values to add to the matrix or rhs vector based on node. More...
 

Static Public Member Functions

static std::string get_info ()
 Gets a string containing the reference information. More...
 
static void print_info (std::ostream &out=libMesh::out)
 Prints the reference information, by default to libMesh::out. More...
 
static unsigned int n_objects ()
 Prints the number of outstanding (created, but not yet destroyed) objects. More...
 
static void enable_print_counter_info ()
 Methods to enable/disable the reference counter output from print_info() More...
 
static void disable_print_counter_info ()
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 Data structure to log the information. More...
 

Protected Member Functions

void increment_constructor_count (const std::string &name)
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name)
 Increments the destruction counter. More...
 

Static Protected Attributes

static Counts _counts
 Actually holds the data. More...
 
static Threads::atomic< unsigned int > _n_objects
 The number of objects. More...
 
static Threads::spin_mutex _mutex
 Mutual exclusion object to enable thread-safe reference counting. More...
 
static bool _enable_print_counter
 Flag to control whether reference count information is printed when print_info is called. More...
 

Detailed Description

Definition at line 127 of file assembly.h.

Member Typedef Documentation

◆ Counts

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information.

The log is identified by the class name.

Definition at line 117 of file reference_counter.h.

Member Function Documentation

◆ boundary_assembly()

virtual void F0::boundary_assembly ( FEMContext )
virtual

Perform the element boundary assembly.

Reimplemented from libMesh::ElemAssembly.

Definition at line 233 of file assembly.h.

References libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_residual(), libMesh::FEMContext::get_side_fe(), libMesh::FEMContext::get_side_qrule(), libMesh::FEMContext::has_side_boundary_id(), and libMesh::QBase::n_points().

234  {
235  if (c.has_side_boundary_id(1)) // Output is calculated on the horn "inlet"
236  {
237  const unsigned int p_var = 0;
238 
239  FEBase * side_fe = nullptr;
240  c.get_side_fe(p_var, side_fe);
241 
242  const std::vector<Real> & JxW_face = side_fe->get_JxW();
243 
244  const std::vector<std::vector<Real>> & phi_face = side_fe->get_phi();
245 
246  // The number of local degrees of freedom in each variable
247  const unsigned int n_p_dofs = c.get_dof_indices(p_var).size();
248 
249  // Now we will build the affine operator
250  unsigned int n_sidepoints = c.get_side_qrule().n_points();
251 
252  for (unsigned int qp=0; qp != n_sidepoints; qp++)
253  for (unsigned int i=0; i != n_p_dofs; i++)
254  c.get_elem_residual()(i) += JxW_face[qp] * phi_face[i][qp];
255  }
256  }
FEGenericBase< Real > FEBase

◆ disable_print_counter_info()

static void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

◆ enable_print_counter_info()

static void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

◆ get_info()

static std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

◆ get_nodal_values()

virtual void libMesh::ElemAssembly::get_nodal_values ( std::vector< dof_id_type > &  ,
DenseMatrix< Number > &  ,
DenseVector< Number > &  ,
const System ,
const Node  
)
virtualinherited

Get values to add to the matrix or rhs vector based on node.

This allows one to impose point loads or springs, for example.

Definition at line 69 of file elem_assembly.h.

74  {
75  // Do nothing by default
76  }

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
protectedinherited

Increments the construction counter.

Should be called in the constructor of any derived class that will be reference counted.

Definition at line 181 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }
std::string name(const ElemQuality q)
static Counts _counts
Actually holds the data.
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
protectedinherited

Increments the destruction counter.

Should be called in the destructor of any derived class that will be reference counted.

Definition at line 194 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }
std::string name(const ElemQuality q)
static Counts _counts
Actually holds the data.
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.

◆ interior_assembly() [1/3]

virtual void F0::interior_assembly ( FEMContext )
virtual

Perform the element interior assembly.

Reimplemented from libMesh::ElemAssembly.

Definition at line 130 of file assembly.h.

References libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_residual(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), and libMesh::QBase::n_points().

131  {
132  const unsigned int u_var = 0;
133 
134  FEBase * elem_fe = nullptr;
135  c.get_element_fe(u_var, elem_fe);
136 
137  const std::vector<Real> & JxW = elem_fe->get_JxW();
138 
139  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
140 
141  // The number of local degrees of freedom in each variable
142  const unsigned int n_u_dofs = c.get_dof_indices(u_var).size();
143 
144  // Now we will build the affine operator
145  unsigned int n_qpoints = c.get_element_qrule().n_points();
146 
147  for (unsigned int qp=0; qp != n_qpoints; qp++)
148  for (unsigned int i=0; i != n_u_dofs; i++)
149  c.get_elem_residual()(i) += JxW[qp] * (1.*phi[i][qp]);
150  }
FEGenericBase< Real > FEBase

◆ interior_assembly() [2/3]

virtual void F0::interior_assembly ( FEMContext )
virtual

Perform the element interior assembly.

Reimplemented from libMesh::ElemAssembly.

Definition at line 158 of file assembly.h.

References libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_residual(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), and libMesh::QBase::n_points().

159  {
160  const unsigned int u_var = 0;
161 
162  FEBase * elem_fe = nullptr;
163  c.get_element_fe(u_var, elem_fe);
164 
165  const std::vector<Real> & JxW = elem_fe->get_JxW();
166 
167  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
168 
169  // The number of local degrees of freedom in each variable
170  const unsigned int n_u_dofs = c.get_dof_indices(u_var).size();
171 
172  // Now we will build the affine operator
173  unsigned int n_qpoints = c.get_element_qrule().n_points();
174 
175  for (unsigned int qp=0; qp != n_qpoints; qp++)
176  for (unsigned int i=0; i != n_u_dofs; i++)
177  c.get_elem_residual()(i) += JxW[qp] * (1.*phi[i][qp]);
178  }
FEGenericBase< Real > FEBase

◆ interior_assembly() [3/3]

virtual void F0::interior_assembly ( FEMContext )
virtual

Perform the element interior assembly.

Reimplemented from libMesh::ElemAssembly.

Definition at line 174 of file assembly.h.

References libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_residual(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), and libMesh::QBase::n_points().

175  {
176  const unsigned int u_var = 0;
177 
178  FEBase * elem_fe = nullptr;
179  c.get_element_fe(u_var, elem_fe);
180 
181  const std::vector<Real> & JxW = elem_fe->get_JxW();
182 
183  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
184 
185  // The number of local degrees of freedom in each variable
186  const unsigned int n_u_dofs = c.get_dof_indices(u_var).size();
187 
188  // Now we will build the affine operator
189  unsigned int n_qpoints = c.get_element_qrule().n_points();
190 
191  for (unsigned int qp=0; qp != n_qpoints; qp++)
192  for (unsigned int i=0; i != n_u_dofs; i++)
193  c.get_elem_residual()(i) += JxW[qp] * (1.*phi[i][qp]);
194  }
FEGenericBase< Real > FEBase

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
staticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects
The number of objects.

◆ print_info()

static void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Member Data Documentation

◆ _counts

Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 141 of file reference_counter.h.

◆ _mutex

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

Threads::atomic<unsigned int> libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects.

Print the reference count information when the number returns to 0.

Definition at line 130 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().


The documentation for this struct was generated from the following file: