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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26