refineable_spherical_advection_diffusion_elements.cc
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
27 
28 namespace oomph
29 {
30  //==========================================================================
31  /// Add the element's contribution to the elemental residual vector
32  /// and/or elemental jacobian matrix.
33  /// This function overloads the standard version so that the possible
34  /// presence of hanging nodes is taken into account.
35  //=========================================================================
38  Vector<double>& residuals,
39  DenseMatrix<double>& jacobian,
40  DenseMatrix<double>& mass_matrix,
41  unsigned flag)
42  {
43  // Find out how many nodes there are in the element
44  const unsigned n_node = nnode();
45 
46  // Get the nodal index at which the unknown is stored
47  const unsigned u_nodal_index = this->u_index_spherical_adv_diff();
48 
49  // Set up memory for the shape and test functions
50  Shape psi(n_node), test(n_node);
51  DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
52 
53  // Set the value of n_intpt
54  const unsigned n_intpt = integral_pt()->nweight();
55 
56  // Set the Vector to hold local coordinates
57  Vector<double> s(2);
58 
59  // Get Peclet number
60  const double scaled_peclet = this->pe();
61 
62  // Get the Peclet number multiplied by the Strouhal number
63  const double scaled_peclet_st = this->pe_st();
64 
65  // Integers used to store the local equation number and local unknown
66  // indices for the residuals and jacobians
67  int local_eqn = 0, local_unknown = 0;
68 
69  // Local storage for pointers to hang_info objects
70  HangInfo *hang_info_pt = 0, *hang_info2_pt = 0;
71 
72  // Local variable to determine the ALE stuff
73  // bool ALE_is_disabled_flag = this->ALE_is_disabled;
74 
75  // Loop over the integration points
76  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
77  {
78  // Assign values of s
79  for (unsigned i = 0; i < 2; i++) s[i] = integral_pt()->knot(ipt, i);
80 
81  // Get the integral weight
82  double w = integral_pt()->weight(ipt);
83 
84  // Call the derivatives of the shape and test functions
86  ipt, psi, dpsidx, test, dtestdx);
87 
88  // Premultiply the weights and the Jacobian
89  double W = w * J;
90 
91  // Calculate local values of the function, initialise to zero
92  double dudt = 0.0;
93  double interpolated_u = 0.0;
94 
95  // These need to be a Vector to be ANSI C++, initialise to zero
97  Vector<double> interpolated_dudx(2, 0.0);
98  // Vector<double> mesh_velocity(DIM,0.0);
99 
100  // Calculate function value and derivatives:
101  //-----------------------------------------
102 
103  // Loop over nodes
104  for (unsigned l = 0; l < n_node; l++)
105  {
106  // Get the value at the node
107  double u_value = this->nodal_value(l, u_nodal_index);
108  interpolated_u += u_value * psi(l);
109  dudt += this->du_dt_spherical_adv_diff(l) * psi(l);
110  // Loop over directions
111  for (unsigned j = 0; j < 2; j++)
112  {
113  interpolated_x[j] += nodal_position(l, j) * psi(l);
114  interpolated_dudx[j] += u_value * dpsidx(l, j);
115  }
116  }
117 
118  // Get the mesh velocity, if required
119  /* if (!ALE_is_disabled_flag)
120  {
121  for(unsigned l=0;l<n_node;l++)
122  {
123  // Loop over directions
124  for(unsigned j=0;j<2;j++)
125  {
126  mesh_velocity[j] += dnodal_position_dt(l,j)*psi(l);
127  }
128  }
129  }*/
130 
131 
132  // Get body force
133  double source;
134  this->get_source_spherical_adv_diff(ipt, interpolated_x, source);
135 
136 
137  // Get wind
138  //--------
139  Vector<double> wind(3);
140  this->get_wind_spherical_adv_diff(ipt, s, interpolated_x, wind);
141 
142  // r is the first position component
143  double r = interpolated_x[0];
144  // theta is the second position component
145  double sin_th = sin(interpolated_x[1]);
146  // dS is the area weighting
147  double dS = r * r * sin_th;
148 
149 
150  // Assemble residuals and Jacobian
151  //================================
152 
153  // Loop over the nodes for the test functions
154  for (unsigned l = 0; l < n_node; l++)
155  {
156  // Local variables to store the number of master nodes and
157  // the weight associated with the shape function if the node is hanging
158  unsigned n_master = 1;
159  double hang_weight = 1.0;
160  // Local bool (is the node hanging)
161  bool is_node_hanging = this->node_pt(l)->is_hanging();
162 
163  // If the node is hanging, get the number of master nodes
164  if (is_node_hanging)
165  {
166  hang_info_pt = this->node_pt(l)->hanging_pt();
167  n_master = hang_info_pt->nmaster();
168  }
169  // Otherwise there is just one master node, the node itself
170  else
171  {
172  n_master = 1;
173  }
174 
175  // Loop over the number of master nodes
176  for (unsigned m = 0; m < n_master; m++)
177  {
178  // Get the local equation number and hang_weight
179  // If the node is hanging
180  if (is_node_hanging)
181  {
182  // Read out the local equation from the master node
183  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
184  u_nodal_index);
185  // Read out the weight from the master node
186  hang_weight = hang_info_pt->master_weight(m);
187  }
188  // If the node is not hanging
189  else
190  {
191  // The local equation number comes from the node itself
192  local_eqn = this->nodal_local_eqn(l, u_nodal_index);
193  // The hang weight is one
194  hang_weight = 1.0;
195  }
196 
197  // If the nodal equation is not a boundary conditino
198  if (local_eqn >= 0)
199  {
200  // Add du/dt and body force/source term here
201  residuals[local_eqn] -= (scaled_peclet_st * dudt + source) * dS *
202  test(l) * W * hang_weight;
203 
204  // The Advection Diffusion bit itself
205  residuals[local_eqn] -=
206  // radial terms
207  (dS * interpolated_dudx[0] *
208  (scaled_peclet * wind[0] * test(l) + dtestdx(l, 0)) +
209  // azimuthal terms
210  (sin_th * interpolated_dudx[1] *
211  (r * scaled_peclet * wind[1] * test(l) + dtestdx(l, 1)))) *
212  W * hang_weight;
213 
214  // Calculate the Jacobian
215  if (flag)
216  {
217  // Local variables to store the number of master nodes
218  // and the weights associated with each hanging node
219  unsigned n_master2 = 1;
220  double hang_weight2 = 1.0;
221  // Loop over the nodes for the variables
222  for (unsigned l2 = 0; l2 < n_node; l2++)
223  {
224  // Local bool (is the node hanging)
225  bool is_node2_hanging = this->node_pt(l2)->is_hanging();
226  // If the node is hanging, get the number of master nodes
227  if (is_node2_hanging)
228  {
229  hang_info2_pt = this->node_pt(l2)->hanging_pt();
230  n_master2 = hang_info2_pt->nmaster();
231  }
232  // Otherwise there is one master node, the node itself
233  else
234  {
235  n_master2 = 1;
236  }
237 
238  // Loop over the master nodes
239  for (unsigned m2 = 0; m2 < n_master2; m2++)
240  {
241  // Get the local unknown and weight
242  // If the node is hanging
243  if (is_node2_hanging)
244  {
245  // Read out the local unknown from the master node
246  local_unknown = this->local_hang_eqn(
247  hang_info2_pt->master_node_pt(m2), u_nodal_index);
248  // Read out the hanging weight from the master node
249  hang_weight2 = hang_info2_pt->master_weight(m2);
250  }
251  // If the node is not hanging
252  else
253  {
254  // The local unknown number comes from the node
255  local_unknown = this->nodal_local_eqn(l2, u_nodal_index);
256  // The hang weight is one
257  hang_weight2 = 1.0;
258  }
259 
260  // If the unknown is not pinned
261  if (local_unknown >= 0)
262  {
263  // Add contribution to Elemental Matrix
264  // Mass matrix du/dt term
265  jacobian(local_eqn, local_unknown) -=
266  scaled_peclet_st * test(l) * psi(l2) *
267  this->node_pt(l2)->time_stepper_pt()->weight(1, 0) * dS *
268  W * hang_weight * hang_weight2;
269 
270  // Add contribution to mass matrix
271  if (flag == 2)
272  {
273  mass_matrix(local_eqn, local_unknown) +=
274  scaled_peclet_st * test(l) * psi(l2) * dS * W *
275  hang_weight * hang_weight2;
276  }
277 
278  // Add contribution to Elemental Matrix
279  // Assemble Jacobian term
280  jacobian(local_eqn, local_unknown) -=
281  // radial terms
282  (dS * dpsidx(l2, 0) *
283  (scaled_peclet * wind[0] * test(l) + dtestdx(l, 0)) +
284  // azimuthal terms
285  (sin_th * dpsidx(l2, 1) *
286  (r * scaled_peclet * wind[1] * test(l) +
287  dtestdx(l, 1)))) *
288  W * hang_weight * hang_weight2;
289  }
290  } // End of loop over master nodes
291  } // End of loop over nodes
292  } // End of Jacobian calculation
293 
294  } // End of non-zero equation
295 
296  } // End of loop over the master nodes for residual
297  } // End of loop over nodes
298 
299  } // End of loop over integration points
300  }
301 
302 
303  //====================================================================
304  // Force build of templates
305  //====================================================================
309 
310 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3992
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition: elements.h:1436
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1967
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Definition: elements.h:2321
Class that contains data for hanging nodes.
Definition: nodes.h:742
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
Definition: nodes.h:1228
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Access function that returns the local equation number for the hanging node variables (values stored ...
Refineable version of QSphericalAdvectionDiffusionElement. Inherit from the standard QSphericalAdvect...
void fill_in_generic_residual_contribution_spherical_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to the elemental residual vector and/or Jacobian matrix flag=1: comput...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
virtual void get_source_spherical_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
virtual double dshape_and_dtest_eulerian_at_knot_spherical_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
double du_dt_spherical_adv_diff(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
const double & pe_st() const
Peclet number multiplied by Strouhal number.
virtual unsigned u_index_spherical_adv_diff() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
virtual void get_wind_spherical_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Get wind at (Eulerian) position x and/or local coordinate s. This function is virtual to allow overlo...
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
//////////////////////////////////////////////////////////////////// ////////////////////////////////...