/[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 2271 - (show annotations)
Mon Feb 16 05:08:29 2009 UTC (10 years, 7 months ago) by jfenwick
File MIME type: text/plain
File size: 4755 byte(s)
Merging version 2269 to trunk

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26