biharmonic_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-2025 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//====================================================================
26// Config header
27#ifdef HAVE_CONFIG_H
28#include <oomph-lib-config.h>
29#endif
30
31// oomph-lib includes
32#include "biharmonic_elements.h"
33
34
35namespace oomph
36{
37 /// Compute element residual Vector only (if JFLAG=and/or element
38 /// Jacobian matrix
39 template<unsigned DIM>
43 {
44 // Find out how many nodes there are
45 unsigned n_node = this->nnode();
46
47 // Find out how many values there
48 unsigned n_value = this->node_pt(0)->nvalue();
49
50 // Find out how many degrees of freedom are associated with 2nd derivative
51 unsigned d2_dof = this->get_d2_dof();
52
53 // set up memory for shape and test functions
56 DShape d2psidx(n_node, n_value, d2_dof);
57 // Shape psi(n_node,n_value);
59 DShape d2psids(n_node, n_value, d2_dof);
60
61 // storage of single local coordinate
63
64 // number of integration points
65 unsigned n_intpt = this->integral_pt()->nweight();
66
67 // loop over integration points
68 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
69 {
70 // set local coordinates to be intergration scheme knot
71 for (unsigned i = 0; i < DIM; i++)
72 s[i] = this->integral_pt()->knot(ipt, i);
73
74 // find weight at knot
75 double w = this->integral_pt()->weight(ipt);
76
77 // find shape fn and derivate fns at knot point,
78 // and Jacobian of local to global mapping
79 double J = this->d2shape_eulerian(s, psi, dpsidx, d2psidx);
80 this->d2shape_local(s, psi, dpsids, d2psids);
81
82 // premultiply weight by jacobian
83 double W = w * J;
84
85 // initialise storage for d2u_interpolated
87 for (unsigned i = 0; i < DIM; i++)
88 {
89 d2u_interpolated[i] = 0.0;
90 }
91
92 // loop over nodes, degrees of freedom, and dimension to calculate
93 // interpolated_d2u
94 for (unsigned n = 0; n < n_node; n++)
95 {
96 for (unsigned k = 0; k < n_value; k++)
97 {
98 for (unsigned d = 0; d < DIM; d++)
99 {
100 d2u_interpolated[d] +=
101 this->node_pt(n)->value(k) * d2psidx(n, k, d);
102 }
103 }
104 }
105
106 // create vector for interpolated position and calculate
109
110 // evaluate source function at knot position
111 double source = 0.0;
112 get_source(ipt, interpolated_position, source);
113
114 // loop over nodes
115 for (unsigned n1 = 0; n1 < n_node; n1++)
116 {
117 // loop over types of degrees of freedom
118 for (unsigned k1 = 0; k1 < n_value; k1++)
119 {
120 // determine the local equation number
122
123 // if local equation number equal to -1 then its a boundary
124 // node(pinned)
125 if (local_eqn_number >= 0)
126 {
127 // add source contribution to residual
128 residuals[local_eqn_number] -= source * psi(n1, k1) * W;
129
130 // loop over dimension of d2u interpolated
131 for (unsigned d1 = 0; d1 < DIM; d1++)
132 {
133 // loop over derivatives of d2psidx
134 for (unsigned d2 = 0; d2 < DIM; d2++)
135 {
136 // add residual contribution
138 d2u_interpolated[d1] * d2psidx(n1, k1, d2) * W;
139 }
140
141 // calculate the jacobian
142 if (JFLAG)
143 {
144 // loop over nodes
145 for (unsigned n2 = 0; n2 < n_node; n2++)
146 {
147 // loop over types of degrees of freedom
148 for (unsigned k2 = 0; k2 < n_value; k2++)
149 {
150 // get local dof, -if -1 then node pinned
152
153 // if local dof # = -1, then its pinned
154 if (local_dof_number >= 0)
155 {
156 // loop over derivatives
157 for (unsigned d2 = 0; d2 < DIM; d2++)
158 {
159 // add contribution to jacobian
161 d2psidx(n1, k1, d1) * d2psidx(n2, k2, d2) * W;
162 }
163 }
164 }
165 }
166 }
167 }
168 }
169 }
170 }
171 }
172 }
173
174 template class BiharmonicElement<2>;
175 template class BiharmonicElement<1>;
176
177} // namespace oomph
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
virtual void fill_in_generic_residual_contribution_biharmonic(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned JFLAG)
Compute element residual Vector only (if JFLAG=and/or element Jacobian matrix.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition shape.h:278
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition nodes.h:483
double d2shape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
Compute the geometric shape functions and also first and second derivatives w.r.t....
Definition elements.cc:3478
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
virtual void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Function to compute the geometric shape functions and also first and second derivatives w....
Definition elements.h:2020
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
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
int local_eqn_number(const unsigned long &ieqn_global) const
Return the local equation number corresponding to the ieqn_global-th global equation number....
Definition elements.h:713
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.
double value(const unsigned &i) const
Return i-th value (dofs or pinned) at this node either directly or via hanging node representation....
Definition nodes.cc:2408
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).