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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1812 - (show annotations)
Fri Sep 26 00:19:18 2008 UTC (11 years, 1 month ago) by ksteube
File MIME type: text/plain
File size: 15558 byte(s)
Re-ordered the methods in paso/src/SolverFCT_solve.c and resolved several compiler warnings

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 /* Paso: Flux correction transport solver
18 *
19 * solves Mu_t=Ku+q
20 * u(0) >=0
21 *
22 * Warning: the program assums sum_{j} k_{ij}=0!!!
23 *
24 */
25 /**************************************************************/
26
27 /* Author: l.gross@uq.edu.au */
28
29 /**************************************************************/
30
31 #include "Paso.h"
32 #include "Solver.h"
33 #include "SolverFCT.h"
34 #include "PasoUtil.h"
35 #include "escript/blocktimer.h"
36 #ifdef _OPENMP
37 #include <omp.h>
38 #endif
39 #ifdef PASO_MPI
40 #include <mpi.h>
41 #endif
42
43 /*
44 * inserts the source term into the problem
45 */
46 void Paso_FCT_setSource(const dim_t n,const double *source, double* sourceN, double* sourceP)
47 {
48 dim_t i;
49 register double rtmp;
50 /*
51 * seperate source into positive and negative part:
52 */
53 #pragma omp parallel for private(i,rtmp)
54 for (i = 0; i < n; ++i) {
55 rtmp=source[i];
56 if (rtmp <0) {
57 sourceN[i]=-rtmp;
58 sourceP[i]=0;
59 } else {
60 sourceN[i]= 0;
61 sourceP[i]= rtmp;
62 }
63 }
64 }
65
66 err_t Paso_FCT_setUpRightHandSide(Paso_FCTransportProblem* fctp, const double dt, const double *u_m, Paso_Coupler* u_m_coupler, double * z_m,
67 Paso_SystemMatrix* flux_matrix, Paso_Coupler* uTilde_coupler, const double *b,
68 Paso_Coupler* QN_coupler, Paso_Coupler* QP_coupler,
69 double *RN_m, Paso_Coupler* RN_m_coupler, double* RP_m, Paso_Coupler* RP_m_coupler, const double *sourceN,
70 Paso_Performance* pp)
71 {
72 dim_t i;
73 const dim_t n=Paso_SystemMatrix_getTotalNumRows(fctp->transport_matrix);
74 register double m, rtmp;
75 /* distribute u */
76 Paso_Coupler_startCollect(u_m_coupler,u_m);
77 Paso_Coupler_finishCollect(u_m_coupler);
78 /*
79 * set the ant diffusion fluxes:
80 *
81 */
82 Paso_FCTransportProblem_setAntiDiffusionFlux(dt,fctp,flux_matrix,u_m_coupler);
83 /*
84 * apply pre flux-correction: f_{ij}:=0 if f_{ij}*(\tilde{u}[i]- \tilde{u}[j])<=0
85 *
86 */
87 Paso_FCTransportProblem_applyPreAntiDiffusionCorrection(flux_matrix,uTilde_coupler);
88 /*
89 * set flux limiters RN_m,RP_m
90 *
91 */
92 Paso_FCTransportProblem_setRs(flux_matrix,fctp->lumped_mass_matrix,QN_coupler,QP_coupler,RN_m,RP_m);
93 Paso_Coupler_startCollect(RN_m_coupler,RN_m);
94 Paso_Coupler_startCollect(RP_m_coupler,RP_m);
95 /*
96 * z_m[i]=b[i] - (m_i*u_m[i] - dt*theta*sum_{j<>i} l_{ij} (u_m[j]-u_m[i]) + dt q^-[i])
97 *
98 * note that iteration_matrix stores the negative values of the
99 * low order transport matrix l therefore a=dt*fctp->theta is used.
100 */
101 Paso_SolverFCT_setMuPaLuPbQ(z_m,fctp->lumped_mass_matrix, u_m_coupler,dt*fctp->theta,fctp->iteration_matrix,dt,sourceN);
102 /* z_m=b-z_m */
103 Paso_Update(n,-1.,z_m,1.,b);
104
105 Paso_Coupler_finishCollect(RN_m_coupler);
106 Paso_Coupler_finishCollect(RP_m_coupler);
107 /* add corrected fluxes into z_m */
108 Paso_FCTransportProblem_addCorrectedFluxes(z_m,flux_matrix,RN_m_coupler,RP_m_coupler);
109 return NO_ERROR;
110 }
111
112 void Paso_SolverFCT_solve(Paso_FCTransportProblem* fctp, double* u, double dt, double* source, Paso_Options* options) {
113 const dim_t FAILURES_MAX=5;
114 err_t error_code;
115 dim_t m,n_substeps, i_substeps, Failed, i, iter;
116 err_t errorCode;
117 double *z_m=NULL, *b_n=NULL, *sourceP=NULL, *sourceN=NULL, *uTilde_n=NULL, *QN_n=NULL, *QP_n=NULL, *RN_m=NULL, *RP_m=NULL, *du_m=NULL;
118 Paso_Coupler *QN_n_coupler=NULL, *QP_n_coupler=NULL, *RN_m_coupler=NULL, *RP_m_coupler=NULL, *uTilde_n_coupler=NULL, *u_m_coupler=NULL;
119 Paso_SystemMatrix *flux_matrix_m=NULL;
120 double dt_max, dt2,t, norm_u_m, omega, norm_du_m, tol;
121 register double mass, rtmp;
122 const dim_t n=Paso_SystemMatrix_getTotalNumRows(fctp->transport_matrix);
123 dim_t blockSize=Paso_FCTransportProblem_getBlockSize(fctp);
124 const double atol=options->absolute_tolerance;
125 const double rtol=options->tolerance;
126 const dim_t max_m=options->iter_max;
127 Paso_Performance pp;
128 bool_t converged=FALSE, max_m_reached=FALSE;
129 if (dt<=0.) {
130 Paso_setError(TYPE_ERROR,"Paso_SolverFCT_solve: dt must be positive.");
131 }
132 dt_max=Paso_FCTransportProblem_getSafeTimeStepSize(fctp);
133 /*
134 * allocate memory
135 *
136 */
137 z_m=TMPMEMALLOC(n,double);
138 Paso_checkPtr(z_m);
139 du_m=TMPMEMALLOC(n,double);
140 Paso_checkPtr(du_m);
141 b_n=TMPMEMALLOC(n,double);
142 Paso_checkPtr(b_n);
143 sourceP=TMPMEMALLOC(n,double);
144 Paso_checkPtr(sourceP);
145 sourceN=TMPMEMALLOC(n,double);
146 Paso_checkPtr(sourceN);
147 uTilde_n=TMPMEMALLOC(n,double);
148 Paso_checkPtr(uTilde_n);
149 QN_n=TMPMEMALLOC(n,double);
150 Paso_checkPtr(QN_n);
151 QP_n=TMPMEMALLOC(n,double);
152 Paso_checkPtr(QP_n);
153 RN_m=TMPMEMALLOC(n,double);
154 Paso_checkPtr(RN_m);
155 RP_m=TMPMEMALLOC(n,double);
156 Paso_checkPtr(RP_m);
157 QN_n_coupler=Paso_Coupler_alloc(Paso_FCTransportProblem_borrowConnector(fctp),blockSize);
158 QP_n_coupler=Paso_Coupler_alloc(Paso_FCTransportProblem_borrowConnector(fctp),blockSize);
159 RN_m_coupler=Paso_Coupler_alloc(Paso_FCTransportProblem_borrowConnector(fctp),blockSize);
160 RP_m_coupler=Paso_Coupler_alloc(Paso_FCTransportProblem_borrowConnector(fctp),blockSize);
161 uTilde_n_coupler=Paso_Coupler_alloc(Paso_FCTransportProblem_borrowConnector(fctp),blockSize);
162 u_m_coupler=Paso_Coupler_alloc(Paso_FCTransportProblem_borrowConnector(fctp),blockSize);
163 flux_matrix_m=Paso_SystemMatrix_alloc(fctp->transport_matrix->type,
164 fctp->transport_matrix->pattern,
165 fctp->transport_matrix->row_block_size,
166 fctp->transport_matrix->col_block_size);
167
168 if (Paso_noError()) {
169 /*
170 * Preparation:
171 *
172 */
173 Paso_FCT_setSource(n, source, sourceN, sourceP);
174 /*
175 * let the show begin!!!!
176 *
177 */
178 t=0;
179 i_substeps=0;
180 Paso_Copy(n,u,fctp->u);
181 norm_u_m=Paso_lsup(n,u, fctp->mpi_info);
182 /* while(i_substeps<n_substeps && Paso_noError()) { */
183 if (fctp->dt_max < LARGE_POSITIVE_FLOAT) {
184 dt2=MIN(dt_max,dt);
185 } else {
186 dt2=dt;
187 }
188 while(t<dt && Paso_noError()) {
189 printf("substep step %d at t=%e (step size= %e)\n",i_substeps+1,t+dt2,dt2);
190 Paso_FCT_setUp(fctp,dt2,sourceN,sourceP,b_n,uTilde_n,uTilde_n_coupler,QN_n,QN_n_coupler,QP_n,QP_n_coupler,
191 options,&pp);
192 /* now the iteration starts */
193 m=0;
194 converged=FALSE;
195 max_m_reached=FALSE;
196 /* tolerance? */
197 while ( (!converged) && (! max_m_reached) && Paso_noError()) {
198 /* set z_m */
199 Paso_FCT_setUpRightHandSide(fctp,dt2,u,u_m_coupler,z_m,flux_matrix_m,uTilde_n_coupler,b_n,
200 QN_n_coupler,QP_n_coupler,RN_m,RN_m_coupler,RP_m,RP_m_coupler,sourceN,&pp);
201 /*
202 * now we solve the linear system to get the correction dt:
203 *
204 */
205 if (fctp->theta > 0) {
206 omega=1./(dt2*fctp->theta);
207 Paso_Solver_solvePreconditioner(fctp->iteration_matrix,du_m,z_m);
208 Paso_Update(n,1.,u,omega,du_m);
209 } else {
210 omega=1;
211 #pragma omp parallel for private(i,mass,rtmp)
212 for (i = 0; i < n; ++i) {
213 mass=fctp->lumped_mass_matrix[i];
214 if (ABS(mass)>0.) {
215 rtmp=z_m[i]/mass;
216 } else {
217 rtmp=0;
218 }
219 du_m[i]=rtmp;
220 u[i]+=rtmp;
221 }
222 }
223 norm_u_m=Paso_lsup(n,u, fctp->mpi_info);
224 norm_du_m=Paso_lsup(n,du_m, fctp->mpi_info)*omega;
225 printf("iteration step %d completed: norm increment= %e (tolerance = %e)\n",m+1, norm_du_m, rtol * norm_u_m + atol);
226
227 max_m_reached=(m>max_m);
228 converged=(norm_du_m <= rtol * norm_u_m + atol);
229 m++;
230 }
231 if (converged) {
232 Failed=0;
233 /* #pragma omp parallel for schedule(static) private(i) */
234 Paso_Copy(n,fctp->u,u);
235 i_substeps++;
236 t+=dt2;
237 if (fctp->dt_max < LARGE_POSITIVE_FLOAT) {
238 dt2=MIN3(dt_max,dt2*1.5,dt-t);
239 } else {
240 dt2=MIN(dt2*1.5,dt-t);
241 }
242 } else if (max_m_reached) {
243 /* if FAILURES_MAX failures in a row: give up */
244 if (Failed > FAILURES_MAX) {
245 Paso_setError(VALUE_ERROR,"Paso_SolverFCT_solve: no convergence after time step reduction.");
246 } else {
247 printf("no convergence in Paso_Solver_NewtonGMRES: Trying smaller time step size.");
248 dt2=dt*0.5;
249 Failed++;
250 }
251 }
252
253 }
254 /*
255 * clean-up:
256 *
257 */
258 MEMFREE(z_m);
259 MEMFREE(du_m);
260 TMPMEMFREE(b_n);
261 Paso_SystemMatrix_free(flux_matrix_m);
262 TMPMEMFREE(sourceP);
263 TMPMEMFREE(sourceN);
264 TMPMEMFREE(uTilde_n);
265 TMPMEMFREE(QN_n);
266 TMPMEMFREE(QP_n);
267 TMPMEMFREE(RN_m);
268 TMPMEMFREE(RP_m);
269 Paso_Coupler_free(QN_n_coupler);
270 Paso_Coupler_free(QP_n_coupler);
271 Paso_Coupler_free(RN_m_coupler);
272 Paso_Coupler_free(RP_m_coupler);
273 Paso_Coupler_free(uTilde_n_coupler);
274 Paso_Coupler_free(u_m_coupler);
275 options->absolute_tolerance=atol;
276 options->tolerance=rtol;
277 }
278 }
279 double Paso_FCTransportProblem_getSafeTimeStepSize(Paso_FCTransportProblem* fctp)
280 {
281 dim_t i, n;
282 double dt_max, dt_max_loc;
283 register double l_ii,m;
284 n=Paso_SystemMatrix_getTotalNumRows(fctp->transport_matrix);
285 if (! fctp->valid_matrices) {
286 fctp->dt_max=LARGE_POSITIVE_FLOAT;
287 /* extract the row sum of the advective part */
288 Paso_SystemMatrix_rowSum(fctp->mass_matrix,fctp->lumped_mass_matrix);
289
290 /* set low order transport operator */
291 Paso_FCTransportProblem_setLowOrderOperator(fctp);
292 /*
293 * calculate time step size:
294 */
295 dt_max=LARGE_POSITIVE_FLOAT;
296 if (fctp->theta < 1.) {
297 #pragma omp parallel private(dt_max_loc)
298 {
299 dt_max_loc=LARGE_POSITIVE_FLOAT;
300 #pragma omp for schedule(static) private(i,l_ii,m)
301 for (i=0;i<n;++i) {
302 l_ii=fctp->main_diagonal_low_order_transport_matrix[i];
303 m=fctp->lumped_mass_matrix[i];
304 if ( (l_ii<0 && m>0.) || (l_ii>0 && m<0) ) {
305 dt_max_loc=MIN(dt_max_loc,-m/l_ii);
306 }
307 }
308 #pragma omp critical
309 {
310 dt_max=MIN(dt_max,dt_max_loc);
311 }
312 }
313 #ifdef PASO_MPI
314 dt_max_loc = dt_max;
315 MPI_Allreduce(&dt_max_loc, &dt_max, 1, MPI_DOUBLE, MPI_MIN, fctp->mpi_info->comm);
316 #endif
317 if (dt_max<LARGE_POSITIVE_FLOAT) dt_max*=1./(1.-fctp->theta);
318 }
319 if (dt_max <= 0.) {
320 Paso_setError(TYPE_ERROR,"Paso_SolverFCT_solve: dt must be positive.");
321 } else {
322 if (dt_max<LARGE_POSITIVE_FLOAT) printf("maximum time step size is %e (theta = %e).\n",dt_max,fctp->theta);
323 }
324 fctp->dt_max=dt_max;
325 fctp->valid_matrices=Paso_noError();
326 }
327 return fctp->dt_max;
328 }
329 void Paso_FCT_setUp(Paso_FCTransportProblem* fctp, const double dt, const double *sourceN, const double *sourceP, double* b, double* uTilde,
330 Paso_Coupler* uTilde_coupler, double *QN, Paso_Coupler* QN_coupler, double *QP, Paso_Coupler* QP_coupler,
331 Paso_Options* options, Paso_Performance* pp)
332 {
333 dim_t i;
334 const dim_t n=Paso_SystemMatrix_getTotalNumRows(fctp->transport_matrix);
335 double omega, factor;
336 register double m, u_tilde_i, rtmp4;
337 /* distribute u */
338 Paso_Coupler_startCollect(fctp->u_coupler,fctp->u);
339 Paso_Coupler_finishCollect(fctp->u_coupler);
340 /*
341 * b^n[i]=m u^n[i] + dt*(1-theta) sum_{j <> i} l_{ij}*(u^n[j]-u^n[i]) + dt*sourceP[i]
342 *
343 * note that iteration_matrix stores the negative values of the
344 * low order transport matrix l therefore a=-dt*(1-fctp->theta) is used.
345 *
346 */
347 Paso_SolverFCT_setMuPaLuPbQ(b,fctp->lumped_mass_matrix,fctp->u_coupler,
348 -dt*(1-fctp->theta),fctp->iteration_matrix,dt,sourceP);
349 /*
350 * uTilde[i]=b[i]/m[i]
351 *
352 * fctp->iteration_matrix[i,i]=m[i]/(dt theta) + \frac{1}{\theta} \frac{q^-[i]}-l[i,i]
353 *
354 */
355 if (fctp->theta > 0) {
356 Paso_solve_free(fctp->iteration_matrix);
357 omega=1./(dt*fctp->theta);
358 factor=dt*omega;
359 #pragma omp parallel for private(i,m,u_tilde_i,rtmp4)
360 for (i = 0; i < n; ++i) {
361 m=fctp->lumped_mass_matrix[i];
362 if (ABS(m)>0.) {
363 u_tilde_i=b[i]/m;
364 } else {
365 u_tilde_i=fctp->u[i];
366 }
367 rtmp4=m*omega-fctp->main_diagonal_low_order_transport_matrix[i];
368 if (ABS(u_tilde_i)>0) rtmp4+=sourceN[i]*factor/u_tilde_i;
369 fctp->iteration_matrix->mainBlock->val[fctp->main_iptr[i]]=rtmp4;
370 uTilde[i]=u_tilde_i;
371 }
372 Performance_startMonitor(pp,PERFORMANCE_PRECONDITIONER_INIT);
373 Paso_Solver_setPreconditioner(fctp->iteration_matrix,options);
374 Performance_stopMonitor(pp,PERFORMANCE_PRECONDITIONER_INIT);
375 } else {
376 #pragma omp parallel for private(i,m,u_tilde_i)
377 for (i = 0; i < n; ++i) {
378 m=fctp->lumped_mass_matrix[i];
379 if (ABS(m)>0.) {
380 u_tilde_i=b[i]/m;
381 } else {
382 u_tilde_i=fctp->u[i];
383 }
384 uTilde[i]=u_tilde_i;
385 }
386 /* no update of iteration_matrix required! */
387 } /* end (fctp->theta > 0) */
388
389 /* distribute uTilde: */
390 Paso_Coupler_startCollect(uTilde_coupler,uTilde);
391 Paso_Coupler_finishCollect(uTilde_coupler);
392 /*
393 * calculate QP[i] max_{j} (\tilde{u}[j]- \tilde{u}[i] )
394 * QN[i] min_{j} (\tilde{u}[j]- \tilde{u}[i] )
395 *
396 */
397 Paso_SolverFCT_setQs(uTilde_coupler,QN,QP,fctp->iteration_matrix);
398 Paso_Coupler_startCollect(QN_coupler,QN);
399 Paso_Coupler_startCollect(QP_coupler,QP);
400 Paso_Coupler_finishCollect(QN_coupler);
401 Paso_Coupler_finishCollect(QP_coupler);
402 }

  ViewVC Help
Powered by ViewVC 1.1.26