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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 123 - (hide annotations)
Fri Jul 8 04:08:13 2005 UTC (14 years, 4 months ago) by jgs
Original Path: trunk/esys2/finley/src/finleyC/Assemble_setNormal.c
File MIME type: text/plain
File size: 4175 byte(s)
Merge of development branch back to main trunk on 2005-07-08

1 jgs 82 /* $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 jgs 123 index_t sign,node_offset;
31     dim_t e,q;
32 jgs 82 if (nodes==NULL || elements==NULL) return;
33 jgs 123 dim_t NN=elements->ReferenceElement->Type->numNodes;
34     dim_t NS=elements->ReferenceElement->Type->numShapes;
35     dim_t numDim=nodes->numDim;
36     dim_t numQuad=elements->ReferenceElement->numQuadNodes;
37     dim_t numDim_local=elements->ReferenceElement->Type->numDim;
38 jgs 82
39     /* set some parameter */
40    
41     if (getFunctionSpaceType(normal)==FINLEY_CONTACT_ELEMENTS_2) {
42     node_offset=NN-NS;
43     sign=-1;
44     } else {
45     node_offset=0;
46     sign=1;
47     }
48     /* check the dimensions of normal */
49 jgs 102 if (! (numDim==numDim_local || numDim-1==numDim_local)) {
50 jgs 82 Finley_ErrorCode=TYPE_ERROR;
51     sprintf(Finley_ErrorMsg,"Cannot calculate normal vector");
52     } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
53     Finley_ErrorCode=TYPE_ERROR;
54     sprintf(Finley_ErrorMsg,"illegal number of samples of normal Data object");
55     } else if (! numSamplesEqual(normal,numQuad,elements->numElements)) {
56     Finley_ErrorCode=TYPE_ERROR;
57     sprintf(Finley_ErrorMsg,"illegal number of samples of normal Data object");
58     } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
59     Finley_ErrorCode=TYPE_ERROR;
60     sprintf(Finley_ErrorMsg,"illegal point data shape of normal Data object");
61     } else if (!isExpanded(normal)) {
62     Finley_ErrorCode=TYPE_ERROR;
63     sprintf(Finley_ErrorMsg,"expanded Data object is expected for normal.");
64     }
65    
66     /* now we can start */
67     if (Finley_ErrorCode==NO_ERROR) {
68     #pragma omp parallel private(local_X,dVdv)
69     {
70     local_X=dVdv=NULL;
71     /* allocation of work arrays */
72 jgs 102 local_X=THREAD_MEMALLOC(NS*numDim,double);
73     dVdv=THREAD_MEMALLOC(numQuad*numDim*numDim_local,double);
74 jgs 82 if (!(Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) ) ) {
75     /* open the element loop */
76     #pragma omp for private(e,q,normal_array) schedule(static)
77     for(e=0;e<elements->numElements;e++) {
78     /* gather local coordinates of nodes into local_X: */
79     Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
80     /* calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q) */
81     Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
82     /* get normalized vector: */
83     normal_array=getSampleData(normal,e);
84     Finley_NormalVector(numQuad,numDim,numDim_local,dVdv,normal_array);
85     for (q=0;q<numQuad*numDim;q++) normal_array[q]*=sign;
86     }
87     }
88     THREAD_MEMFREE(local_X);
89     THREAD_MEMFREE(dVdv);
90     }
91     }
92     }
93     /*
94     * $Log$
95 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
96     * Merge of development branch back to main trunk on 2005-07-08
97     *
98 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
99 jgs 97 * *** empty log message ***
100 jgs 123 * Revision 1.1.1.1.2.3 2005/06/29 02:34:48 gross
101     * some changes towards 64 integers in finley
102 jgs 82 *
103 jgs 123 * Revision 1.1.1.1.2.2 2004/11/24 01:37:13 gross
104     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
105 jgs 97 *
106 jgs 82 *
107 jgs 123 *
108 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26