/[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 777 - (hide annotations)
Wed Jul 12 08:54:45 2006 UTC (13 years, 3 months ago) by gross
Original Path: trunk/finley/src/Assemble_setNormal.c
File MIME type: text/plain
File size: 4798 byte(s)
integration with persistent jacobeans is running now
1 jgs 150 /*
2 elspeth 616 ************************************************************
3     * Copyright 2006 by ACcESS MNRF *
4     * *
5     * http://www.access.edu.au *
6     * Primary Business: Queensland, Australia *
7     * Licensed under the Open Software License version 3.0 *
8     * http://www.opensource.org/licenses/osl-3.0.php *
9     * *
10     ************************************************************
11 jgs 150 */
12 jgs 82
13     /**************************************************************/
14    
15     /* assemblage routines: calculates the normal vector at quadrature points on face elements */
16    
17     /**************************************************************/
18    
19 jgs 150 /* Author: gross@access.edu.au */
20     /* Version: $Id$ */
21 jgs 82
22     /**************************************************************/
23    
24     #include "Assemble.h"
25     #include "Util.h"
26     #ifdef _OPENMP
27     #include <omp.h>
28     #endif
29    
30     /**************************************************************/
31    
32    
33     void Finley_Assemble_setNormal(Finley_NodeFile* nodes, Finley_ElementFile* elements, escriptDataC* normal) {
34     double *local_X=NULL, *dVdv=NULL,*normal_array;
35 jgs 123 index_t sign,node_offset;
36     dim_t e,q;
37 jgs 82 if (nodes==NULL || elements==NULL) return;
38 jgs 123 dim_t NN=elements->ReferenceElement->Type->numNodes;
39     dim_t NS=elements->ReferenceElement->Type->numShapes;
40     dim_t numDim=nodes->numDim;
41     dim_t numQuad=elements->ReferenceElement->numQuadNodes;
42     dim_t numDim_local=elements->ReferenceElement->Type->numDim;
43 jgs 150 Finley_resetError();
44 jgs 82
45     /* set some parameter */
46    
47     if (getFunctionSpaceType(normal)==FINLEY_CONTACT_ELEMENTS_2) {
48     node_offset=NN-NS;
49     sign=-1;
50     } else {
51     node_offset=0;
52     sign=1;
53     }
54     /* check the dimensions of normal */
55 jgs 102 if (! (numDim==numDim_local || numDim-1==numDim_local)) {
56 gross 777 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: Cannot calculate normal vector");
57 jgs 82 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
58 gross 777 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: illegal number of samples of normal Data object");
59 jgs 82 } else if (! numSamplesEqual(normal,numQuad,elements->numElements)) {
60 gross 777 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: illegal number of samples of normal Data object");
61 jgs 82 } else if (! isDataPointShapeEqual(normal,1,&(numDim))) {
62 gross 777 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: illegal point data shape of normal Data object");
63 jgs 82 } else if (!isExpanded(normal)) {
64 gross 777 Finley_setError(TYPE_ERROR,"Finley_Assemble_setNormal: expanded Data object is expected for normal.");
65 jgs 82 }
66    
67     /* now we can start */
68 jgs 150 if (Finley_noError()) {
69 jgs 82 #pragma omp parallel private(local_X,dVdv)
70     {
71     local_X=dVdv=NULL;
72     /* allocation of work arrays */
73 jgs 102 local_X=THREAD_MEMALLOC(NS*numDim,double);
74     dVdv=THREAD_MEMALLOC(numQuad*numDim*numDim_local,double);
75 jgs 82 if (!(Finley_checkPtr(local_X) || Finley_checkPtr(dVdv) ) ) {
76     /* open the element loop */
77     #pragma omp for private(e,q,normal_array) schedule(static)
78     for(e=0;e<elements->numElements;e++) {
79     /* gather local coordinates of nodes into local_X: */
80     Finley_Util_Gather_double(NS,&(elements->Nodes[INDEX2(node_offset,e,NN)]),numDim,nodes->Coordinates,local_X);
81     /* calculate dVdv(i,j,q)=local_X(i,n)*DSDv(n,j,q) */
82     Finley_Util_SmallMatMult(numDim,numDim_local*numQuad,dVdv,NS,local_X,elements->ReferenceElement->dSdv);
83     /* get normalized vector: */
84     normal_array=getSampleData(normal,e);
85     Finley_NormalVector(numQuad,numDim,numDim_local,dVdv,normal_array);
86     for (q=0;q<numQuad*numDim;q++) normal_array[q]*=sign;
87     }
88     }
89     THREAD_MEMFREE(local_X);
90     THREAD_MEMFREE(dVdv);
91     }
92     }
93     }
94     /*
95     * $Log$
96 jgs 150 * Revision 1.6 2005/09/15 03:44:21 jgs
97     * Merge of development branch dev-02 back to main trunk on 2005-09-15
98     *
99     * Revision 1.5.2.1 2005/09/07 06:26:18 gross
100     * the solver from finley are put into the standalone package paso now
101     *
102 jgs 123 * Revision 1.5 2005/07/08 04:07:48 jgs
103     * Merge of development branch back to main trunk on 2005-07-08
104     *
105 jgs 102 * Revision 1.4 2004/12/15 07:08:32 jgs
106 jgs 97 * *** empty log message ***
107 jgs 123 * Revision 1.1.1.1.2.3 2005/06/29 02:34:48 gross
108     * some changes towards 64 integers in finley
109 jgs 82 *
110 jgs 123 * Revision 1.1.1.1.2.2 2004/11/24 01:37:13 gross
111     * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
112 jgs 97 *
113 jgs 82 *
114 jgs 123 *
115 jgs 82 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26