void straightenMeshUpdateEdgesOnBoundaryIsolated( ElementSpaceType & straightener, mpl::int_<3> /**/ ) { typedef typename ElementSpaceType::functionspace_type space_type; typedef typename space_type::mesh_type mesh_type; typedef typename space_type::dof_type::fe_type fe_type; auto const ncdof = space_type::dof_type::nComponents; auto const dofshift = fe_type::nDofPerVertex*mesh_type::element_type::numVertices; auto mesh = straightener.mesh(); auto const myrank = mesh->worldComm().localRank(); std::set<size_type> edgeIdFoundToUpdate; // search if there are some edges on boundary which are the interface of two ghost boundary faces // in this case this edges must not be straigthen auto itedge = mesh->beginEdgeOnBoundary(); auto const enedge = mesh->endEdgeOnBoundary(); for ( ; itedge!=enedge ; ++itedge ) { if ( itedge->numberOfProcGhost()==0 ) continue; auto const theedgeid = itedge->id(); std::set<size_type> ghostFaceIdFoundOnBoundary; auto itprocghost=itedge->elementsGhost().begin(); auto const enprocghost=itedge->elementsGhost().end(); for ( ; itprocghost!=enprocghost ; ++itprocghost) { auto iteltghost = itprocghost->second.begin(); auto const eneltghost = itprocghost->second.end(); for ( ; iteltghost!=eneltghost ; ++iteltghost ) { auto const& eltGhost = mesh->element(*iteltghost,itprocghost->first); for ( uint16_type f = 0 ; f < mesh_type::element_type::numTopologicalFaces ; ++f ) { auto const& theface = eltGhost.face(f); if ( theface.isOnBoundary() ) { bool findEdge=false; for ( uint16_type e = 0; e < mesh_type::face_type::numEdges && !findEdge ; ++e ) { if ( theface.edge(e).id() == theedgeid) { findEdge=true; ghostFaceIdFoundOnBoundary.insert(theface.id());} } } } } } // for ( ; itprocghost!=enprocghost ; ++itprocghost) // if 2 faces are find then the edge must not be straigten if (ghostFaceIdFoundOnBoundary.size()==2) edgeIdFoundToUpdate.insert(theedgeid); } // for ( ; itedge!=enedge ; ++itedge ) // apply a null displacement on egdes founded (not straighten) if (edgeIdFoundToUpdate.size() > 0) { auto range = elements(mesh,EntityProcessType::ALL); auto iteltactif = range.template get<1>(); auto const eneltactif = range.template get<2>(); for ( ; iteltactif!=eneltactif ; ++iteltactif ) { auto const& curElt = boost::unwrap_ref(*iteltactif); if ( !curElt.isOnBoundary() ) continue; for ( uint16_type e = 0; e < mesh_type::element_type::numEdges ; ++e ) { if ( edgeIdFoundToUpdate.find(curElt.edge(e).id()) != edgeIdFoundToUpdate.end()) { // find edge auto const idEltFind = curElt.id(); for ( uint16_type locdof = 0 ; locdof<fe_type::nDofPerEdge ; ++locdof ) { auto const local_id = dofshift + e*fe_type::nDofPerEdge + locdof; for ( uint16_type comp = 0; comp < ncdof; ++comp ) { const size_type globdof = straightener.functionSpace()->dof()->localToGlobal( idEltFind, local_id, comp ).template get<0>(); straightener( globdof ) = 0; } } } } } // for ( ; iteltactif!=eneltactif ; ++iteltactif ) } // if (edgeIdFoundToUpdate.size() > 0) } // straightenMeshUpdateEdgesOnBoundaryIsolated
typename CRBTrilinear<TruthModelType>::convergence_type CRBTrilinear<TruthModelType>::offline() { int proc_number = this->worldComm().globalRank(); bool rebuild_database = boption(_name="crb.rebuild-database") ; bool orthonormalize_primal = boption(_name="crb.orthonormalize-primal") ; boost::timer ti; if( this->worldComm().isMasterRank() ) std::cout << "Offline CRBTrilinear starts, this may take a while until Database is computed..."<<std::endl; LOG(INFO) << "[CRBTrilinear::offline] Starting offline for output " << this->M_output_index << "\n"; LOG(INFO) << "[CRBTrilinear::offline] initialize underlying finite element model\n"; //M_model->initModel(); LOG( INFO )<< " -- model init done in " << ti.elapsed() << "s"; parameter_type mu( this->M_Dmu ); double delta_pr; double delta_du; size_type index; //if M_N == 0 then there is not an already existing database if ( rebuild_database || this->M_N == 0) { ti.restart(); LOG(INFO) << "[CRBTrilinear::offline] compute random sampling\n"; int total_proc = this->worldComm().globalSize(); std::string sampling_mode = soption("crb.sampling-mode"); bool all_proc_same_sampling=boption("crb.all-procs-have-same-sampling"); int sampling_size = ioption("crb.sampling-size"); std::string file_name = ( boost::format("M_Xi_%1%_"+sampling_mode+"-proc%2%on%3%") % sampling_size %proc_number %total_proc ).str(); if( all_proc_same_sampling ) file_name+="-all-proc-have-same-sampling"; std::ifstream file ( file_name ); if( ! file ) { // random sampling std::string supersamplingname =(boost::format("Dmu-%1%-generated-by-master-proc") %sampling_size ).str(); if( sampling_mode == "log-random" ) this->M_Xi->randomize( sampling_size , all_proc_same_sampling , supersamplingname ); else if( sampling_mode == "log-equidistribute" ) this->M_Xi->logEquidistribute( sampling_size , all_proc_same_sampling , supersamplingname ); else if( sampling_mode == "equidistribute" ) this->M_Xi->equidistribute( sampling_size , all_proc_same_sampling , supersamplingname ); else throw std::logic_error( "[CRBTrilinear::offline] ERROR invalid option crb.sampling-mode, please select between log-random, log-equidistribute or equidistribute" ); //M_Xi->equidistribute( this->vm()["crb.sampling-size"].template as<int>() ); this->M_Xi->writeOnFile(file_name); } else { this->M_Xi->clear(); this->M_Xi->readFromFile(file_name); } this->M_WNmu->setSuperSampling( this->M_Xi ); LOG( INFO )<<"[CRBTrilinear offline] M_error_type = "<<this->M_error_type<<std::endl; LOG(INFO) << " -- sampling init done in " << ti.elapsed() << "s"; ti.restart(); // empty sets this->M_WNmu->clear(); if( this->M_error_type == CRB_NO_RESIDUAL ) mu = this->M_Dmu->element(); else { // start with M_C = { arg min mu, mu \in Xi } boost::tie( mu, index ) = this->M_Xi->min(); } int size = mu.size(); //std::cout << " -- WN size : " << M_WNmu->size() << "\n"; // dimension of reduced basis space this->M_N = 0; this->M_maxerror = 1e10; delta_pr = 0; delta_du = 0; //boost::tie( M_maxerror, mu, index ) = maxErrorBounds( N ); LOG(INFO) << "[CRBTrilinear::offline] allocate reduced basis data structures\n"; this->M_Aqm_pr.resize( this->M_model->Qa() ); for(int q=0; q<this->M_model->Qa(); q++) { this->M_Aqm_pr[q].resize( 1 ); } M_Aqm_tril_pr.resize( this->M_model->QaTri() ); //for(int q=0; q<this->M_model->QaTri(); q++) this->M_Fqm_pr.resize( this->M_model->Ql( 0 ) ); for(int q=0; q<this->M_model->Ql( 0 ); q++) { this->M_Fqm_pr[q].resize( 1 ); } this->M_Lqm_pr.resize( this->M_model->Ql( this->M_output_index ) ); for(int q=0; q<this->M_model->Ql( this->M_output_index ); q++) this->M_Lqm_pr[q].resize( 1 ); }//end of if( rebuild_database ) #if 1 else { mu = this->M_current_mu; if( proc_number == 0 ) { std::cout<<"we are going to enrich the reduced basis"<<std::endl; std::cout<<"there are "<<this->M_N<<" elements in the database"<<std::endl; } LOG(INFO) <<"we are going to enrich the reduced basis"<<std::endl; LOG(INFO) <<"there are "<<this->M_N<<" elements in the database"<<std::endl; }//end of else associated to if ( rebuild_databse ) #endif LOG(INFO) << "[CRBTrilinear::offline] compute affine decomposition\n"; std::vector< std::vector<sparse_matrix_ptrtype> > Aqm; std::vector< std::vector<sparse_matrix_ptrtype> > Aqm_tril; std::vector< std::vector<std::vector<vector_ptrtype> > > Fqm; boost::tie( boost::tuples::ignore, Aqm, Fqm ) = this->M_model->computeAffineDecomposition(); element_ptrtype u( new element_type( this->M_model->functionSpace() ) ); LOG(INFO) << "[CRBTrilinear::offline] starting offline adaptive loop\n"; bool reuse_prec = this->vm()["crb.reuse-prec"].template as<bool>() ; bool use_predefined_WNmu = this->vm()["crb.use-predefined-WNmu"].template as<bool>() ; int N_log_equi = this->vm()["crb.use-logEquidistributed-WNmu"].template as<int>() ; int N_equi = this->vm()["crb.use-equidistributed-WNmu"].template as<int>() ; int N_random = ioption( "crb.use-random-WNmu" ); /* if( N_log_equi > 0 || N_equi > 0 ) use_predefined_WNmu = true;*/ // file where the sampling is savec std::string file_name = ( boost::format("SamplingWNmu") ).str(); std::ifstream file ( file_name ); this->M_WNmu->clear(); if ( use_predefined_WNmu ) // In this case we want to read the sampling { if( ! file ) // The user forgot to give the sampling file throw std::logic_error( "[CRB::offline] ERROR the file SamplingWNmu doesn't exist so it's impossible to known which parameters you want to use to build the database" ); else { int sampling_size = this->M_WNmu->readFromFile(file_name); if( Environment::isMasterRank() ) std::cout<<"[CRB::offline] Read WNmu ( sampling size : " << sampling_size <<" )"<<std::endl; LOG( INFO )<<"[CRB::offline] Read WNmu ( sampling size : " << sampling_size <<" )"; } } else // We generate the sampling with choosen strategy { if ( N_log_equi>0 ) { this->M_WNmu->logEquidistribute( N_log_equi , true ); if( Environment::isMasterRank() ) std::cout<<"[CRB::offline] Log-Equidistribute WNmu ( sampling size : " <<N_log_equi<<" )"<<std::endl; LOG( INFO )<<"[CRB::offline] Log-Equidistribute WNmu ( sampling size : " <<N_log_equi<<" )"; } else if ( N_equi>0 ) { this->M_WNmu->equidistribute( N_equi , true ); if( Environment::isMasterRank() ) std::cout<<"[CRB::offline] Equidistribute WNmu ( sampling size : " <<N_equi<<" )"<<std::endl; LOG( INFO )<<"[CRB::offline] Equidistribute WNmu ( sampling size : " <<N_equi<<" )"; } else if ( N_random>0 ) { this->M_WNmu->randomize( N_random , true ); if( Environment::isMasterRank() ) std::cout<<"[CRB::offline] Randomize WNmu ( sampling size : " <<N_random<<" )"<<std::endl; LOG( INFO )<<"[CRB::offline] Randomize WNmu ( sampling size : " <<N_random<<" )"; } else // In this case we don't know what sampling to use throw std::logic_error( "[CRB::offline] ERROR : You have to choose an appropriate strategy for the offline sampling : random, equi, logequi or predefined" ); this->M_WNmu->writeOnFile(file_name); /* if( ! file ) { this->M_WNmu->clear(); std::vector< parameter_type > V; parameter_type __mu; __mu = this->M_Dmu->element(); __mu(0)= 1 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 111112 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 222223 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 333334 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 444445 , __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 555556 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 666667 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 777778 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 888889 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 1e+06 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 8123 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)= 9123 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=1.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=2.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=4.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=912 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=1.123e3 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=4.123e3 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=7.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=2123 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=6.123e3 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=3.123e3 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=3.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=5.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=9.123e4 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=812 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=5.111e3 ; __mu(1)= 1 ; V.push_back( __mu ); __mu(0)=5.124e2 ; __mu(1)= 1 ; V.push_back( __mu ); this->M_WNmu->setElements( V ); this->M_iter_max = this->M_WNmu->size(); this->M_WNmu->writeOnFile(file_name); }*/ use_predefined_WNmu=true; } //build sampling this->M_iter_max = this->M_WNmu->size(); mu = this->M_WNmu->at( this->M_N ); // first element if( this->M_error_type == CRB_NO_RESIDUAL || use_predefined_WNmu ) { //in this case it makes no sens to check the estimated error this->M_maxerror = 1e10; } LOG(INFO) << "[CRBTrilinear::offline] strategy "<< this->M_error_type <<"\n"; while ( this->M_maxerror > this->M_tolerance && this->M_N < this->M_iter_max ) { boost::timer timer, timer2; LOG(INFO) <<"========================================"<<"\n"; if( proc_number == this->worldComm().masterRank() ) std::cout<<"construction of "<<this->M_N<<"/"<<this->M_iter_max<<" basis "<<std::endl; LOG(INFO) << "N=" << this->M_N << "/" << this->M_iter_max << "( nb proc : "<<worldComm().globalSize()<<")"; // for a given parameter \p mu assemble the left and right hand side u->setName( ( boost::format( "fem-primal-N%1%-proc%2%" ) % (this->M_N) % proc_number ).str() ); mu.check(); u->zero(); timer2.restart(); LOG(INFO) << "[CRB::offline] solving primal" << "\n"; *u = this->M_model->solve( mu ); //if( proc_number == this->worldComm().masterRank() ) std::cout << " -- primal problem solved in " << timer2.elapsed() << "s\n"; timer2.restart(); if( ! use_predefined_WNmu ) this->M_WNmu->push_back( mu, index ); this->M_WNmu_complement = this->M_WNmu->complement(); this->M_model->rBFunctionSpace()->addPrimalBasisElement( *u ); //WARNING : the dual element is not the real dual solution ! //no dual problem was solved this->M_model->rBFunctionSpace()->addDualBasisElement( *u ); int number_of_added_elements=1; this->M_N+=number_of_added_elements; if ( orthonormalize_primal ) { this->orthonormalize( this->M_N, this->M_model->rBFunctionSpace()->primalRB() , number_of_added_elements ); this->orthonormalize( this->M_N, this->M_model->rBFunctionSpace()->primalRB() , number_of_added_elements ); this->orthonormalize( this->M_N, this->M_model->rBFunctionSpace()->primalRB() , number_of_added_elements ); } LOG(INFO) << "[CRB::offline] compute Aq_pr, Aq_du, Aq_pr_du" << "\n"; for (size_type q = 0; q < this->M_model->Qa(); ++q ) { this->M_Aqm_pr[q][0].conservativeResize( this->M_N, this->M_N ); // only compute the last line and last column of reduced matrices for ( size_type i = this->M_N-number_of_added_elements; i < this->M_N; i++ ) { for ( size_type j = 0; j < this->M_N; ++j ) { this->M_Aqm_pr[q][0]( i, j ) = Aqm[q][0]->energy( this->M_model->rBFunctionSpace()->primalBasisElement(i) , this->M_model->rBFunctionSpace()->primalBasisElement(j) ); } } for ( size_type j=this->M_N-number_of_added_elements; j < this->M_N; j++ ) { for ( size_type i = 0; i < this->M_N; ++i ) { this->M_Aqm_pr[q][0]( i, j ) = Aqm[q][0]->energy( this->M_model->rBFunctionSpace()->primalBasisElement(i), this->M_model->rBFunctionSpace()->primalBasisElement(j) ); } } }//loop over q LOG(INFO) << "[CRBTrilinear::offline] compute Fq_pr" << "\n"; for ( size_type q = 0; q < this->M_model->Ql( 0 ); ++q ) { this->M_Fqm_pr[q][0].conservativeResize( this->M_N ); for ( size_type l = 1; l <= number_of_added_elements; ++l ) { int index = this->M_N-l; this->M_Fqm_pr[q][0]( index ) = this->M_model->Fqm( 0, q, 0, this->M_model->rBFunctionSpace()->primalBasisElement(index) ); } }//loop over q LOG(INFO) << "[CRB::offline] compute Lq_pr" << "\n"; for ( size_type q = 0; q < this->M_model->Ql( this->M_output_index ); ++q ) { this->M_Lqm_pr[q][0].conservativeResize( this->M_N ); for ( size_type l = 1; l <= number_of_added_elements; ++l ) { int index = this->M_N-l; this->M_Lqm_pr[q][0]( index ) = this->M_model->Fqm( this->M_output_index, q, 0, this->M_model->rBFunctionSpace()->primalBasisElement(index) ); } }//loop over q sparse_matrix_ptrtype trilinear_form; for (size_type q = 0; q < this->M_model->QaTri(); ++q ) { M_Aqm_tril_pr[q].resize( this->M_N ); for (int k=0 ; k<this->M_N; k++) { //bring back the matrix associated to the trilinear form for a given basis function //we do this here to use only one matrix trilinear_form = this->M_model->computeTrilinearForm( this->M_model->rBFunctionSpace()->primalBasisElement(k) ); M_Aqm_tril_pr[q][k].conservativeResize( this->M_N, this->M_N ); for ( int i = 0; i < this->M_N; ++i ) { for ( int j = 0; j < this->M_N; ++j ) { M_Aqm_tril_pr[q][k]( i, j ) = trilinear_form->energy( this->M_model->rBFunctionSpace()->primalBasisElement(j), this->M_model->rBFunctionSpace()->primalBasisElement(i) ); }//j }//i }//k }// q timer2.restart(); if ( ! use_predefined_WNmu ) { bool already_exist; do { //initialization already_exist=false; //pick randomly an element mu = this->M_Dmu->element(); //make sure that the new mu is not already is M_WNmu BOOST_FOREACH( auto _mu, *this->M_WNmu ) { if( mu == _mu ) already_exist=true; } } while( already_exist ); this->M_current_mu = mu; } else { //remmber that in this case M_iter_max = sampling size if( this->M_N < this->M_iter_max )