Line data Source code
1 : /********************************************************************/ 2 : /* SOFTWARE COPYRIGHT NOTIFICATION */ 3 : /* Cardinal */ 4 : /* */ 5 : /* (c) 2021 UChicago Argonne, LLC */ 6 : /* ALL RIGHTS RESERVED */ 7 : /* */ 8 : /* Prepared by UChicago Argonne, LLC */ 9 : /* Under Contract No. DE-AC02-06CH11357 */ 10 : /* With the U. S. Department of Energy */ 11 : /* */ 12 : /* Prepared by Battelle Energy Alliance, LLC */ 13 : /* Under Contract No. DE-AC07-05ID14517 */ 14 : /* With the U. S. Department of Energy */ 15 : /* */ 16 : /* See LICENSE for full restrictions */ 17 : /********************************************************************/ 18 : 19 : #ifdef ENABLE_NEK_COUPLING 20 : 21 : #include "NekVolumetricSource.h" 22 : 23 : registerMooseObject("CardinalApp", NekVolumetricSource); 24 : 25 : InputParameters 26 198 : NekVolumetricSource::validParams() 27 : { 28 198 : auto params = ConservativeFieldTransfer::validParams(); 29 396 : params.addParam<Real>( 30 : "initial_source_integral", 31 396 : 0, 32 : "Initial value to use for the 'postprocessor_to_conserve'; this initial value will be " 33 : "overridden once the coupled app executes its transfer of the volumetric source term " 34 : "integral into the 'postprocessor_to_conserve'. You may want to use this parameter if NekRS " 35 : "runs first, or if you are running NekRS in isolation but still want to apply a source term " 36 : "via Cardinal. Remember that this parameter is only used to normalize the source term " 37 : "'source_variable', so you will need to populate an initial shape."); 38 198 : params.addClassDescription("Reads/writes volumetric source data between NekRS and MOOSE."); 39 198 : return params; 40 0 : } 41 : 42 102 : NekVolumetricSource::NekVolumetricSource(const InputParameters & parameters) 43 : : ConservativeFieldTransfer(parameters), 44 204 : _initial_source_integral(getParam<Real>("initial_source_integral")) 45 : { 46 102 : if (_direction == "to_nek") 47 : { 48 203 : addExternalVariable(_usrwrk_slot[0], _variable); 49 101 : indices.heat_source = _usrwrk_slot[0] * nekrs::fieldOffset(); 50 : 51 101 : if (_usrwrk_slot.size() > 1) 52 2 : paramError("usrwrk_slot", 53 : "'usrwrk_slot' must be of length 1 for volumetric source transfers; you have " 54 1 : "entered a vector of length " + 55 0 : Moose::stringify(_usrwrk_slot.size())); 56 : } 57 : 58 100 : if (!_nek_mesh->volume()) 59 1 : mooseError("The NekVolumetricSource object can only be used when there is volumetric coupling " 60 : "of NekRS with MOOSE, i.e. when 'volume = true' in NekRSMesh."); 61 : 62 99 : if (_direction == "from_nek") 63 0 : paramError("direction", 64 : "The NekVolumetricSource currently only supports transfers 'to_nek'; contact the " 65 : "Cardinal developer team if you require reading of NekRS volumetric source terms."); 66 : 67 : // Check that there is a udf function providing the source for the passive scalar 68 : // equations. NOTE: This check is imperfect, because even if there is a source kernel, 69 : // we cannot tell _which_ passive scalar equation that it is applied to (we have 70 : // source kernels for the RANS passive scalar equations, for instance). 71 99 : if (nekrs::hasTemperatureSolve() && !nekrs::hasHeatSourceKernel()) 72 1 : mooseError("In order to send a volumetric heat source to NekRS, you must have an OCCA source " 73 : "kernel in the passive scalar equations!"); 74 : 75 98 : if (!nekrs::hasTemperatureVariable()) 76 1 : mooseError("In order to send a volumetric heat source to NekRS, your case files must have a " 77 1 : "[TEMPERATURE] block. Note that you can set 'solver = none' in '" + 78 1 : _nek_problem.casename() + ".par' if you don't want to solve for temperature."); 79 : 80 97 : if (!nekrs::hasTemperatureSolve()) 81 14 : mooseWarning("By setting 'solver = none' for temperature in '" + _nek_problem.casename() + 82 : ".par', NekRS will not solve for temperature. The volumetric heat source sent by " 83 : "this object will be unused."); 84 : 85 96 : addExternalPostprocessor(_postprocessor_name, _initial_source_integral); 86 96 : _source_integral = &getPostprocessorValueByName(_postprocessor_name); 87 : 88 96 : _source_elem = (double *)calloc(_n_per_vol, sizeof(double)); 89 96 : } 90 : 91 192 : NekVolumetricSource::~NekVolumetricSource() { freePointer(_source_elem); } 92 : 93 : bool 94 1132 : NekVolumetricSource::normalizeVolumetricSource(const double moose, 95 : double nek, 96 : double & normalized_nek) 97 : { 98 : auto dimension_multiplier = 99 1132 : nekrs::referenceVolume() * nekrs::nondimensionalDivisor(field::heat_source); 100 : 101 : // scale the nek source to dimensional form for the sake of normalizing against 102 : // a dimensional MOOSE source 103 1132 : nek *= dimension_multiplier; 104 : 105 : // avoid divide-by-zero 106 1132 : if (std::abs(nek) < _abs_tol) 107 : return true; 108 : 109 1108 : nekrs::scaleUsrwrk(indices.heat_source, moose / nek); 110 : 111 : // check that the normalization worked properly 112 1108 : normalized_nek = 113 1108 : nekrs::usrwrkVolumeIntegral(indices.heat_source, nek_mesh::all) * dimension_multiplier; 114 1108 : bool low_rel_err = std::abs(normalized_nek - moose) / moose < _rel_tol; 115 1108 : bool low_abs_err = std::abs(normalized_nek - moose) < _abs_tol; 116 : 117 1108 : return low_rel_err && low_abs_err; 118 : } 119 : 120 : void 121 1132 : NekVolumetricSource::sendDataToNek() 122 : { 123 1132 : _console << "Sending volumetric source to NekRS..." << std::endl; 124 : 125 1132 : auto d = nekrs::nondimensionalDivisor(field::heat_source); 126 1132 : auto a = nekrs::nondimensionalAdditive(field::heat_source); 127 3157788 : for (unsigned int e = 0; e < _nek_mesh->numVolumeElems(); e++) 128 : { 129 : // We can only write into the nekRS scratch space if that face is "owned" by the current process 130 3156656 : if (nekrs::commRank() != _nek_mesh->volumeCoupling().processor_id(e)) 131 2759908 : continue; 132 : 133 396748 : _nek_problem.mapVolumeDataToNekVolume(e, _variable_number[_variable], d, a, &_source_elem); 134 396748 : _nek_problem.writeVolumeSolution(e, field::heat_source, _source_elem); 135 : } 136 : 137 : // Because the NekRSMesh may be quite different from that used in the app solving for 138 : // the heat source, we will need to normalize the total source on the nekRS side by the 139 : // total source computed by the coupled MOOSE app. 140 1132 : const Real scale_cubed = _nek_mesh->scaling() * _nek_mesh->scaling() * _nek_mesh->scaling(); 141 1132 : const double nek_source = nekrs::usrwrkVolumeIntegral(indices.heat_source, nek_mesh::all); 142 1132 : const double moose_source = *_source_integral; 143 : 144 : // For the sake of printing diagnostics to the screen regarding source normalization, 145 : // we first scale the nek source by any unit changes and then by the reference source 146 : const double nek_source_print_mult = 147 1132 : scale_cubed * nekrs::nondimensionalDivisor(field::heat_source); 148 1132 : double normalized_nek_source = 0.0; 149 : bool successful_normalization; 150 : 151 1132 : _console << "[volume]: Normalizing total NekRS heat source of " 152 1132 : << Moose::stringify(nek_source * nek_source_print_mult) 153 2264 : << " to the conserved MOOSE value of " + Moose::stringify(moose_source) << std::endl; 154 : 155 : // Any unit changes (for DIMENSIONAL nekRS runs) are automatically accounted for 156 : // here because moose_source is an integral on the MOOSE mesh, while nek_source is 157 : // an integral on the nek mesh 158 : successful_normalization = 159 1132 : normalizeVolumetricSource(moose_source, nek_source, normalized_nek_source); 160 : 161 : // If before normalization, there is a large difference between the nekRS imposed source 162 : // and the MOOSE source, this could mean that there is a poor match between the domains, 163 : // even if neither value is zero. For instance, if you forgot that the nekRS mesh is in 164 : // units of centimeters, but you're coupling to an app based in meters, the sources will 165 : // be very different from one another. 166 1132 : if (moose_source && 167 1132 : (std::abs(nek_source * nek_source_print_mult - moose_source) / moose_source) > 0.25) 168 24 : mooseDoOnce( 169 : mooseWarning("NekRS source differs from MOOSE source by more than 25\%! This is NOT " 170 : "necessarily a problem - but it could indicate that your geometries don't " 171 : "line up properly or something is amiss with your transfer. We recommend " 172 : "opening the output files to visually inspect the volumetric source in both " 173 : "the main and sub applications to check that the fields look correct.")); 174 : 175 1132 : if (!successful_normalization) 176 0 : mooseError("Volumetric source normalization process failed! NekRS integrated source: ", 177 : normalized_nek_source, 178 : " MOOSE integrated source: ", 179 : moose_source, 180 : ".\n\n", 181 0 : normalizationHint()); 182 1132 : } 183 : 184 : #endif