/[escript]/branches/RW_WIN32/finley/src/finleyC/Assemble_gradient.c
ViewVC logotype

Contents of /branches/RW_WIN32/finley/src/finleyC/Assemble_gradient.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 210 - (show annotations)
Wed Nov 23 09:54:02 2005 UTC (14 years, 7 months ago) by robwdcock
File MIME type: text/plain
File size: 8752 byte(s)
PARTIAL WIN32 BUILD SYSTEM AND PORT
+ All libraries now build under new build system. No unit test routines yet
+ Reversed some incorrect refactorings in Bruce.cpp

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26