/[escript]/trunk/finley/src/Assemble_setNormal.c
ViewVC logotype

Contents of /trunk/finley/src/Assemble_setNormal.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2748 - (show annotations)
Tue Nov 17 07:32:59 2009 UTC (10 years ago) by gross
File MIME type: text/plain
File size: 3814 byte(s)
Macro elements are implemented now. VTK writer for macro elements still needs testing.
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2009 by University of Queensland
5 * Earth Systems Science Computational Center (ESSCC)
6 * http://www.uq.edu.au/esscc
7 *
8 * Primary Business: Queensland, Australia
9 * Licensed under the Open Software License version 3.0
10 * http://www.opensource.org/licenses/osl-3.0.php
11 *
12 *******************************************************/
13
14
15 /**************************************************************/
16
17 /* assemblage routines: calculates the normal vector at quadrature points on face elements */
18
19 /**************************************************************/
20
21 #include "Assemble.h"
22 #include "Util.h"
23 #ifdef _OPENMP
24 #include <omp.h>
25 #endif
26
27 /**************************************************************/
28
29
30 void Finley_Assemble_setNormal(Finley_NodeFile* nodes, Finley_ElementFile* elements, escriptDataC* normal) {
31 double *local_X=NULL, *dVdv=NULL,*normal_array;
32 index_t sign,node_offset;
33 Finley_ReferenceElement* reference_element=NULL;
34 dim_t e,q, NN, NS, numDim, numQuad, numDim_local;
35 if (nodes==NULL || elements==NULL) return;
36 Finley_resetError();
37 reference_element=Finley_ReferenceElementSet_borrowReferenceElement(elements->referenceElementSet, Finley_Assemble_reducedIntegrationOrder(normal));
38 NN=elements->numNodes;
39 numDim=nodes->numDim;
40
41 numQuad=reference_element->Parametrization->numQuadNodes;
42 numDim_local=reference_element->Parametrization->Type->numDim;
43 NS=reference_element->Parametrization->Type->numShapes;
44
45
46
47 /* set some parameter */
48
49 if (getFunctionSpaceType(normal)==FINLEY_CONTACT_ELEMENTS_2) {
50 node_offset=reference_element->Type->offsets[1];
51 sign=-1;
52 } else {
53 node_offset=reference_element->Type->offsets[0];
54 sign=1;
55 }
56 /* check the dimensions of normal */
57 if (! (numDim==numDim_local || numDim-1==numDim_local)) {
58 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: Cannot calculate normal vector");
59 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
60 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: illegal number of samples of normal Data object");
61 } else if (! numSamplesEqual(normal,numQuad,elements->numElements)) {
62 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: illegal number of samples of normal Data object");
63 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
64 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: illegal point data shape of normal Data object");
65 } else if (!isExpanded(normal)) {
66 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: expanded Data object is expected for normal.");
67 }
68
69 /* now we can start */
70 if (Finley_noError()) {
71 requireWrite(normal);
72 #pragma omp parallel private(local_X,dVdv)
73 {
74 local_X=dVdv=NULL;
75 /* allocation of work arrays */
76 local_X=THREAD_MEMALLOC(NS*numDim,double);
77 dVdv=THREAD_MEMALLOC(numQuad*numDim*numDim_local,double);
78 if (!(Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) ) ) {
79 /* open the element loop */
80 #pragma omp for private(e,q,normal_array) schedule(static)
81 for(e=0;e<elements->numElements;e++) {
82 /* gather local coordinates of nodes into local_X: */
83 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
84 /* calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q) */
85 Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,reference_element->Parametrization->dSdv);
86 /* get normalized vector: */
87 normal_array=getSampleDataRW(normal,e);
88 Finley_NormalVector(numQuad,numDim,numDim_local,dVdv,normal_array);
89 for (q=0;q<numQuad*numDim;q++) normal_array[q]*=sign;
90 }
91 }
92 THREAD_MEMFREE(local_X);
93 THREAD_MEMFREE(dVdv);
94 }
95 }
96 }
97

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26