/[escript]/branches/trilinos_from_5897/dudley/src/Assemble_getSize.cpp
ViewVC logotype

Contents of /branches/trilinos_from_5897/dudley/src/Assemble_getSize.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 6009 - (show annotations)
Wed Mar 2 04:13:26 2016 UTC (3 years, 1 month ago) by caltinay
File size: 3404 byte(s)
Much needed sync with trunk...

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 /****************************************************************************
18
19 Assemblage routines: calculates the minimum distance between two vertices
20 of elements and assigns the value to each quadrature point in out
21
22 *****************************************************************************/
23
24 #include "Assemble.h"
25 #include "Util.h"
26
27 namespace dudley {
28
29 void Assemble_getSize(Dudley_NodeFile* nodes, Dudley_ElementFile* elements, escript::Data* out)
30 {
31 if (!nodes || !elements)
32 return;
33
34 const int numDim = nodes->numDim;
35
36 // now we look up what type of elements we need based on the function space
37 // of out. If it is DUDLEY_REDUCED_ELEMENTS or
38 // DUDLEY_REDUCED_FACE_ELEMENTS then we have single quad point
39 int numQuad;
40 if (Assemble_reducedIntegrationOrder(out)) {
41 numQuad = 1;
42 } else {
43 numQuad = elements->numDim + 1;
44 }
45 const int NN = elements->numNodes;
46 const int NS = elements->numDim + 1;
47 const int NVertices = elements->numDim + 1;
48
49 // check the dimensions of out
50 if (!out->numSamplesEqual(numQuad, elements->numElements)) {
51 throw DudleyException("Assemble_getSize: illegal number of samples of element size Data object");
52 } else if (!out->isDataPointShapeEqual(0, &numDim)) {
53 throw DudleyException("Assemble_getSize: illegal data point shape of element size Data object");
54 } else if (!out->actsExpanded()) {
55 throw DudleyException("Assemble_getSize: expanded Data object is expected for element size.");
56 }
57
58 // now we can start
59 out->requireWrite();
60 #pragma omp parallel
61 {
62 std::vector<double> local_X(NN * numDim);
63 #pragma omp for
64 for (index_t e = 0; e < elements->numElements; e++) {
65 // gather local coordinates of nodes into local_X(numDim,NN)
66 Dudley_Util_Gather_double(NS, &elements->Nodes[INDEX2(0, e, NN)],
67 numDim, nodes->Coordinates, &local_X[0]);
68 // calculate minimal differences
69 double max_diff = 0;
70 for (int n0 = 0; n0 < NVertices; n0++) {
71 for (int n1 = n0 + 1; n1 < NVertices; n1++) {
72 double diff = 0;
73 for (int i = 0; i < numDim; i++) {
74 const double d = local_X[INDEX2(i, n0, numDim)] - local_X[INDEX2(i, n1, numDim)];
75 diff += d * d;
76 }
77
78 max_diff = std::max(max_diff, diff);
79 }
80 }
81 max_diff = sqrt(max_diff);
82 // set all values to max_diff
83 double* out_array = out->getSampleDataRW(e);
84 for (int q = 0; q < numQuad; q++)
85 out_array[q] = max_diff;
86 }
87 } // end of parallel region
88 }
89
90 } // namespace dudley
91

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision
svn:mergeinfo /branches/4.0fordebian/dudley/src/Assemble_getSize.cpp:5567-5588 /branches/lapack2681/finley/src/Assemble_getSize.cpp:2682-2741 /branches/pasowrap/dudley/src/Assemble_getSize.cpp:3661-3674 /branches/py3_attempt2/dudley/src/Assemble_getSize.cpp:3871-3891 /branches/restext/finley/src/Assemble_getSize.cpp:2610-2624 /branches/ripleygmg_from_3668/dudley/src/Assemble_getSize.cpp:3669-3791 /branches/stage3.0/finley/src/Assemble_getSize.cpp:2569-2590 /branches/symbolic_from_3470/dudley/src/Assemble_getSize.cpp:3471-3974 /branches/symbolic_from_3470/ripley/test/python/dudley/src/Assemble_getSize.cpp:3517-3974 /release/3.0/finley/src/Assemble_getSize.cpp:2591-2601 /release/4.0/dudley/src/Assemble_getSize.cpp:5380-5406 /trunk/dudley/src/Assemble_getSize.cpp:4257-4344,5898-6007 /trunk/ripley/test/python/dudley/src/Assemble_getSize.cpp:3480-3515

  ViewVC Help
Powered by ViewVC 1.1.26