1 /* Copyright (C) 2005 The Scalable Software Infrastructure Project. All rights reserved.
2 
3    Redistribution and use in source and binary forms, with or without
4    modification, are permitted provided that the following conditions are met:
5    1. Redistributions of source code must retain the above copyright
6       notice, this list of conditions and the following disclaimer.
7    2. Redistributions in binary form must reproduce the above copyright
8       notice, this list of conditions and the following disclaimer in the
9       documentation and/or other materials provided with the distribution.
10    3. Neither the name of the project nor the names of its contributors
11       may be used to endorse or promote products derived from this software
12       without specific prior written permission.
13 
14    THIS SOFTWARE IS PROVIDED BY THE SCALABLE SOFTWARE INFRASTRUCTURE PROJECT
15    ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
16    TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
17    PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE SCALABLE SOFTWARE INFRASTRUCTURE
18    PROJECT BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
19    OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
20    SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
21    INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
22    CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
23    ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
24    POSSIBILITY OF SUCH DAMAGE.
25 */
26 
27 #ifdef HAVE_CONFIG_H
28 	#include "lis_config.h"
29 #else
30 #ifdef HAVE_CONFIG_WIN_H
31 	#include "lis_config_win.h"
32 #endif
33 #endif
34 
35 #include <stdio.h>
36 #include <stdlib.h>
37 #ifdef HAVE_MALLOC_H
38         #include <malloc.h>
39 #endif
40 #include <string.h>
41 #include <stdarg.h>
42 #ifdef _OPENMP
43 	#include <omp.h>
44 #endif
45 #ifdef USE_MPI
46 	#include <mpi.h>
47 #endif
48 #include "lislib.h"
49 
50 /************************************************
51  * lis_precon_create
52  * lis_psolve
53  * lis_psolveh
54  ************************************************/
55 
56 #undef __FUNC__
57 #define __FUNC__ "lis_precon_create_ssor"
lis_precon_create_ssor(LIS_SOLVER solver,LIS_PRECON precon)58 LIS_INT lis_precon_create_ssor(LIS_SOLVER solver, LIS_PRECON precon)
59 {
60 	LIS_INT	err;
61 	LIS_SCALAR w;
62 	LIS_MATRIX A;
63 
64 	LIS_DEBUG_FUNC_IN;
65 
66 	A   = solver->A;
67 	w   = solver->params[LIS_PARAMS_SSOR_OMEGA-LIS_OPTIONS_LEN];
68 
69 	err = lis_matrix_convert_self(solver);
70 	if( err ) return err;
71 
72 	err = lis_matrix_split(A);
73 	if( err )
74 	{
75 		return err;
76 	}
77 	if( A->use_wd!=LIS_SOLVER_SOR )
78 	{
79 		if( !A->WD )
80 		{
81 			err = lis_matrix_diag_duplicate(A->D,&A->WD);
82 			if( err ) return err;
83 		}
84 		lis_matrix_diag_copy(A->D,A->WD);
85 		lis_matrix_diag_scale(w,A->WD);
86 		lis_matrix_diag_inverse(A->WD);
87 		A->use_wd = LIS_SOLVER_SOR;
88 	}
89 
90 	precon->A       = A;
91 	precon->is_copy = LIS_FALSE;
92 
93 	LIS_DEBUG_FUNC_OUT;
94     return LIS_SUCCESS;
95 }
96 
97 #undef __FUNC__
98 #define __FUNC__ "lis_psolve_ssor"
lis_psolve_ssor(LIS_SOLVER solver,LIS_VECTOR B,LIS_VECTOR X)99 LIS_INT lis_psolve_ssor(LIS_SOLVER solver, LIS_VECTOR B, LIS_VECTOR X)
100 {
101 	LIS_MATRIX A;
102 
103 	/*
104 	 *  Mx = b
105 	 *  M  = (D/w + L) * (I + w*D^-1 * U)
106 	 */
107 
108 	LIS_DEBUG_FUNC_IN;
109 
110 	A = solver->precon->A;
111 
112 	lis_matrix_solve(A,B,X,LIS_MATRIX_SSOR);
113 
114 
115 	LIS_DEBUG_FUNC_OUT;
116 	return LIS_SUCCESS;
117 }
118 
119 #undef __FUNC__
120 #define __FUNC__ "lis_psolveh_ssor"
lis_psolveh_ssor(LIS_SOLVER solver,LIS_VECTOR B,LIS_VECTOR X)121 LIS_INT lis_psolveh_ssor(LIS_SOLVER solver, LIS_VECTOR B, LIS_VECTOR X)
122 {
123 	LIS_MATRIX A;
124 
125 	/*
126 	 *  M'x = b
127 	 *  M'  = (I + U' * w*D^-1) * (D/w + L')
128 	 */
129 
130 	LIS_DEBUG_FUNC_IN;
131 
132 	A = solver->precon->A;
133 
134 	lis_matrix_solveh(A,B,X,LIS_MATRIX_SSOR);
135 
136 	LIS_DEBUG_FUNC_OUT;
137 	return LIS_SUCCESS;
138 }
139