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

Annotation of /trunk/finley/src/Assemble_interpolate.c

Parent Directory Parent Directory | Revision Log Revision Log


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

1 jgs 150 /*
2 elspeth 616 ************************************************************
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 jgs 150 */
12 jgs 82
13     /**************************************************************/
14    
15     /* assemblage routines: interpolates nodal data in a data array onto elements (=integration points) */
16    
17     /**************************************************************/
18    
19 jgs 150 /* Author: gross@access.edu.au */
20     /* Version: $Id$ */
21 jgs 82
22     /**************************************************************/
23    
24     #include "Assemble.h"
25     #include "Util.h"
26     #ifdef _OPENMP
27     #include <omp.h>
28     #endif
29    
30     /**************************************************************/
31    
32    
33     void Finley_Assemble_interpolate(Finley_NodeFile *nodes, Finley_ElementFile* elements,escriptDataC* data,escriptDataC* interpolated_data) {
34     double* local_data=NULL,*S=NULL,*data_array;
35 jgs 123 index_t dof_offset,*resort_nodes;
36     dim_t q,i,NS_DOF,NN_DOF,numNodes,e;
37     type_t type;
38 jgs 82 #define NODES 0
39     #define DOF 1
40     #define REDUCED_DOF 2
41     if (nodes==NULL || elements==NULL) return;
42 jgs 123 dim_t NN=elements->ReferenceElement->Type->numNodes;
43     dim_t NS=elements->ReferenceElement->Type->numShapes;
44     dim_t numComps=getDataPointSize(data);
45     type_t data_type=getFunctionSpaceType(data);
46     dim_t numQuad=elements->ReferenceElement->numQuadNodes;
47     index_t id[NN];
48 jgs 82 for (i=0;i<NN;i++) id[i]=i;
49 jgs 150 Finley_resetError();
50 jgs 82
51     /* set some parameter */
52    
53     if (data_type==FINLEY_NODES) {
54     type=NODES;
55     resort_nodes=id;
56     NN_DOF=elements->ReferenceElement->Type->numNodes;
57     NS_DOF=elements->ReferenceElement->Type->numShapes;
58     S=elements->ReferenceElement->S;
59     numNodes=nodes->numNodes;
60     } else if (data_type==FINLEY_DEGREES_OF_FREEDOM) {
61     type=DOF;
62     resort_nodes=id;
63     NN_DOF=elements->ReferenceElement->Type->numNodes;
64     NS_DOF=elements->ReferenceElement->Type->numShapes;
65     S=elements->ReferenceElement->S;
66     numNodes=nodes->numDegreesOfFreedom;
67     } else if (data_type==FINLEY_REDUCED_DEGREES_OF_FREEDOM) {
68     type=REDUCED_DOF;
69     resort_nodes=elements->ReferenceElement->Type->linearNodes;
70     NN_DOF=elements->LinearReferenceElement->Type->numNodes;
71     NS_DOF=elements->LinearReferenceElement->Type->numShapes;
72     S=elements->LinearReferenceElement->S;
73     numNodes=nodes->reducedNumDegreesOfFreedom;
74     } else {
75 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: Cannot interpolate data");
76 jgs 82 }
77    
78     if (getFunctionSpaceType(interpolated_data)==FINLEY_CONTACT_ELEMENTS_2) {
79     dof_offset=NN_DOF-NS_DOF;
80     } else {
81     dof_offset=0;
82     }
83    
84     /* check the dimensions of interpolated_data and data */
85    
86     if (! numSamplesEqual(interpolated_data,numQuad,elements->numElements)) {
87 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: illegal number of samples of output Data object");
88 jgs 82 } else if (! numSamplesEqual(data,1,numNodes)) {
89 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: illegal number of samples of input Data object");
90 jgs 82 } else if (numComps!=getDataPointSize(interpolated_data)) {
91 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: number of components of input and interpolated Data do not match.");
92 jgs 82 } else if (!isExpanded(interpolated_data)) {
93 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: expanded Data object is expected for output data.");
94 jgs 82 }
95    
96     /* now we can start */
97    
98 jgs 150 if (Finley_noError()) {
99 jgs 82 #pragma omp parallel private(local_data)
100     {
101     local_data=NULL;
102     /* allocation of work arrays */
103 jgs 102 local_data=THREAD_MEMALLOC(NS*numComps,double);
104 jgs 82 if (! Finley_checkPtr(local_data)) {
105    
106     /* open the element loop */
107    
108     #pragma omp for private(e,q,i,data_array) schedule(static)
109     for(e=0;e<elements->numElements;e++) {
110     /* gather local data into local_data(numComps,NS_DOF): */
111     switch (type) {
112     case NODES:
113     for (q=0;q<NS_DOF;q++) {
114     i=elements->Nodes[INDEX2(resort_nodes[dof_offset+q],e,NN)];
115     data_array=getSampleData(data,i);
116     Finley_copyDouble(numComps,data_array,local_data+q*numComps);
117 jgs 102 }
118 jgs 82 break;
119     case DOF:
120     for (q=0;q<NS_DOF;q++) {
121     i=elements->Nodes[INDEX2(resort_nodes[dof_offset+q],e,NN)];
122     data_array=getSampleData(data,nodes->degreeOfFreedom[i]);
123     Finley_copyDouble(numComps,data_array,local_data+q*numComps);
124     }
125     break;
126     case REDUCED_DOF:
127     for (q=0;q<NS_DOF;q++) {
128     i=elements->Nodes[INDEX2(resort_nodes[dof_offset+q],e,NN)];
129     data_array=getSampleData(data,nodes->reducedDegreeOfFreedom[i]);
130     Finley_copyDouble(numComps,data_array,local_data+q*numComps);
131     }
132     break;
133     }
134    
135     /* calculate interpolated_data=local_data*S */
136    
137     Finley_Util_SmallMatMult(numComps,numQuad,getSampleData(interpolated_data,e),NS_DOF,local_data,S);
138    
139     } /* end of element loop */
140    
141     }
142     THREAD_MEMFREE(local_data);
143 jgs 102 } /* end of parallel region */
144 jgs 82 }
145     #undef NODES
146     #undef DOF
147     #undef REDUCED_DOF
148     }
149     /*
150     * $Log$
151 jgs 150 * Revision 1.6 2005/09/15 03:44:21 jgs
152     * Merge of development branch dev-02 back to main trunk on 2005-09-15
153     *
154     * Revision 1.5.2.1 2005/09/07 06:26:18 gross
155     * the solver from finley are put into the standalone package paso now
156     *
157 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
158     * Merge of development branch back to main trunk on 2005-07-08
159     *
160 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
161 jgs 97 * *** empty log message ***
162 jgs 123 * Revision 1.1.1.1.2.3 2005/06/29 02:34:48 gross
163     * some changes towards 64 integers in finley
164 jgs 82 *
165 jgs 123 * Revision 1.1.1.1.2.2 2004/11/24 01:37:12 gross
166     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
167 jgs 97 *
168 jgs 82 *
169 jgs 123 *
170 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26