escript  Revision_
dudley/src/Assemble.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2018 by The University of Queensland
5 * http://www.uq.edu.au
6 *
7 * Primary Business: Queensland, Australia
8 * Licensed under the Apache License, version 2.0
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 * Development 2012-2013 by School of Earth Sciences
13 * Development from 2014 by Centre for Geoscience Computing (GeoComp)
14 *
15 *****************************************************************************/
16 
17 /****************************************************************************
18 
19  Assemblage routines: header file
20 
21 *****************************************************************************/
22 
23 #ifndef __DUDLEY_ASSEMBLE_H__
24 #define __DUDLEY_ASSEMBLE_H__
25 
26 #include "Dudley.h"
27 #include "ElementFile.h"
28 #include "NodeFile.h"
29 #include <escript/AbstractSystemMatrix.h>
30 
31 namespace dudley {
32 
33 struct AssembleParameters
34 {
35  AssembleParameters(const NodeFile* nodes, const ElementFile* ef,
37  bool reducedOrder);
38 
40  const ElementFile* elements;
46  int numQuad;
48  int numDim;
50  int NN;
52  int numEqu;
54  const index_t* DOF;
59  int numShapes;
60  const double* shapeFns;
61 };
62 
63 void Assemble_PDE(const NodeFile* nodes, const ElementFile* elements,
65  const escript::Data& A, const escript::Data& B,
66  const escript::Data& C, const escript::Data& D,
67  const escript::Data& X, const escript::Data& Y);
68 
69 template<typename Scalar = double>
71  const escript::Data& d_dirac,
72  const escript::Data& y_dirac);
73 
74 template<typename Scalar = double>
76  const escript::Data& A, const escript::Data& B,
77  const escript::Data& C, const escript::Data& D,
78  const escript::Data& X, const escript::Data& Y);
79 
80 template<typename Scalar = double>
82  const escript::Data& A, const escript::Data& B,
83  const escript::Data& C, const escript::Data& D,
84  const escript::Data& X, const escript::Data& Y);
85 
86 template<typename Scalar = double>
88  const escript::Data& A, const escript::Data& B,
89  const escript::Data& C, const escript::Data& D,
90  const escript::Data& X, const escript::Data& Y);
91 
92 template<typename Scalar = double>
94  const escript::Data& A, const escript::Data& B,
95  const escript::Data& C, const escript::Data& D,
96  const escript::Data& X, const escript::Data& Y);
97 
98 
104 template<typename Scalar>
106  const std::vector<index_t>& Nodes, int numEq,
107  const std::vector<Scalar>& array);
108 
112 void Assemble_LumpedSystem(const NodeFile* nodes, const ElementFile* elements,
113  escript::Data& lumpedMat, const escript::Data& D,
114  bool useHRZ);
115 
117 template<typename Scalar>
118 void Assemble_AverageElementData(const ElementFile* elements,
119  escript::Data& out, const escript::Data& in);
120 
122 template<typename Scalar>
123 void Assemble_CopyElementData(const ElementFile* elements, escript::Data& out,
124  const escript::Data& in);
125 
127 template<typename Scalar>
128 void Assemble_CopyNodalData(const NodeFile* nodes, escript::Data& out,
129  const escript::Data& in);
130 
132 void Assemble_NodeCoordinates(const NodeFile* nodes, escript::Data& x);
133 
135 void Assemble_getNormal(const NodeFile* nodes, const ElementFile* elements,
136  escript::Data& normals);
137 
140 void Assemble_getSize(const NodeFile* nodes, const ElementFile* elements,
141  escript::Data& size);
142 
145 template<typename Scalar>
146 void Assemble_gradient(const NodeFile* nodes, const ElementFile* elements,
147  escript::Data& gradient, const escript::Data& data);
148 
150 template<typename Scalar>
151 void Assemble_integrate(const NodeFile* nodes, const ElementFile* elements,
152  const escript::Data& data, std::vector<Scalar>& integrals);
153 
155 template<typename Scalar>
156 void Assemble_interpolate(const NodeFile* nodes, const ElementFile* elements,
157  const escript::Data& data, escript::Data& output);
158 
159 void Assemble_jacobians_2D(const double* coordinates, int numQuad,
160  dim_t numElements, int numNodes,
161  const index_t* nodes, double* dTdX, double* absD,
162  double* quadWeight, const index_t* elementId);
163 
164 void Assemble_jacobians_2D_M1D_E1D(const double* coordinates, int numQuad,
165  dim_t numElements, int numNodes,
166  const index_t* nodes, double* dTdX, double* absD,
167  double* quadWeight, const index_t* elementId);
168 
169 void Assemble_jacobians_3D(const double* coordinates, int numQuad,
170  dim_t numElements, int numNodes,
171  const index_t* nodes, double* dTdX, double* abs_D,
172  double* quadWeight, const index_t* elementId);
173 
174 void Assemble_jacobians_3D_M2D_E2D(const double* coordinates, int numQuad,
175  dim_t numElements, int numNodes,
176  const index_t* nodes, double* dTdX, double* absD,
177  double* quadWeight, const index_t* elementId);
178 
179 
180 } // namespace dudley
181 
182 #endif // __DUDLEY_ASSEMBLE_H__
183 
dudley::ElementFile_Jacobians
Definition: dudley/src/ElementFile.h:38
dudley::Assemble_jacobians_3D_M2D_E2D
void Assemble_jacobians_3D_M2D_E2D(const double *coordinates, int numQuad, dim_t numElements, int numNodes, const index_t *nodes, double *dTdX, double *absD, double *quadWeight, const index_t *elementId)
Definition: dudley/src/Assemble_jacobians.cpp:253
dudley::Assemble_LumpedSystem
void Assemble_LumpedSystem(const NodeFile *nodes, const ElementFile *elements, escript::Data &lumpedMat, const escript::Data &D, bool useHRZ)
Definition: dudley/src/Assemble_LumpedSystem.cpp:36
dudley::AssembleParameters::DOF_UpperBound
dim_t DOF_UpperBound
number of local degrees of freedom
Definition: dudley/src/Assemble.h:86
dudley::Assemble_getNormal
void Assemble_getNormal(const NodeFile *nodes, const ElementFile *elements, escript::Data &normals)
calculates the normal vector at quadrature points on face elements
Definition: dudley/src/Assemble_getNormal.cpp:36
dudley::NodeFile
Definition: dudley/src/NodeFile.h:37
dudley::Nodes
Definition: Dudley.h:63
dudley::Assemble_CopyElementData
void Assemble_CopyElementData(const ElementFile *elements, escript::Data &out, const escript::Data &in)
copies data between different types of elements
Definition: dudley/src/Assemble_CopyElementData.cpp:37
dudley::AssembleParameters::numShapes
int numShapes
Definition: dudley/src/Assemble.h:89
dudley::Assemble_interpolate
void Assemble_interpolate(const NodeFile *nodes, const ElementFile *elements, const escript::Data &data, escript::Data &output)
interpolates nodal data in a data array onto elements (=integration points)
Definition: dudley/src/Assemble_interpolate.cpp:37
dudley::AssembleParameters::S
escript::AbstractSystemMatrix * S
system matrix to be updated
Definition: dudley/src/Assemble.h:72
dudley::Assemble_AverageElementData
void Assemble_AverageElementData(const ElementFile *elements, escript::Data &out, const escript::Data &in)
averages data
Definition: dudley/src/Assemble_AverageElementData.cpp:40
dudley::Assemble_PDE_Single_3D
void Assemble_PDE_Single_3D(const AssembleParameters &p, const escript::Data &A, const escript::Data &B, const escript::Data &C, const escript::Data &D, const escript::Data &X, const escript::Data &Y)
Definition: dudley/src/Assemble_PDE_Single_3D.cpp:72
dudley::AssembleParameters
Definition: dudley/src/Assemble.h:48
dudley::Assemble_NodeCoordinates
void Assemble_NodeCoordinates(const NodeFile *nodes, escript::Data &x)
copies node coordinates into expanded Data object x
Definition: dudley/src/Assemble_NodeCoordinates.cpp:37
dudley::Assemble_addToSystemMatrix
void Assemble_addToSystemMatrix(escript::AbstractSystemMatrix *S, const std::vector< index_t > &Nodes, int numEq, const std::vector< Scalar > &array)
Dudley.h
dudley::AssembleParameters::elements
const ElementFile * elements
element file these parameters apply to
Definition: dudley/src/Assemble.h:70
dudley::Assemble_PDE_System_3D
void Assemble_PDE_System_3D(const AssembleParameters &p, const escript::Data &A, const escript::Data &B, const escript::Data &C, const escript::Data &D, const escript::Data &X, const escript::Data &Y)
Definition: dudley/src/Assemble_PDE_System_3D.cpp:78
dudley
A suite of factory methods for creating 2D and 3D dudley domains.
Definition: dudley/src/Assemble.h:31
dudley::Assemble_getSize
void Assemble_getSize(const NodeFile *nodes, const ElementFile *elements, escript::Data &size)
Definition: dudley/src/Assemble_getSize.cpp:35
escript::Data
Data represents a collection of datapoints.
Definition: Data.h:62
dudley::Assemble_jacobians_2D
void Assemble_jacobians_2D(const double *coordinates, int numQuad, dim_t numElements, int numNodes, const index_t *nodes, double *dTdX, double *absD, double *quadWeight, const index_t *elementId)
Definition: dudley/src/Assemble_jacobians.cpp:54
escript::DataTypes::dim_t
index_t dim_t
Definition: DataTypes.h:87
dudley::Assemble_jacobians_3D
void Assemble_jacobians_3D(const double *coordinates, int numQuad, dim_t numElements, int numNodes, const index_t *nodes, double *dTdX, double *abs_D, double *quadWeight, const index_t *elementId)
Definition: dudley/src/Assemble_jacobians.cpp:178
escript::AbstractSystemMatrix
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:53
dudley::ElementFile
Definition: dudley/src/ElementFile.h:62
ElementFile.h
dudley::AssembleParameters::NN
int NN
leading dimension of element node table
Definition: dudley/src/Assemble.h:80
dudley::Assemble_PDE_System_2D
void Assemble_PDE_System_2D(const AssembleParameters &p, const escript::Data &A, const escript::Data &B, const escript::Data &C, const escript::Data &D, const escript::Data &X, const escript::Data &Y)
Definition: dudley/src/Assemble_PDE_System_2D.cpp:78
dudley::AssembleParameters::DOF
const index_t * DOF
row and column degrees of freedom
Definition: dudley/src/Assemble.h:84
dudley::Assemble_CopyNodalData
void Assemble_CopyNodalData(const NodeFile *nodes, escript::Data &out, const escript::Data &in)
copies data between different types of nodal representations
Definition: dudley/src/Assemble_CopyNodalData.cpp:37
NodeFile.h
dudley::Assemble_jacobians_2D_M1D_E1D
void Assemble_jacobians_2D_M1D_E1D(const double *coordinates, int numQuad, dim_t numElements, int numNodes, const index_t *nodes, double *dTdX, double *absD, double *quadWeight, const index_t *elementId)
Definition: dudley/src/Assemble_jacobians.cpp:132
dudley::Assemble_integrate
void Assemble_integrate(const NodeFile *nodes, const ElementFile *elements, const escript::Data &data, std::vector< Scalar > &integrals)
integrates data on quadrature points
Definition: dudley/src/Assemble_integrate.cpp:39
dudley::AssembleParameters::AssembleParameters
AssembleParameters(const NodeFile *nodes, const ElementFile *ef, escript::ASM_ptr sm, escript::Data &rhs, bool reducedOrder)
Definition: dudley/src/Assemble_getAssembleParameters.cpp:32
dudley::Assemble_PDE_Points
void Assemble_PDE_Points(const AssembleParameters &p, const escript::Data &d_dirac, const escript::Data &y_dirac)
Definition: dudley/src/Assemble_PDE_Points.cpp:71
dudley::AssembleParameters::numDim
int numDim
number of spatial dimensions
Definition: dudley/src/Assemble.h:78
escript::DataTypes::index_t
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:82
dudley::AssembleParameters::numQuad
int numQuad
number of quadrature nodes
Definition: dudley/src/Assemble.h:76
dudley::AssembleParameters::jac
const ElementFile_Jacobians * jac
reference to jacobians
Definition: dudley/src/Assemble.h:88
dudley::Assemble_gradient
void Assemble_gradient(const NodeFile *nodes, const ElementFile *elements, escript::Data &gradient, const escript::Data &data)
Definition: dudley/src/Assemble_gradient.cpp:39
dudley::AssembleParameters::F
escript::Data & F
right-hand side to be updated
Definition: dudley/src/Assemble.h:74
dudley::Assemble_PDE
void Assemble_PDE(const NodeFile *nodes, const ElementFile *elements, escript::ASM_ptr S, escript::Data &F, const escript::Data &A, const escript::Data &B, const escript::Data &C, const escript::Data &D, const escript::Data &X, const escript::Data &Y)
Definition: dudley/src/Assemble_PDE.cpp:120
S
#define S(_J_, _I_)
Definition: ShapeFunctions.cpp:134
escript::ASM_ptr
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:43
dudley::AssembleParameters::shapeFns
const double * shapeFns
Definition: dudley/src/Assemble.h:90
dudley::Assemble_PDE_Single_2D
void Assemble_PDE_Single_2D(const AssembleParameters &p, const escript::Data &A, const escript::Data &B, const escript::Data &C, const escript::Data &D, const escript::Data &X, const escript::Data &Y)
Definition: dudley/src/Assemble_PDE_Single_2D.cpp:72
dudley::AssembleParameters::numEqu
int numEqu
number of equations (= matrix row/column block size)
Definition: dudley/src/Assemble.h:82