refineable_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 //=========================================================================
42 template<unsigned DIM>
45  DenseMatrix<double> &jacobian,
47  &mass_matrix,
48  unsigned flag)
49 {
50 
51 //Find out how many nodes there are in the element
52 unsigned n_node = nnode();
53 
54 //Get the nodal index at which the unknown is stored
55  unsigned u_nodal_index = this->u_index_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,DIM), dtestdx(n_node,DIM);
60 
61 //Set the value of n_intpt
62 unsigned n_intpt = integral_pt()->nweight();
63 
64 //Set the Vector to hold local coordinates
65 Vector<double> s(DIM);
66 
67 //Get Peclet number
68 double peclet = this->pe();
69 
70 //Get the Peclet multiplied by the Strouhal number
71 double 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<DIM;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 =
95  this->dshape_and_dtest_eulerian_at_knot_adv_diff(ipt,psi,dpsidx,
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
106  Vector<double> interpolated_x(DIM,0.0);
107  Vector<double> interpolated_dudx(DIM,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_adv_diff(l)*psi(l);
120  // Loop over directions
121  for(unsigned j=0;j<DIM;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<DIM;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_adv_diff(ipt,interpolated_x,source);
145 
146 
147  //Get wind
148  //--------
149  Vector<double> wind(DIM);
150  this->get_wind_adv_diff(ipt,s,interpolated_x,wind);
151 
152  // Assemble residuals and Jacobian
153  //================================
154 
155  // Loop over the nodes for the test functions
156  for(unsigned l=0;l<n_node;l++)
157  {
158  //Local variables to store the number of master nodes and
159  //the weight associated with the shape function if the node is hanging
160  unsigned n_master=1; double hang_weight=1.0;
161  //Local bool (is the node hanging)
162  bool is_node_hanging = this->node_pt(l)->is_hanging();
163 
164  //If the node is hanging, get the number of master nodes
165  if(is_node_hanging)
166  {
167  hang_info_pt = this->node_pt(l)->hanging_pt();
168  n_master = hang_info_pt->nmaster();
169  }
170  //Otherwise there is just one master node, the node itself
171  else
172  {
173  n_master = 1;
174  }
175 
176  //Loop over the number of master nodes
177  for(unsigned m=0;m<n_master;m++)
178  {
179  //Get the local equation number and hang_weight
180  //If the node is hanging
181  if(is_node_hanging)
182  {
183  //Read out the local equation from the master node
184  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
185  u_nodal_index);
186  //Read out the weight from the master node
187  hang_weight = hang_info_pt->master_weight(m);
188  }
189  //If the node is not hanging
190  else
191  {
192  //The local equation number comes from the node itself
193  local_eqn = this->nodal_local_eqn(l,u_nodal_index);
194  //The hang weight is one
195  hang_weight = 1.0;
196  }
197 
198  //If the nodal equation is not a boundary conditino
199  if(local_eqn >= 0)
200  {
201  //Add du/dt and body force/source term here
202  residuals[local_eqn] -=
203  (peclet_st*dudt + source)*test(l)*W*hang_weight;
204 
205  //The Mesh velocity and Advection--Diffusion bit
206  for(unsigned k=0;k<DIM;k++)
207  {
208  //Terms that multiply the test function
209  double tmp = peclet*wind[k];
210  //If the mesh is moving need to subtract the mesh velocity
211  if(!ALE_is_disabled_flag) tmp -= peclet_st*mesh_velocity[k];
212  //Now combine the terms into the residual
213  residuals[local_eqn] -= interpolated_dudx[k]*
214  (tmp*test(l) + dtestdx(l,k))*W*hang_weight;
215  }
216 
217  // Calculate the Jacobian
218  if(flag)
219  {
220  //Local variables to store the number of master nodes
221  //and the weights associated with each hanging node
222  unsigned n_master2=1; double hang_weight2=1.0;
223  //Loop over the nodes for the variables
224  for(unsigned l2=0;l2<n_node;l2++)
225  {
226  //Local bool (is the node hanging)
227  bool is_node2_hanging = this->node_pt(l2)->is_hanging();
228  //If the node is hanging, get the number of master nodes
229  if(is_node2_hanging)
230  {
231  hang_info2_pt = this->node_pt(l2)->hanging_pt();
232  n_master2 = hang_info2_pt->nmaster();
233  }
234  //Otherwise there is one master node, the node itself
235  else
236  {
237  n_master2 = 1;
238  }
239 
240  //Loop over the master nodes
241  for(unsigned m2=0;m2<n_master2;m2++)
242  {
243  //Get the local unknown and weight
244  //If the node is hanging
245  if(is_node2_hanging)
246  {
247  //Read out the local unknown from the master node
248  local_unknown =
249  this->local_hang_eqn(hang_info2_pt->master_node_pt(m2),
250  u_nodal_index);
251  //Read out the hanging weight from the master node
252  hang_weight2 = hang_info2_pt->master_weight(m2);
253  }
254  //If the node is not hanging
255  else
256  {
257  //The local unknown number comes from the node
258  local_unknown = this->nodal_local_eqn(l2,u_nodal_index);
259  //The hang weight is one
260  hang_weight2 = 1.0;
261  }
262 
263  //If the unknown is not pinned
264  if(local_unknown >= 0)
265  {
266  //Add contribution to Elemental Matrix
267  // Mass matrix du/dt term
268  jacobian(local_eqn,local_unknown)
269  -= peclet_st*test(l)*psi(l2)*
270  this->node_pt(l2)->time_stepper_pt()->weight(1,0)
271  *W*hang_weight*hang_weight2;
272 
273  //Add contribution to mass matrix
274  if(flag==2)
275  {
276  mass_matrix(local_eqn,local_unknown) +=
277  peclet_st*test(l)*psi(l2)*W*hang_weight*hang_weight2;
278  }
279 
280  //Add contribution to Elemental Matrix
281  for(unsigned k=0;k<DIM;k++)
282  {
283  //Terms that multiply the test function
284  double tmp = peclet*wind[k];
285  //If the mesh is moving, subtract the mesh velocity
286  if(!ALE_is_disabled_flag)
287  {tmp -= peclet_st*mesh_velocity[k];}
288  //Construct the jacobian term
289  jacobian(local_eqn,local_unknown) -=
290  dpsidx(l2,k)*(tmp*test(l) + dtestdx(l,k))*W*
291  hang_weight*hang_weight2;
292  }
293  }
294  } //End of loop over master nodes
295  } //End of loop over nodes
296  } //End of Jacobian calculation
297 
298  } //End of non-zero equation
299 
300  } //End of loop over the master nodes for residual
301  } //End of loop over nodes
302 
303  } // End of loop over integration points
304 }
305 
306 
307 
308 //====================================================================
309 // Force build of templates
310 //====================================================================
314 
318 
319 }
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed...
const double & pe() const
Peclet number.
cstr elem_len * i
Definition: cfortran.h:607
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:753
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:733
void fill_in_generic_residual_contribution_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...
const double & pe_st() const
Peclet number multiplied by Strouhal number.
static char t char * s
Definition: cfortran.h:572
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:736
Class that contains data for hanging nodes.
Definition: nodes.h:684
Refineable version of QAdvectionDiffusionElement. Inherit from the standard QAdvectionDiffusionElemen...