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

Contents of /trunk/finley/src/Assemble_integrate.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1811 - (show annotations)
Thu Sep 25 23:11:13 2008 UTC (10 years, 11 months ago) by ksteube
File MIME type: text/plain
File size: 4365 byte(s)
Copyright updated in all files

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: integrates data on quadrature points */
18
19 /**************************************************************/
20
21 #include "Assemble.h"
22 #include "Util.h"
23 #ifdef _OPENMP
24 #include <omp.h>
25 #endif
26
27 /**************************************************************/
28
29 void Finley_Assemble_integrate(Finley_NodeFile* nodes, Finley_ElementFile* elements,escriptDataC* data,double* out) {
30 type_t data_type=getFunctionSpaceType(data);
31 dim_t numComps=getDataPointSize(data);
32 Finley_ElementFile_Jacobeans* jac=NULL;
33 Paso_MPI_rank my_mpi_rank;
34
35 Finley_resetError();
36 if (nodes==NULL || elements==NULL) return;
37 my_mpi_rank = nodes->MPIInfo->rank;
38 /* set some parameter */
39 jac=Finley_ElementFile_borrowJacobeans(elements,nodes,FALSE,Finley_Assemble_reducedIntegrationOrder(data));
40 if (Finley_noError()) {
41
42 /* check the shape of the data */
43 if (! numSamplesEqual(data,jac->ReferenceElement->numQuadNodes,elements->numElements)) {
44 Finley_setError(TYPE_ERROR,"Finley_Assemble_integrate: illegal number of samples of integrant kernel Data object");
45 }
46 /* now we can start */
47
48 if (Finley_noError()) {
49 dim_t q,e,i;
50 double *out_local=NULL, rtmp,*data_array=NULL;
51 for (q=0;q<numComps;q++) out[q]=0;
52 #pragma omp parallel private(q,i,rtmp,data_array,out_local)
53 {
54 out_local=THREAD_MEMALLOC(numComps,double);
55 if (! Finley_checkPtr(out_local) ) {
56 /* initialize local result */
57
58 for (i=0;i<numComps;i++) out_local[i]=0;
59
60 /* open the element loop */
61
62 if (isExpanded(data)) {
63 #pragma omp for private(e) schedule(static)
64 for(e=0;e<elements->numElements;e++) {
65 if (elements->Owner[e] == my_mpi_rank) {
66 data_array=getSampleData(data,e);
67 for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) {
68 for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
69 }
70 }
71 }
72 } else {
73 #pragma omp for private(e) schedule(static)
74 for(e=0;e<elements->numElements;e++) {
75 if (elements->Owner[e] == my_mpi_rank) {
76 data_array=getSampleData(data,e);
77 rtmp=0.;
78 for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) rtmp+=jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
79 for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
80 }
81 }
82 }
83 /* add local results to global result */
84 #pragma omp critical
85 for (i=0;i<numComps;i++) out[i]+=out_local[i];
86 }
87 THREAD_MEMFREE(out_local);
88 }
89 }
90 }
91 }
92
93 /*
94 * $Log$
95 * Revision 1.6 2005/09/15 03:44:21 jgs
96 * Merge of development branch dev-02 back to main trunk on 2005-09-15
97 *
98 * Revision 1.5.2.1 2005/09/07 06:26:18 gross
99 * the solver from finley are put into the standalone package paso now
100 *
101 * Revision 1.5 2005/07/08 04:07:48 jgs
102 * Merge of development branch back to main trunk on 2005-07-08
103 *
104 * Revision 1.4 2004/12/15 07:08:32 jgs
105 * *** empty log message ***
106 * Revision 1.1.1.1.2.2 2005/06/29 02:34:48 gross
107 * some changes towards 64 integers in finley
108 *
109 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 gross
110 * some changes dealing with the integer overflow in memory allocation. Finley solves 4M unknowns now
111 *
112 *
113 *
114 */

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26