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 MODEL::LastNode() { GRID* rg = region; // allocate temporary memory for counter vector and initialize ------------------------- int* counter = (int*) MEMORY::memo.Array_nd( rg->Getnp() ); for( int i=0; i<rg->Getnp(); i++ ) counter[i] = 0; // initialize the elem->isLast array --------------------------------------------------- for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; int nnd = el->Getnnd(); for( int i=0; i<nnd; i++ ) el->isLast[0][i] = 0; } // determine the number of elements associated to nodes ---------------------------------- for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; int nnd = el->Getnnd(); for( int i=0; i<nnd; i++ ) { int no = el->nd[i]->Getno(); counter[no]++; } } // determine the last occurence of nodes ----------------------------------------------- for( int e=0; e<ne; e++ ) { ELEM* el = elem[e]; int nnd = el->Getnnd(); for( int i=0; i<nnd; i++ ) { int no = el->nd[i]->Getno(); counter[no]--; if( counter[no] < 0 ) REPORT::rpt.Error( kUnexpectedFault, "wrong count of nodes (MODEL::lastNode - 2)" ); if( counter[no] == 0 ) el->isLast[0][i] = 1; } } MEMORY::memo.Detach( counter ); char text[100]; sprintf( text, "\n%-25s%s\n", " (MODEL::LastNode)", "last occurence of nodes determined" ); REPORT::rpt.Output( text, 2 ); }