void EQS_KD2D::Execute( PROJECT* project, int steadyFlow, int shape ) { switch( shape ) { case 0: linearShape = false; quarterShape = false; break; case 1: linearShape = true; quarterShape = false; break; case 2: linearShape = false; quarterShape = true; break; } MODEL* model = project->M2D; GRID* rg = project->M2D->region; NODE** node = model->node; ELEM** elem = model->elem; int np = model->np; int ne = model->ne; int diverged_cg = 0; double* B = NULL; double* xKD = NULL; double* xKDo = NULL; double* cxKo = NULL; double* cxDo = NULL; model->Incinit(); // print information on actual iteration ----------------------------------------------- project->PrintTheCycle( 1 ); REPORT::rpt.PrintTime( 1 ); // set parameters according to time integration and relaxation ------------------------- double th = project->timeint.thetaTurb; double dt = project->timeint.incTime.Getsec(); if( fabs(th) < 1.0e-10 && fabs(dt) < 1.0e-10 ) { REPORT::rpt.Error( kParameterFault, "theta and timeInterval too small (EQS_KD2D::execute - 1)" ); } double thdt = 1.0 / dt / th; double dt_KD; double relaxDt_KD = project->timeint.relaxTimeTurb.Getsec(); int relaxMethod = project->relaxMethod; if( steadyFlow ) { if( relaxMethod >= 3 ) { dt_KD = relaxDt_KD; relaxThdt_KD = 1.0 / dt_KD; } else { dt_KD = dt; relaxThdt_KD = 0.0; } for( int i=0; i<np; i++ ) { node[i]->v.dKdt = 0.0; node[i]->v.dDdt = 0.0; } } else { if( relaxMethod >= 3 ) { dt_KD = relaxDt_KD; } else { dt_KD = dt; } relaxThdt_KD = 1.0 / dt_KD / th; // time prediction ------------------------------------------------------------------- for( int i=0; i<np; i++ ) { VARS* v = &(node[i]->v); //VARS* vo = &(node[i]->vo); //v->K = vo->K + vo->dKdt * dt; //v->D = vo->D + vo->dDdt * dt; //v->dKdt = (1.0 - 1.0/th)*vo->dKdt + thdt*(v->K - vo->K); //v->dDdt = (1.0 - 1.0/th)*vo->dDdt + thdt*(v->D - vo->D); v->dKdt = 0.0; v->dDdt = 0.0; } } // check KD-values for validity (>= 0) ------------------------------------------------- Validate( project, np, node ); // determine friction coefficients ----------------------------------------------------- model->DoFriction( project ); // initialize Reynolds stresses and eddy viscosity ------------------------------------- rg->Turbulence( project ); // set KD boundary conditions ---------------------------------------------------------- model->SetBoundKD( project ); // ------------------------------------------------------------------------------------- // in case of quartered elements: // compute averaged values of U, V, S, K and D for the virtual center node in quads int nq = 0; double* cxK = NULL; double* cxD = NULL; if( quarterShape ) { // count number of quadrilaterals for( int e=0; e<ne; e++ ) { if( elem[e]->Getncn() == 4 ) nq++; } // allocate memory for nq center nodes if( !cbuf ) cbuf = new NODE [nq]; if( !cbuf ) { REPORT::rpt.Error( kMemoryFault, "cannot allocate memory (EQS_KD2D::execute - 2)" ); } cent = (NODE**) MEMORY::memo.Array_el( ne ); cxK = (double*) MEMORY::memo.Array_el( ne ); cxD = (double*) MEMORY::memo.Array_el( ne ); nq = 0; for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { int no = el->Getno(); int ncn = el->Getncn(); cent[no] = NULL; if( ncn == 4 ) { cent[no] = &cbuf[nq++]; cent[no]->Setno( no ); cent[no]->x = 0.0; cent[no]->y = 0.0; cent[no]->z = 0.0; cent[no]->cf = 0.0; cent[no]->v.U = 0.0; cent[no]->v.V = 0.0; cent[no]->v.S = 0.0; cent[no]->v.K = 0.0; cent[no]->v.D = 0.0; cent[no]->v.dKdt = 0.0; cent[no]->v.dDdt = 0.0; for( int i=0; i<ncn; i++ ) { cent[no]->x += el->nd[i]->x; cent[no]->y += el->nd[i]->y; cent[no]->z += el->nd[i]->z; cent[no]->cf += el->nd[i]->cf; cent[no]->v.U += el->nd[i]->v.U; cent[no]->v.V += el->nd[i]->v.V; cent[no]->v.S += el->nd[i]->v.S; cent[no]->v.K += el->nd[i]->v.K; cent[no]->v.D += el->nd[i]->v.D; cent[no]->v.dKdt += el->nd[i]->v.dKdt; cent[no]->v.dDdt += el->nd[i]->v.dDdt; } cent[no]->x /= ncn; cent[no]->y /= ncn; cent[no]->z /= ncn; cent[no]->cf /= ncn; cent[no]->v.U /= ncn; cent[no]->v.V /= ncn; cent[no]->v.S /= ncn; cent[no]->v.K /= ncn; cent[no]->v.D /= ncn; cent[no]->v.dKdt /= ncn; cent[no]->v.dDdt /= ncn; } } } } // ------------------------------------------------------------------------------------- // iteration loop double dt_KDo = 0.0; double relaxo = 1.0; double maxKDo = -1.0; int conv; for( int it=0; it<project->actualCycit; it++ ) { conv = true; // print information on actual iteration --------------------------------------------- if( it > 0 ) { project->PrintTheCycle( it+1 ); REPORT::rpt.PrintTime( 1 ); } // initialize Reynolds stresses and eddy viscosity ----------------------------------- if( it > 0 && isFS(project->actualTurb, BCONSET::kVtIterat) && !isFS(project->actualTurb, BCONSET::kVtPrandtlKol) ) { rg->Turbulence( project ); } // set up equation numbers ----------------------------------------------------------- if( model->Getinit() != modelInit ) { initStructure = true; modelInit = model->Getinit(); project->fix[0] = BCON::kFixK; project->fix[1] = BCON::kFixD; project->elemKind = ELEM::kRegion; if( linearShape ) { SetEqno( model, 2, 0, 0, project->fix, project->elemKind ); } else { SetEqno( model, 2, 2, 0, project->fix, project->elemKind ); } if( B ) MEMORY::memo.Detach( B ); if( xKD ) MEMORY::memo.Detach( xKD ); B = (double*) MEMORY::memo.Array_eq( neq ); xKD = (double*) MEMORY::memo.Array_eq( neq ); // allocate memory for relaxed Newton-Rahpson if( relaxMethod >= 3 ) { xKDo = (double*) MEMORY::memo.Array_eq( neq ); if( quarterShape ) { cxKo = (double*) MEMORY::memo.Array_el( ne ); cxDo = (double*) MEMORY::memo.Array_el( ne ); } } } // solve equations with frontal solving algorithm ------------------------------------ for( int i=0; i<neq; i++ ) xKD[i] = 0.0; diverged_cg = Solve( model, neq, B, xKD, project ); // ----------------------------------------------------------------------------------- // determine new values for K and D at virtual center nodes if( quarterShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { int no = el->Getno(); int ncn = el->Getncn(); if( ncn == 4 ) { Coefs( el, project, estifm, force ); cxK[no] = force[16]; cxD[no] = force[17]; for( int i=0; i<8; i++ ) { int eqK = GetEqno( el->nd[i], 0 ); int eqD = GetEqno( el->nd[i], 1 ); cxK[no] -= estifm[16][i] * xKD[eqK] + estifm[16][i+8] * xKD[eqD]; cxD[no] -= estifm[17][i] * xKD[eqK] + estifm[17][i+8] * xKD[eqD]; } cxK[no] /= estifm[16][16]; cxD[no] -= estifm[17][16] * cxK[no]; cxD[no] /= estifm[17][17]; } } } } // ----------------------------------------------------------------------------------- // statistics REPORT::rpt.Message( 2, "\n\n%-25s%s\n\n %s\n\n", " (EQS_KD2D::Execute)", "convergence parameters ...", " variable node average maximum" ); // compute averaged changes and standard deviation of K and D ------------------------ double stdevK = 0.0; double stdevD = 0.0; double aveAbsK = 0.0; double aveAbsD = 0.0; for( int i=0; i<np; i++ ) { int eqno; double dK = 0.0; double dD = 0.0; eqno = GetEqno( node[i], 0 ); if( eqno >= 0 ) dK = xKD[eqno]; eqno = GetEqno( node[i], 1 ); if( eqno >= 0 ) dD = xKD[eqno]; aveAbsK += dK; stdevK += dK*dK; aveAbsD += dD; stdevD += dD*dD; } if( quarterShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { int no = el->Getno(); int ncn = el->Getncn(); if( ncn == 4 ) { double dK = cxK[no]; double dD = cxD[no]; aveAbsK += dK; stdevK += dK*dK; aveAbsD += dD; stdevD += dD*dD; } } } } // ----------------------------------------------------------------------------------- int nptot = 0; for( int i=0; i<np; i++ ) { NODE* nd = rg->Getnode(i); if( !isFS(nd->flag, NODE::kInface_DN) ) nptot++; } ////////////////////////////////////////////////////////////////////////////////////// // MPI: broadcast statistic # ifdef _MPI_ aveAbsK = project->subdom.Mpi_sum( aveAbsK ); stdevK = project->subdom.Mpi_sum( stdevK ); aveAbsD = project->subdom.Mpi_sum( aveAbsD ); stdevD = project->subdom.Mpi_sum( stdevD ); nptot = project->subdom.Mpi_sum( nptot ); # endif ////////////////////////////////////////////////////////////////////////////////////// aveAbsK /= nptot; stdevK = sqrt( stdevK / nptot ); aveAbsD /= nptot; stdevD = sqrt( stdevD / nptot ); double norm = stdevK + stdevD; double fractK = 2.0 * sqrt( stdevK ); double fractD = 2.0 * sqrt( stdevD ); // compute maximum changes of K and D limited to ~95% fractile ----------------------- int cntK = 0; int cntD = 0; int noPerK = 0; double avePerK = 0.0; double maxPerK = 0.0; int noAbsK = 0; double maxAbsK = 0.0; int noPerD = 0; double avePerD = 0.0; double maxPerD = 0.0; int noAbsD = 0; double maxAbsD = 0.0; for( int i=0; i<np; i++ ) { int eqno; double dK = 0.0; double dD = 0.0; eqno = GetEqno( node[i], 0 ); if( eqno >= 0 ) { dK = xKD[eqno]; if( fabs(dK) > project->maxDeltaKD && fabs(dK) > fractK ) { xKD[eqno] = dK/fabs(dK) * fractK; } } eqno = GetEqno( node[i], 1 ); if( eqno >= 0 ) { dD = xKD[eqno]; if( fabs(dD) > project->maxDeltaKD && fabs(dD) > fractD ) { xKD[eqno] = dD/fabs(dD) * fractD; } } // maximum changes and percentage if( node[i]->v.K + dK > 0.0 ) { if( fabs(dK) > fabs(maxAbsK) ) { maxAbsK = dK; noAbsK = node[i]->Getname(); } if( node[i]->v.K > 0.0 ) { double per = dK / node[i]->v.K; avePerK += fabs(per); cntK++; if( fabs(per) > fabs(maxPerK) ) { maxPerK = per; noPerK = node[i]->Getname(); } } } if( node[i]->v.D + dD > 0.0 ) { if( fabs(dD) > fabs(maxAbsD) ) { maxAbsD = dD; noAbsD = node[i]->Getname(); } if( node[i]->v.D > 0.0 ) { double per = dD / node[i]->v.D; avePerD += fabs(per); cntD++; if( fabs(per) > fabs(maxPerD) ) { maxPerD = per; noPerD = node[i]->Getname(); } } } } if( quarterShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { int no = el->Getno(); int ncn = el->Getncn(); if( ncn == 4 ) { double dK = cxK[no]; double dD = cxD[no]; if( fabs(dK) > project->maxDeltaKD && fabs(dK) > fractK ) { cxK[no] = dK/fabs(dK) * fractK; } if( fabs(dD) > project->maxDeltaKD && fabs(dD) > fractD ) { cxD[no] = dD/fabs(dD) * fractD; } // maximum changes and percentage if( cent[no]->v.K + dK > 0.0 ) { if( fabs(dK) > fabs(maxAbsK) ) { maxAbsK = dK; noAbsK = -(no+1); } if( cent[no]->v.K > 0.0 ) { double per = dK / cent[no]->v.K; avePerK += fabs(per); cntK++; if( fabs(per) > fabs(maxPerK) ) { maxPerK = per; noPerK = -(no+1); } } } if( cent[no]->v.D + dD > 0.0 ) { if( fabs(dD) > fabs(maxAbsD) ) { maxAbsD = dD; noAbsD = -(no+1); } if( cent[no]->v.D > 0.0 ) { double per = dD / cent[no]->v.D; avePerD += fabs(per); cntD++; if( fabs(per) > fabs(maxPerD) ) { maxPerD = per; noPerD = -(no+1); } } } } } } } ////////////////////////////////////////////////////////////////////////////////////// // MPI: broadcast statistic # ifdef _MPI_ cntK = project->subdom.Mpi_sum( cntK ); maxAbsK = project->subdom.Mpi_maxabs( maxAbsK ); avePerK = project->subdom.Mpi_sum( avePerK ); maxPerK = project->subdom.Mpi_maxabs( maxPerK ); cntD = project->subdom.Mpi_sum( cntD ); maxAbsD = project->subdom.Mpi_maxabs( maxAbsD ); avePerD = project->subdom.Mpi_sum( avePerD ); maxPerD = project->subdom.Mpi_maxabs( maxPerD ); # endif ////////////////////////////////////////////////////////////////////////////////////// avePerK /= cntK; avePerD /= cntD; REPORT::rpt.Message( 2, " %1c %5d %12.5le %12.5le %s\n", 'K', noAbsK+1, aveAbsK, maxAbsK, " (abs)" ); REPORT::rpt.Message( 2, " %1c %5d %12.5lf %12.5lf %s\n\n", ' ', noPerK+1, avePerK, maxPerK, " ( % )" ); REPORT::rpt.Message( 2, " %1c %5d %12.5le %12.5le %s\n", 'D', noAbsD+1, aveAbsD, maxAbsD, " (abs)" ); REPORT::rpt.Message( 2, " %1c %5d %12.5lf %12.5lf %s\n\n", ' ', noPerD+1, avePerD, maxPerD, " ( % )" ); if( fabs(maxAbsK) > project->convKD ) conv = false; if( fabs(maxAbsD) > project->convKD ) conv = false; // determine relaxation parameter for NEWTON-RAPHSON --------------------------------- double relax; double maxKD = fabs(maxAbsK); if( fabs(maxAbsD) > maxKD ) maxKD = fabs(maxAbsD); switch( relaxMethod ) { default: REPORT::rpt.Warning( kParameterFault, "relaxation method %d not supported", relaxMethod ); case 0: relax = 1.0; REPORT::rpt.Message( 2, "\n%-25s%s %12.4le\n%-25s%s %12.4le\n\n", " ", "relaxation (0): norm =", norm, " ", " relax =", relax ); break; case 2: relax = project->maxDeltaKD / maxKD; if( relax > 1.0 ) relax = 1.0; REPORT::rpt.Message( 2, "\n%-25s%s %12.4le\n%-25s%s %12.4le\n\n", " ", "relaxation (2): norm =", norm, " ", " relax =", relax ); break; case 3: REPORT::rpt.Message( 2, "\n%-25s%s %12.4le\n", " ", "relaxed time: dt_KD =", dt_KD ); if( dt_KD < dt ) conv = false; relax = project->maxDeltaKD / maxKD; dt_KD *= relax; if( dt_KD > dt ) dt_KD = dt; if( dt_KD < relaxDt_KD ) dt_KD = relaxDt_KD; if( steadyFlow ) relaxThdt_KD = 1.0 / dt_KD; else relaxThdt_KD = 1.0 / dt_KD / th; if( relax > 1.0 ) relax = 1.0; REPORT::rpt.Message( 2, "\n%-25s%s %12.4le\n%-25s%s %12.4le\n\n", " ", "relaxation (3): norm =", norm, " ", " relax =", relax ); break; case 4: REPORT::rpt.Message( 2, "\n%-25s%s %12.4le\n", " ", "relaxed time: dt_KD =", dt_KD ); if( dt_KD < dt ) conv = false; relax = project->maxDeltaKD / maxKD; if( relax < project->relaxMax ) { if( maxKDo < 0.0 || maxKD < maxKDo ) // initialisation: maxKDo = -1 { if( relax < project->relaxMin ) relax = project->relaxMin; dt_KDo = dt_KD; maxKDo = maxKD; } else { if( relaxo > 1.01 * project->relaxMin ) { // relaxed Newton-Raphson: restore K and D from previous iteration for( int i=0; i<np; i++ ) { int n; n = GetEqno( node[i], 0 ); if( n >= 0 ) { node[i]->v.K -= relaxo * xKDo[n]; xKD[n] = xKDo[n]; } n = GetEqno( node[i], 1 ); if( n >= 0 ) { node[i]->v.D -= relaxo * xKDo[n]; xKD[n] = xKDo[n]; } } if( quarterShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { if( el->Getncn() == 4 ) { int no = el->Getno(); cent[no]->v.K -= relaxo * cxKo[no]; cent[no]->v.D -= relaxo * cxDo[no]; cxK[no] = cxKo[no]; cxD[no] = cxDo[no]; } } } } } if( relax < project->relaxMin ) relax = project->relaxMin; dt_KDo = dt_KD; dt_KD *= relax; // decrease dt_KD if( dt_KD < relaxDt_KD ) dt_KD = relaxDt_KD; maxKDo = maxKD; } } else { dt_KDo = dt_KD; if( relax > 1.0 ) dt_KD *= relax; // increase dt_KD if( dt_KD > dt ) dt_KD = dt; relax = 1.0; maxKDo = maxKD; // if( relax > relaxo ) // { // relax = 0.1 * relaxo; // if( relax >= project->relaxMin ) // { // // relaxed Newton-Raphson: restore K and D from previous iteration // for( int i=0; i<np; i++ ) // { // int n; // n = GetEqno( node[i], 0 ); // if( n >= 0 ) // { // node[i]->v.K -= relaxo * xKDo[n]; // xKD[n] = xKDo[n]; // } // n = GetEqno( node[i], 1 ); // if( n >= 0 ) // { // node[i]->v.D -= relaxo * xKDo[n]; // xKD[n] = xKDo[n]; // } // } // if( quarterShape ) // { // for( int e=0; e<ne; e++ ) // { // ELEM* el = elem[e]; // if( isFS(el->flag, ELEM::kRegion) ) // { // if( el->Getncn() == 4 ) // { // int no = el->Getno(); // cent[no]->v.K -= relaxo * cxKo[no]; // cent[no]->v.D -= relaxo * cxDo[no]; // cxK[no] = cxKo[no]; // cxD[no] = cxDo[no]; // } // } // } // } // maxKDo = maxKD; // } // else // { // maxKDo = -1.0; // } // } // else // { // relax = project->relaxMax; // maxKDo = maxKD; // } // if( relax > project->relaxMax ) relax = project->relaxMax; // if( relax < project->relaxMin ) relax = project->relaxMin; // if( relax < 0.999 * relaxo ) conv = false; } if( steadyFlow ) relaxThdt_KD = 1.0 / dt_KD; else relaxThdt_KD = 1.0 / dt_KD / th; REPORT::rpt.Message( 2, "\n%-25s%s %12.4le\n%-25s%s %12.4le\n\n", " ", "relaxation (4): norm =", norm, " ", " relax =", relax ); // relaxed Newton-Raphson: store xKD, cxK and cxD for( int i=0; i<np; i++ ) { int n; n = GetEqno( node[i], 0 ); if( n >= 0 ) xKDo[n] = xKD[n]; n = GetEqno( node[i], 1 ); if( n >= 0 ) xKDo[n] = xKD[n]; } if( quarterShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { if( el->Getncn() == 4 ) { int no = el->Getno(); cxKo[no] = cxK[no]; cxDo[no] = cxD[no]; } } } } break; } // update ---------------------------------------------------------------------------- relaxo = relax; for( int i=0; i<np; i++ ) { int n; n = GetEqno( node[i], 0 ); if( n >= 0 ) node[i]->v.K += relax * xKD[n]; n = GetEqno( node[i], 1 ); if( n >= 0 ) node[i]->v.D += relax * xKD[n]; } if( quarterShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; if( isFS(el->flag, ELEM::kRegion) ) { if( el->Getncn() == 4 ) { int no = el->Getno(); cent[no]->v.K += relax * cxK[no]; cent[no]->v.D += relax * cxD[no]; } } } } // check for range of values --------------------------------------------------------- double minK = 0.0; double minD = 0.0; double maxK = 0.0; double maxD = 0.0; int first = true; int jk = 0; int jd = 0; for( int i=0; i<np; i++ ) { if( first ) { first = false; minK = maxK = node[i]->v.K; minD = maxD = node[i]->v.D; } else { if( node[i]->v.K < minK ) minK = node[i]->v.K; if( node[i]->v.K > maxK ) maxK = node[i]->v.K; if( node[i]->v.D < minD ) minD = node[i]->v.D; if( node[i]->v.D > maxD ) maxD = node[i]->v.D; } if( node[i]->v.K <= 0.0 ) { if( GetEqno(node[i],0) >= 0 ) jk++; } if( node[i]->v.D <= 0.0 ) { if( GetEqno(node[i],1) >= 0 ) jd++; } } ////////////////////////////////////////////////////////////////////////////////////// // MPI: broadcast statistic # ifdef _MPI_ jk = project->subdom.Mpi_sum( jk ); minK = project->subdom.Mpi_min( minK ); maxK = project->subdom.Mpi_max( maxK ); jd = project->subdom.Mpi_sum( jd ); minD = project->subdom.Mpi_min( minD ); maxD = project->subdom.Mpi_max( maxD ); # endif ////////////////////////////////////////////////////////////////////////////////////// // compute useful values for K and D where they are negative ------------------------- if( jk || jd ) { REPORT::rpt.Message( 3, "%-25s%d %s\n\n", " ", jk, "nodes with K out of range" ); REPORT::rpt.Message( 3, "%-25s%d %s\n\n", " ", jd, "nodes with D out of range" ); } REPORT::rpt.Message( 3, "%-25sminimum of K: %le\n", " ", minK ); REPORT::rpt.Message( 3, "%-25smaximum of K: %le\n", " ", maxK ); REPORT::rpt.Message( 3, "%-25sminimum of D: %le\n", " ", minD ); REPORT::rpt.Message( 3, "%-25smaximum of D: %le\n", " ", maxD ); Validate( project, np, node, ne, cent, elem ); // compute midside values (linear interpolation) ------------------------------------- if( linearShape ) { for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; int ncn = el->Getncn(); int nnd = el->Getnnd(); for( int i=ncn; i<nnd; i++ ) { // get left and right corner node to midside node i int il, ir; double left, rght; el->GetQShape()->getCornerNodes( i, &il, &ir ); left = el->nd[il]->v.K; rght = el->nd[ir]->v.K; el->nd[i]->v.K = 0.5 * (left + rght); left = el->nd[il]->v.D; rght = el->nd[ir]->v.D; el->nd[i]->v.D = 0.5 * (left + rght); } } } // compute time derivatives ---------------------------------------------------------- rg->ReportCuPe( dt, project->vk ); if( !steadyFlow ) { double iTheta = 1.0 - 1.0 / project->timeint.thetaTurb; for( int i=0; i<np; i++ ) { if( isFS(node[i]->flag, NODE::kDry) || isFS(node[i]->flag, NODE::kMarsh) ) { node[i]->v.dKdt = node[i]->v.dDdt = 0.0; } else { double K, pK, D, pD, pdKdt, pdDdt; VARS *v, *vo; v = &(node[i]->v); vo = &(node[i]->vo); K = v->K; pK = vo->K; pdKdt = vo->dKdt; D = v->D; pD = vo->D; pdDdt = vo->dDdt; // compute derivatives at actual time step node[i]->v.dKdt = iTheta*pdKdt + thdt*(K - pK); node[i]->v.dDdt = iTheta*pdDdt + thdt*(D - pD); } } } // ----------------------------------------------------------------------------------- if( conv || it == project->actualCycit-1 ) { if( REPORT::rpt.level == 1 && it == project->actualCycit-1 ) { REPORT::rpt.Message( 1, "\n\n%-25s%s: %d\n\n", " (EQS_KD2D::Execute)", "finished in iteration step", it+1 ); REPORT::rpt.Message( 1, "\n\n%-25s%s\n\n %s\n\n", " (EQS_KD2D::Execute)", "convergence parameters ...", " variable node average maximum" ); REPORT::rpt.Message( 1, " %1c %5d %12.5le %12.5le %s\n", 'K', noAbsK+1, aveAbsK, maxAbsK, " (abs)" ); REPORT::rpt.Message( 1, " %1c %5d %12.5lf %12.5lf %s\n\n", ' ', noPerK+1, avePerK, maxPerK, " ( % )" ); REPORT::rpt.Message( 1, " %1c %5d %12.5le %12.5le %s\n", 'D', noAbsD+1, aveAbsD, maxAbsD, " (abs)" ); REPORT::rpt.Message( 1, " %1c %5d %12.5lf %12.5lf %s\n\n", ' ', noPerD+1, avePerD, maxPerD, " ( % )" ); } break; } } // ------------------------------------------------------------------------------------- // finally: // compute eddy viscosity from revised turbulence parameters rg->Turbulence( project ); // ------------------------------------------------------------------------------------- MEMORY::memo.Detach( B ); MEMORY::memo.Detach( xKD ); if( cxK ) MEMORY::memo.Detach( cxK ); if( cxD ) MEMORY::memo.Detach( cxD ); if( cent ) MEMORY::memo.Detach( cent ); if( xKDo ) MEMORY::memo.Detach( xKDo ); if( cxKo ) MEMORY::memo.Detach( cxKo ); if( cxDo ) MEMORY::memo.Detach( cxDo ); // ------------------------------------------------------------------------------------- if( !conv ) project->errLevel |= kErr_some_errors | kErr_no_conv_nr; if( diverged_cg ) project->errLevel |= diverged_cg | kErr_no_conv_cg; }
void PRECO_ILUT::Solve( PROJECT* project, EQS* eqs, double* B, double* X ) { int neq = crsi->m_neq; int neq_dn = crsi->m_neq_dn; int neq_up = crsi->m_neq_up; int* width = crsi->m_width; int** index = crsi->m_index; REALPR** ILU = crsi->m_A; //////////////////////////////////////////////////////////////////////////////////////// // 1. forward solve: manipulate vector B with lower part of matrix ILU // // 1.1 MPI: forward solve for interior subdomain nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 1.1\n" ); # endif for( int i=0; i<neq_up; i++ ) { X[i] = B[i]; for( int j=1; j<width[i]; j++ ) { int eq = index[i][j]; if( eq < i ) // L-matrix { X[i] += ILU[i][j] * X[eq]; } } } //////////////////////////////////////////////////////////////////////////////////////// // 1.2 MPI communication: // receive from subdomains with smaller pid (2 <- 1) // s < pid ==> upstream interface nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 1.2\n" ); # endif # ifdef _MPI_ if( project->subdom.npr > 1 ) { SUBDOM* subdom = &project->subdom; INFACE* inface = subdom->inface; MPI_Status status; int df = eqs->dfcn; for( int s=subdom->npr-1; s>=0; s-- ) { int np = inface[s].np; // number of interface nodes if( np > 0 && s < subdom->pid ) // upstream interface nodes { int cnt = 0; # ifdef _MPI_DBG { char text[200]; sprintf( text, " ### receiving from %d:", s+1 ); REPORT::rpt.Output( text ); } # endif MPI_Recv( &cnt, 1, MPI_INT, s, 1, MPI_COMM_WORLD, &status ); # ifdef _MPI_DBG { char text[200]; sprintf( text, " %d values\n", cnt ); REPORT::rpt.Output( text ); } # endif if( cnt ) { MPI_Recv( inface[s].recv, cnt, MPI_DOUBLE, s, 2, MPI_COMM_WORLD, &status ); cnt = 0; for( int n=0; n<np; n++ ) { NODE* nd = inface[s].node[n]; if( !isFS(nd->flag, NODE::kDry) ) // nothing received if node is dry... { SUB* sub = nd->sub; while( sub ) { if( sub->no == s ) break; sub = sub->next; } if( !sub->dry ) // ...or if the node is dry in { // any adjacent subdomain for( int e=0; e<df; e++ ) { int eqno = eqs->GetEqno( nd, e ); if( eqno >= 0 ) { X[eqno] = inface[s].recv[cnt]; cnt++; } } } } } } } } } //////////////////////////////////////////////////////////////////////////////////////// // 1.3 MPI: forward solve for upstream interface nodes // Note: X[i] is not initialized with B[i] since // this was performed in prior subdomain # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 1.3\n" ); # endif for( int i=neq_up; i<neq_dn; i++ ) { for( int j=1; j<width[i]; j++ ) { int eq = index[i][j]; if( eq < i ) { X[i] += ILU[i][j] * X[eq]; } } } //////////////////////////////////////////////////////////////////////////////////////// // 1.4 MPI: faktorisation of downstream interface nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 1.4\n" ); # endif for( int i=neq_dn; i<neq; i++ ) { X[i] = B[i]; for( int j=1; j<width[i]; j++ ) { int eq = index[i][j]; if( eq < neq_dn ) { X[i] += ILU[i][j] * X[eq]; } } } //////////////////////////////////////////////////////////////////////////////////////// // 1.5 MPI communication: // send to subdomains with larger pid (2 -> 3) // s > pid ==> downstream interface nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 1.5\n" ); # endif if( project->subdom.npr > 1 ) { SUBDOM* subdom = &project->subdom; INFACE* inface = subdom->inface; int df = eqs->dfcn; for( int s=0; s<subdom->npr; s++ ) { int np = inface[s].np; // number of interface nodes if( np > 0 && s > subdom->pid ) // downstream interface nodes { int cnt = 0; for( int n=0; n<np; n++ ) { NODE* nd = inface[s].node[n]; if( !isFS(nd->flag, NODE::kDry) ) // nothing to send if the node is dry... { SUB* sub = nd->sub; while( sub ) { if( sub->no == s ) break; sub = sub->next; } if( !sub->dry ) // ...or if the node is dry in { // any adjacent subdomain for( int e=0; e<df; e++ ) { int eqno = eqs->GetEqno( nd, e ); if( eqno >= 0 ) { inface[s].send[cnt] = X[eqno]; cnt++; } } } } } # ifdef _MPI_DBG { char text[200]; sprintf( text, " ### sending %d values to %d\n", cnt, s+1 ); REPORT::rpt.Output( text ); } # endif MPI_Send( &cnt, 1, MPI_INT, s, 1, MPI_COMM_WORLD ); if( cnt ) MPI_Send( inface[s].send, cnt, MPI_DOUBLE, s, 2, MPI_COMM_WORLD ); } } } # endif //////////////////////////////////////////////////////////////////////////////////////// # ifdef kDebug { FILE* dbg; char filename[80]; MODEL* model = project->M2D; GRID* region = model->region; if( project->subdom.npr == 1 ) sprintf( filename, "forwX.dbg" ); else sprintf( filename, "forwX_%02d.dbg", project->subdom.pid ); dbg = fopen( filename, "w" ); fprintf( dbg, "%d\n", region->Getnp() ); for( int n=0; n<region->Getnp(); n++ ) { NODE* nd = region->Getnode( n ); for( int e=0; e<eqs->dfcn; e++ ) { int eqno = eqs->GetEqno( nd, e ); fprintf( dbg, "%5d %1d ", nd->Getname(), e ); if( eqno >= 0 ) fprintf( dbg, "%14.6le\n", X[eqno] ); else fprintf( dbg, "%14.6le\n", 0.0 ); } } fclose( dbg ); } # endif # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) forward solve finished\n" ); # endif //////////////////////////////////////////////////////////////////////////////////////// // 2. solve for X with upper part of matrix ILU (backward substitution) // // 2.1 MPI communication: // receive from subdomains with larger pid (2 <- 3) // s > pid ==> downstream interface nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 2.1\n" ); # endif # ifdef _MPI_ if( project->subdom.npr > 1 ) { SUBDOM* subdom = &project->subdom; INFACE* inface = subdom->inface; MPI_Status status; int df = eqs->dfcn; for( int s=subdom->npr-1; s>=0; s-- ) { int np = inface[s].np; // number of interface nodes if( np > 0 && s > subdom->pid ) // downstream interface nodes { int cnt = 0; # ifdef _MPI_DBG { char text[200]; sprintf( text, " ### receiving from %d:", s+1 ); REPORT::rpt.Output( text ); } # endif MPI_Recv( &cnt, 1, MPI_INT, s, 1, MPI_COMM_WORLD, &status ); # ifdef _MPI_DBG { char text[200]; sprintf( text, " %d values\n", cnt ); REPORT::rpt.Output( text ); } # endif if( cnt ) { MPI_Recv( inface[s].recv, cnt, MPI_DOUBLE, s, 2, MPI_COMM_WORLD, &status ); cnt = 0; for( int n=0; n<np; n++ ) { NODE* nd = inface[s].node[n]; if( !isFS(nd->flag, NODE::kDry) ) // nothing received if node is dry... { SUB* sub = nd->sub; while( sub ) { if( sub->no == s ) break; sub = sub->next; } if( !sub->dry ) // ...or if the node is dry in { // any adjacent subdomain for( int e=0; e<df; e++ ) { int req = eqs->GetEqno( nd, e ); if( req >= 0 ) { X[req] = inface[s].recv[cnt]; cnt++; } } } } } } } } } //////////////////////////////////////////////////////////////////////////////////////// // 2.2 MPI: solve for upstream interface nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 2.2\n" ); # endif for( int i=neq_dn-1; i>=neq_up; i-- ) { for( int j=1; j<width[i]; j++ ) { int eq = index[i][j]; if( eq > i ) { X[i] -= ILU[i][j] * X[eq]; } } if( fabs(ILU[i][0]) < kZero ) REPORT::rpt.Error( kParameterFault, "division by zero - EQS::ILU_solver(1)" ); X[i] /= ILU[i][0]; } //////////////////////////////////////////////////////////////////////////////////////// // 2.3 MPI communication: // send to subdomains with smaller pid (2 -> 1) // s < pid ==> upstream interface nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 2.3\n" ); # endif if( project->subdom.npr > 1 ) { SUBDOM* subdom = &project->subdom; INFACE* inface = subdom->inface; int df = eqs->dfcn; for( int s=0; s<subdom->npr; s++ ) { int np = inface[s].np; // number of interface nodes if( np > 0 && s < subdom->pid ) // upstream interface nodes { int cnt = 0; for( int n=0; n<np; n++ ) { NODE* nd = inface[s].node[n]; if( !isFS(nd->flag, NODE::kDry) ) // nothing to send if node is dry... { SUB* sub = nd->sub; while( sub ) { if( sub->no == s ) break; sub = sub->next; } if( !sub->dry ) // ...or if the node is dry in { // any adjacent subdomain for( int e=0; e<df; e++ ) { int eqno = eqs->GetEqno( nd, e ); if( eqno >= 0 ) { inface[s].send[cnt] = X[eqno]; cnt++; } } } } } # ifdef _MPI_DBG { char text[200]; sprintf( text, " ### sending %d values to %d\n", cnt, s+1 ); REPORT::rpt.Output( text ); } # endif MPI_Send( &cnt, 1, MPI_INT, s, 1, MPI_COMM_WORLD ); if( cnt ) MPI_Send( inface[s].send, cnt, MPI_DOUBLE, s, 2, MPI_COMM_WORLD ); } } } # endif //////////////////////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////////////////////// // 2.4 MPI: backward solve for interior nodes # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) starting with 2.4\n" ); # endif for( int i=neq_up-1; i>=0; i-- ) { for( int j=1; j<width[i]; j++ ) { int eq = index[i][j]; if( eq > i ) { X[i] -= ILU[i][j] * X[eq]; } } if( fabs(ILU[i][0]) < kZero ) REPORT::rpt.Error( kParameterFault, "division by zero - EQS::ILU_solver(2)" ); X[i] /= ILU[i][0]; } # ifdef kDebug { FILE* dbg; char filename[80]; MODEL* model = project->M2D; GRID* region = model->region; if( project->subdom.npr == 1 ) sprintf( filename, "backX.dbg" ); else sprintf( filename, "backX_%02d.dbg", project->subdom.pid ); dbg = fopen( filename, "w" ); fprintf( dbg, "%d\n", region->Getnp() ); for( int n=0; n<region->Getnp(); n++ ) { NODE* nd = region->Getnode( n ); for( int e=0; e<eqs->dfcn; e++ ) { fprintf( dbg, "%5d %1d ", nd->Getname(), e ); int eqno = eqs->GetEqno( nd, e ); if( eqno >= 0 ) fprintf( dbg, "%14.6le\n", X[eqno] ); else fprintf( dbg, "%14.6le\n", 0.0 ); } } fclose( dbg ); } MPI_Barrier( MPI_COMM_WORLD ); # endif # ifdef _MPI_DBG REPORT::rpt.Output( " (PRECO_ILUT::Solve) backward solve finished\n" ); # endif }
void EQS_KL2D::Dissipation( PROJECT* project ) { double cm = project->KD.cm; double cd = project->KD.cd; MODEL* model = project->M2D; GRID* rg = model->region; // initialization ---------------------------------------------------------------------- for( int n=0; n<rg->Getnp(); n++ ) { NODE* nd = rg->Getnode(n); nd->v.D = 0.0; } double* ndarea = (double*) MEMORY::memo.Array_nd( rg->Getnp() ); for( int n=0; n<rg->Getnp(); n++ ) ndarea[n] = 0.0; // ------------------------------------------------------------------------------------- // loop on elements // ------------------------------------------------------------------------------------- for( int e=0; e<rg->Getne(); e++ ) { ELEM* el = rg->Getelem(e); if( !isFS(el->flag, ELEM::kDry) ) { TYPE* type = TYPE::Getid( el->type ); int nnd = el->Getnnd(); NODE** nd = el->nd; // --------------------------------------------------------------------------------- // compute coordinates relative to first node // --------------------------------------------------------------------------------- double xmin, xmax, ymin, ymax; double x[kMaxNodes2D], y[kMaxNodes2D]; xmin = xmax = x[0] = nd[0]->x; ymin = ymax = y[0] = nd[0]->y; for( int i=1; i<nnd; i++ ) { if( nd[i]->x < xmin ) xmin = nd[i]->x; if( nd[i]->x > xmax ) xmax = nd[i]->x; if( nd[i]->y < ymin ) ymin = nd[i]->y; if( nd[i]->y > ymax ) ymax = nd[i]->y; x[i] = nd[i]->x - x[0]; y[i] = nd[i]->y - y[0]; } x[0] = y[0] = 0.0; double lm = type->lm; // grid depending lm --------------------------------------------------------------- if( isFS(project->actualTurb, BCONSET::kVtLES) ) { lm *= sqrt( (xmax-xmin)*(ymax-ymin) ); } double area = el->area(); for( int i=0; i<nnd; i++ ) { double K = nd[i]->v.K; if( K <= 0.0 ) K = project->minK; double ls = lm / sqrt( sqrt(cm*cm*cm/cd) ); nd[i]->v.D += area * cd * sqrt(K*K*K) / ls; ndarea[nd[i]->Getno()] += area; } } } for( int n=0; n<rg->Getnp(); n++ ) { NODE* nd = rg->Getnode(n); if( ndarea[n] > 0.0 ) nd->v.D /= ndarea[n]; } MEMORY::memo.Detach( ndarea ); }
void EQS_D2D::Execute( PROJECT* project, int steadyFlow, int linearShape ) { MODEL* model = project->M2D; GRID* rg = project->M2D->region; NODE** node = model->node; ELEM** elem = model->elem; int np = model->np; int ne = model->ne; int div_cg = 0; double* B = NULL; double* X = NULL; model->Incinit(); // print information on actual iteration ----------------------------------------------- project->PrintTheCycle( 1 ); REPORT::rpt.PrintTime( 1 ); // set parameters according to time integration and relaxation ------------------------- double th = project->timeint.thetaTurb; double dt = project->timeint.incTime.Getsec(); if( fabs(th) < 1.0e-10 && fabs(dt) < 1.0e-10 ) { REPORT::rpt.Error( kParameterFault, "theta and timeInterval too small (EQS_D2D::execute - 1)" ); } double thdt = 1.0 / dt / th; double dt_KD; double relaxDt_KD = project->timeint.relaxTimeTurb.Getsec(); int relaxMethod = project->relaxMethod; if( steadyFlow ) { if( relaxMethod >= 3 ) { dt_KD = relaxDt_KD; relaxThdt_KD = 1.0 / dt_KD; } else { dt_KD = dt; relaxThdt_KD = 0.0; } for( int i=0; i<np; i++ ) { node[i]->v.dDdt = 0.0; } } else { if( relaxMethod >= 3 ) { dt_KD = relaxDt_KD; } else { dt_KD = dt; } relaxThdt_KD = 1.0 / dt_KD / th; // time prediction ------------------------------------------------------------------- for( int i=0; i<np; i++ ) { VARS* v = &node[i]->v; //VARS* vo = &(node[i]->vo); //v->D = vo->D + vo->dDdt * dt; //v->dDdt = (1.0 - 1.0/th)*vo->dDdt + thdt*(v->D - vo->D); v->dDdt = 0.0; } } // check KD-values for validity (>= 0) ------------------------------------------------- Validate( np, node, project ); // determine friction coefficients ----------------------------------------------------- model->DoFriction( project ); // initialize Reynolds stresses and eddy viscosity ------------------------------------- rg->Turbulence( project ); // set KD boundary conditions ---------------------------------------------------------- model->SetBoundKD( project ); // ------------------------------------------------------------------------------------- int conv = true; for( int it=0; it<project->actualCycit; it++ ) { conv = true; // print information on actual iteration --------------------------------------------- if( it > 0 ) { project->PrintTheCycle( it+1 ); REPORT::rpt.PrintTime( 1 ); } // initialize Reynolds stresses and eddy viscosity ----------------------------------- if( it > 0 && isFS(project->actualTurb, BCONSET::kVtIterat) && !isFS(project->actualTurb, BCONSET::kVtPrandtlKol) ) { rg->Turbulence( project ); } // set up equation numbers ----------------------------------------------------------- if( model->Getinit() != modelInit ) { initStructure = true; modelInit = model->Getinit(); project->fix[0] = BCON::kFixD; project->elemKind = ELEM::kRegion; SetEqno( model, 1, 1, 0, project->fix, project->elemKind ); if( B ) MEMORY::memo.Detach( B ); if( X ) MEMORY::memo.Detach( X ); B = (double*) MEMORY::memo.Array_eq( neq ); X = (double*) MEMORY::memo.Array_eq( neq ); } // solve equations with frontal solving algorithm ------------------------------------ for( int i=0; i<neq; i++ ) X[i] = 0.0; div_cg = Solve( model, neq, B, X, project ); // ----------------------------------------------------------------------------------- // statistics REPORT::rpt.Message( 1, "\n\n%-25s%s\n\n %s\n\n", " (EQS_D2D::Execute)","convergence parameters ...", " variable node average maximum" ); // compute averaged changes and standard deviation of D ------------------------------ double stdevD = 0.0; double aveAbsD = 0.0; for( int i=0; i<np; i++ ) { int eqno; double dD = 0.0; eqno = GetEqno( node[i], 0 ); if( eqno >= 0 ) dD = X[eqno]; aveAbsD += dD; stdevD += dD*dD; } int nptot = 0; for( int i=0; i<np; i++ ) { NODE* nd = rg->Getnode(i); if( !isFS(nd->flag, NODE::kInface_DN) ) nptot++; } ////////////////////////////////////////////////////////////////////////////////////// // MPI: broadcast statistic # ifdef _MPI_ aveAbsD = project->subdom.Mpi_sum( aveAbsD ); stdevD = project->subdom.Mpi_sum( stdevD ); nptot = project->subdom.Mpi_sum( nptot ); # endif ////////////////////////////////////////////////////////////////////////////////////// aveAbsD /= nptot; stdevD = sqrt( stdevD / nptot ); double norm = stdevD; double fractD = 2.0 * sqrt( stdevD ); // compute maximum changes of K and D limited to ~95% fractile ----------------------- int cntD = 0; int noPerD = 0; double avePerD = 0.0; double maxPerD = 0.0; int noAbsD = 0; double maxAbsD = 0.0; for( int i=0; i<np; i++ ) { int eqno; double dK = 0.0; double dD = 0.0; eqno = GetEqno( node[i], 0 ); if( eqno >= 0 ) { dD = X[eqno]; if( fabs(dD) > project->maxDeltaKD && fabs(dD) > fractD ) { X[eqno] = dD/fabs(dD) * fractD; } } if( node[i]->v.D + dD > 0.0 ) { if( fabs(dD) > fabs(maxAbsD) ) { maxAbsD = dD; noAbsD = node[i]->Getname(); } if( node[i]->v.D > 0.0 ) { double per = dD / node[i]->v.D; avePerD += fabs(per); cntD++; if( fabs(per) > fabs(maxPerD) ) { maxPerD = per; noPerD = node[i]->Getname(); } } } } ////////////////////////////////////////////////////////////////////////////////////// // MPI: broadcast statistic # ifdef _MPI_ cntD = project->subdom.Mpi_sum( cntD ); maxAbsD = project->subdom.Mpi_max( maxAbsD ); avePerD = project->subdom.Mpi_sum( avePerD ); maxPerD = project->subdom.Mpi_max( maxPerD ); # endif ////////////////////////////////////////////////////////////////////////////////////// avePerD /= cntD; REPORT::rpt.Message( 1, " %1c %5d %12.5le %12.5le %s\n", 'D', noAbsD+1, aveAbsD, maxAbsD, " (abs)" ); REPORT::rpt.Message( 1, " %1c %5d %12.5lf %12.5lf %s\n\n", ' ', noPerD+1, avePerD, maxPerD, " ( % )" ); if( fabs(maxAbsD) > project->convKD ) conv = false; // determine relaxation parameter for NEWTON-RAPHSON --------------------------------- double relax; double maxKD = fabs(maxAbsD); switch( relaxMethod ) { default: REPORT::rpt.Warning( kParameterFault, "relaxation method %d not supported", relaxMethod ); case 0: relax = 1.0; REPORT::rpt.Message( 1, "\n%-25s%s %12.4le\n %s %12.4le\n\n", " ", "relaxation: norm =", norm, " ", " relax =", relax ); break; case 2: relax = project->maxDeltaKD / maxKD; if( relax > 1.0 ) relax = 1.0; REPORT::rpt.Message( 1, "\n%-25s%s %12.4le\n%-25s%s %12.4le\n\n", " ", "relaxation: norm =", norm, " ", " relax =", relax ); break; case 3: case 4: REPORT::rpt.Message( 1, "\n%-25s%s %12.4le\n", " ", "relaxed time: dt_KD =", dt_KD ); if( dt_KD < dt ) conv = false; relax = project->maxDeltaKD / maxKD; dt_KD *= relax; if( dt_KD > dt ) dt_KD = dt; if( dt_KD < relaxDt_KD ) dt_KD = relaxDt_KD; if( steadyFlow ) relaxThdt_KD = 1.0 / dt_KD; else relaxThdt_KD = 1.0 / dt_KD / th; if( relax > 1.0 ) relax = 1.0; REPORT::rpt.Message( 1, "\n%-25s%s %12.4le\n%-25s%s %12.4le\n\n", " ", "relaxation: norm =", norm, " ", " relax =", relax ); break; } // update ---------------------------------------------------------------------------- for( int i=0; i<np; i++ ) { int n = GetEqno( node[i], 0 ); if( n >= 0 ) node[i]->v.D += relax * X[n]; } // check for range of values --------------------------------------------------------- double minD = 0.0; double maxD = 0.0; int first = true; int jd = 0; for( int i=0; i<np; i++ ) { if( first ) { first = false; minD = maxD = node[i]->v.D; } else { if( node[i]->v.D < minD ) minD = node[i]->v.D; if( node[i]->v.D > maxD ) maxD = node[i]->v.D; } if( node[i]->v.D <= 0.0 ) { // node[i]->v.D = project->minD; if( GetEqno(node[i],1) >= 0 ) jd++; } } // compute useful values for K and D where they are negative ------------------------- if( jd ) { REPORT::rpt.Message( 1, "%-25s%d %s\n\n", " ", jd, "nodes with D out of range" ); Validate( np, node, project ); } REPORT::rpt.Message( 1, "%-25sminimum of D: %le\n", " ", minD ); REPORT::rpt.Message( 1, "%-25smaximum of D: %le\n", " ", maxD ); // compute time derivatives ---------------------------------------------------------- rg->ReportCuPe( dt, project->vk ); if( !steadyFlow ) { double iTheta; iTheta = 1.0 - 1.0 / project->timeint.thetaTurb; for( int i=0; i<np; i++ ) { if( isFS(node[i]->flag, NODE::kDry) || isFS(node[i]->flag, NODE::kMarsh) ) { node[i]->v.dDdt = 0.0; } else { double D, pD, pdDdt; VARS *v, *vo; v = &(node[i]->v); vo = &(node[i]->vo); D = v->D; pD = vo->D; pdDdt = vo->dDdt; // compute derivatives at actual time step node[i]->v.dDdt = iTheta*pdDdt + thdt*(D - pD); } } } if( conv ) break; } // finally: compute eddy viscosity from revised turbulence parameters K,D -------------- rg->Turbulence( project ); // ------------------------------------------------------------------------------------- MEMORY::memo.Detach( B ); MEMORY::memo.Detach( X ); if( !conv ) project->errLevel |= kErr_some_errors | kErr_no_conv_nr; if( div_cg ) project->errLevel |= div_cg | kErr_no_conv_cg; }