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