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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 6079 - (show annotations)
Mon Mar 21 12:22:38 2016 UTC (2 years, 11 months ago) by caltinay
File size: 2904 byte(s)
Big commit - making dudley much more like finley to make it more
managable. Fixed quite a few issues that had been fixed in finley.
Disposed of all ReducedNode/ReducedDOF entities that dudley never supported.
Compiles and passes tests.

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_getNormal(const NodeFile* nodes, const ElementFile* elements,
24 escript::Data& normal)
25 {
26 if (!nodes || !elements)
27 return;
28
29 const int NN = elements->numNodes;
30 const int numDim = nodes->numDim;
31 const int numQuad = (hasReducedIntegrationOrder(normal) ? 1 : NN);
32 const int numDim_local = elements->numLocalDim;
33 const int NS = elements->numDim + 1;
34
35 const double *dSdv = NULL;
36 switch (elements->numDim) {
37 case 2:
38 dSdv = &DTDV_2D[0][0];
39 break;
40 case 3:
41 dSdv = &DTDV_3D[0][0];
42 break;
43 default:
44 dSdv = &DTDV_1D[0][0];
45 break;
46 }
47
48 // check the dimensions of normal
49 if (!(numDim == numDim_local || numDim - 1 == numDim_local)) {
50 throw DudleyException("Assemble_setNormal: Cannot calculate normal vector");
51 } else if (!normal.isDataPointShapeEqual(1, &numDim)) {
52 throw DudleyException("Assemble_setNormal: illegal point data shape of normal Data object");
53 } else if (!normal.numSamplesEqual(numQuad, elements->numElements)) {
54 throw DudleyException("Assemble_setNormal: illegal number of samples of normal Data object");
55 } else if (!normal.actsExpanded()) {
56 throw DudleyException("Assemble_setNormal: expanded Data object is expected for normal.");
57 }
58
59 normal.requireWrite();
60 #pragma omp parallel
61 {
62 std::vector<double> local_X(NS * numDim);
63 std::vector<double> dVdv(numQuad * numDim * numDim_local);
64 #pragma omp for
65 for (index_t e = 0; e < elements->numElements; e++) {
66 // gather local coordinates of nodes into local_X
67 util::gather(NS, &elements->Nodes[INDEX2(0, e, NN)], numDim,
68 nodes->Coordinates, &local_X[0]);
69
70 // calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q)
71 util::smallMatMult(numDim, numDim_local * numQuad,
72 &dVdv[0], NS, &local_X[0], dSdv);
73 // get normalized vector
74 double* normal_array = normal.getSampleDataRW(e);
75 util::normalVector(numQuad, numDim, numDim_local, &dVdv[0], normal_array);
76 }
77 }
78 }
79
80 } // namespace dudley
81

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26