SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT > Class Template Reference

#include <poisson_elements_with_singularity.h>

+ Inheritance diagram for SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >:

Public Types

typedef double(* PoissonSingularFctPt) (const Vector< double > &x)
 Function pointer to the singular function: More...
 
typedef Vector< double >(* PoissonGradSingularFctPt) (const Vector< double > &x)
 Function pointer to the gradient of the singular function: More...
 

Public Member Functions

 SingularPoissonSolutionElement ()
 Constructor. More...
 
boolsingular_function_satisfies_laplace_equation ()
 
double c () const
 Find the value of the unknown C. More...
 
void set_wrapped_poisson_element_pt (WRAPPED_POISSON_ELEMENT *wrapped_poisson_el_pt, const Vector< double > &s, unsigned *direction_pt)
 
PoissonSingularFctPtsingular_fct_pt ()
 Access function to pointer to singular function. More...
 
double singular_function (const Vector< double > &x) const
 Evaluate singular function at Eulerian position x. More...
 
double u_bar (const Vector< double > &x)
 
PoissonGradSingularFctPtgrad_singular_fct_pt ()
 Access function to pointer to the gradient of sing function. More...
 
Vector< doublegrad_singular_function (const Vector< double > &x) const
 Evaluate the gradient of the sing function at Eulerian position x. More...
 
Vector< doublegrad_u_bar (const Vector< double > &x)
 Gradient of ubar (including the constant!) More...
 
void fill_in_contribution_to_residuals (Vector< double > &residual)
 
void fill_in_contribution_to_jacobian (Vector< double > &residual, DenseMatrix< double > &jacobian)
 

Private Member Functions

void fill_in_generic_contribution_to_residuals (Vector< double > &residual, DenseMatrix< double > &jacobian, const unsigned &flag)
 Compute local residual, and, if flag=1, local jacobian matrix. More...
 

Private Attributes

WRAPPED_POISSON_ELEMENT * Wrapped_poisson_el_pt
 Pointer to Poisson element that contains the singularity. More...
 
PoissonSingularFctPt Singular_fct_pt
 Pointer to singular function. More...
 
PoissonGradSingularFctPt Grad_singular_fct_pt
 Pointer to the gradient of the singular function. More...
 
Vector< doubleS_in_wrapped_poisson_element
 Local coordinates of singularity in the wrapped Poisson element. More...
 
unsignedDirection_pt
 Direction of the derivative used for the residual of the element. More...
 
bool Singular_function_satisfies_laplace_equation
 Does singular fct satisfy Laplace's eqn? More...
 

Detailed Description

template<class WRAPPED_POISSON_ELEMENT>
class SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// We consider a singularity in the solution at the point O.

This class defines the functions u_sing and u_bar = C*u_sing and their gradients.

The class also defines the function that computes the residual associated with C. R_C = \frac{\partial u_FE}{\partial x_d}(O) where u_FE is : \sum_{i} U_i \Psi_i. This sets the derivative of the finite-element part of the solution to zero in the specified direction (d) and thus imposes regularity on that part.

Member Typedef Documentation

◆ PoissonGradSingularFctPt

template<class WRAPPED_POISSON_ELEMENT >
typedef Vector<double>(* SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::PoissonGradSingularFctPt) (const Vector< double > &x)

Function pointer to the gradient of the singular function:

◆ PoissonSingularFctPt

template<class WRAPPED_POISSON_ELEMENT >
typedef double(* SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::PoissonSingularFctPt) (const Vector< double > &x)

Function pointer to the singular function:

Constructor & Destructor Documentation

◆ SingularPoissonSolutionElement()

template<class WRAPPED_POISSON_ELEMENT >
SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::SingularPoissonSolutionElement ( )
inline

Constructor.

63  {
64  // Initialise Function pointer to singular function
66 
67  // Initialise Function pointer to gradient of the singular function
69 
70  // Initalise pointer to the wrapped Poisson element which will be used to
71  // compute the residual and which includes the point O
73 
74  // Initialise the pointer to the direction of the derivative used
75  // for the residual of this element
76  Direction_pt=0;
77 
78  // Create a single item of internal Data, storing one unknown which
79  // represents the unknown C.
80  add_internal_data(new Data(1));
81 
82  // Safe assumption: Singular fct does not satisfy Laplace's eqn
84 
85  } // End of constructor
WRAPPED_POISSON_ELEMENT * Wrapped_poisson_el_pt
Pointer to Poisson element that contains the singularity.
Definition: poisson_elements_with_singularity.h:288
PoissonGradSingularFctPt Grad_singular_fct_pt
Pointer to the gradient of the singular function.
Definition: poisson_elements_with_singularity.h:294
bool Singular_function_satisfies_laplace_equation
Does singular fct satisfy Laplace's eqn?
Definition: poisson_elements_with_singularity.h:303
unsigned * Direction_pt
Direction of the derivative used for the residual of the element.
Definition: poisson_elements_with_singularity.h:300
PoissonSingularFctPt Singular_fct_pt
Pointer to singular function.
Definition: poisson_elements_with_singularity.h:291

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Direction_pt, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Grad_singular_fct_pt, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Singular_fct_pt, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Singular_function_satisfies_laplace_equation, and SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Wrapped_poisson_el_pt.

