24 #ifndef __FINLEY_ASSEMBLE_H__ 25 #define __FINLEY_ASSEMBLE_H__ 30 #include "paso/SystemMatrix.h" 126 const index_t* Nodes_Equa,
const int num_Equa,
const int NN_Sol,
127 const index_t* Nodes_Sol,
const int num_Sol,
const double* array);
160 const double* QuadWeights,
int numShape,
162 const double* DSDv,
int numTest,
const double* DTDv,
163 double* dTdX,
double* volume,
const index_t* elementId);
165 const double* QuadWeights,
int numShape,
167 const double* DSDv,
int numTest,
const double* DTDv,
168 double* dTdX,
double* volume,
const index_t* elementId);
170 const double* QuadWeights,
int numShape,
172 const double* DSDv,
int numTest,
const double* DTDv,
173 double* dTdX,
double* volume,
const index_t* elementId);
175 const double* QuadWeights,
int numShape,
177 const double* DSDv,
int numTest,
const double* DTDv,
178 double* dTdX,
double* volume,
const index_t* elementId);
180 const double* QuadWeights,
int numShape,
182 const double* DSDv,
int numTest,
const double* DTDv,
183 double* dTdX,
double* volume,
const index_t* elementId);
185 const double* QuadWeights,
int numShape,
187 const double* DSDv,
int numTest,
const double* DTDv,
188 double* dTdX,
double* volume,
const index_t* elementId);
190 const double* QuadWeights,
int numShape,
192 const double* DSDv,
int numTest,
const double* DTDv,
193 double* dTdX,
double* volume,
const index_t* elementId);
195 const double* QuadWeights,
int numShape,
197 const double* DSDv,
int numTest,
const double* DTDv,
198 double* dTdX,
double* volume,
const index_t* elementId);
200 const double* QuadWeights,
int numShape,
202 const double* DSDv,
int numTest,
const double* DTDv,
203 double* dTdX,
double* volume,
const index_t* elementId);
205 const double* QuadWeights,
int numShape,
207 const double* DSDv,
int numTest,
const double* DTDv,
208 double* dTdX,
double* volume,
const index_t* elementId);
210 const double* QuadWeights,
int numShape,
212 const double* DSDv,
int numTest,
const double* DTDv,
213 double* dTdX,
double* volume,
const index_t* elementId);
217 #endif // __FINLEY_ASSEMBLE_H__ 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:46
void Assemble_PDE_Points(const AssembleParameters &p, const escript::Data &d_dirac, const escript::Data &y_dirac)
Definition: finley/src/Assemble_PDE_Points.cpp:44
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: Assemble_PDE_System_3D.cpp:49
void Assemble_jacobians_2D_M1D_E2D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:248
void Assemble_PDE_Single_C(const AssembleParameters &p, const escript::Data &D, const escript::Data &Y)
Definition: Assemble_PDE_Single_C.cpp:42
void Assemble_AverageElementData(const ElementFile *elements, escript::Data &out, const escript::Data &in)
Definition: finley/src/Assemble_AverageElementData.cpp:33
void Assemble_integrate(const NodeFile *nodes, const ElementFile *elements, const escript::Data &data, double *integrals)
Definition: finley/src/Assemble_integrate.cpp:33
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:50
void Assemble_LumpedSystem(const NodeFile *nodes, const ElementFile *elements, escript::Data &lumpedMat, const escript::Data &D, bool useHRZ)
Definition: finley/src/Assemble_LumpedSystem.cpp:38
void Assemble_getSize(const NodeFile *nodes, const ElementFile *elements, escript::Data &size)
Definition: finley/src/Assemble_getSize.cpp:32
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:33
void Assemble_getNormal(const NodeFile *nodes, const ElementFile *elements, escript::Data &normals)
Definition: Assemble_getNormal.cpp:33
void Assemble_addToSystemMatrix(paso::SystemMatrix_ptr S, const int NN_Equa, const index_t *Nodes_Equa, const int num_Equa, const int NN_Sol, const index_t *Nodes_Sol, const int num_Sol, const double *array)
Definition: finley/src/Assemble_addToSystemMatrix.cpp:56
void Assemble_jacobians_2D_M1D_E2D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:302
Definition: finley/src/ElementFile.h:27
void Assemble_jacobians_1D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:51
void Assemble_jacobians_3D_M2D_E2D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:735
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)
Definition: finley/src/Assemble_CopyNodalData.cpp:34
void Assemble_NodeCoordinates(const NodeFile *nodes, escript::Data &out)
Definition: finley/src/Assemble_NodeCoordinates.cpp:35
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
AssembleParameters(const NodeFile *nodes, const ElementFile *ef, paso::SystemMatrix_ptr sm, escript::Data &rhs, bool reducedOrder)
Definition: finley/src/Assemble_getAssembleParameters.cpp:31
Definition: finley/src/Assemble.h:32
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)
Definition: finley/src/Assemble_CopyElementData.cpp:33
index_t row_DOF_UpperBound
Definition: finley/src/Assemble.h:62
boost::shared_ptr< SystemMatrix > SystemMatrix_ptr
Definition: SystemMatrix.h:38
ElementFile_Jacobians * col_jac
Definition: finley/src/Assemble.h:70
paso::SystemMatrix_ptr S
system matrix to be updated
Definition: finley/src/Assemble.h:42
void Assemble_jacobians_3D_M2D_E3D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:535
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: Assemble_PDE_Single_2D.cpp:46
void Assemble_PDE_System_C(const AssembleParameters &p, const escript::Data &D, const escript::Data &Y)
Definition: Assemble_PDE_System_C.cpp:46
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: Assemble_PDE_Single_3D.cpp:46
void Assemble_jacobians_2D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:88
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)
Definition: finley/src/Assemble_interpolate.cpp:34
void Assemble_PDE(const NodeFile *nodes, const ElementFile *elements, paso::SystemMatrix_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:88
Data represents a collection of datapoints.
Definition: Data.h:68
Definition: finley/src/NodeFile.h:30
int numSides
number of sides
Definition: finley/src/Assemble.h:50
void Assemble_jacobians_2D_M1D_E1D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:142
int numQuadTotal
total number of quadrature nodes = numQuadSub * numQuadSub
Definition: finley/src/Assemble.h:46
const ElementFile * elements
element file these parameters apply to
Definition: finley/src/Assemble.h:40
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: Assemble_PDE_System_2D.cpp:49
Definition: finley/src/Assemble.h:34
int numSub
number of sub-elements
Definition: finley/src/Assemble.h:52
int index_t
Definition: types.h:24
ElementFile_Jacobians * row_jac
Definition: finley/src/Assemble.h:63
void Assemble_jacobians_3D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:378
index_t dim_t
Definition: types.h:27
void Assemble_jacobians_3D_M2D_E3D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:454
int col_numShapesTotal
Definition: finley/src/Assemble.h:72
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
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_jacobians_2D_M1D_E1D_C(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:188
escript::Data & F
right-hand side to be updated
Definition: finley/src/Assemble.h:44
Definition: finley/src/ElementFile.h:60
void Assemble_jacobians_3D_M2D_E2D(const double *coordinates, int numQuad, const double *QuadWeights, int numShape, dim_t numElements, dim_t numNodes, const index_t *nodes, const double *DSDv, int numTest, const double *DTDv, double *dTdX, double *volume, const index_t *elementId)
Definition: Assemble_jacobians.cpp:668