Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
AuxArray.c
Go to the documentation of this file.
1
14#include <math.h>
15
16#ifdef _OPENMP
17#include <omp.h>
18#endif
19
20#include "fasp.h"
21#include "fasp_functs.h"
22
23/*---------------------------------*/
24/*-- Public Functions --*/
25/*---------------------------------*/
26
41void fasp_darray_set(const INT n, REAL* x, const REAL val)
42{
43 SHORT use_openmp = FALSE;
44
45#ifdef _OPENMP
46 INT nthreads = 1;
47
48 if (n > OPENMP_HOLDS) {
49 use_openmp = TRUE;
50 nthreads = fasp_get_num_threads();
51 }
52#endif
53
54 if (val == 0.0) {
55 if (use_openmp) {
56#ifdef _OPENMP
57 INT mybegin, myend, myid;
58#pragma omp parallel for private(myid, mybegin, myend)
59 for (myid = 0; myid < nthreads; myid++) {
60 fasp_get_start_end(myid, nthreads, n, &mybegin, &myend);
61 memset(&x[mybegin], 0x0, sizeof(REAL) * (myend - mybegin));
62 }
63#endif
64 } else
65 memset(x, 0x0, sizeof(REAL) * n);
66 } else {
67 INT i;
68
69 if (use_openmp) {
70#ifdef _OPENMP
71 INT mybegin, myend, myid;
72#pragma omp parallel for private(myid, mybegin, myend, i)
73 for (myid = 0; myid < nthreads; myid++) {
74 fasp_get_start_end(myid, nthreads, n, &mybegin, &myend);
75 for (i = mybegin; i < myend; ++i) x[i] = val;
76 }
77#endif
78 } else {
79 for (i = 0; i < n; ++i) x[i] = val;
80 }
81 }
82}
83
98void fasp_iarray_set(const INT n, INT* x, const INT val)
99{
100 SHORT use_openmp = FALSE;
101
102#ifdef _OPENMP
103 INT nthreads = 1;
104
105 if (n > OPENMP_HOLDS) {
106 use_openmp = TRUE;
107 nthreads = fasp_get_num_threads();
108 }
109#endif
110
111 if (val == 0) {
112 if (use_openmp) {
113#ifdef _OPENMP
114 INT mybegin, myend, myid;
115#pragma omp parallel for private(myid, mybegin, myend)
116 for (myid = 0; myid < nthreads; myid++) {
117 fasp_get_start_end(myid, nthreads, n, &mybegin, &myend);
118 memset(&x[mybegin], 0, sizeof(INT) * (myend - mybegin));
119 }
120#endif
121 } else {
122 memset(x, 0, sizeof(INT) * n);
123 }
124 } else {
125 INT i;
126
127 if (use_openmp) {
128#ifdef _OPENMP
129 INT mybegin, myend, myid;
130#pragma omp parallel for private(myid, mybegin, myend, i)
131 for (myid = 0; myid < nthreads; myid++) {
132 fasp_get_start_end(myid, nthreads, n, &mybegin, &myend);
133 for (i = mybegin; i < myend; ++i) x[i] = val;
134 }
135#endif
136 } else {
137 for (i = 0; i < n; ++i) x[i] = val;
138 }
139 }
140}
141
154void fasp_ldarray_set(const INT n, LONGREAL* x, const LONGREAL val)
155{
156 SHORT use_openmp = FALSE;
157
158#ifdef _OPENMP
159 INT nthreads = 1;
160
161 if (n > OPENMP_HOLDS) {
162 use_openmp = TRUE;
163 nthreads = fasp_get_num_threads();
164 }
165#endif
166
167 if (val == 0) {
168 if (use_openmp) {
169#ifdef _OPENMP
170 INT mybegin, myend, myid;
171#pragma omp parallel for private(myid, mybegin, myend)
172 for (myid = 0; myid < nthreads; myid++) {
173 fasp_get_start_end(myid, nthreads, n, &mybegin, &myend);
174 memset(&x[mybegin], 0, sizeof(LONGREAL) * (myend - mybegin));
175 }
176#endif
177 } else {
178 memset(x, 0, sizeof(LONGREAL) * n);
179 }
180 } else {
181 INT i;
182
183 if (use_openmp) {
184#ifdef _OPENMP
185 INT mybegin, myend, myid;
186#pragma omp parallel for private(myid, mybegin, myend, i)
187 for (myid = 0; myid < nthreads; myid++) {
188 fasp_get_start_end(myid, nthreads, n, &mybegin, &myend);
189 for (i = mybegin; i < myend; ++i) x[i] = val;
190 }
191#endif
192 } else {
193 for (i = 0; i < n; ++i) x[i] = val;
194 }
195 }
196}
197
210void fasp_darray_cp(const INT n, const REAL* x, REAL* y)
211{
212 memcpy(y, x, n * sizeof(REAL));
213}
214
227void fasp_iarray_cp(const INT n, const INT* x, INT* y)
228{
229 memcpy(y, x, n * sizeof(INT));
230}
231
232/*---------------------------------*/
233/*-- End of File --*/
234/*---------------------------------*/
void fasp_darray_set(const INT n, REAL *x, const REAL val)
Set initial value for an array to be x=val.
Definition: AuxArray.c:41
void fasp_iarray_set(const INT n, INT *x, const INT val)
Set initial value for an array to be x=val.
Definition: AuxArray.c:98
void fasp_ldarray_set(const INT n, LONGREAL *x, const LONGREAL val)
Set initial value for an array to be x=val.
Definition: AuxArray.c:154
void fasp_darray_cp(const INT n, const REAL *x, REAL *y)
Copy an array to the other y=x.
Definition: AuxArray.c:210
void fasp_iarray_cp(const INT n, const INT *x, INT *y)
Copy an array to the other y=x.
Definition: AuxArray.c:227
void fasp_get_start_end(const INT procid, const INT nprocs, const INT n, INT *start, INT *end)
Assign Load to each thread.
Definition: AuxThreads.c:93
Main header file for the FASP project.
#define REAL
Definition: fasp.h:75
#define SHORT
FASP integer and floating point numbers.
Definition: fasp.h:71
#define LONGREAL
Definition: fasp.h:76
#define INT
Definition: fasp.h:72
#define OPENMP_HOLDS
Definition: fasp_const.h:269
#define TRUE
Definition of logic type.
Definition: fasp_const.h:61
#define FALSE
Definition: fasp_const.h:62