refineable_gen_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  //=========================================================================
36  template<unsigned DIM>
39  Vector<double>& residuals,
40  DenseMatrix<double>& jacobian,
41  DenseMatrix<double>& mass_matrix,
42  unsigned flag)
43  {
44  // Find out how many nodes there are in the element
45  unsigned n_node = nnode();
46 
47  // Get the nodal index at which the unknown is stored
48  unsigned u_nodal_index = this->u_index_cons_adv_diff();
49 
50  // Set up memory for the shape and test functions
51  Shape psi(n_node), test(n_node);
52  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
53 
54  // Set the value of n_intpt
55  unsigned n_intpt = integral_pt()->nweight();
56 
57  // Set the Vector to hold local coordinates
58  Vector<double> s(DIM);
59 
60  // Get Peclet number
61  double peclet = this->pe();
62 
63  // Get the Peclet multiplied by the Strouhal number
64  double peclet_st = this->pe_st();
65 
66  // Integers used to store the local equation number and local unknown
67  // indices for the residuals and jacobians
68  int local_eqn = 0, local_unknown = 0;
69 
70  // Local storage for pointers to hang_info objects
71  HangInfo *hang_info_pt = 0, *hang_info2_pt = 0;
72 
73  // Local variable to determine the ALE stuff
74  bool ALE_is_disabled_flag = this->ALE_is_disabled;
75 
76  // Loop over the integration points
77  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
78  {
79  // Assign values of s
80  for (unsigned i = 0; i < DIM; i++) s[i] = integral_pt()->knot(ipt, i);
81 
82  // Get the integral weight
83  double w = integral_pt()->weight(ipt);
84 
85  // Call the derivatives of the shape and test functions
86  double J = this->dshape_and_dtest_eulerian_at_knot_cons_adv_diff(
87  ipt, psi, dpsidx, test, dtestdx);
88 
89  // Premultiply the weights and the Jacobian
90  double W = w * J;
91 
92  // Calculate local values of the function, initialise to zero
93  double dudt = 0.0;
94  double interpolated_u = 0.0;
95 
96  // These need to be a Vector to be ANSI C++, initialise to zero
97  Vector<double> interpolated_x(DIM, 0.0);
98  Vector<double> interpolated_dudx(DIM, 0.0);
99  Vector<double> mesh_velocity(DIM, 0.0);
100 
101  // Calculate function value and derivatives:
102  //-----------------------------------------
103 
104  // Loop over nodes
105  for (unsigned l = 0; l < n_node; l++)
106  {
107  // Get the value at the node
108  double u_value = this->nodal_value(l, u_nodal_index);
109  interpolated_u += u_value * psi(l);
110  dudt += this->du_dt_cons_adv_diff(l) * psi(l);
111  // Loop over directions
112  for (unsigned j = 0; j < DIM; j++)
113  {
114  interpolated_x[j] += nodal_position(l, j) * psi(l);
115  interpolated_dudx[j] += u_value * dpsidx(l, j);
116  }
117  }
118 
119  // Get the mesh velocity, if required
120  if (!ALE_is_disabled_flag)
121  {
122  for (unsigned l = 0; l < n_node; l++)
123  {
124  // Loop over directions
125  for (unsigned j = 0; j < DIM; j++)
126  {
127  mesh_velocity[j] += dnodal_position_dt(l, j) * psi(l);
128  }
129  }
130  }
131 
132 
133  // Get body force
134  double source;
135  this->get_source_cons_adv_diff(ipt, interpolated_x, source);
136 
137 
138  // Get wind
139  //--------
140  Vector<double> wind(DIM);
141  this->get_wind_cons_adv_diff(ipt, s, interpolated_x, wind);
142 
143  // Get the conserved wind (non-divergence free)
144  Vector<double> conserved_wind(DIM);
145  this->get_conserved_wind_cons_adv_diff(
146  ipt, s, interpolated_x, conserved_wind);
147 
148  // Get diffusivity tensor
149  DenseMatrix<double> D(DIM, DIM);
150  this->get_diff_cons_adv_diff(ipt, s, interpolated_x, D);
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;
161  double hang_weight = 1.0;
162  // Local bool (is the node hanging)
163  bool is_node_hanging = this->node_pt(l)->is_hanging();
164 
165  // If the node is hanging, get the number of master nodes
166  if (is_node_hanging)
167  {
168  hang_info_pt = this->node_pt(l)->hanging_pt();
169  n_master = hang_info_pt->nmaster();
170  }
171  // Otherwise there is just one master node, the node itself
172  else
173  {
174  n_master = 1;
175  }
176 
177  // Loop over the number of master nodes
178  for (unsigned m = 0; m < n_master; m++)
179  {
180  // Get the local equation number and hang_weight
181  // If the node is hanging
182  if (is_node_hanging)
183  {
184  // Read out the local equation from the master node
185  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
186  u_nodal_index);
187  // Read out the weight from the master node
188  hang_weight = hang_info_pt->master_weight(m);
189  }
190  // If the node is not hanging
191  else
192  {
193  // The local equation number comes from the node itself
194  local_eqn = this->nodal_local_eqn(l, u_nodal_index);
195  // The hang weight is one
196  hang_weight = 1.0;
197  }
198 
199  // If the nodal equation is not a boundary conditino
200  if (local_eqn >= 0)
201  {
202  // Add du/dt and body force/source term here
203  residuals[local_eqn] -=
204  (peclet_st * dudt + source) * test(l) * W * hang_weight;
205 
206  // The Mesh velocity and GeneralisedAdvection--Diffusion bit
207  for (unsigned k = 0; k < DIM; k++)
208  {
209  // Terms that multiply the test function
210  double tmp = peclet * wind[k];
211  // If the mesh is moving need to subtract the mesh velocity
212  if (!ALE_is_disabled_flag)
213  {
214  tmp -= peclet_st * mesh_velocity[k];
215  }
216  tmp *= interpolated_dudx[k];
217 
218  // Terms that multiply the derivative of the test function
219  // Advective term
220  double tmp2 = -conserved_wind[k] * interpolated_u;
221  // Now the diuffusive term
222  for (unsigned j = 0; j < DIM; j++)
223  {
224  tmp2 += D(k, j) * interpolated_dudx[j];
225  }
226  // Now construct the contribution to the residuals
227  residuals[local_eqn] -=
228  (tmp * test(l) + tmp2 * dtestdx(l, k)) * W * hang_weight;
229  }
230 
231  // Calculate the Jacobian
232  if (flag)
233  {
234  // Local variables to store the number of master nodes
235  // and the weights associated with each hanging node
236  unsigned n_master2 = 1;
237  double hang_weight2 = 1.0;
238  // Loop over the nodes for the variables
239  for (unsigned l2 = 0; l2 < n_node; l2++)
240  {
241  // Local bool (is the node hanging)
242  bool is_node2_hanging = this->node_pt(l2)->is_hanging();
243  // If the node is hanging, get the number of master nodes
244  if (is_node2_hanging)
245  {
246  hang_info2_pt = this->node_pt(l2)->hanging_pt();
247  n_master2 = hang_info2_pt->nmaster();
248  }
249  // Otherwise there is one master node, the node itself
250  else
251  {
252  n_master2 = 1;
253  }
254 
255  // Loop over the master nodes
256  for (unsigned m2 = 0; m2 < n_master2; m2++)
257  {
258  // Get the local unknown and weight
259  // If the node is hanging
260  if (is_node2_hanging)
261  {
262  // Read out the local unknown from the master node
263  local_unknown = this->local_hang_eqn(
264  hang_info2_pt->master_node_pt(m2), u_nodal_index);
265  // Read out the hanging weight from the master node
266  hang_weight2 = hang_info2_pt->master_weight(m2);
267  }
268  // If the node is not hanging
269  else
270  {
271  // The local unknown number comes from the node
272  local_unknown = this->nodal_local_eqn(l2, u_nodal_index);
273  // The hang weight is one
274  hang_weight2 = 1.0;
275  }
276 
277  // If the unknown is not pinned
278  if (local_unknown >= 0)
279  {
280  // Add contribution to Elemental Matrix
281  // Mass matrix du/dt term
282  jacobian(local_eqn, local_unknown) -=
283  peclet_st * test(l) * psi(l2) *
284  this->node_pt(l2)->time_stepper_pt()->weight(1, 0) * W *
285  hang_weight * hang_weight2;
286 
287  // Add contribution to mass matrix
288  if (flag == 2)
289  {
290  mass_matrix(local_eqn, local_unknown) +=
291  peclet_st * test(l) * psi(l2) * W * hang_weight *
292  hang_weight2;
293  }
294 
295  // Add contribution to Elemental Matrix
296  for (unsigned k = 0; k < DIM; k++)
297  {
298  // Temporary term used in assembly
299  double tmp = peclet * wind[k];
300  if (!ALE_is_disabled_flag)
301  {
302  tmp -= peclet_st * mesh_velocity[k];
303  }
304  tmp *= dpsidx(l2, k);
305 
306  double tmp2 = -conserved_wind[k] * psi(l2);
307  // Now the diffusive term
308  for (unsigned j = 0; j < DIM; j++)
309  {
310  tmp2 += D(k, j) * dpsidx(l2, j);
311  }
312 
313  // Now assemble Jacobian term
314  jacobian(local_eqn, local_unknown) -=
315  (tmp * test(l) + tmp2 * dtestdx(l, k)) * W *
316  hang_weight * hang_weight2;
317  }
318  }
319  } // End of loop over master nodes
320  } // End of loop over nodes
321  } // End of Jacobian calculation
322 
323  } // End of non-zero equation
324 
325  } // End of loop over the master nodes for residual
326  } // End of loop over nodes
327 
328  } // End of loop over integration points
329  }
330 
331 
332  //====================================================================
333  // Force build of templates
334  //====================================================================
338 
342 
343 } // 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
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
void fill_in_generic_residual_contribution_cons_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...
Refineable version of QGeneralisedAdvectionDiffusionElement. Inherit from the standard QGeneralisedAd...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
//////////////////////////////////////////////////////////////////// ////////////////////////////////...