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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 6009 - (show annotations)
Wed Mar 2 04:13:26 2016 UTC (2 years, 11 months ago) by caltinay
File size: 3253 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 normal vector at quadrature points on
20 face elements
21
22 *****************************************************************************/
23
24 #include "Assemble.h"
25 #include "ShapeTable.h"
26 #include "Util.h"
27
28 namespace dudley {
29
30 void Assemble_setNormal(Dudley_NodeFile* nodes, Dudley_ElementFile* elements, escript::Data* normal)
31 {
32 if (!nodes || !elements)
33 return;
34
35 const int NN = elements->numNodes;
36 const int numDim = nodes->numDim;
37 bool reduced_integration = Assemble_reducedIntegrationOrder(normal);
38 const int numQuad = (!reduced_integration) ? (elements->numDim + 1) : 1;
39 const int numDim_local = elements->numLocalDim;
40 const int NS = elements->numDim + 1;
41
42 const double *dSdv = NULL;
43 switch (elements->numDim) {
44 case 2:
45 dSdv = &DTDV_2D[0][0];
46 break;
47 case 3:
48 dSdv = &DTDV_3D[0][0];
49 break;
50 default:
51 dSdv = &DTDV_1D[0][0];
52 break;
53 }
54
55 // check the dimensions of normal
56 if (!(numDim == numDim_local || numDim - 1 == numDim_local)) {
57 throw DudleyException("Assemble_setNormal: Cannot calculate normal vector");
58 } else if (!normal->isDataPointShapeEqual(1, &numDim)) {
59 throw DudleyException("Assemble_setNormal: illegal point data shape of normal Data object");
60 } else if (!normal->numSamplesEqual(numQuad, elements->numElements)) {
61 throw DudleyException("Assemble_setNormal: illegal number of samples of normal Data object");
62 } else if (!normal->actsExpanded()) {
63 throw DudleyException("Assemble_setNormal: expanded Data object is expected for normal.");
64 }
65
66 normal->requireWrite();
67 #pragma omp parallel
68 {
69 std::vector<double> local_X(NS * numDim);
70 std::vector<double> dVdv(numQuad * numDim * numDim_local);
71 #pragma omp for
72 for (index_t e = 0; e < elements->numElements; e++) {
73 // gather local coordinates of nodes into local_X
74 Dudley_Util_Gather_double(NS,
75 &elements->Nodes[INDEX2(0, e, NN)], numDim,
76 nodes->Coordinates, &local_X[0]);
77
78 // calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q)
79 Dudley_Util_SmallMatMult(numDim, numDim_local * numQuad,
80 &dVdv[0], NS, &local_X[0], dSdv);
81 // get normalized vector
82 double* normal_array = normal->getSampleDataRW(e);
83 Dudley_NormalVector(numQuad, numDim, numDim_local, &dVdv[0], normal_array);
84 }
85 }
86 }
87
88 } // namespace dudley
89

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