escript  Revision_
finley/src/Assemble.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2016 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 __FINLEY_ASSEMBLE_H__
24 #define __FINLEY_ASSEMBLE_H__
25 
26 #include "Finley.h"
27 #include "ElementFile.h"
28 #include "NodeFile.h"
29 #include <escript/AbstractSystemMatrix.h>
30 
31 namespace finley {
32 
34 {
35  AssembleParameters(const NodeFile* nodes, const ElementFile* ef,
37  bool reducedOrder);
38 
50  int numSides;
52  int numSub;
54  int numDim;
56  int NN;
59 
60  int numEqu;
61  const index_t* row_DOF;
64  const int* row_node;
67  int numComp;
68  const index_t* col_DOF;
71  const int* col_node;
74 };
75 
76 
80 void Assemble_PDE(const NodeFile* nodes, const ElementFile* elements,
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>
88  const escript::Data& d_dirac,
89  const escript::Data& y_dirac);
90 
92  const escript::Data& A, const escript::Data& B,
93  const escript::Data& C, const escript::Data& D,
94  const escript::Data& X, const escript::Data& Y);
95 
96 template<typename Scalar>
98  const escript::Data& A, const escript::Data& B,
99  const escript::Data& C, const escript::Data& D,
100  const escript::Data& X, const escript::Data& Y);
101 
102 template<typename Scalar>
104  const escript::Data& A, const escript::Data& B,
105  const escript::Data& C, const escript::Data& D,
106  const escript::Data& X, const escript::Data& Y);
107 
108 template<typename Scalar>
110  const escript::Data& Y);
111 
113  const escript::Data& A, const escript::Data& B,
114  const escript::Data& C, const escript::Data& D,
115  const escript::Data& X, const escript::Data& Y);
116 
117 template<typename Scalar>
119  const escript::Data& A, const escript::Data& B,
120  const escript::Data& C, const escript::Data& D,
121  const escript::Data& X, const escript::Data& Y);
122 
123 template<typename Scalar>
125  const escript::Data& A, const escript::Data& B,
126  const escript::Data& C, const escript::Data& D,
127  const escript::Data& X, const escript::Data& Y);
128 
129 template<typename Scalar>
131  const escript::Data& Y);
132 
133 template<typename Scalar = double>
135  const index_t* Nodes_Equa, int num_Equa, int NN_Sol,
136  const index_t* Nodes_Sol, int num_Sol, const Scalar* array);
137 
138 void Assemble_LumpedSystem(const NodeFile* nodes, const ElementFile* elements,
139  escript::Data& lumpedMat, const escript::Data& D,
140  bool useHRZ);
141 
144  escript::Data& out, const escript::Data& in);
145 
148  const escript::Data& in);
149 
151 void Assemble_CopyNodalData(const NodeFile* nodes, escript::Data& out,
152  const escript::Data& in);
153 
155 void Assemble_NodeCoordinates(const NodeFile* nodes, escript::Data& x);
156 
158 void Assemble_getNormal(const NodeFile* nodes, const ElementFile* elements,
159  escript::Data& normals);
160 
163 void Assemble_getSize(const NodeFile* nodes, const ElementFile* elements,
164  escript::Data& size);
165 
168 void Assemble_gradient(const NodeFile* nodes, const ElementFile* elements,
169  escript::Data& gradient, const escript::Data& data);
170 
172 void Assemble_integrate(const NodeFile* nodes, const ElementFile* elements,
173  const escript::Data& data, double* integrals);
174 
176 void Assemble_interpolate(const NodeFile* nodes, const ElementFile* elements,
177  const escript::Data& data, escript::Data& output);
178 
179 void Assemble_jacobians_1D(const double* coordinates, int numQuad,
180  const double* QuadWeights, int numShape,
181  dim_t numElements, int numNodes, const index_t* nodes,
182  const double* DSDv, int numTest, const double* DTDv,
183  double* dTdX, double* volume, const index_t* elementId);
184 
185 void Assemble_jacobians_2D(const double* coordinates, int numQuad,
186  const double* QuadWeights, int numShape,
187  dim_t numElements, int numNodes, const index_t* nodes,
188  const double* DSDv, int numTest, const double* DTDv,
189  double* dTdX, double* volume, const index_t* elementId);
190 
191 void Assemble_jacobians_2D_M1D_E1D(const double* coordinates, int numQuad,
192  const double* QuadWeights, int numShape,
193  dim_t numElements, int numNodes, const index_t* nodes,
194  const double* DSDv, int numTest, const double* DTDv,
195  double* dTdX, double* volume, const index_t* elementId);
196 
197 void Assemble_jacobians_2D_M1D_E1D_C(const double* coordinates, int numQuad,
198  const double* QuadWeights, int numShape,
199  dim_t numElements, int numNodes, const index_t* nodes,
200  const double* DSDv, int numTest, const double* DTDv,
201  double* dTdX, double* volume, const index_t* elementId);
202 
203 void Assemble_jacobians_2D_M1D_E2D(const double* coordinates, int numQuad,
204  const double* QuadWeights, int numShape,
205  dim_t numElements, int numNodes, const index_t* nodes,
206  const double* DSDv, int numTest, const double* DTDv,
207  double* dTdX, double* volume, const index_t* elementId);
208 
209 void Assemble_jacobians_2D_M1D_E2D_C(const double* coordinates, int numQuad,
210  const double* QuadWeights, int numShape,
211  dim_t numElements, int numNodes, const index_t* nodes,
212  const double* DSDv, int numTest, const double* DTDv,
213  double* dTdX, double* volume, const index_t* elementId);
214 
215 void Assemble_jacobians_3D(const double* coordinates, int numQuad,
216  const double* QuadWeights, int numShape,
217  dim_t numElements, int numNodes, const index_t* nodes,
218  const double* DSDv, int numTest, const double* DTDv,
219  double* dTdX, double* volume, const index_t* elementId);
220 
221 void Assemble_jacobians_3D_M2D_E2D(const double* coordinates, int numQuad,
222  const double* QuadWeights, int numShape,
223  dim_t numElements, int numNodes, const index_t* nodes,
224  const double* DSDv, int numTest, const double* DTDv,
225  double* dTdX, double* volume, const index_t* elementId);
226 
227 void Assemble_jacobians_3D_M2D_E2D_C(const double* coordinates, int numQuad,
228  const double* QuadWeights, int numShape,
229  dim_t numElements, int numNodes, const index_t* nodes,
230  const double* DSDv, int numTest, const double* DTDv,
231  double* dTdX, double* volume, const index_t* elementId);
232 
233 void Assemble_jacobians_3D_M2D_E3D(const double* coordinates, int numQuad,
234  const double* QuadWeights, int numShape,
235  dim_t numElements, int numNodes, const index_t* nodes,
236  const double* DSDv, int numTest, const double* DTDv,
237  double* dTdX, double* volume, const index_t* elementId);
238 
239 void Assemble_jacobians_3D_M2D_E3D_C(const double* coordinates, int numQuad,
240  const double* QuadWeights, int numShape,
241  dim_t numElements, int numNodes, const index_t* nodes,
242  const double* DSDv, int numTest, const double* DTDv,
243  double* dTdX, double* volume, const index_t* elementId);
244 
245 } // namespace finley
246 
247 #endif // __FINLEY_ASSEMBLE_H__
248 
void Assemble_PDE_Points(const AssembleParameters &p, const escript::Data &d_dirac, const escript::Data &y_dirac)
Definition: finley/src/Assemble_PDE_Points.cpp:43
void Assemble_PDE_Single_1D(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: Assemble_PDE_Single_1D.cpp:45
void Assemble_AverageElementData(const ElementFile *elements, escript::Data &out, const escript::Data &in)
averages data
Definition: finley/src/Assemble_AverageElementData.cpp:31
void Assemble_integrate(const NodeFile *nodes, const ElementFile *elements, const escript::Data &data, double *integrals)
integrates data on quadrature points
Definition: finley/src/Assemble_integrate.cpp:31
int col_numShapes
Definition: finley/src/Assemble.h:73
void Assemble_PDE_System_1D(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: Assemble_PDE_System_1D.cpp:48
void Assemble_LumpedSystem(const NodeFile *nodes, const ElementFile *elements, escript::Data &lumpedMat, const escript::Data &D, bool useHRZ)
Definition: finley/src/Assemble_LumpedSystem.cpp:36
void Assemble_jacobians_2D_M1D_E1D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:139
void Assemble_getSize(const NodeFile *nodes, const ElementFile *elements, escript::Data &size)
Definition: finley/src/Assemble_getSize.cpp:32
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: finley/src/Assemble_PDE_Single_2D.cpp:45
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: finley/src/Assemble_PDE_Single_3D.cpp:45
void Assemble_PDE_System_C(const AssembleParameters &p, const escript::Data &D, const escript::Data &Y)
Definition: Assemble_PDE_System_C.cpp:42
int numDim
number of spatial dimensions
Definition: finley/src/Assemble.h:54
void Assemble_gradient(const NodeFile *nodes, const ElementFile *elements, escript::Data &gradient, const escript::Data &data)
Definition: finley/src/Assemble_gradient.cpp:32
void Assemble_getNormal(const NodeFile *nodes, const ElementFile *elements, escript::Data &normals)
calculates the normal vector at quadrature points on face elements
Definition: finley/src/Assemble_getNormal.cpp:32
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: finley/src/Assemble_PDE_System_2D.cpp:48
Definition: finley/src/ElementFile.h:27
void Assemble_NodeCoordinates(const NodeFile *nodes, escript::Data &x)
copies node coordinates into expanded Data object x
Definition: finley/src/Assemble_NodeCoordinates.cpp:33
void Assemble_jacobians_2D_M1D_E2D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:245
const int * col_node
Definition: finley/src/Assemble.h:71
const index_t * col_DOF
Definition: finley/src/Assemble.h:68
void Assemble_CopyNodalData(const NodeFile *nodes, escript::Data &out, const escript::Data &in)
copies data between different types of nodal representations
Definition: finley/src/Assemble_CopyNodalData.cpp:22
index_t col_DOF_UpperBound
Definition: finley/src/Assemble.h:69
int row_numShapes
Definition: finley/src/Assemble.h:66
int numComp
Definition: finley/src/Assemble.h:67
escript::ASM_ptr S
system matrix to be updated
Definition: finley/src/Assemble.h:42
A suite of factory methods for creating various finley domains.
Definition: finley/src/Assemble.h:31
dim_t numElements
number of elements
Definition: finley/src/Assemble.h:58
void Assemble_CopyElementData(const ElementFile *elements, escript::Data &out, const escript::Data &in)
copies data between different types of elements
Definition: finley/src/Assemble_CopyElementData.cpp:29
void Assemble_jacobians_3D_M2D_E3D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:532
index_t row_DOF_UpperBound
Definition: finley/src/Assemble.h:62
void Assemble_jacobians_1D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:49
ElementFile_Jacobians * col_jac
Definition: finley/src/Assemble.h:70
void Assemble_addToSystemMatrix(escript::ASM_ptr S, int NN_Equa, const index_t *Nodes_Equa, int num_Equa, int NN_Sol, const index_t *Nodes_Sol, int num_Sol, const Scalar *array)
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:59
int row_numShapesTotal
Definition: finley/src/Assemble.h:65
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: finley/src/Assemble_interpolate.cpp:32
Data represents a collection of datapoints.
Definition: Data.h:63
Definition: finley/src/NodeFile.h:40
int numSides
number of sides
Definition: finley/src/Assemble.h:50
int numQuadTotal
total number of quadrature nodes = numQuadSub * numQuadSub
Definition: finley/src/Assemble.h:46
void Assemble_jacobians_3D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:375
const ElementFile * elements
element file these parameters apply to
Definition: finley/src/Assemble.h:40
void Assemble_jacobians_2D_M1D_E1D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:185
Definition: finley/src/Assemble.h:33
int numSub
number of sub-elements
Definition: finley/src/Assemble.h:52
void Assemble_jacobians_3D_M2D_E3D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:451
void Assemble_jacobians_2D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:86
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: finley/src/Assemble_PDE.cpp:85
void Assemble_jacobians_2D_M1D_E2D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:299
ElementFile_Jacobians * row_jac
Definition: finley/src/Assemble.h:63
Data Scalar(double value, const FunctionSpace &what, bool expanded)
A collection of factory functions for creating Data objects which contain data points of various shap...
Definition: DataFactory.cpp:32
void Assemble_jacobians_3D_M2D_E2D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:732
AssembleParameters(const NodeFile *nodes, const ElementFile *ef, escript::ASM_ptr sm, escript::Data &rhs, bool reducedOrder)
Definition: finley/src/Assemble_getAssembleParameters.cpp:33
int col_numShapesTotal
Definition: finley/src/Assemble.h:72
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:32
int numQuadSub
number of quadrature nodes per subelements
Definition: finley/src/Assemble.h:48
const index_t * row_DOF
Definition: finley/src/Assemble.h:61
void Assemble_jacobians_3D_M2D_E2D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, int numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: finley/src/Assemble_jacobians.cpp:665
int NN
leading dimension of element node table
Definition: finley/src/Assemble.h:56
const int * row_node
Definition: finley/src/Assemble.h:64
int numEqu
Definition: finley/src/Assemble.h:60
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: finley/src/Assemble_PDE_System_3D.cpp:48
escript::Data & F
right-hand side to be updated
Definition: finley/src/Assemble.h:44
index_t dim_t
Definition: DataTypes.h:64
Definition: finley/src/ElementFile.h:61
void Assemble_PDE_Single_C(const AssembleParameters &p, const escript::Data &D, const escript::Data &Y)
Definition: Assemble_PDE_Single_C.cpp:41