/[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 1223 - (show annotations)
Fri Aug 3 02:40:39 2007 UTC (11 years, 8 months ago) by gross
File MIME type: text/plain
File size: 6225 byte(s)
first attemt towards an improved MPI version.  

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 index_t isPrepared; /* UNKNOWN, UNPREPARED, PREPARED to indicate that the element table has been pertpared for calculation (maybe not optimized) */
42 Finley_RefElement* ReferenceElement; /* the reference element, see Reference element.c */
43 Finley_RefElement* ReferenceElementReducedOrder; /* the reference element with reduced integration order, see Reference element.c */
44 Finley_RefElement* LinearReferenceElement; /* the reference element for the linear mesh. it is vital that it is using the same quadrature scheme like ReferenceElement*/
45 Finley_RefElement* LinearReferenceElementReducedOrder; /* the reference element for the linear mesh. it is vital that it is using the same quadrature
46 \ scheme like LinearReferenceElementReducedIntegration*/
47
48 dim_t numElements; /* number of elements. */
49
50 index_t *Id; /* Id[i] is the id nmber of
51 node i. this number is not
52 used but useful when
53 elements are resorted. in
54 the entire code the term
55 'element id' refers to i
56 but nor to Id[i] if not
57 explicitly stated
58 otherwise. */
59
60 index_t *Tag; /* Tag[i] is the tag of
61 element i. */
62
63 dim_t numNodes; /* number of nodes per element = ReferenceElement.Type.numNodes */
64 index_t *Nodes; /* Nodes[INDEX(k, i, numNodes)]
65 is the k-the node in the
66 i-the element. note that
67 in the way the nodes are
68 ordered Nodes[INDEX(k, i,
69 LinearReferenceElement.Type.numNodes)
70 is k-the node of element i
71 when refering to the
72 linear version of the
73 mesh. */
74 index_t minColor; /* minimum color */
75 index_t maxColor; /* maximum color */
76 index_t *Color; /* assigns each element a color. elements with the same color */
77 /* are don't share a node so they can be processed simultaneously */
78 /* at anytime Color must provide a valid value. In any case one can set */
79 /* Color[e]=e for all e */
80 index_t order; /* order of the element integration scheme*/
81 index_t reduced_order; /* order of the reduced element integration scheme*/
82
83 Finley_ElementFile_Jacobeans* jacobeans; /* element jacobeans */
84 Finley_ElementFile_Jacobeans* jacobeans_reducedS; /* element jacobeans for reduced order of shape function*/
85 Finley_ElementFile_Jacobeans* jacobeans_reducedQ; /* element jacobeans for reduced integration order*/
86 Finley_ElementFile_Jacobeans* jacobeans_reducedS_reducedQ; /* element jacobeans for reduced integration order and reduced order of shape function*/
87
88 };
89
90 typedef struct Finley_ElementFile Finley_ElementFile;
91 Finley_ElementFile* Finley_ElementFile_alloc( ElementTypeId, index_t, index_t, Paso_MPIInfo* );
92 void Finley_ElementFile_free(Finley_ElementFile*);
93 void Finley_ElementFile_allocTable(Finley_ElementFile*,dim_t);
94 void Finley_ElementFile_freeTable(Finley_ElementFile*);
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 /*=================== */
101 void Finley_ElementFile_improveColoring(Finley_ElementFile* in,dim_t numNodes,dim_t* degreeOfFreedom);
102 void Finley_ElementFile_optimizeDistribution(Finley_ElementFile** in);
103 /*=================== */
104 void Finley_ElementFile_setNodeRange(dim_t*,dim_t*,Finley_ElementFile*);
105 void Finley_ElementFile_relableNodes(dim_t*,dim_t,Finley_ElementFile*);
106 void Finley_ElementFile_markNodes(dim_t*,dim_t,Finley_ElementFile*,dim_t);
107 void Finley_ElementFile_scatter(dim_t*,Finley_ElementFile*,Finley_ElementFile*);
108 void Finley_ElementFile_gather(dim_t*,Finley_ElementFile*,Finley_ElementFile*);
109 void Finley_ElementFile_copyTable(dim_t,Finley_ElementFile*,dim_t,dim_t,Finley_ElementFile*);
110
111 void Finley_ElementFile_setTags(Finley_ElementFile*,const int,escriptDataC*);
112 Finley_ElementFile_Jacobeans* Finley_ElementFile_Jacobeans_alloc(Finley_RefElement*);
113 void Finley_ElementFile_Jacobeans_dealloc(Finley_ElementFile_Jacobeans*);
114 Finley_ElementFile_Jacobeans* Finley_ElementFile_borrowJacobeans(Finley_ElementFile*, Finley_NodeFile*, bool_t, bool_t);
115
116
117 #endif /* #ifndef INC_FINLEY_ELEMENTFILE */
118

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26