4646 *
4747 */
4848
49- #define FINE_NODE_EXTRAPOLATED_PROLONGATION () \
50- do { \
51- if (i_r & 1 ) { \
52- if (i_theta & 1 ) { \
53- /* (odd, odd) -> node in center of coarse cell */ \
54- double value = \
55- 0.5 * (coarse_values[coarse_grid.index (i_r_coarse + 1 , i_theta_coarse)] + /* Bottom right */ \
56- coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse + 1 )] /* Top left */ \
57- ); \
58- fine_result[fine_grid.index (i_r, i_theta)] = value; \
59- } \
60- else { \
61- /* (odd, even) -> between coarse nodes in radial direction */ \
62- double value = 0.5 * (coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse)] + /* Left */ \
63- coarse_values[coarse_grid.index (i_r_coarse + 1 , i_theta_coarse)] /* Right */ \
64- ); \
65- fine_result[fine_grid.index (i_r, i_theta)] = value; \
66- } \
67- } \
68- else { \
69- if (i_theta & 1 ) { \
70- /* (even, odd) -> between coarse nodes in angular direction */ \
71- double value = 0.5 * (coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse)] + /* Bottom */ \
72- coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse + 1 )] /* Top */ \
73- ); \
74- fine_result[fine_grid.index (i_r, i_theta)] = value; \
75- } \
76- else { \
77- /* (even, even) -> node lies exactly on coarse grid */ \
78- fine_result[fine_grid.index (i_r, i_theta)] = \
79- coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse)]; /* Center */ \
80- } \
81- } \
82- } while (0 )
49+ inline void fineNodeExtrapolatedProlongation (int i_r, int i_theta, int i_r_coarse, int i_theta_coarse,
50+ const PolarGrid& coarse_grid, const PolarGrid& fine_grid,
51+ Vector<double >& fine_result, ConstVector<double >& coarse_values)
52+ {
53+ if (i_r & 1 ) {
54+ if (i_theta & 1 ) { /* (odd, odd) -> node in center of coarse cell */
55+ double value = 0.5 * (coarse_values[coarse_grid.index (i_r_coarse + 1 , i_theta_coarse)] + /* Bottom right */
56+ coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse + 1 )] /* Top left */
57+ );
58+ fine_result[fine_grid.index (i_r, i_theta)] = value;
59+ }
60+ else { /* (odd, even) -> between coarse nodes in radial direction */
61+ double value = 0.5 * (coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse)] + /* Left */
62+ coarse_values[coarse_grid.index (i_r_coarse + 1 , i_theta_coarse)] /* Right */
63+ );
64+ fine_result[fine_grid.index (i_r, i_theta)] = value;
65+ }
66+ }
67+ else {
68+ if (i_theta & 1 ) { /* (even, odd) -> between coarse nodes in angular direction */
69+ double value = 0.5 * (coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse)] + /* Bottom */
70+ coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse + 1 )] /* Top */
71+ );
72+ fine_result[fine_grid.index (i_r, i_theta)] = value;
73+ }
74+ else { /* (even, even) -> node lies exactly on coarse grid */
75+ fine_result[fine_grid.index (i_r, i_theta)] =
76+ coarse_values[coarse_grid.index (i_r_coarse, i_theta_coarse)]; /* Center */
77+ }
78+ }
79+ }
8380
8481void Interpolation::applyExtrapolatedProlongation (const PolarGrid& coarse_grid, const PolarGrid& fine_grid,
8582 Vector<double > fine_result, ConstVector<double > coarse_values) const
@@ -98,7 +95,8 @@ void Interpolation::applyExtrapolatedProlongation(const PolarGrid& coarse_grid,
9895 int i_r_coarse = i_r / 2 ;
9996 for (int i_theta = 0 ; i_theta < fine_grid.ntheta (); i_theta++) {
10097 int i_theta_coarse = i_theta / 2 ;
101- FINE_NODE_EXTRAPOLATED_PROLONGATION ();
98+ fineNodeExtrapolatedProlongation (i_r, i_theta, i_r_coarse, i_theta_coarse, coarse_grid, fine_grid,
99+ fine_result, coarse_values);
102100 }
103101 }
104102
@@ -108,7 +106,8 @@ void Interpolation::applyExtrapolatedProlongation(const PolarGrid& coarse_grid,
108106 int i_theta_coarse = i_theta / 2 ;
109107 for (int i_r = fine_grid.numberSmootherCircles (); i_r < fine_grid.nr (); i_r++) {
110108 int i_r_coarse = i_r / 2 ;
111- FINE_NODE_EXTRAPOLATED_PROLONGATION ();
109+ fineNodeExtrapolatedProlongation (i_r, i_theta, i_r_coarse, i_theta_coarse, coarse_grid, fine_grid,
110+ fine_result, coarse_values);
112111 }
113112 }
114113 }
0 commit comments