escript  Revision_
dudley/src/Assemble.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2020 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-2017 by Centre for Geoscience Computing (GeoComp)
14 * Development from 2019 by School of Earth and Environmental Sciences
15 **
16 *****************************************************************************/
17 
18 /****************************************************************************
19 
20  Assemblage routines: header file
21 
22 *****************************************************************************/
23 
24 #ifndef __DUDLEY_ASSEMBLE_H__
25 #define __DUDLEY_ASSEMBLE_H__
26 
27 #include "Dudley.h"
28 #include "ElementFile.h"
29 #include "NodeFile.h"
30 #include <escript/AbstractSystemMatrix.h>
31 
32 namespace dudley {
33 
35 {
36  AssembleParameters(const NodeFile* nodes, const ElementFile* ef,
38  bool reducedOrder);
39 
47  int numQuad;
49  int numDim;
51  int NN;
53  int numEqu;
55  const index_t* DOF;
60  int numShapes;
61  const double* shapeFns;
62 };
63 
64 void Assemble_PDE(const NodeFile* nodes, const ElementFile* elements,
66  const escript::Data& A, const escript::Data& B,
67  const escript::Data& C, const escript::Data& D,
68  const escript::Data& X, const escript::Data& Y);
69 
70 template<typename Scalar = double>
72  const escript::Data& d_dirac,
73  const escript::Data& y_dirac);
74 
75 template<typename Scalar = double>
77  const escript::Data& A, const escript::Data& B,
78  const escript::Data& C, const escript::Data& D,
79  const escript::Data& X, const escript::Data& Y);
80 
81 template<typename Scalar = double>
83  const escript::Data& A, const escript::Data& B,
84  const escript::Data& C, const escript::Data& D,
85  const escript::Data& X, const escript::Data& Y);
86 
87 template<typename Scalar = double>
89  const escript::Data& A, const escript::Data& B,
90  const escript::Data& C, const escript::Data& D,
91  const escript::Data& X, const escript::Data& Y);
92 
93 template<typename Scalar = double>
95  const escript::Data& A, const escript::Data& B,
96  const escript::Data& C, const escript::Data& D,
97  const escript::Data& X, const escript::Data& Y);
98 
99 
105 template<typename Scalar>
107  const std::vector<index_t>& Nodes, int numEq,
108  const std::vector<Scalar>& array);
109 
113 void Assemble_LumpedSystem(const NodeFile* nodes, const ElementFile* elements,
114  escript::Data& lumpedMat, const escript::Data& D,
115  bool useHRZ);
116 
118 template<typename Scalar>
119 void Assemble_AverageElementData(const ElementFile* elements,
120  escript::Data& out, const escript::Data& in);
121 
123 template<typename Scalar>
124 void Assemble_CopyElementData(const ElementFile* elements, escript::Data& out,
125  const escript::Data& in);
126 
128 template<typename Scalar>
129 void Assemble_CopyNodalData(const NodeFile* nodes, escript::Data& out,
130  const escript::Data& in);
131 
133 void Assemble_NodeCoordinates(const NodeFile* nodes, escript::Data& x);
134 
136 void Assemble_getNormal(const NodeFile* nodes, const ElementFile* elements,
137  escript::Data& normals);
138 
141 void Assemble_getSize(const NodeFile* nodes, const ElementFile* elements,
142  escript::Data& size);
143 
146 template<typename Scalar>
147 void Assemble_gradient(const NodeFile* nodes, const ElementFile* elements,
148  escript::Data& gradient, const escript::Data& data);
149 
151 template<typename Scalar>
152 void Assemble_integrate(const NodeFile* nodes, const ElementFile* elements,
153  const escript::Data& data, std::vector<Scalar>& integrals);
154 
156 template<typename Scalar>
157 void Assemble_interpolate(const NodeFile* nodes, const ElementFile* elements,
158  const escript::Data& data, escript::Data& output);
159 
160 void Assemble_jacobians_2D(const double* coordinates, int numQuad,
161  dim_t numElements, int numNodes,
162  const index_t* nodes, double* dTdX, double* absD,
163  double* quadWeight, const index_t* elementId);
164 
165 void Assemble_jacobians_2D_M1D_E1D(const double* coordinates, int numQuad,
166  dim_t numElements, int numNodes,
167  const index_t* nodes, double* dTdX, double* absD,
168  double* quadWeight, const index_t* elementId);
169 
170 void Assemble_jacobians_3D(const double* coordinates, int numQuad,
171  dim_t numElements, int numNodes,
172  const index_t* nodes, double* dTdX, double* abs_D,
173  double* quadWeight, const index_t* elementId);
174 
175 void Assemble_jacobians_3D_M2D_E2D(const double* coordinates, int numQuad,
176  dim_t numElements, int numNodes,
177  const index_t* nodes, double* dTdX, double* absD,
178  double* quadWeight, const index_t* elementId);
179 
180 
181 } // namespace dudley
182 
183 #endif // __DUDLEY_ASSEMBLE_H__
184 
dudley::ElementFile_Jacobians
Definition: dudley/src/ElementFile.h:29
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:255
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:26
dudley::AssembleParameters::DOF_UpperBound
dim_t DOF_UpperBound
number of local degrees of freedom
Definition: dudley/src/Assemble.h:57
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:26
dudley::NodeFile
Definition: dudley/src/NodeFile.h:40
dudley::Nodes
@ Nodes
Definition: Dudley.h:53
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:27
dudley::AssembleParameters::numShapes
int numShapes
Definition: dudley/src/Assemble.h:60
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:27
NodeFile.h
dudley::AssembleParameters::S
escript::AbstractSystemMatrix * S
system matrix to be updated
Definition: dudley/src/Assemble.h:43
dudley::Assemble_AverageElementData
void Assemble_AverageElementData(const ElementFile *elements, escript::Data &out, const escript::Data &in)
averages data
Definition: dudley/src/Assemble_AverageElementData.cpp:29
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:46
dudley::AssembleParameters
Definition: dudley/src/Assemble.h:35
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:27
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:41
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:49
dudley
A suite of factory methods for creating 2D and 3D dudley domains.
Definition: dudley/src/Assemble.h:32
dudley::Assemble_getSize
void Assemble_getSize(const NodeFile *nodes, const ElementFile *elements, escript::Data &size)
Definition: dudley/src/Assemble_getSize.cpp:25
escript::Data
Data represents a collection of datapoints.
Definition: Data.h:64
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:56
escript::DataTypes::dim_t
index_t dim_t
Definition: DataTypes.h:65
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:180
escript::AbstractSystemMatrix
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:44
dudley::ElementFile
Definition: dudley/src/ElementFile.h:53
ElementFile.h
dudley::AssembleParameters::NN
int NN
leading dimension of element node table
Definition: dudley/src/Assemble.h:51
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:49
dudley::AssembleParameters::DOF
const index_t * DOF
row and column degrees of freedom
Definition: dudley/src/Assemble.h:55
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:27
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:134
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:29
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:35
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:46
dudley::AssembleParameters::numDim
int numDim
number of spatial dimensions
Definition: dudley/src/Assemble.h:49
escript::DataTypes::index_t
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:60
dudley::AssembleParameters::numQuad
int numQuad
number of quadrature nodes
Definition: dudley/src/Assemble.h:47
dudley::AssembleParameters::jac
const ElementFile_Jacobians * jac
reference to jacobians
Definition: dudley/src/Assemble.h:59
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:29
dudley::AssembleParameters::F
escript::Data & F
right-hand side to be updated
Definition: dudley/src/Assemble.h:45
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:84
S
#define S(_J_, _I_)
Definition: ShapeFunctions.cpp:122
escript::ASM_ptr
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:33
dudley::AssembleParameters::shapeFns
const double * shapeFns
Definition: dudley/src/Assemble.h:61
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:46
dudley::AssembleParameters::numEqu
int numEqu
number of equations (= matrix row/column block size)
Definition: dudley/src/Assemble.h:53