Member Function Documentation

◆ c()

template<class WRAPPED_POISSON_ELEMENT >
double SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::c ( ) const
inline

Find the value of the unknown C.

98  {
99  return internal_data_pt(0)->value(0);
100  }

Referenced by SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::grad_u_bar(), and SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::u_bar().

◆ fill_in_contribution_to_jacobian()

template<class WRAPPED_POISSON_ELEMENT >
void SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residual,
DenseMatrix< double > &  jacobian 
)
inline
232  {
233  fill_in_generic_contribution_to_residuals(residual,jacobian,1);
234  }
void fill_in_generic_contribution_to_residuals(Vector< double > &residual, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute local residual, and, if flag=1, local jacobian matrix.
Definition: poisson_elements_with_singularity.h:240

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::fill_in_generic_contribution_to_residuals().

◆ fill_in_contribution_to_residuals()

template<class WRAPPED_POISSON_ELEMENT >
void SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residual)
inline

◆ fill_in_generic_contribution_to_residuals()

template<class WRAPPED_POISSON_ELEMENT >
void SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::fill_in_generic_contribution_to_residuals ( Vector< double > &  residual,
DenseMatrix< double > &  jacobian,
const unsigned flag 
)
inlineprivate

Compute local residual, and, if flag=1, local jacobian matrix.

243  {
244  // Get the local eqn number of our one-and-only
245  // unknown
246  int eqn_number=internal_local_eqn(0,0);
247 
248  // Get the derivative
249  unsigned cached_dim=Wrapped_poisson_el_pt->dim();
250  Vector<double> flux(cached_dim);
252  double derivative = flux[*Direction_pt];
253 
254  // fill in the contribution to the residual
255  residual[eqn_number] = derivative;
256 
257  if (flag)
258  {
259  // Find the number of nodes in the Poisson element associated to the
260  // SingularPoissonSolutionElement class
261  unsigned n_node = Wrapped_poisson_el_pt->nnode();
262 
263  // Compute the derivatives of the shape functions of the poisson
264  // element associated to the SingularPoissonSolutionElement at the
265  // local coordinate S_in_wrapped_poisson_element_pt
266  Shape psi(n_node);
267  DShape dpsidx(n_node,cached_dim);
268  Wrapped_poisson_el_pt->dshape_eulerian(
269  S_in_wrapped_poisson_element,psi,dpsidx);
270 
271  // Loop over the nodes
272  for (unsigned j=0;j<n_node;j++)
273  {
274  // Find the local equation number of the node
275  int node_eqn_number = external_local_eqn(j,0);
276  if (node_eqn_number>=0)
277  {
278  // Add the contribution of the node to the local jacobian
279  jacobian(eqn_number,node_eqn_number) = dpsidx(j,*Direction_pt);
280  }
281 
282  }
283  }
284 
285  }
Vector< double > S_in_wrapped_poisson_element
Local coordinates of singularity in the wrapped Poisson element.
Definition: poisson_elements_with_singularity.h:297
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Direction_pt, ProblemParameters::flux(), j, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::S_in_wrapped_poisson_element, and SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Wrapped_poisson_el_pt.

Referenced by SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::fill_in_contribution_to_jacobian(), and SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::fill_in_contribution_to_residuals().

◆ grad_singular_fct_pt()

template<class WRAPPED_POISSON_ELEMENT >
PoissonGradSingularFctPt& SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::grad_singular_fct_pt ( )
inline

Access function to pointer to the gradient of sing function.

172  {return Grad_singular_fct_pt;}

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Grad_singular_fct_pt.

◆ grad_singular_function()

template<class WRAPPED_POISSON_ELEMENT >
Vector<double> SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::grad_singular_function ( const Vector< double > &  x) const
inline

Evaluate the gradient of the sing function at Eulerian position x.

