/[escript]/trunk-mpi-branch/finley/src/ElementFile.h
ViewVC logotype

Contents of /trunk-mpi-branch/finley/src/ElementFile.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1267 - (show annotations)
Tue Aug 21 22:01:21 2007 UTC (12 years, 9 months ago) by gross
File MIME type: text/plain
File size: 6274 byte(s)
check for preperation status removed. is not really needed.
1 /*
2 ************************************************************
3 * Copyright 2006 by ACcESS MNRF *
4 * *
5 * http://www.access.edu.au *
6 * Primary Business: Queensland, Australia *
7 * Licensed under the Open Software License version 3.0 *
8 * http://www.opensource.org/licenses/osl-3.0.php *
9 * *
10 ************************************************************
11 */
12
13 /* Version: $Id$ */
14
15 #ifndef INC_FINLEY_ELEMENTFILE
16 #define INC_FINLEY_ELEMENTFILE
17
18 #include "Finley.h"
19 #include "NodeFile.h"
20 #include "ReferenceElements.h"
21 #include "escript/DataC.h"
22
23 #ifdef PASO_MPI
24 #include "paso/Paso_MPI.h"
25 #endif
26
27 struct Finley_ElementFile_Jacobeans {
28 Finley_Status_t status; /* status of mesh when jacobeans where updated last time */
29 dim_t numDim; /* spatial dimension */
30 Finley_RefElement* ReferenceElement; /* reference elemnt used to calculate jacobeans (this is a borrowd reference) */
31 double* volume; /* local volume */
32 double* DSDX; /* derivatives of shape functions in global coordinates at quadrature points*/
33 };
34
35 typedef struct Finley_ElementFile_Jacobeans Finley_ElementFile_Jacobeans;
36
37 struct Finley_ElementFile {
38 Paso_MPIInfo *MPIInfo;
39 Paso_MPI_rank *Owner;
40
41 Finley_RefElement* ReferenceElement; /* the reference element, see Reference element.c */
42 Finley_RefElement* ReferenceElementReducedOrder; /* the reference element with reduced integration order, see Reference element.c */
43 Finley_RefElement* LinearReferenceElement; /* the reference element for the linear mesh. it is vital that it is using the same quadrature scheme like ReferenceElement*/
44 Finley_RefElement* LinearReferenceElementReducedOrder; /* the reference element for the linear mesh. it is vital that it is using the same quadrature
45 \ scheme like LinearReferenceElementReducedIntegration*/
46
47 dim_t numElements; /* number of elements. */
48
49 index_t *Id; /* Id[i] is the id nmber of
50 node i. this number is not
51 used but useful when
52 elements are resorted. in
53 the entire code the term
54 'element id' refers to i
55 but nor to Id[i] if not
56 explicitly stated
57 otherwise. */
58
59 index_t *Tag; /* Tag[i] is the tag of
60 element i. */
61
62 dim_t numNodes; /* number of nodes per element = ReferenceElement.Type.numNodes */
63 index_t *Nodes; /* Nodes[INDEX(k, i, numNodes)]
64 is the k-the node in the
65 i-the element. note that
66 in the way the nodes are
67 ordered Nodes[INDEX(k, i,
68 LinearReferenceElement.Type.numNodes)
69 is k-the node of element i
70 when refering to the
71 linear version of the
72 mesh. */
73 index_t minColor; /* minimum color */
74 index_t maxColor; /* maximum color */
75 index_t *Color; /* assigns each element a color. elements with the same color */
76 /* are don't share a node so they can be processed simultaneously */
77 /* at anytime Color must provide a valid value. In any case one can set */
78 /* Color[e]=e for all e */
79 index_t order; /* order of the element integration scheme*/
80 index_t reduced_order; /* order of the reduced element integration scheme*/
81
82 Finley_ElementFile_Jacobeans* jacobeans; /* element jacobeans */
83 Finley_ElementFile_Jacobeans* jacobeans_reducedS; /* element jacobeans for reduced order of shape function*/
84 Finley_ElementFile_Jacobeans* jacobeans_reducedQ; /* element jacobeans for reduced integration order*/
85 Finley_ElementFile_Jacobeans* jacobeans_reducedS_reducedQ; /* element jacobeans for reduced integration order and reduced order of shape function*/
86
87 };
88
89 typedef struct Finley_ElementFile Finley_ElementFile;
90 Finley_ElementFile* Finley_ElementFile_alloc( ElementTypeId, index_t, index_t, Paso_MPIInfo* );
91 void Finley_ElementFile_free(Finley_ElementFile*);
92 void Finley_ElementFile_allocTable(Finley_ElementFile*,dim_t);
93 void Finley_ElementFile_freeTable(Finley_ElementFile*);
94 void Finley_ElementFile_setElementDistribution(Finley_ElementFile* in, dim_t* distribution);
95 dim_t Finley_ElementFile_getGlobalNumElements(Finley_ElementFile* in);
96 dim_t Finley_ElementFile_getMyNumElements(Finley_ElementFile* in);
97 index_t Finley_ElementFile_getFirstElement(Finley_ElementFile* in);
98 void Finley_ElementFile_distributeByRankOfDOF(Finley_ElementFile* self, Paso_MPI_rank* mpiRankOfDOF, index_t *Id);
99
100 void Finley_ElementFile_createColoring(Finley_ElementFile* in,dim_t numNodes,dim_t* degreeOfFreedom);
101 void Finley_ElementFile_optimizeOrdering(Finley_ElementFile** in);
102 void Finley_ElementFile_setNodeRange(dim_t*,dim_t*,Finley_ElementFile*);
103 void Finley_ElementFile_relableNodes(dim_t*,dim_t,Finley_ElementFile*);
104 void Finley_ElementFile_markNodes(dim_t*,dim_t,Finley_ElementFile*,dim_t);
105 void Finley_ElementFile_scatter(dim_t*,Finley_ElementFile*,Finley_ElementFile*);
106 void Finley_ElementFile_gather(dim_t*,Finley_ElementFile*,Finley_ElementFile*);
107 void Finley_ElementFile_copyTable(dim_t,Finley_ElementFile*,dim_t,dim_t,Finley_ElementFile*);
108 void Finley_ElementFile_markDOFsConnectedToRange(index_t* mask,index_t offset,index_t marker,index_t firstDOF,index_t lastDOF,index_t *dofIndex,Finley_ElementFile*in ,bool_t useLinear);
109
110 void Finley_ElementFile_setTags(Finley_ElementFile*,const int,escriptDataC*);
111 Finley_ElementFile_Jacobeans* Finley_ElementFile_Jacobeans_alloc(Finley_RefElement*);
112 void Finley_ElementFile_Jacobeans_dealloc(Finley_ElementFile_Jacobeans*);
113 Finley_ElementFile_Jacobeans* Finley_ElementFile_borrowJacobeans(Finley_ElementFile*, Finley_NodeFile*, bool_t, bool_t);
114
115
116 #endif /* #ifndef INC_FINLEY_ELEMENTFILE */
117

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26