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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 155 - (hide annotations)
Wed Nov 9 02:02:19 2005 UTC (13 years, 11 months ago) by jgs
Original Path: trunk/finley/src/finleyC/Assemble_setNormal.c
File MIME type: text/plain
File size: 5035 byte(s)
move all directories from trunk/esys2 into trunk and remove esys2

1 jgs 150 /*
2     ******************************************************************************
3     * *
4     * COPYRIGHT ACcESS 2003,2004,2005 - All Rights Reserved *
5     * *
6     * This software is the property of ACcESS. No part of this code *
7     * may be copied in any form or by any means without the expressed written *
8     * consent of ACcESS. Copying, use or modification of this software *
9     * by any unauthorised person is illegal unless that person has a software *
10     * license agreement with ACcESS. *
11     * *
12     ******************************************************************************
13     */
14 jgs 82
15     /**************************************************************/
16    
17     /* assemblage routines: calculates the normal vector at quadrature points on face elements */
18    
19     /**************************************************************/
20    
21 jgs 150 /* Author: gross@access.edu.au */
22     /* Version: $Id$ */
23 jgs 82
24     /**************************************************************/
25    
26     #include "Assemble.h"
27     #include "Util.h"
28     #ifdef _OPENMP
29     #include <omp.h>
30     #endif
31    
32     /**************************************************************/
33    
34    
35     void Finley_Assemble_setNormal(Finley_NodeFile* nodes, Finley_ElementFile* elements, escriptDataC* normal) {
36     double *local_X=NULL, *dVdv=NULL,*normal_array;
37 jgs 123 index_t sign,node_offset;
38     dim_t e,q;
39 jgs 82 if (nodes==NULL || elements==NULL) return;
40 jgs 123 dim_t NN=elements->ReferenceElement->Type->numNodes;
41     dim_t NS=elements->ReferenceElement->Type->numShapes;
42     dim_t numDim=nodes->numDim;
43     dim_t numQuad=elements->ReferenceElement->numQuadNodes;
44     dim_t numDim_local=elements->ReferenceElement->Type->numDim;
45 jgs 150 Finley_resetError();
46 jgs 82
47     /* set some parameter */
48    
49     if (getFunctionSpaceType(normal)==FINLEY_CONTACT_ELEMENTS_2) {
50     node_offset=NN-NS;
51     sign=-1;
52     } else {
53     node_offset=0;
54     sign=1;
55     }
56     /* check the dimensions of normal */
57 jgs 102 if (! (numDim==numDim_local || numDim-1==numDim_local)) {
58 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: Cannot calculate normal vector");
59 jgs 82 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
60 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: illegal number of samples of normal Data object");
61 jgs 82 } else if (! numSamplesEqual(normal,numQuad,elements->numElements)) {
62 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: illegal number of samples of normal Data object");
63 jgs 82 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
64 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: illegal point data shape of normal Data object");
65 jgs 82 } else if (!isExpanded(normal)) {
66 jgs 150 Finley_setError(TYPE_ERROR,"__FILE__: expanded Data object is expected for normal.");
67 jgs 82 }
68    
69     /* now we can start */
70 jgs 150 if (Finley_noError()) {
71 jgs 82 #pragma omp parallel private(local_X,dVdv)
72     {
73     local_X=dVdv=NULL;
74     /* allocation of work arrays */
75 jgs 102 local_X=THREAD_MEMALLOC(NS*numDim,double);
76     dVdv=THREAD_MEMALLOC(numQuad*numDim*numDim_local,double);
77 jgs 82 if (!(Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) ) ) {
78     /* open the element loop */
79     #pragma omp for private(e,q,normal_array) schedule(static)
80     for(e=0;e<elements->numElements;e++) {
81     /* gather local coordinates of nodes into local_X: */
82     Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
83     /* calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q) */
84     Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
85     /* get normalized vector: */
86     normal_array=getSampleData(normal,e);
87     Finley_NormalVector(numQuad,numDim,numDim_local,dVdv,normal_array);
88     for (q=0;q<numQuad*numDim;q++) normal_array[q]*=sign;
89     }
90     }
91     THREAD_MEMFREE(local_X);
92     THREAD_MEMFREE(dVdv);
93     }
94     }
95     }
96     /*
97     * $Log$
98 jgs 150 * Revision 1.6 2005/09/15 03:44:21 jgs
99     * Merge of development branch dev-02 back to main trunk on 2005-09-15
100     *
101     * Revision 1.5.2.1 2005/09/07 06:26:18 gross
102     * the solver from finley are put into the standalone package paso now
103     *
104 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
105     * Merge of development branch back to main trunk on 2005-07-08
106     *
107 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
108 jgs 97 * *** empty log message ***
109 jgs 123 * Revision 1.1.1.1.2.3 2005/06/29 02:34:48 gross
110     * some changes towards 64 integers in finley
111 jgs 82 *
112 jgs 123 * Revision 1.1.1.1.2.2 2004/11/24 01:37:13 gross
113     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
114 jgs 97 *
115 jgs 82 *
116 jgs 123 *
117 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26