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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 616 - (show annotations)
Wed Mar 22 02:46:56 2006 UTC (13 years, 8 months ago) by elspeth
File MIME type: text/plain
File size: 8617 byte(s)
Copyright added to more source files.

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 /**************************************************************/
14
15 /* assemblage routines: calculate the gradient of nodal data at quadrature points */
16
17 /**************************************************************/
18
19 /* Copyrights by ACcESS Australia, 2003,2004,2005 */
20 /* author: gross@access.edu.au */
21 /* version: $Id$ */
22
23 /**************************************************************/
24
25 #include "Assemble.h"
26 #include "Util.h"
27 #ifdef _OPENMP
28 #include <omp.h>
29 #endif
30 /*****************************************************************/
31
32
33 #define NODES 0
34 #define DOF 1
35 #define REDUCED_DOF 2
36
37 void Finley_Assemble_gradient(Finley_NodeFile* nodes, Finley_ElementFile* elements,
38 escriptDataC* grad_data,escriptDataC* data) {
39
40 double *local_X=NULL, *local_data=NULL, *dVdv=NULL, *dvdV=NULL, *Vol=NULL, *d_datadv=NULL, *gradS=NULL,*data_array;
41 index_t node_offset,*resort_nodes=FALSE,dof_offset;
42 dim_t numNodes=0,e,i,q,NS_DOF=0,NN_DOF=0;
43 type_t type=DOF;
44 if (nodes==NULL || elements==NULL) return;
45 dim_t NN=elements->ReferenceElement->Type->numNodes;
46 dim_t NS=elements->ReferenceElement->Type->numShapes;
47 index_t id[NN];
48 dim_t numDim=nodes->numDim;
49 type_t data_type=getFunctionSpaceType(data);
50 dim_t numComps=getDataPointSize(data);
51 dim_t numQuad=elements->ReferenceElement->numQuadNodes;
52 for (i=0;i<NN;i++) id[i]=i;
53 Finley_resetError();
54
55 /* set some parameter */
56
57 if (data_type==FINLEY_NODES) {
58 type=NODES;
59 resort_nodes=id;
60 NN_DOF=elements->ReferenceElement->Type->numNodes;
61 NS_DOF=elements->ReferenceElement->Type->numShapes;
62 gradS=elements->ReferenceElement->dSdv;
63 numNodes=nodes->numNodes;
64 } else if (data_type==FINLEY_DEGREES_OF_FREEDOM) {
65 type=DOF;
66 resort_nodes=id;
67 NN_DOF=elements->ReferenceElement->Type->numNodes;
68 NS_DOF=elements->ReferenceElement->Type->numShapes;
69 gradS=elements->ReferenceElement->dSdv;
70 numNodes=nodes->numDegreesOfFreedom;
71 } else if (data_type==FINLEY_REDUCED_DEGREES_OF_FREEDOM) {
72 type=REDUCED_DOF;
73 resort_nodes=elements->ReferenceElement->Type->linearNodes;
74 NN_DOF=elements->LinearReferenceElement->Type->numNodes;
75 NS_DOF=elements->LinearReferenceElement->Type->numShapes;
76 gradS=elements->LinearReferenceElement->dSdv;
77 numNodes=nodes->reducedNumDegreesOfFreedom;
78 } else {
79 Finley_setError(TYPE_ERROR,"Finley_Assemble_gradient: Cannot calculate gradient of data");
80 }
81 if (getFunctionSpaceType(grad_data)==FINLEY_CONTACT_ELEMENTS_2) {
82 node_offset=NN-NS;
83 dof_offset=NN_DOF-NS_DOF;
84 } else {
85 node_offset=0;
86 dof_offset=0;
87 }
88
89 /* check the dimensions of interpolated_data and data */
90
91 if (numDim!=elements->ReferenceElement->Type->numDim) {
92 Finley_setError(TYPE_ERROR,"Finley_Assemble_gradient: Spatial and element dimension must match.");
93 } else if (! numSamplesEqual(grad_data,numQuad,elements->numElements)) {
94 Finley_setError(TYPE_ERROR,"Finley_Assemble_gradient: illegal number of samples in gradient Data object");
95 } else if (! numSamplesEqual(data,1,numNodes)) {
96 Finley_setError(TYPE_ERROR,"Finley_Assemble_gradient: illegal number of samples of input Data object");
97 } else if (numDim*numComps!=getDataPointSize(grad_data)) {
98 Finley_setError(TYPE_ERROR,"Finley_Assemble_gradient: illegal number of components in gradient data object.");
99 } else if (!isExpanded(grad_data)) {
100 Finley_setError(TYPE_ERROR,"Finley_Assemble_gradient: expanded Data object is expected for output data.");
101 }
102
103 /* now we can start */
104
105 if (Finley_noError()) {
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 grad_data(l,i,q)=local_data(l,n)* DSDV(n,i,q) */
152 // Finley_Util_SmallMatMult(numQuad,numComps,numQuad*numDim,getSampleData(grad_data,e),NS_DOF,local_data,dSdV);
153
154 /* calculate d_datadv(l,i,q)=local_data(l,n)*DSDv(n,i,q) */
155 Finley_Util_SmallMatMult(numComps,numDim*numQuad,d_datadv,NS_DOF,local_data,gradS);
156 /* calculate grad_data(l,i,q)=d_datadv(l,k,q)*dvdV(k,i,q) */
157 Finley_Util_SmallMatSetMult(numQuad,numComps,numDim,getSampleData(grad_data,e),numDim,d_datadv,dvdV);
158 } /* for */
159 }
160 THREAD_MEMFREE(local_X);
161 THREAD_MEMFREE(dVdv);
162 THREAD_MEMFREE(dvdV);
163 THREAD_MEMFREE(Vol);
164 THREAD_MEMFREE(local_data);
165 THREAD_MEMFREE(d_datadv);
166 }
167 }
168 }
169 #undef NODES
170 #undef DOF
171 #undef REDUCED_DOF
172 /*
173 * $Log$
174 * Revision 1.6 2005/09/15 03:44:21 jgs
175 * Merge of development branch dev-02 back to main trunk on 2005-09-15
176 *
177 * Revision 1.5.2.1 2005/09/07 06:26:17 gross
178 * the solver from finley are put into the standalone package paso now
179 *
180 * Revision 1.5 2005/07/08 04:07:47 jgs
181 * Merge of development branch back to main trunk on 2005-07-08
182 *
183 * Revision 1.4 2004/12/15 07:08:32 jgs
184 * *** empty log message ***
185 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
186 * some changes towards 64 integers in finley
187 *
188 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
189 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
190 *
191 *
192 *
193 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26