176  {
177  // Find the dimension of the problem
178  unsigned cached_dim=Wrapped_poisson_el_pt->dim();
179 
180  // Declare the gradient of sing
181  Vector<double> dsingular(cached_dim);
182 
183 #ifdef PARANOID
184  if (Grad_singular_fct_pt==0)
185  {
186  std::stringstream error_stream;
187  error_stream
188  << "Pointer to gradient of singular function hasn't been defined!"
189  << std::endl;
190  throw OomphLibError(
191  error_stream.str(),
194  }
195 #endif
196 
197  // Evaluate gradient of the sing function
198  return (*Grad_singular_fct_pt)(x);
199  }
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Grad_singular_fct_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Wrapped_poisson_el_pt, and plotDoE::x.

Referenced by SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::grad_u_bar().

◆ grad_u_bar()

template<class WRAPPED_POISSON_ELEMENT >
Vector<double> SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::grad_u_bar ( const Vector< double > &  x)
inline

Gradient of ubar (including the constant!)

203  {
204  // Compute the gradient of sing
205  Vector<double> dsingular = grad_singular_function(x);
206 
207  // Initialise the gradient of ubar
208  unsigned cached_dim=Wrapped_poisson_el_pt->dim();
209  Vector<double> dubar(cached_dim);
210 
211  // Find the value of C
212  double c=internal_data_pt(0)->value(0);
213  for (unsigned i=0;i<cached_dim;i++)
214  {
215  dubar[i] = c*dsingular[i];
216  }
217 
218  return dubar;
219  } // End of function
int i
Definition: BiCGSTAB_step_by_step.cpp:9
double c() const
Find the value of the unknown C.
Definition: poisson_elements_with_singularity.h:97
Vector< double > grad_singular_function(const Vector< double > &x) const
Evaluate the gradient of the sing function at Eulerian position x.
Definition: poisson_elements_with_singularity.h:175

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::c(), SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::grad_singular_function(), i, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Wrapped_poisson_el_pt, and plotDoE::x.

◆ set_wrapped_poisson_element_pt()

template<class WRAPPED_POISSON_ELEMENT >
void SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::set_wrapped_poisson_element_pt ( WRAPPED_POISSON_ELEMENT *  wrapped_poisson_el_pt,
const Vector< double > &  s,
unsigned direction_pt 
)
inline

Set pointer to associated wrapped Poisson element which contains the singularity (at local coordinate s). Also specify the direction in which the slope of the FE part of the solution is set to zero

109  {
110  // Assign the pointer to the variable Wrapped_poisson_el_pt
111  Wrapped_poisson_el_pt = wrapped_poisson_el_pt;
112 
113  // Find number of nodes in the element
114  unsigned nnod=wrapped_poisson_el_pt->nnode();
115 
116  // Loop over the nodes of the element
117  for (unsigned j=0;j<nnod;j++)
118  {
119  // Add the node as external data in the SingularPoissonSolutionElement class
120  // because they affect the slope that we set to zero to
121  // determine the value of the amplitude C
122  add_external_data(Wrapped_poisson_el_pt->node_pt(j));
123  }
124 
125  // Assign the pointer to the local coordinate at which the residual
126  // will be computed
128 
129  // Assign the pointer to the direction at which the derivative used
130  // in the residual will be computed
131  Direction_pt = direction_pt;
132  }
RealScalar s
Definition: level1_cplx_impl.h:130

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Direction_pt, j, s, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::S_in_wrapped_poisson_element, and SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Wrapped_poisson_el_pt.

◆ singular_fct_pt()

template<class WRAPPED_POISSON_ELEMENT >
PoissonSingularFctPt& SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::singular_fct_pt ( )
inline

Access function to pointer to singular function.

136 {return Singular_fct_pt;}

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Singular_fct_pt.

◆ singular_function()

template<class WRAPPED_POISSON_ELEMENT >
double SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::singular_function ( const Vector< double > &  x) const
inline

Evaluate singular function at Eulerian position x.

140  {
141 #ifdef PARANOID
142  if (Singular_fct_pt==0)
143  {
144  std::stringstream error_stream;
145  error_stream
146  << "Pointer to singular function hasn't been defined!"
147  << std::endl;
148  throw OomphLibError(
149  error_stream.str(),
152  }
153 #endif
154 
155  // Evaluate singular function
156  return (*Singular_fct_pt)(x);
157  }

References OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Singular_fct_pt, and plotDoE::x.

Referenced by SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::u_bar().

◆ singular_function_satisfies_laplace_equation()

template<class WRAPPED_POISSON_ELEMENT >
bool& SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::singular_function_satisfies_laplace_equation ( )
inline

Does the singular fct satisfy Laplace's eqn? Default assumption is that it doesn't; user can overwrite this here to avoid unnecessary computation of integrals in residual

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Singular_function_satisfies_laplace_equation.

◆ u_bar()

template<class WRAPPED_POISSON_ELEMENT >
double SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::u_bar ( const Vector< double > &  x)
inline

The singular function (including its amplitude): ubar = C * sing where sing is defined via the function pointer.

162  {
163  // Find the value of C
164  double c=internal_data_pt(0)->value(0);
165 
166  // Value of ubar at the position x
167  return c*singular_function(x);
168  } // End of function
double singular_function(const Vector< double > &x) const
Evaluate singular function at Eulerian position x.
Definition: poisson_elements_with_singularity.h:139

References SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::c(), SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::singular_function(), and plotDoE::x.

Member Data Documentation

◆ Direction_pt

◆ Grad_singular_fct_pt

◆ S_in_wrapped_poisson_element

template<class WRAPPED_POISSON_ELEMENT >
Vector<double> SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::S_in_wrapped_poisson_element
private

◆ Singular_fct_pt

◆ Singular_function_satisfies_laplace_equation

template<class WRAPPED_POISSON_ELEMENT >
bool SingularPoissonSolutionElement< WRAPPED_POISSON_ELEMENT >::Singular_function_satisfies_laplace_equation
private

◆ Wrapped_poisson_el_pt


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