1 |
|
2 |
/* $Id$ */ |
3 |
|
4 |
/******************************************************* |
5 |
* |
6 |
* Copyright 2003-2007 by ACceSS MNRF |
7 |
* Copyright 2007 by University of Queensland |
8 |
* |
9 |
* http://esscc.uq.edu.au |
10 |
* Primary Business: Queensland, Australia |
11 |
* Licensed under the Open Software License version 3.0 |
12 |
* http://www.opensource.org/licenses/osl-3.0.php |
13 |
* |
14 |
*******************************************************/ |
15 |
|
16 |
/**************************************************************/ |
17 |
|
18 |
/* assemblage routines: interpolates nodal data in a data array onto elements (=integration points) */ |
19 |
|
20 |
/**************************************************************/ |
21 |
|
22 |
#include "Assemble.h" |
23 |
#include "Util.h" |
24 |
#ifdef _OPENMP |
25 |
#include <omp.h> |
26 |
#endif |
27 |
|
28 |
/**************************************************************/ |
29 |
|
30 |
|
31 |
void Finley_Assemble_interpolate(Finley_NodeFile *nodes, Finley_ElementFile* elements,escriptDataC* data,escriptDataC* interpolated_data) { |
32 |
double* local_data=NULL,*S=NULL,*data_array; |
33 |
index_t dof_offset, NN, NS; |
34 |
bool_t reduced_integration=FALSE; |
35 |
dim_t q,i,NS_DOF,NN_DOF,numNodes,e, numQuad; |
36 |
Finley_RefElement* reference_element=NULL; |
37 |
dim_t numComps=getDataPointSize(data); |
38 |
index_t id[MAX_numNodes], *resort_nodes, *map; |
39 |
type_t data_type=getFunctionSpaceType(data); |
40 |
type_t type; |
41 |
size_t numComps_size; |
42 |
Finley_resetError(); |
43 |
#define NODES 0 |
44 |
#define REDUCED_NODES 3 |
45 |
#define DOF 1 |
46 |
#define REDUCED_DOF 2 |
47 |
if (nodes==NULL || elements==NULL) return; |
48 |
|
49 |
NN=elements->numNodes; |
50 |
NS=elements->ReferenceElement->Type->numShapes; |
51 |
reduced_integration = Finley_Assemble_reducedIntegrationOrder(interpolated_data); |
52 |
for (i=0;i<NN;i++) id[i]=i; |
53 |
|
54 |
/* set some parameter */ |
55 |
|
56 |
if (data_type==FINLEY_NODES) { |
57 |
type=NODES; |
58 |
resort_nodes=id; |
59 |
if (reduced_integration) { |
60 |
reference_element=elements->ReferenceElementReducedOrder; |
61 |
} else { |
62 |
reference_element=elements->ReferenceElement; |
63 |
} |
64 |
numNodes=Finley_NodeFile_getNumNodes(nodes); |
65 |
map=Finley_NodeFile_borrowTargetNodes(nodes); |
66 |
} else if (data_type==FINLEY_REDUCED_NODES) { |
67 |
type=REDUCED_NODES; |
68 |
resort_nodes=elements->ReferenceElement->Type->linearNodes; |
69 |
if (reduced_integration) { |
70 |
reference_element=elements->LinearReferenceElementReducedOrder; |
71 |
} else { |
72 |
reference_element=elements->LinearReferenceElement; |
73 |
} |
74 |
numNodes=Finley_NodeFile_getNumReducedNodes(nodes); |
75 |
map=Finley_NodeFile_borrowTargetReducedNodes(nodes); |
76 |
} else if (data_type==FINLEY_DEGREES_OF_FREEDOM) { |
77 |
if (elements->MPIInfo->size>1) { |
78 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: for more than one processor DEGREES_OF_FREEDOM data are not accepted as input."); |
79 |
return; |
80 |
} |
81 |
type=DOF; |
82 |
resort_nodes=id; |
83 |
if (reduced_integration) { |
84 |
reference_element=elements->ReferenceElementReducedOrder; |
85 |
} else { |
86 |
reference_element=elements->ReferenceElement; |
87 |
} |
88 |
numNodes=Finley_NodeFile_getNumDegreesOfFreedom(nodes); |
89 |
map=Finley_NodeFile_borrowTargetDegreesOfFreedom(nodes); |
90 |
} else if (data_type==FINLEY_REDUCED_DEGREES_OF_FREEDOM) { |
91 |
if (elements->MPIInfo->size>1) { |
92 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: for more than one processor REDUCED_DEGREES_OF_FREEDOM data are not accepted as input."); |
93 |
return; |
94 |
} |
95 |
type=REDUCED_DOF; |
96 |
resort_nodes=elements->ReferenceElement->Type->linearNodes; |
97 |
if (reduced_integration) { |
98 |
reference_element=elements->LinearReferenceElementReducedOrder; |
99 |
} else { |
100 |
reference_element=elements->LinearReferenceElement; |
101 |
} |
102 |
numNodes=Finley_NodeFile_getNumReducedDegreesOfFreedom(nodes); |
103 |
map=Finley_NodeFile_borrowTargetReducedDegreesOfFreedom(nodes); |
104 |
} else { |
105 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: Cannot interpolate data"); |
106 |
} |
107 |
NN_DOF=reference_element->Type->numNodes; |
108 |
NS_DOF=reference_element->Type->numShapes; |
109 |
S=reference_element->S; |
110 |
numQuad=reference_element->numQuadNodes; |
111 |
if (getFunctionSpaceType(interpolated_data)==FINLEY_CONTACT_ELEMENTS_2) { |
112 |
dof_offset=NN_DOF-NS_DOF; |
113 |
} else { |
114 |
dof_offset=0; |
115 |
} |
116 |
|
117 |
/* check the dimensions of interpolated_data and data */ |
118 |
|
119 |
if (! numSamplesEqual(interpolated_data,numQuad,elements->numElements)) { |
120 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: illegal number of samples of output Data object"); |
121 |
} else if (! numSamplesEqual(data,1,numNodes)) { |
122 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: illegal number of samples of input Data object"); |
123 |
} else if (numComps!=getDataPointSize(interpolated_data)) { |
124 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: number of components of input and interpolated Data do not match."); |
125 |
} else if (!isExpanded(interpolated_data)) { |
126 |
Finley_setError(TYPE_ERROR,"Finley_Assemble_interpolate: expanded Data object is expected for output data."); |
127 |
} |
128 |
|
129 |
/* now we can start */ |
130 |
|
131 |
if (Finley_noError()) { |
132 |
#pragma omp parallel private(local_data, numComps_size) |
133 |
{ |
134 |
local_data=NULL; |
135 |
/* allocation of work arrays */ |
136 |
local_data=THREAD_MEMALLOC(NS*numComps,double); |
137 |
if (! Finley_checkPtr(local_data)) { |
138 |
|
139 |
numComps_size=(size_t)numComps*sizeof(double); |
140 |
/* open the element loop */ |
141 |
|
142 |
#pragma omp for private(e,q,i,data_array) schedule(static) |
143 |
for(e=0;e<elements->numElements;e++) { |
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,map[i]); |
147 |
memcpy(local_data+q*numComps, data_array, numComps_size); |
148 |
} |
149 |
/* calculate interpolated_data=local_data*S */ |
150 |
|
151 |
Finley_Util_SmallMatMult(numComps,numQuad,getSampleData(interpolated_data,e),NS_DOF,local_data,S); |
152 |
|
153 |
} /* end of element loop */ |
154 |
|
155 |
} |
156 |
THREAD_MEMFREE(local_data); |
157 |
} /* end of parallel region */ |
158 |
} |
159 |
#undef NODES |
160 |
#undef REDUCED_NODES |
161 |
#undef DOF |
162 |
#undef REDUCED_DOF |
163 |
} |