1 |
|
2 |
/***************************************************************************** |
3 |
* |
4 |
* Copyright (c) 2003-2016 by The University of Queensland |
5 |
* http://www.uq.edu.au |
6 |
* |
7 |
* Primary Business: Queensland, Australia |
8 |
* Licensed under the Open Software License version 3.0 |
9 |
* http://www.opensource.org/licenses/osl-3.0.php |
10 |
* |
11 |
* Development until 2012 by Earth Systems Science Computational Center (ESSCC) |
12 |
* Development 2012-2013 by School of Earth Sciences |
13 |
* Development from 2014 by Centre for Geoscience Computing (GeoComp) |
14 |
* |
15 |
*****************************************************************************/ |
16 |
|
17 |
#include "Assemble.h" |
18 |
#include "ShapeTable.h" |
19 |
#include "Util.h" |
20 |
|
21 |
namespace dudley { |
22 |
|
23 |
void Assemble_AverageElementData(const ElementFile* elements, |
24 |
escript::Data& out, const escript::Data& in) |
25 |
{ |
26 |
if (!elements) |
27 |
return; |
28 |
|
29 |
double wq; |
30 |
int numQuad_in, numQuad_out; |
31 |
if (hasReducedIntegrationOrder(in)) { |
32 |
numQuad_in = QuadNums[elements->numDim][0]; |
33 |
wq = QuadWeight[elements->numDim][0]; |
34 |
} else { |
35 |
numQuad_in = QuadNums[elements->numDim][1]; |
36 |
wq = QuadWeight[elements->numDim][1]; |
37 |
} |
38 |
if (hasReducedIntegrationOrder(out)) { |
39 |
numQuad_out = QuadNums[elements->numDim][0]; |
40 |
} else { |
41 |
numQuad_out = QuadNums[elements->numDim][1]; |
42 |
} |
43 |
|
44 |
// check out and in |
45 |
const dim_t numElements = elements->numElements; |
46 |
const int numComps = out.getDataPointSize(); |
47 |
|
48 |
if (numComps != in.getDataPointSize()) { |
49 |
throw DudleyException("Assemble_AverageElementData: number of components of input and output Data do not match."); |
50 |
} else if (!in.numSamplesEqual(numQuad_in, numElements)) { |
51 |
throw DudleyException("Assemble_AverageElementData: illegal number of samples of input Data object"); |
52 |
} else if (!out.numSamplesEqual(numQuad_out, numElements)) { |
53 |
throw DudleyException("Assemble_AverageElementData: illegal number of samples of output Data object"); |
54 |
} else if (!out.actsExpanded()) { |
55 |
throw DudleyException("Assemble_AverageElementData: expanded Data object is expected for output data."); |
56 |
} else { |
57 |
out.requireWrite(); |
58 |
if (in.actsExpanded()) { |
59 |
const double vol = wq * numQuad_in; |
60 |
const double volinv = 1. / vol; |
61 |
#pragma omp parallel for |
62 |
for (index_t n = 0; n < numElements; n++) { |
63 |
const double* in_array = in.getSampleDataRO(n); |
64 |
double* out_array = out.getSampleDataRW(n); |
65 |
for (int i = 0; i < numComps; ++i) { |
66 |
double rtmp = 0.; |
67 |
for (int q = 0; q < numQuad_in; ++q) |
68 |
rtmp += in_array[INDEX2(i, q, numComps)] * wq; |
69 |
rtmp *= volinv; |
70 |
for (int q = 0; q < numQuad_out; ++q) |
71 |
out_array[INDEX2(i, q, numComps)] = rtmp; |
72 |
} |
73 |
} |
74 |
} else { // constant data |
75 |
const size_t numComps_size = numComps * sizeof(double); |
76 |
#pragma omp parallel for |
77 |
for (index_t n = 0; n < numElements; n++) { |
78 |
const double* in_array = in.getSampleDataRO(n); |
79 |
double* out_array = out.getSampleDataRW(n); |
80 |
for (int q = 0; q < numQuad_out; q++) |
81 |
memcpy(out_array + q * numComps, in_array, numComps_size); |
82 |
} |
83 |
} |
84 |
} |
85 |
} |
86 |
|
87 |
} // namespace dudley |
88 |
|