/// <summary> /// writes the dammed result of the integration to the sparse matrix /// </summary> protected override void SaveIntegrationResults(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { int GAMMA = this.m_CodomainMap.BasisS.Count; // GAMMA: number of codom/row/test variables LECQuadratureLevelSet <IMutableMatrixEx, double[]> .CompOffsets(i0, Length, out int[] offsetRow, out int[] RowNonxN, m_CodomainMap); SpeciesId[] spcS = new[] { this.SpeciesA, this.SpeciesB }; int[] _i0aff = new int[2]; int[] _iEaff = new int[2]; // loop over cells... for (int i = 0; i < Length; i++) { int jCell = i + i0; #if DEBUG Debug.Assert(RowNonxN.ListEquals(m_CodomainMap.GetNonXBasisLengths(jCell))); #endif _i0aff[0] = i; _iEaff[0] = -1; for (int gamma = 0; gamma < GAMMA; gamma++) // loop over rows... { for (int cr = 0; cr < 2; cr++) // loop over neg/pos species row... { SpeciesId rowSpc = spcS[cr]; int Row0 = m_CodomainMap.LocalUniqueCoordinateIndex(m_lsTrk, gamma, jCell, rowSpc, 0); _i0aff[1] = offsetRow[gamma] + RowNonxN[gamma] * cr; // the 'RowXbSw' is 0 for non-xdg, so both species will be added _iEaff[1] = _i0aff[1] + RowNonxN[gamma] - 1; var BlockRes = ResultsOfIntegration.ExtractSubArrayShallow(_i0aff, _iEaff); for (int r = BlockRes.GetLength(0) - 1; r >= 0; r--) { ResultVector[Row0 + r] += BlockRes[r]; } } } } }
protected override void Evaluate(int i0, int Len, QuadRule QR, MultidimensionalArray EvalResult) { NodeSet QuadNodes = QR.Nodes; int D = gridData.SpatialDimension; int NoOfNodes = QuadNodes.NoOfNodes; int GAMMA = m_CodomainMap.NoOfVariables; // GAMMA: number of codom variables // Evaluate Domain & Parameter fields // -------------------------------- Field_Eval.Start(); for (int i = 0; i < m_DomainAndParamFields.Length; i++) { if (m_ValueRequired[i]) { DGField _Field = m_DomainAndParamFields[i]; if (_Field != null) { if (_Field is XDGField) { // jump in parameter i at level-set: separate evaluation for both sides var _xField = _Field as XDGField; _xField.GetSpeciesShadowField(this.SpeciesA).Evaluate(i0, Len, QuadNodes, m_FieldValuesNeg[i]); _xField.GetSpeciesShadowField(this.SpeciesB).Evaluate(i0, Len, QuadNodes, m_FieldValuesPos[i]); } else { // no jump at level set: positive and negative limit of parameter i are equal _Field.Evaluate(i0, Len, QuadNodes, m_FieldValuesPos[i]); m_FieldValuesNeg[i].Set(m_FieldValuesPos[i]); } } else { m_FieldValuesPos[i].Clear(); m_FieldValuesNeg[i].Clear(); } } if (m_GradientRequired[i]) { DGField _Field = m_DomainAndParamFields[i]; if (_Field != null) { if (_Field is XDGField) { // jump in parameter i at level-set: separate evaluation for both sides var _xField = _Field as XDGField; _xField.GetSpeciesShadowField(this.SpeciesA).EvaluateGradient(i0, Len, QuadNodes, m_FieldGradientValuesNeg[i]); _xField.GetSpeciesShadowField(this.SpeciesB).EvaluateGradient(i0, Len, QuadNodes, m_FieldGradientValuesPos[i]); } else { // no jump at level set: positive and negative limit of parameter i are equal _Field.EvaluateGradient(i0, Len, QuadNodes, m_FieldGradientValuesPos[i]); m_FieldGradientValuesNeg[i].Set(m_FieldGradientValuesPos[i]); } } else { m_FieldGradientValuesPos[i].Clear(); m_FieldGradientValuesNeg[i].Clear(); } } } Field_Eval.Stop(); // Evaluate level sets and normals // ------------------------------- ParametersAndNormals.Start(); var NoOfLevSets = m_lsTrk.LevelSets.Count; MultidimensionalArray Normals = m_lsTrk.DataHistories[m_LevSetIdx].Current.GetLevelSetNormals(QuadNodes, i0, Len); MultidimensionalArray NodesGlobal = gridData.GlobalNodes.GetValue_Cell(QuadNodes, i0, Len); ParametersAndNormals.Stop(); // Evaluate basis and test functions // --------------------------------- bool[] ReqV = new bool[GAMMA]; bool[] ReqGradV = new bool[GAMMA]; for (int gamma = 0; gamma < GAMMA; gamma++) { if (Koeff_V[gamma] != null) { ReqV[gamma] = true; } if (Koeff_GradV[gamma] != null) { ReqGradV[gamma] = true; } } MultidimensionalArray[] TestValues; // index: codom variable/test function MultidimensionalArray[] TestGradientValues; // index: codom variable/test function int[,] sectionsTest; Basis_Eval.Start(); LECQuadratureLevelSet <IMutableMatrixEx, double[]> .EvalBasis(i0, Len, this.m_CodomainMap.BasisS, ReqV, ReqGradV, out TestValues, out TestGradientValues, out sectionsTest, QuadNodes); Basis_Eval.Stop(); // Evaluate Integral components // ---------------------------- // loop over codomain variables ... for (int gamma = 0; gamma < GAMMA; gamma++) { // prepare parameters // - - - - - - - - - // set Normal's LevSetIntParams _inParams = new LevSetIntParams(); _inParams.Normal = Normals; // set Nodes Global _inParams.X = NodesGlobal; _inParams.time = this.time; _inParams.LsTrk = this.m_lsTrk; _inParams.i0 = i0; Debug.Assert(_inParams.Len == Len); // clear summation buffers // - - - - - - - - - - - - if (Koeff_V[gamma] != null) { Koeff_V[gamma].Clear(); } if (Koeff_GradV[gamma] != null) { Koeff_GradV[gamma].Clear(); } // Evaluate Bilin. forms // - - - - - - - - - - - { EvalComponent(_inParams, gamma, this.m_NonlinLsForm_V[gamma], this.m_NonlinLsForm_V_Watches[gamma], Koeff_V[gamma], m_FieldValuesPos, m_FieldValuesNeg, m_FieldGradientValuesPos, m_FieldGradientValuesNeg, DELTA, Flux_Eval, delegate(INonlinLevelSetForm_V _comp, LevSetIntParams inp, MultidimensionalArray[] uA, MultidimensionalArray[] uB, MultidimensionalArray[] Grad_uA, MultidimensionalArray[] Grad_uB, MultidimensionalArray SumBuf) { _comp.LevelSetForm_V(_inParams, uA, uB, Grad_uA, Grad_uB, SumBuf); }); } { EvalComponent(_inParams, gamma, this.m_NonlinLsForm_GradV[gamma], this.m_NonlinLsForm_GradV_Watches[gamma], Koeff_GradV[gamma], m_FieldValuesPos, m_FieldValuesNeg, m_FieldGradientValuesPos, m_FieldGradientValuesNeg, DELTA, Flux_Eval, delegate(INonlinLevelSetForm_GradV _comp, LevSetIntParams inp, MultidimensionalArray[] uA, MultidimensionalArray[] uB, MultidimensionalArray[] Grad_uA, MultidimensionalArray[] Grad_uB, MultidimensionalArray SumBuf) { _comp.LevelSetForm_GradV(_inParams, uA, uB, Grad_uA, Grad_uB, SumBuf); }); } } // Summation Loops: multiply with test and trial functions // ------------------------------------------------------- int[] offsetCod = new int[GAMMA]; LECQuadratureLevelSet <IMutableMatrixEx, double[]> . CompOffsets(i0, Len, offsetCod, m_CodomainMap); for (int gamma = 0; gamma < GAMMA; gamma++) { // Evaluate Integrand // - - - - - - - - - var TestVal = TestValues[gamma]; var TestGradVal = TestGradientValues[gamma]; int N; if (TestVal != null) { N = TestVal.GetLength(2); } else if (TestGradVal != null) { N = TestGradVal.GetLength(2); } else { N = 0; } Loops.Start(); // affine offset for (int cr = 0; cr < 2; cr++) // loop over negative/positive species { int[] extr0 = new int[] { 0, 0, sectionsTest[gamma, cr] * N + offsetCod[gamma] }; int[] extrE = new int[] { Len - 1, NoOfNodes - 1, extr0[2] + N - 1 }; var SubRes = EvalResult.ExtractSubArrayShallow(extr0, extrE); if (Koeff_V[gamma] != null) { var Sum_Koeff_V_Cr = Koeff_V[gamma].ExtractSubArrayShallow(-1, -1, cr); SubRes.Multiply(1.0, Sum_Koeff_V_Cr, TestVal, 1.0, "jkn", "jk", "jkn"); } if (Koeff_GradV[gamma] != null) { var Sum_Koeff_NablaV_Cr = Koeff_GradV[gamma].ExtractSubArrayShallow(-1, -1, cr, -1); SubRes.Multiply(1.0, Sum_Koeff_NablaV_Cr, TestGradVal, 1.0, "jkn", "jkd", "jknd"); } } Loops.Stop(); } }
/// <summary> /// ctor. /// </summary> internal NECQuadratureLevelSet(IGridData context, XSpatialOperatorMk2 DiffOp, V __ResultVector, IList <DGField> __DomainFields, IList <DGField> __Parameters, UnsetteledCoordinateMapping CodomainMap, LevelSetTracker lsTrk, int _iLevSet, Tuple <SpeciesId, SpeciesId> SpeciesPair, ICompositeQuadRule <QuadRule> domAndRule) // : base(new int[] { CodomainMap.GetNonXBasisLengths(0).Sum() * 2 }, // we always integrate over species in pairs (neg + pos), so we need to alloc mem only 2 species context, domAndRule) // { MPICollectiveWatchDog.Watch(); // ----------------------------------- // set members / check ctor parameters // ----------------------------------- m_lsTrk = lsTrk; this.m_LevSetIdx = _iLevSet; this.m_SpeciesPair = SpeciesPair; this.ResultVector = __ResultVector; m_CodomainMap = CodomainMap; var _Parameters = (__Parameters != null) ? __Parameters.ToArray() : new DGField[0]; if (__DomainFields.Count != DiffOp.DomainVar.Count) { throw new ArgumentException("mismatch between number of domain variables in spatial operator and given domain variables"); } if (_Parameters.Length != DiffOp.ParameterVar.Count) { throw new ArgumentException("mismatch between number of parameter variables in spatial operator and given parameters"); } if (m_CodomainMap.NoOfVariables != DiffOp.CodomainVar.Count) { throw new ArgumentException("mismatch between number of codomain variables in spatial operator and given codomain mapping"); } var _DomainAndParamFields = ArrayTools.Cat(__DomainFields, _Parameters); this.DELTA = __DomainFields.Count; m_DomainAndParamFieldsA = new ConventionalDGField[_DomainAndParamFields.Length]; m_DomainAndParamFieldsB = new ConventionalDGField[_DomainAndParamFields.Length]; for (int i = 0; i < m_DomainAndParamFieldsA.Length; i++) { var f = _DomainAndParamFields[i]; if (f == null) { m_DomainAndParamFieldsA[i] = null; m_DomainAndParamFieldsB[i] = null; } else if (f is XDGField xf) { m_DomainAndParamFieldsA[i] = xf.GetSpeciesShadowField(this.SpeciesA); m_DomainAndParamFieldsB[i] = xf.GetSpeciesShadowField(this.SpeciesB); } else if (f is ConventionalDGField cf) { m_DomainAndParamFieldsA[i] = cf; m_DomainAndParamFieldsB[i] = null; } else { throw new NotImplementedException("missing implementation for " + f.GetType().Name); } } LECQuadratureLevelSet <IMutableMatrix, double[]> .TestNegativeAndPositiveSpecies(domAndRule, m_lsTrk, SpeciesA, SpeciesB, m_LevSetIdx); // ------------------------ // sort equation components // ------------------------ int Gamma = m_CodomainMap.NoOfVariables; m_NonlinLsForm_V = EquationComponentArgMapping <INonlinLevelSetForm_V> .GetArgMapping(DiffOp, true, eq => ((eq.LevelSetTerms & (TermActivationFlags.V | TermActivationFlags.UxV | TermActivationFlags.GradUxV)) != 0) && Compfilter(eq), eq => (eq is ILevelSetForm)?new NonlinearLevelSetFormVectorizer((ILevelSetForm)eq, lsTrk) : null); m_NonlinLsForm_GradV = EquationComponentArgMapping <INonlinLevelSetForm_GradV> .GetArgMapping(DiffOp, true, eq => ((eq.LevelSetTerms & (TermActivationFlags.GradV | TermActivationFlags.UxGradV | TermActivationFlags.GradUxGradV)) != 0) && Compfilter(eq), eq => (eq is ILevelSetForm)?new NonlinearLevelSetFormVectorizer((ILevelSetForm)eq, lsTrk) : null); m_ValueRequired = new bool[m_DomainAndParamFieldsA.Length]; m_GradientRequired = new bool[m_DomainAndParamFieldsA.Length]; m_NonlinLsForm_V.DetermineReqFields(m_GradientRequired, comp => ((comp.LevelSetTerms & (TermActivationFlags.GradUxGradV | TermActivationFlags.GradUxV)) != 0)); m_NonlinLsForm_GradV.DetermineReqFields(m_GradientRequired, comp => ((comp.LevelSetTerms & (TermActivationFlags.GradUxGradV | TermActivationFlags.GradUxV)) != 0)); m_NonlinLsForm_V.DetermineReqFields(m_ValueRequired, comp => ((comp.LevelSetTerms & (TermActivationFlags.UxGradV | TermActivationFlags.UxV)) != 0)); m_NonlinLsForm_GradV.DetermineReqFields(m_ValueRequired, comp => ((comp.LevelSetTerms & (TermActivationFlags.UxGradV | TermActivationFlags.UxV)) != 0)); for (int i = __DomainFields.Count; i < m_DomainAndParamFieldsA.Length; i++) { m_ValueRequired[i] = true; // parameters are always required! } // ----- // alloc // ----- Koeff_V = new MultidimensionalArray[Gamma]; Koeff_GradV = new MultidimensionalArray[Gamma]; for (int gamma = 0; gamma < Gamma; gamma++) { Koeff_V[gamma] = new MultidimensionalArray(3); Koeff_GradV[gamma] = new MultidimensionalArray(4); } Debug.Assert(m_DomainAndParamFieldsA.Length == m_DomainAndParamFieldsB.Length); int L = m_DomainAndParamFieldsA.Length; m_FieldValuesPos = new MultidimensionalArray[L]; m_FieldValuesNeg = new MultidimensionalArray[L]; m_FieldGradientValuesPos = new MultidimensionalArray[L]; m_FieldGradientValuesNeg = new MultidimensionalArray[L]; for (int l = 0; l < L; l++) { if (m_ValueRequired[l]) { m_FieldValuesPos[l] = new MultidimensionalArray(2); m_FieldValuesNeg[l] = new MultidimensionalArray(2); } if (m_GradientRequired[l]) { m_FieldGradientValuesPos[l] = new MultidimensionalArray(3); m_FieldGradientValuesNeg[l] = new MultidimensionalArray(3); } } // ------------------ // init custom timers // ------------------ base.CustomTimers = new Stopwatch[] { new Stopwatch(), new Stopwatch(), new Stopwatch(), new Stopwatch(), new Stopwatch() }; base.CustomTimers_Names = new string[] { "Flux-Eval", "Basis-Eval", "Loops", "ParametersAndNormals", "Field-Eval" }; base.CustomTimers_RootPointer = new int[5]; ArrayTools.SetAll(base.CustomTimers_RootPointer, -1); this.m_NonlinLsForm_V_Watches = this.m_NonlinLsForm_V.InitStopWatches(0, this); this.m_NonlinLsForm_GradV_Watches = this.m_NonlinLsForm_GradV.InitStopWatches(0, this); Flux_Eval = base.CustomTimers[0]; Basis_Eval = base.CustomTimers[1]; Loops = base.CustomTimers[2]; ParametersAndNormals = base.CustomTimers[3]; Field_Eval = base.CustomTimers[4]; }