/** * @brief Set the source for each FSR and energy group to some value. * @param value the value to assign to each FSR source */ void CPUSolver::flattenFSRSources(FP_PRECISION value) { #pragma omp parallel for schedule(guided) for (int r=0; r < _num_FSRs; r++) { for (int e=0; e < _num_groups; e++) { _source(r,e) = value; _old_source(r,e) = value; } } return; }
/** * @brief Computes the total source (fission and scattering) in each FSR. * @details This method computes the total source in each FSR based on * this iteration's current approximation to the scalar flux. A * residual for the source with respect to the source compute on * the previous iteration is computed and returned. The residual * is determined as follows: * /f$ res = \sqrt{\frac{\displaystyle\sum \displaystyle\sum * \left(\frac{Q^i - Q^{i-1}{Q^i}\right)^2}{\# FSRs}}} \f$ * * @return the residual between this source and the previous source */ FP_PRECISION CPUSolver::computeFSRSources() { int tid; Material* material; FP_PRECISION scatter_source; FP_PRECISION fission_source; FP_PRECISION* nu_sigma_f; FP_PRECISION* sigma_s; FP_PRECISION* sigma_t; FP_PRECISION* chi; FP_PRECISION source_residual = 0.0; FP_PRECISION inverse_k_eff = 1.0 / _k_eff; /* For all FSRs, find the source */ #pragma omp parallel for private(tid, material, nu_sigma_f, chi, \ sigma_s, sigma_t, fission_source, scatter_source) schedule(guided) for (int r=0; r < _num_FSRs; r++) { tid = omp_get_thread_num(); material = _FSR_materials[r]; nu_sigma_f = material->getNuSigmaF(); chi = material->getChi(); sigma_s = material->getSigmaS(); sigma_t = material->getSigmaT(); /* Initialize the source residual to zero */ _source_residuals[r] = 0.; /* Compute fission source for each group */ if (material->isFissionable()) { for (int e=0; e < _num_groups; e++) _fission_sources(r,e) = _scalar_flux(r,e) * nu_sigma_f[e]; fission_source = pairwise_sum<FP_PRECISION>(&_fission_sources(r,0), _num_groups); fission_source *= inverse_k_eff; } else fission_source = 0.0; /* Compute total scattering source for group G */ for (int G=0; G < _num_groups; G++) { scatter_source = 0; for (int g=0; g < _num_groups; g++) _scatter_sources(tid,g) = sigma_s[G*_num_groups+g] * _scalar_flux(r,g); scatter_source=pairwise_sum<FP_PRECISION>(&_scatter_sources(tid,0), _num_groups); /* Set the total source for FSR r in group G */ _source(r,G) = (fission_source * chi[G] + scatter_source) * ONE_OVER_FOUR_PI; _reduced_source(r,G) = _source(r,G) / sigma_t[G]; /* Compute the norm of residual of the source in the FSR */ if (fabs(_source(r,G)) > 1E-10) _source_residuals[r] += pow((_source(r,G) - _old_source(r,G)) / _source(r,G), 2); /* Update the old source */ _old_source(r,G) = _source(r,G); } } /* Sum up the residuals from each FSR */ source_residual = pairwise_sum<FP_PRECISION>(_source_residuals, _num_FSRs); source_residual = sqrt(source_residual / (_num_FSRs * _num_groups)); return source_residual; }
/** * @brief Computes the total source (fission and scattering) in each FSR. * @details This method computes the total source in each FSR based on * this iteration's current approximation to the scalar flux. A * residual for the source with respect to the source compute on * the previous iteration is computed and returned. The residual * is determined as follows: * /f$ res = \sqrt{\frac{\displaystyle\sum \displaystyle\sum * \left(\frac{Q^i - Q^{i-1}{Q^i}\right)^2}{# FSRs}}} \f$ * * @return the residual between this source and the previous source */ FP_PRECISION VectorizedSolver::computeFSRSources() { int tid; FP_PRECISION scatter_source; FP_PRECISION fission_source; FP_PRECISION* nu_sigma_f; FP_PRECISION* sigma_s; FP_PRECISION* sigma_t; FP_PRECISION* chi; Material* material; FP_PRECISION source_residual = 0.0; FP_PRECISION inverse_k_eff = 1.0 / _k_eff; /* For all FSRs, find the source */ #pragma omp parallel for private(material, nu_sigma_f, chi, \ sigma_s, sigma_t, fission_source, scatter_source) schedule(guided) for (int r=0; r < _num_FSRs; r++) { tid = omp_get_thread_num(); material = _FSR_materials[r]; nu_sigma_f = material->getNuSigmaF(); chi = material->getChi(); sigma_s = material->getSigmaS(); sigma_t = material->getSigmaT(); /* Initialize the source residual to zero */ _source_residuals[r] = 0.; /* Compute fission source for each group */ if (material->isFissionable()) { for (int v=0; v < _num_vector_lengths; v++) { /* Compute fission source for each group */ #pragma simd vectorlength(VEC_LENGTH) for (int e=v*VEC_LENGTH; e < (v+1)*VEC_LENGTH; e++) _fission_sources(r,e) = _scalar_flux(r,e) * nu_sigma_f[e]; } #ifdef SINGLE fission_source = cblas_sasum(_num_groups, &_fission_sources(r,0),1); #else fission_source = cblas_dasum(_num_groups, &_fission_sources(r,0),1); #endif fission_source *= inverse_k_eff; } else fission_source = 0.0; /* Compute total scattering source for group G */ for (int G=0; G < _num_groups; G++) { scatter_source = 0; for (int v=0; v < _num_vector_lengths; v++) { #pragma simd vectorlength(VEC_LENGTH) for (int g=v*VEC_LENGTH; g < (v+1)*VEC_LENGTH; g++) _scatter_sources(tid,g) = sigma_s[G*_num_groups+g] * _scalar_flux(r,g); } #ifdef SINGLE scatter_source=cblas_sasum(_num_groups,&_scatter_sources(tid,0),1); #else scatter_source=cblas_dasum(_num_groups,&_scatter_sources(tid,0),1); #endif /* Set the total source for FSR r in group G */ _source(r,G) = (fission_source * chi[G] + scatter_source) * ONE_OVER_FOUR_PI; _reduced_source(r,G) = _source(r,G) / sigma_t[G]; /* Compute the norm of residual of the source in the FSR */ if (fabs(_source(r,G)) > 1E-10) _source_residuals[r] += pow((_source(r,G) - _old_source(r,G)) / _source(r,G), 2); /* Update the old source */ _old_source(r,G) = _source(r,G); } } /* Sum up the residuals from each group and in each FSR */ #ifdef SINGLE source_residual = cblas_sasum(_num_FSRs,_source_residuals,1); #else source_residual = cblas_dasum(_num_FSRs,_source_residuals,1); #endif source_residual = sqrt(source_residual / (_num_groups * _num_FSRs)); return source_residual; }