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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 102 - (show annotations)
Wed Dec 15 07:08:39 2004 UTC (14 years, 10 months ago) by jgs
Original Path: trunk/esys2/finley/src/finleyC/Assemble_setNormal.c
File MIME type: text/plain
File size: 3788 byte(s)
*** empty log message ***

1 /* $Id$ */
2
3 /**************************************************************/
4
5 /* assemblage routines: calculates the normal vector at quadrature points on face elements */
6
7 /**************************************************************/
8
9 /* Copyrights by ACcESS Australia, 2003,2004 */
10 /* author: gross@access.edu.au */
11 /* Version: $Id$ */
12
13 /**************************************************************/
14
15 #include "escript/Data/DataC.h"
16 #include "Finley.h"
17 #include "Assemble.h"
18 #include "NodeFile.h"
19 #include "ElementFile.h"
20 #include "Util.h"
21 #ifdef _OPENMP
22 #include <omp.h>
23 #endif
24
25 /**************************************************************/
26
27
28 void Finley_Assemble_setNormal(Finley_NodeFile* nodes, Finley_ElementFile* elements, escriptDataC* normal) {
29 double *local_X=NULL, *dVdv=NULL,*normal_array;
30 int e,sign,q,node_offset;
31 if (nodes==NULL || elements==NULL) return;
32 int NN=elements->ReferenceElement->Type->numNodes;
33 int NS=elements->ReferenceElement->Type->numShapes;
34 int numDim=nodes->numDim;
35 int numQuad=elements->ReferenceElement->numQuadNodes;
36 int numDim_local=elements->ReferenceElement->Type->numDim;
37
38 /* set some parameter */
39
40 if (getFunctionSpaceType(normal)==FINLEY_CONTACT_ELEMENTS_2) {
41 node_offset=NN-NS;
42 sign=-1;
43 } else {
44 node_offset=0;
45 sign=1;
46 }
47 /* check the dimensions of normal */
48 if (! (numDim==numDim_local || numDim-1==numDim_local)) {
49 Finley_ErrorCode=TYPE_ERROR;
50 sprintf(Finley_ErrorMsg,"Cannot calculate normal vector");
51 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
52 Finley_ErrorCode=TYPE_ERROR;
53 sprintf(Finley_ErrorMsg,"illegal number of samples of normal Data object");
54 } else if (! numSamplesEqual(normal,numQuad,elements->numElements)) {
55 Finley_ErrorCode=TYPE_ERROR;
56 sprintf(Finley_ErrorMsg,"illegal number of samples of normal Data object");
57 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
58 Finley_ErrorCode=TYPE_ERROR;
59 sprintf(Finley_ErrorMsg,"illegal point data shape of normal Data object");
60 } else if (!isExpanded(normal)) {
61 Finley_ErrorCode=TYPE_ERROR;
62 sprintf(Finley_ErrorMsg,"expanded Data object is expected for normal.");
63 }
64
65 /* now we can start */
66 if (Finley_ErrorCode==NO_ERROR) {
67 #pragma omp parallel private(local_X,dVdv)
68 {
69 local_X=dVdv=NULL;
70 /* allocation of work arrays */
71 local_X=THREAD_MEMALLOC(NS*numDim,double);
72 dVdv=THREAD_MEMALLOC(numQuad*numDim*numDim_local,double);
73 if (!(Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) ) ) {
74 /* open the element loop */
75 #pragma omp for private(e,q,normal_array) schedule(static)
76 for(e=0;e<elements->numElements;e++) {
77 /* gather local coordinates of nodes into local_X: */
78 Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
79 /* calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q) */
80 Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
81 /* get normalized vector: */
82 normal_array=getSampleData(normal,e);
83 Finley_NormalVector(numQuad,numDim,numDim_local,dVdv,normal_array);
84 for (q=0;q<numQuad*numDim;q++) normal_array[q]*=sign;
85 }
86 }
87 THREAD_MEMFREE(local_X);
88 THREAD_MEMFREE(dVdv);
89 }
90 }
91 }
92 /*
93 * $Log$
94 * Revision 1.4 2004/12/15 07:08:32 jgs
95 * *** empty log message ***
96 *
97 *
98 *
99 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26