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

Contents of /trunk-mpi-branch/finley/src/Assemble_integrate.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1223 - (show annotations)
Fri Aug 3 02:40:39 2007 UTC (11 years, 8 months ago) by gross
File MIME type: text/plain
File size: 4783 byte(s)
first attemt towards an improved MPI version.  

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: integrates data on quadrature points */
16
17 /**************************************************************/
18
19 /* Copyrights by ACcESS Australia 2003,2004,2005 */
20 /* Author: gross@access.edu.au */
21 /* version: $Id$ */
22
23 /**************************************************************/
24
25 #include "Assemble.h"
26 #include "Util.h"
27 #ifdef _OPENMP
28 #include <omp.h>
29 #endif
30
31 /**************************************************************/
32
33 void Finley_Assemble_integrate(Finley_NodeFile* nodes, Finley_ElementFile* elements,escriptDataC* data,double* out) {
34 type_t datacase;
35 index_t node_offset;
36 bool_t reducedIntegrationOrder;
37 type_t data_type=getFunctionSpaceType(data);
38 dim_t numComps=getDataPointSize(data);
39 Finley_ElementFile_Jacobeans* jac=NULL;
40 Paso_MPI_rank my_mpi_rank;
41
42 Finley_resetError();
43 if (nodes==NULL || elements==NULL) return;
44 my_mpi_rank = nodes->MPIInfo->rank;
45 /* set some parameter */
46 jac=Finley_ElementFile_borrowJacobeans(elements,nodes,FALSE,Finley_Assemble_reducedIntegrationOrder(data));
47 if (Finley_noError()) {
48
49 /* check the shape of the data */
50 if (! numSamplesEqual(data,jac->ReferenceElement->numQuadNodes,elements->numElements)) {
51 Finley_setError(TYPE_ERROR,"Finley_Assemble_integrate: illegal number of samples of integrant kernel Data object");
52 }
53 /* now we can start */
54
55 if (Finley_noError()) {
56 dim_t q,e,i;
57 double *out_local=NULL, rtmp,*data_array=NULL;
58 for (q=0;q<numComps;q++) out[q]=0;
59 #pragma omp parallel private(q,i,rtmp,data_array,out_local)
60 {
61 out_local=THREAD_MEMALLOC(numComps,double);
62 if (! Finley_checkPtr(out_local) ) {
63 /* initialize local result */
64
65 for (i=0;i<numComps;i++) out_local[i]=0;
66
67 /* open the element loop */
68
69 if (isExpanded(data)) {
70 #pragma omp for private(e) schedule(static)
71 for(e=0;e<elements->numElements;e++) {
72 if (elements->Owner[e] == my_mpi_rank) {
73 data_array=getSampleData(data,e);
74 for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) {
75 for (i=0;i<numComps;i++) out_local[i]+=data_array[INDEX2(i,q,numComps)]*jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
76 }
77 }
78 }
79 } else {
80 #pragma omp for private(e) schedule(static)
81 for(e=0;e<elements->numElements;e++) {
82 if (elements->Owner[e] == my_mpi_rank) {
83 data_array=getSampleData(data,e);
84 rtmp=0.;
85 for (q=0;q<jac->ReferenceElement->numQuadNodes;q++) rtmp+=jac->volume[INDEX2(q,e,jac->ReferenceElement->numQuadNodes)];
86 for (i=0;i<numComps;i++) out_local[i]+=data_array[i]*rtmp;
87 }
88 }
89 }
90 /* add local results to global result */
91 #pragma omp critical
92 for (i=0;i<numComps;i++) out[i]+=out_local[i];
93 }
94 THREAD_MEMFREE(out_local);
95 }
96 }
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.2 2005/06/29 02:34:48 gross
114 * some changes towards 64 integers in finley
115 *
116 * Revision 1.1.1.1.2.1 2004/11/24 01:37:12 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