/[escript]/trunk/paso/src/SolverFCT.c
ViewVC logotype

Contents of /trunk/paso/src/SolverFCT.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1415 - (show annotations)
Thu Feb 21 04:57:17 2008 UTC (11 years, 7 months ago) by ksteube
File MIME type: text/plain
File size: 5907 byte(s)
Recent fixes have remedied problems with running OpenMP/MPI mixed mode programming.

1 /* $Id:$ */
2
3 /*******************************************************
4 *
5 * Copyright 2007 by University of Queensland
6 *
7 * http://esscc.uq.edu.au
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 /* Paso: FCTransportProblem */
17
18 /**************************************************************/
19
20 /* Author: l.gross@uq.edu.au */
21
22 /**************************************************************/
23
24
25 #include "Paso.h"
26 #include "SolverFCT.h"
27 #include "PasoUtil.h"
28
29
30 /**************************************************************/
31
32 /* free all memory used by */
33 void Paso_FCTransportProblem_free(Paso_FCTransportProblem* in) {
34 if (in!=NULL) {
35 in->reference_counter--;
36 if (in->reference_counter<=0) {
37
38 Paso_SystemMatrix_free(in->transport_matrix);
39 Paso_SystemMatrix_free(in->mass_matrix);
40 Paso_SystemMatrix_free(in->iteration_matrix);
41 Paso_MPIInfo_free(in->mpi_info);
42 MEMFREE(in->u);
43 MEMFREE(in->main_iptr);
44 MEMFREE(in->lumped_mass_matrix);
45 MEMFREE(in->main_diagonal_low_order_transport_matrix);
46 MEMFREE(in);
47 }
48 }
49 }
50
51 Paso_FCTransportProblem* Paso_FCTransportProblem_getReference(Paso_FCTransportProblem* in) {
52 if (in!=NULL) {
53 ++(in->reference_counter);
54 }
55 return in;
56 }
57
58 Paso_SystemMatrix* Paso_FCTransportProblem_borrowTransportMatrix(Paso_FCTransportProblem* in) {
59 return in->transport_matrix;
60 }
61 Paso_SystemMatrix* Paso_FCTransportProblem_borrowMassMatrix(Paso_FCTransportProblem* in) {
62 return in->mass_matrix;
63 }
64
65 double* Paso_FCTransportProblem_borrowLumpedMassMatrix(Paso_FCTransportProblem* in) {
66 return in->lumped_mass_matrix;
67 }
68
69 dim_t Paso_FCTransportProblem_getTotalNumRows(Paso_FCTransportProblem* in) {
70 return Paso_SystemMatrix_getTotalNumRows(in->transport_matrix);
71 }
72
73 Paso_FCTransportProblem* Paso_FCTransportProblem_alloc(double theta, Paso_SystemMatrixPattern *pattern, int block_size
74
75
76 ) {
77 Paso_SystemMatrixType matrix_type=MATRIX_FORMAT_DEFAULT+MATRIX_FORMAT_BLK1; /* at the moment only block size 1 is supported */
78 Paso_FCTransportProblem* out=NULL;
79 dim_t n,i;
80 index_t iptr,iptr_main,k;
81
82 if ((theta<0.) || (theta >1.)) {
83 Paso_setError(TYPE_ERROR,"Paso_FCTransportProblem_alloc: theta needs to be between 0. and. 1.");
84 return NULL;
85 }
86
87 out=MEMALLOC(1,Paso_FCTransportProblem);
88 if (Paso_checkPtr(out)) return NULL;
89
90
91 out->theta=theta;
92 out->dt_max=LARGE_POSITIVE_FLOAT;
93 out->valid_matrices=FALSE;
94 out->transport_matrix=Paso_SystemMatrix_alloc(matrix_type,pattern,block_size,block_size);
95 Paso_SystemMatrix_allocBuffer(out->transport_matrix);
96 out->mass_matrix=Paso_SystemMatrix_alloc(matrix_type,pattern,block_size,block_size);
97 out->iteration_matrix=NULL;
98
99 out->mpi_info=Paso_MPIInfo_getReference(pattern->mpi_info);
100 out->u=NULL;
101 out->u_min=0.;
102 out->main_iptr=NULL;
103 out->lumped_mass_matrix=NULL;
104 out->main_diagonal_low_order_transport_matrix=NULL;
105
106 if (Paso_noError()) {
107 n=Paso_SystemMatrix_getTotalNumRows(out->transport_matrix);
108
109 out->u=MEMALLOC(n,double);
110 out->main_iptr=MEMALLOC(n,index_t);
111 out->lumped_mass_matrix=MEMALLOC(n,double);
112 out->main_diagonal_low_order_transport_matrix=MEMALLOC(n,double);
113
114 if ( ! (Paso_checkPtr(out->u) || Paso_checkPtr(out->main_iptr) ||
115 Paso_checkPtr(out->lumped_mass_matrix) || Paso_checkPtr(out->main_diagonal_low_order_transport_matrix)) ) {
116
117 #pragma omp parallel for schedule(static) private(i)
118 for (i = 0; i < n; ++i) {
119 out->lumped_mass_matrix[i]=0.;
120 out->main_diagonal_low_order_transport_matrix[i]=0.;
121 out->u[i]=0.;
122 }
123 /* identify the main diagonals */
124 #pragma omp parallel for schedule(static) private(i,iptr,iptr_main,k)
125 for (i = 0; i < n; ++i) {
126 iptr_main=pattern->mainPattern->ptr[0]-1;
127 for (iptr=pattern->mainPattern->ptr[i];iptr<pattern->mainPattern->ptr[i+1]; iptr++) {
128 if (pattern->mainPattern->index[iptr]==i) {
129 iptr_main=iptr;
130 break;
131 }
132 }
133 out->main_iptr[i]=iptr_main;
134 if (iptr_main==pattern->mainPattern->ptr[0]-1)
135 Paso_setError(VALUE_ERROR, "Paso_FCTransportProblem_alloc: no main diagonal");
136 }
137
138 }
139
140 }
141 if (Paso_noError()) {
142 return out;
143 } else {
144 Paso_FCTransportProblem_free(out);
145 return NULL;
146 }
147 }
148
149 void Paso_FCTransportProblem_checkinSolution(Paso_FCTransportProblem* in, double* u) {
150 dim_t i, n;
151 double local_u_min,u_min;
152
153 n=Paso_FCTransportProblem_getTotalNumRows(in);
154 u_min=LARGE_POSITIVE_FLOAT;
155 #pragma omp parallel private(local_u_min)
156 {
157 local_u_min=0.;
158 #pragma omp for schedule(static) private(i)
159 for (i=0;i<n;++i) local_u_min=MIN(local_u_min,u[i]);
160 #pragma omp critical
161 {
162 u_min=MIN(u_min,local_u_min);
163 }
164 }
165 #ifdef PASO_MPI
166 local_u_min=u_min;
167 MPI_Allreduce(&local_u_min,&u_min, 1, MPI_DOUBLE, MPI_MIN, in->mpi_info->comm);
168 #endif
169 in->u_min=u_min;
170 #pragma omp parallel for schedule(static) private(i)
171 for (i=0;i<n;++i) {
172 in->u[i]=u[i]-u_min;
173 }
174 }

  ViewVC Help
Powered by ViewVC 1.1.26