/[escript]/trunk/esys2/finley/src/finleyC/Assemble_gradient.c
ViewVC logotype

Contents of /trunk/esys2/finley/src/finleyC/Assemble_gradient.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 123 - (show annotations)
Fri Jul 8 04:08:13 2005 UTC (14 years, 9 months ago) by jgs
File MIME type: text/plain
File size: 7771 byte(s)
Merge of development branch back to main trunk on 2005-07-08

1 /* $Id$ */
2
3 /**************************************************************/
4
5 /* assemblage routines: calculate the gradient of nodal data at quadrature points */
6
7 /**************************************************************/
8
9 /* Copyrights by ACcESS Australia, 2003,2004 */
10 /* author: gross@access.edu.au */
11 /* Version: $Id$ */
12
13 /**************************************************************/
14
15 #include "escript/Data/DataC.h"
16 #include "Common.h"
17 #include "Finley.h"
18 #include "Assemble.h"
19 #include "NodeFile.h"
20 #include "ElementFile.h"
21 #include "Util.h"
22 #ifdef _OPENMP
23 #include <omp.h>
24 #endif
25 /*****************************************************************/
26
27
28 #define NODES 0
29 #define DOF 1
30 #define REDUCED_DOF 2
31
32 void Finley_Assemble_gradient(Finley_NodeFile* nodes, Finley_ElementFile* elements,
33 escriptDataC* grad_data,escriptDataC* data) {
34
35 double *local_X=NULL, *local_data=NULL, *dVdv=NULL, *dvdV=NULL, *Vol=NULL, *d_datadv=NULL, *gradS=NULL,*data_array;
36 index_t node_offset,*resort_nodes=FALSE,dof_offset;
37 dim_t numNodes=0,e,i,q,NS_DOF=0,NN_DOF=0;
38 type_t type=DOF;
39 if (nodes==NULL || elements==NULL) return;
40 dim_t NN=elements->ReferenceElement->Type->numNodes;
41 dim_t NS=elements->ReferenceElement->Type->numShapes;
42 index_t id[NN];
43 dim_t numDim=nodes->numDim;
44 type_t data_type=getFunctionSpaceType(data);
45 dim_t numComps=getDataPointSize(data);
46 dim_t numQuad=elements->ReferenceElement->numQuadNodes;
47 for (i=0;i<NN;i++) id[i]=i;
48
49 /* set some parameter */
50
51 if (data_type==FINLEY_NODES) {
52 type=NODES;
53 resort_nodes=id;
54 NN_DOF=elements->ReferenceElement->Type->numNodes;
55 NS_DOF=elements->ReferenceElement->Type->numShapes;
56 gradS=elements->ReferenceElement->dSdv;
57 numNodes=nodes->numNodes;
58 } else if (data_type==FINLEY_DEGREES_OF_FREEDOM) {
59 type=DOF;
60 resort_nodes=id;
61 NN_DOF=elements->ReferenceElement->Type->numNodes;
62 NS_DOF=elements->ReferenceElement->Type->numShapes;
63 gradS=elements->ReferenceElement->dSdv;
64 numNodes=nodes->numDegreesOfFreedom;
65 } else if (data_type==FINLEY_REDUCED_DEGREES_OF_FREEDOM) {
66 type=REDUCED_DOF;
67 resort_nodes=elements->ReferenceElement->Type->linearNodes;
68 NN_DOF=elements->LinearReferenceElement->Type->numNodes;
69 NS_DOF=elements->LinearReferenceElement->Type->numShapes;
70 gradS=elements->LinearReferenceElement->dSdv;
71 numNodes=nodes->reducedNumDegreesOfFreedom;
72 } else {
73 Finley_ErrorCode=TYPE_ERROR;
74 sprintf(Finley_ErrorMsg,"Cannot calculate gradient of data");
75 }
76 if (getFunctionSpaceType(grad_data)==FINLEY_CONTACT_ELEMENTS_2) {
77 node_offset=NN-NS;
78 dof_offset=NN_DOF-NS_DOF;
79 } else {
80 node_offset=0;
81 dof_offset=0;
82 }
83
84 /* check the dimensions of interpolated_data and data */
85
86 if (numDim!=elements->ReferenceElement->Type->numDim) {
87 Finley_ErrorCode=TYPE_ERROR;
88 sprintf(Finley_ErrorMsg,"Gradient: Spatial and element dimension must match.");
89 } else if (! numSamplesEqual(grad_data,numQuad,elements->numElements)) {
90 Finley_ErrorCode=TYPE_ERROR;
91 sprintf(Finley_ErrorMsg,"illegal number of samples in gradient Data object");
92 } else if (! numSamplesEqual(data,1,numNodes)) {
93 Finley_ErrorCode=TYPE_ERROR;
94 sprintf(Finley_ErrorMsg,"illegal number of samples of input Data object");
95 } else if (numDim*numComps!=getDataPointSize(grad_data)) {
96 Finley_ErrorCode=TYPE_ERROR;
97 sprintf(Finley_ErrorMsg,"illegal number of components in gradient data object.");
98 } else if (!isExpanded(grad_data)) {
99 Finley_ErrorCode=TYPE_ERROR;
100 sprintf(Finley_ErrorMsg,"expanded Data object is expected for output data.");
101 }
102
103 /* now we can start */
104
105 if (Finley_ErrorCode==NO_ERROR) {
106 #pragma omp parallel private(local_X,local_data,dvdV,dVdv,Vol,d_datadv)
107 {
108 local_X=local_data=dVdv=dvdV=Vol=d_datadv=NULL;
109 /* allocation of work arrays */
110 local_X=THREAD_MEMALLOC(NS*numDim,double);
111 local_data=THREAD_MEMALLOC(NS*numComps,double);
112 dVdv=THREAD_MEMALLOC(numQuad*numDim*numDim,double);
113 dvdV=THREAD_MEMALLOC(numQuad*numDim*numDim,double);
114 Vol=THREAD_MEMALLOC(numQuad,double);
115 d_datadv=THREAD_MEMALLOC(numQuad*numComps*numDim,double);
116 if (!(Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) || Finley_checkPtr(dvdV) || Finley_checkPtr(Vol) || Finley_checkPtr(d_datadv) || Finley_checkPtr(local_data) )) {
117 /* open the element loop */
118 #pragma omp for private(e,i,q,data_array) schedule(static)
119 for(e=0;e<elements->numElements;e++) {
120 /* gather local coordinates of nodes into local_X: */
121 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
122 /* calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q) */
123 Finley_Util_SmallMatMult(numDim,numDim*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
124 /* dvdV=invert(dVdv) */
125 Finley_Util_InvertSmallMat(numQuad,numDim,dVdv,dvdV,Vol);
126 /* gather local data into local_data(numComps,NS_DOF): */
127 switch (type) {
128 case NODES:
129 for (q=0;q<NS_DOF;q++) {
130 i=elements->Nodes[INDEX2(resort_nodes[dof_offset+q],e,NN)];
131 data_array=getSampleData(data,i);
132 Finley_copyDouble(numComps,data_array,local_data+q*numComps);
133 }
134 break;
135 case DOF:
136 for (q=0;q<NS_DOF;q++) {
137 i=elements->Nodes[INDEX2(resort_nodes[dof_offset+q],e,NN)];
138 data_array=getSampleData(data,nodes->degreeOfFreedom[i]);
139 Finley_copyDouble(numComps,data_array,local_data+q*numComps);
140
141 }
142 break;
143 case REDUCED_DOF:
144 for (q=0;q<NS_DOF;q++) {
145 i=elements->Nodes[INDEX2(resort_nodes[dof_offset+q],e,NN)];
146 data_array=getSampleData(data,nodes->reducedDegreeOfFreedom[i]);
147 Finley_copyDouble(numComps,data_array,local_data+q*numComps);
148 }
149 break;
150 }
151 /* calculate d_datadv(l,i,q)=local_data(l,n)*DSDv(n,i,q) */
152 Finley_Util_SmallMatMult(numComps,numDim*numQuad,d_datadv,NS_DOF,local_data,gradS);
153 /* calculate grad_data(l,i)=d_datadv(l,i,q)*dvdV(i,k,q) */
154 Finley_Util_SmallMatSetMult(numQuad,numComps,numDim,getSampleData(grad_data,e),numDim,d_datadv,dvdV);
155 } /* for */
156 }
157 THREAD_MEMFREE(local_X);
158 THREAD_MEMFREE(dVdv);
159 THREAD_MEMFREE(dvdV);
160 THREAD_MEMFREE(Vol);
161 THREAD_MEMFREE(local_data);
162 THREAD_MEMFREE(d_datadv);
163 }
164 }
165 }
166 #undef NODES
167 #undef DOF
168 #undef REDUCED_DOF
169 /*
170 * $Log$
171 * Revision 1.5 2005/07/08 04:07:47 jgs
172 * Merge of development branch back to main trunk on 2005-07-08
173 *
174 * Revision 1.4 2004/12/15 07:08:32 jgs
175 * *** empty log message ***
176 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
177 * some changes towards 64 integers in finley
178 *
179 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
180 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
181 *
182 *
183 *
184 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26